Loading...
Searching...
No Matches
waveModel.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015 IH-Cantabria
9 Copyright (C) 2016-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "waveModel.H"
30#include "fvMesh.H"
31#include "polyPatch.H"
32#include "gravityMeshObject.H"
33#include "volFields.H"
34#include "fvPatchFields.H"
36using namespace Foam::constant;
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
46const Foam::word Foam::waveModel::dictName("waveProperties");
47
48
51 return dictName + '.' + patchName;
52}
53
54
55// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
56
58{
59 // Determine local patch coordinate system given by:
60 // - X: streamwise: patch normal
61 // - Y: spanwise: Z^X
62 // - Z: up: (negative) gravity direction
63 vector x = normalised(-gAverage(patch_.faceAreas()));
64 vector z = -g_/mag(g_);
65 vector y = z^x;
66
67 // Rotation from global<->local about global origin
68 Rlg_ = tensor(x, y, z);
69 Rgl_ = Rlg_.T();
70
71 // Global patch extents
72 const vectorField& Cp = patch_.localPoints();
73 const vectorField CpLocal(Rgl_ & Cp);
74 boundBox bb(CpLocal, true);
75 const scalar xMin = bb.min().x();
76 const scalar xMax = bb.max().x();
77 const scalar yMin = bb.min().y();
78 const scalar yMax = bb.max().y();
79 zSpan_ = bb.max().z() - bb.min().z();
80
81 // Global x, y positions of the paddle centres
82 xPaddle_.setSize(nPaddle_, 0);
83 yPaddle_.setSize(nPaddle_, 0);
84 const scalar xMid = xMin + 0.5*(xMax - xMin);
85 const scalar paddleDy = (yMax - yMin)/scalar(nPaddle_);
86 for (label paddlei = 0; paddlei < nPaddle_; ++paddlei)
87 {
88 xPaddle_[paddlei] = xMid;
89 yPaddle_[paddlei] = paddlei*paddleDy + yMin + 0.5*paddleDy;
90 }
91
92 // Local face centres
93 const vectorField CfLocal(Rgl_ & patch_.faceCentres());
94 z_ = CfLocal.component(2);
95
96 // Local face extents in z-direction
97 zMin_.setSize(patch_.size());
98 zMax_.setSize(patch_.size());
99 const faceList& faces = patch_.localFaces();
100 forAll(faces, facei)
101 {
102 const face& f = faces[facei];
103 const label nPoint = f.size();
104 zMin_[facei] = CpLocal[f[0]].z();
105 zMax_[facei] = CpLocal[f[0]].z();
106
107 for (label fpi = 1; fpi < nPoint; ++fpi)
108 {
109 const label pointi = f[fpi];
110 zMin_[facei] = min(zMin_[facei], CpLocal[pointi].z());
111 zMax_[facei] = max(zMax_[facei], CpLocal[pointi].z());
112 }
113 }
114
115 // Set minimum z reference level
116 zMin0_ = gMin(zMin_);
117
118 // Local paddle-to-face addressing
119 faceToPaddle_.setSize(patch_.size(), -1);
121 {
122 faceToPaddle_[facei] = floor((CfLocal[facei].y() - yMin)/paddleDy);
123 }
124}
125
126
128{
129 // Note: initialising as initial depth
130 auto tlevel = tmp<scalarField>::New(nPaddle_, initialDepth_);
131 auto& level = tlevel.ref();
132
133 const volScalarField& alpha =
134 mesh_.lookupObject<volScalarField>(alphaName_);
135 const fvPatchScalarField& alphap = alpha.boundaryField()[patch_.index()];
136 const scalarField alphac(alphap.patchInternalField());
137
138 const scalarField& magSf = alphap.patch().magSf();
139 scalarList paddleMagSf(nPaddle_, Zero);
140 scalarList paddleWettedMagSf(nPaddle_, Zero);
141
142 forAll(alphac, facei)
143 {
144 label paddlei = faceToPaddle_[facei];
145 paddleMagSf[paddlei] += magSf[facei];
146 paddleWettedMagSf[paddlei] += magSf[facei]*alphac[facei];
147 }
148
149 forAll(paddleMagSf, paddlei)
150 {
151 reduce(paddleMagSf[paddlei], sumOp<scalar>());
152 reduce(paddleWettedMagSf[paddlei], sumOp<scalar>());
153 level[paddlei] +=
154 paddleWettedMagSf[paddlei]*zSpan_
155 /(paddleMagSf[paddlei] + ROOTVSMALL);
156 }
157
158 return tlevel;
159}
160
161
162void Foam::waveModel::setAlpha(const scalarField& level)
163{
164 forAll(alpha_, facei)
165 {
166 const label paddlei = faceToPaddle_[facei];
167 const scalar paddleCalc = level[paddlei];
168
169 const scalar zMin0 = zMin_[facei] - zMin0_;
170 const scalar zMax0 = zMax_[facei] - zMin0_;
171
172 if (zMax0 < paddleCalc)
173 {
174 alpha_[facei] = 1.0;
175 }
176 else if (zMin0 > paddleCalc)
177 {
178 alpha_[facei] = 0.0;
179 }
180 else
181 {
182 scalar dz = paddleCalc - zMin0;
183 alpha_[facei] = dz/(zMax0 - zMin0);
184 }
185 }
186}
187
188
190(
191 const scalarField& level,
192 const label facei,
193 scalar& fraction,
194 scalar& z
195) const
196{
197 const label paddlei = faceToPaddle_[facei];
198 const scalar paddleCalc = level[paddlei];
199 const scalar paddleHeight = min(paddleCalc, waterDepthRef_);
200 const scalar zMin = zMin_[facei] - zMin0_;
201 const scalar zMax = zMax_[facei] - zMin0_;
202
203 fraction = 1;
204 z = 0;
205
206 if (zMax < paddleHeight)
207 {
208 z = z_[facei] - zMin0_;
209 }
210 else if (zMin > paddleCalc)
211 {
212 fraction = -1;
213 }
214 else
215 {
216 if (paddleCalc < waterDepthRef_)
217 {
218 if ((zMax > paddleCalc) && (zMin < paddleCalc))
219 {
220 scalar dz = paddleCalc - zMin;
221 fraction = dz/(zMax - zMin);
222 z = z_[facei] - zMin0_;
223 }
224 }
225 else
226 {
227 if (zMax < paddleCalc)
228 {
229 z = waterDepthRef_;
230 }
231 else if ((zMax > paddleCalc) && (zMin < paddleCalc))
232 {
233 scalar dz = paddleCalc - zMin;
234 fraction = dz/(zMax - zMin);
235 z = waterDepthRef_;
236 }
238 }
239}
240
241
242// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
243
245(
246 const dictionary& dict,
247 const fvMesh& mesh,
248 const polyPatch& patch,
249 const bool readFields
250)
251:
253 (
255 (
256 modelName(patch.name()),
257 Time::timeName(mesh.time().startTime().value()),
258 "uniform",
259 mesh,
260 IOobject::NO_READ,
261 IOobject::AUTO_WRITE
262 )
263 ),
264 mesh_(mesh),
265 patch_(patch),
266 g_(meshObjects::gravity::New(mesh.time()).value()),
267 UName_("U"),
268 alphaName_("alpha"),
269 Rgl_(tensor::I),
270 Rlg_(tensor::I),
271 nPaddle_(1),
272 xPaddle_(),
273 yPaddle_(),
274 z_(),
275 zSpan_(0),
276 zMin_(),
277 zMax_(),
278 waterDepthRef_(0),
279 initialDepth_(0),
280 currTimeIndex_(-1),
281 activeAbsorption_(false),
282 U_(patch.size(), Zero),
283 alpha_(patch.size(), Zero)
284{
285 if (readFields)
286 {
288 }
289}
290
291
292// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
293
294bool Foam::waveModel::readDict(const dictionary& overrideDict)
295{
297 if (headerOk())
298 {
299 IOdictionary::regIOobject::read();
300 }
301
302 merge(overrideDict);
303
304 readIfPresent("U", UName_);
305 readIfPresent("alpha", alphaName_);
306
307 readEntry("nPaddle", nPaddle_);
308 if (nPaddle_ < 1)
309 {
311 << "Number of paddles must be greater than zero. Supplied"
312 << " value nPaddles = " << nPaddle_
313 << exit(FatalIOError);
314 }
315
316 readIfPresent("initialDepth", initialDepth_);
317
318 // Need to initialise the geometry before calling waterLevel()
319 initialiseGeometry();
320
321 // Set the reference water depth
322 if (!readIfPresent("waterDepthRef", waterDepthRef_))
323 {
324 scalar waterDepth = 0;
325 if (readIfPresent("waterDepth", waterDepth))
326 {
327 waterDepthRef_ = waterDepth;
328 }
329 else
330 {
331 const scalarField level(waterLevel());
332 if (level.size())
333 {
334 waterDepthRef_ = level.first();
335 }
336 }
337
338 // Avoid potential zero...
339 waterDepthRef_ += SMALL;
340
341 // Insert the reference water depth into [this] to enable restart
342 add("waterDepthRef", waterDepthRef_);
343 }
344
345 return true;
346}
347
348
349void Foam::waveModel::correct(const scalar t)
350{
351 if (mesh_.time().timeIndex() != currTimeIndex_)
352 {
353 Info<< "Updating " << type() << " wave model for patch "
354 << patch_.name() << endl;
355
356 // Time ramp weight
357 const scalar tCoeff = timeCoeff(t);
358
359 // Reset the velocity and phase fraction fields
360 U_ = vector::zero;
361 alpha_ = 0;
362
363 // Update the calculated water level field
364 scalarField calculatedLevel(nPaddle_, Zero);
365
366 if (patch_.size())
367 {
368 // Set wave level
369 setLevel(t, tCoeff, calculatedLevel);
370
371 // Update the velocity field
372 setVelocity(t, tCoeff, calculatedLevel);
373
374 // Update the phase fraction field
375 setAlpha(calculatedLevel);
376 }
377
378 if (activeAbsorption_)
379 {
380 const scalarField activeLevel(this->waterLevel());
381
382 forAll(U_, facei)
383 {
384 const label paddlei = faceToPaddle_[facei];
385
386 if (zMin_[facei] - zMin0_ < activeLevel[paddlei])
387 {
388 scalar UCorr =
389 (calculatedLevel[paddlei] - activeLevel[paddlei])
390 *sqrt(mag(g_)/activeLevel[paddlei]);
391
392 U_[facei].x() += UCorr;
393 }
394 else
395 {
396 U_[facei].x() = 0;
397 }
398 }
399 }
400
401 // Transform velocity into global coordinate system
402 U_ = Rlg_ & U_;
403
404 currTimeIndex_ = mesh_.time().timeIndex();
405 }
406}
407
410{
411 return U_;
412}
413
416{
417 return alpha_;
418}
419
420
422{
423 os << "Wave model: patch " << patch_.name() << nl
424 << " Type : " << type() << nl
425 << " Velocity field name : " << UName_ << nl
426 << " Phase fraction field name : " << alphaName_ << nl
427 << " Transformation from local to global system : " << Rlg_ << nl
428 << " Number of paddles: " << nPaddle_ << nl
429 << " Reference water depth : " << waterDepthRef_ << nl
430 << " Active absorption: " << activeAbsorption_ << nl;
431}
432
433
434// ************************************************************************* //
scalar y
propsDict readIfPresent("fields", acceptFields)
Foam::label startTime
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
readOption readOpt() const noexcept
Get the read option.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
T & first()
Access first element of the list, position [0].
Definition UList.H:957
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form zero
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
bool merge(const dictionary &dict)
Merge entries from the given dictionary.
Definition dictionary.C:817
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
bool headerOk()
Read and check header info. Does not check the headerClassName.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
Base class for waveModels.
Definition waveModel.H:55
virtual void initialiseGeometry()
Initialise.
Definition waveModel.C:50
scalarField zMin_
Minimum z (point) height per patch face / [m].
Definition waveModel.H:128
scalarField z_
Patch face centre z coordinates / [m].
Definition waveModel.H:118
virtual scalar timeCoeff(const scalar t) const =0
Return the time scaling coefficient.
const polyPatch & patch_
Reference to the polyPatch.
Definition waveModel.H:68
tensor Rlg_
Rotation tensor from local to global system.
Definition waveModel.H:93
const fvMesh & mesh_
Reference to the mesh database.
Definition waveModel.H:63
scalar zSpan_
Overall (point) span in z-direction / [m].
Definition waveModel.H:123
word UName_
Name of velocity field, default = "U".
Definition waveModel.H:78
const vector & g_
Gravity.
Definition waveModel.H:73
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
Definition waveModel.C:408
virtual void correct(const scalar t)
Correct the model for time, t[s].
Definition waveModel.C:342
scalar zMin0_
Minimum z reference level / [m].
Definition waveModel.H:138
scalarField alpha_
Wave indicator field.
Definition waveModel.H:171
bool activeAbsorption_
Active wave absorption switch.
Definition waveModel.H:158
virtual const vectorField & U() const
Return the latest wave velocity prediction.
Definition waveModel.C:402
tensor Rgl_
Rotation tensor from global to local system.
Definition waveModel.H:88
vectorField U_
Velocity field.
Definition waveModel.H:166
scalar waterDepthRef_
Reference water depth / [m].
Definition waveModel.H:143
virtual void setPaddlePropeties(const scalarField &level, const label facei, scalar &fraction, scalar &z) const
Set the paddle coverage fraction and reference height.
Definition waveModel.C:183
scalarField yPaddle_
Paddle y coordinates / [m].
Definition waveModel.H:108
static const word dictName
Dictionary name.
Definition waveModel.H:294
scalarField xPaddle_
Paddle x coordinates / [m].
Definition waveModel.H:103
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)=0
Calculate the wave model velocity.
waveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition waveModel.C:238
static autoPtr< waveModel > New(const word &dictName, const fvMesh &mesh, const polyPatch &patch)
Return a reference to the selected wave model.
scalar initialDepth_
Initial depth / [m].
Definition waveModel.H:148
scalarField zMax_
Maximum z (point) height per patch face / [m].
Definition waveModel.H:133
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const =0
Set the water level.
static word modelName(const word &patchName)
Utility function to construct the model name.
Definition waveModel.C:42
label currTimeIndex_
Time index used for updating.
Definition waveModel.H:153
labelList faceToPaddle_
Addressing from patch face index to paddle index.
Definition waveModel.H:113
virtual tmp< scalarField > waterLevel() const
Water level.
Definition waveModel.C:120
word alphaName_
Name of phase fraction field, default = "alpha".
Definition waveModel.H:83
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
Definition waveModel.C:155
label nPaddle_
Number of paddles.
Definition waveModel.H:98
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition waveModel.C:287
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & Cp
Definition EEqn.H:7
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
word timeName
Definition getTimeIndex.H:3
Different types of constants.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
static const Identity< scalar > I
Definition Identity.H:100
Tensor< scalar > tensor
Definition symmTensor.H:57
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
volScalarField & alpha
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const scalar xMin
const scalar xMax