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"
35
36using namespace Foam::constant;
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
44}
45
46const Foam::word Foam::waveModel::dictName("waveProperties");
47
48
50{
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& Cf = patch_.faceCentres();
94 const vectorField CfLocal(Rgl_ & Cf);
95 z_ = CfLocal.component(2);
96
97 // Local face extents in z-direction
98 zMin_.setSize(patch_.size());
99 zMax_.setSize(patch_.size());
100 const faceList& faces = patch_.localFaces();
101 forAll(faces, facei)
102 {
103 const face& f = faces[facei];
104 const label nPoint = f.size();
105 zMin_[facei] = CpLocal[f[0]].z();
106 zMax_[facei] = CpLocal[f[0]].z();
107
108 for (label fpi = 1; fpi < nPoint; ++fpi)
109 {
110 const label pointi = f[fpi];
111 zMin_[facei] = min(zMin_[facei], CpLocal[pointi].z());
112 zMax_[facei] = max(zMax_[facei], CpLocal[pointi].z());
113 }
114 }
115
116 // Set minimum z reference level
117 zMin0_ = gMin(zMin_);
118
119 // Local paddle-to-face addressing
120 faceToPaddle_.setSize(patch_.size(), -1);
121 forAll(faceToPaddle_, facei)
122 {
123 faceToPaddle_[facei] = floor((CfLocal[facei].y() - yMin)/paddleDy);
124 }
125}
126
127
129{
130 // Note: initialising as initial depth
131 auto tlevel = tmp<scalarField>::New(nPaddle_, initialDepth_);
132 auto& level = tlevel.ref();
133
134 const volScalarField& alpha =
135 mesh_.lookupObject<volScalarField>(alphaName_);
136 const fvPatchScalarField& alphap = alpha.boundaryField()[patch_.index()];
137 const scalarField alphac(alphap.patchInternalField());
138
139 const scalarField& magSf = alphap.patch().magSf();
140 scalarList paddleMagSf(nPaddle_, Zero);
141 scalarList paddleWettedMagSf(nPaddle_, Zero);
142
143 forAll(alphac, facei)
144 {
145 label paddlei = faceToPaddle_[facei];
146 paddleMagSf[paddlei] += magSf[facei];
147 paddleWettedMagSf[paddlei] += magSf[facei]*alphac[facei];
148 }
149
150 forAll(paddleMagSf, paddlei)
151 {
152 reduce(paddleMagSf[paddlei], sumOp<scalar>());
153 reduce(paddleWettedMagSf[paddlei], sumOp<scalar>());
154 level[paddlei] +=
155 paddleWettedMagSf[paddlei]*zSpan_
156 /(paddleMagSf[paddlei] + ROOTVSMALL);
157 }
158
159 return tlevel;
160}
161
162
164{
165 forAll(alpha_, facei)
166 {
167 const label paddlei = faceToPaddle_[facei];
168 const scalar paddleCalc = level[paddlei];
169
170 const scalar zMin0 = zMin_[facei] - zMin0_;
171 const scalar zMax0 = zMax_[facei] - zMin0_;
172
173 if (zMax0 < paddleCalc)
174 {
175 alpha_[facei] = 1.0;
176 }
177 else if (zMin0 > paddleCalc)
178 {
179 alpha_[facei] = 0.0;
180 }
181 else
182 {
183 scalar dz = paddleCalc - zMin0;
184 alpha_[facei] = dz/(zMax0 - zMin0);
185 }
186 }
187}
188
189
191(
192 const scalarField& level,
193 const label facei,
194 scalar& fraction,
195 scalar& z
196) const
197{
198 const label paddlei = faceToPaddle_[facei];
199 const scalar paddleCalc = level[paddlei];
200 const scalar paddleHeight = min(paddleCalc, waterDepthRef_);
201 const scalar zMin = zMin_[facei] - zMin0_;
202 const scalar zMax = zMax_[facei] - zMin0_;
203
204 fraction = 1;
205 z = 0;
206
207 if (zMax < paddleHeight)
208 {
209 z = z_[facei] - zMin0_;
210 }
211 else if (zMin > paddleCalc)
212 {
213 fraction = -1;
214 }
215 else
216 {
217 if (paddleCalc < waterDepthRef_)
218 {
219 if ((zMax > paddleCalc) && (zMin < paddleCalc))
220 {
221 scalar dz = paddleCalc - zMin;
222 fraction = dz/(zMax - zMin);
223 z = z_[facei] - zMin0_;
224 }
225 }
226 else
227 {
228 if (zMax < paddleCalc)
229 {
230 z = waterDepthRef_;
231 }
232 else if ((zMax > paddleCalc) && (zMin < paddleCalc))
233 {
234 scalar dz = paddleCalc - zMin;
235 fraction = dz/(zMax - zMin);
236 z = waterDepthRef_;
237 }
238 }
239 }
240}
241
242
243// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244
246(
247 const dictionary& dict,
248 const fvMesh& mesh,
249 const polyPatch& patch,
250 const bool readFields
251)
252:
254 (
256 (
257 modelName(patch.name()),
258 Time::timeName(mesh.time().startTime().value()),
259 "uniform",
260 mesh,
261 IOobject::NO_READ,
262 IOobject::AUTO_WRITE
263 )
264 ),
265 mesh_(mesh),
266 patch_(patch),
267 g_(meshObjects::gravity::New(mesh.time()).value()),
268 UName_("U"),
269 alphaName_("alpha"),
270 Rgl_(tensor::I),
271 Rlg_(tensor::I),
272 nPaddle_(1),
273 xPaddle_(),
274 yPaddle_(),
275 z_(),
276 zSpan_(0),
277 zMin_(),
278 zMax_(),
279 waterDepthRef_(0),
280 initialDepth_(0),
281 currTimeIndex_(-1),
282 activeAbsorption_(false),
283 U_(patch.size(), Zero),
284 alpha_(patch.size(), Zero)
285{
286 if (readFields)
287 {
288 readDict(dict);
289 }
290}
291
292
293// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
294
295bool Foam::waveModel::readDict(const dictionary& overrideDict)
296{
298 if (headerOk())
299 {
300 IOdictionary::regIOobject::read();
301 }
302
303 merge(overrideDict);
304
305 readIfPresent("U", UName_);
306 readIfPresent("alpha", alphaName_);
307
308 readEntry("nPaddle", nPaddle_);
309 if (nPaddle_ < 1)
310 {
312 << "Number of paddles must be greater than zero. Supplied"
313 << " value nPaddles = " << nPaddle_
314 << exit(FatalIOError);
315 }
316
317 readIfPresent("initialDepth", initialDepth_);
318
319 // Need to initialise the geometry before calling waterLevel()
320 initialiseGeometry();
321
322 // Set the reference water depth
323 if (!readIfPresent("waterDepthRef", waterDepthRef_))
324 {
325 scalar waterDepth = 0;
326 if (readIfPresent("waterDepth", waterDepth))
327 {
328 waterDepthRef_ = waterDepth;
329 }
330 else
331 {
332 const scalarField level(waterLevel());
333 if (level.size())
334 {
335 waterDepthRef_ = level.first();
336 }
337 }
338
339 // Avoid potential zero...
340 waterDepthRef_ += SMALL;
341
342 // Insert the reference water depth into [this] to enable restart
343 add("waterDepthRef", waterDepthRef_);
344 }
345
346 return true;
347}
348
349
350void Foam::waveModel::correct(const scalar t)
351{
352 if (mesh_.time().timeIndex() != currTimeIndex_)
353 {
354 Info<< "Updating " << type() << " wave model for patch "
355 << patch_.name() << endl;
356
357 // Time ramp weight
358 const scalar tCoeff = timeCoeff(t);
359
360 // Reset the velocity and phase fraction fields
361 U_ = vector::zero;
362 alpha_ = 0;
363
364 // Update the calculated water level field
365 scalarField calculatedLevel(nPaddle_, Zero);
366
367 if (patch_.size())
368 {
369 // Set wave level
370 setLevel(t, tCoeff, calculatedLevel);
371
372 // Update the velocity field
373 setVelocity(t, tCoeff, calculatedLevel);
374
375 // Update the phase fraction field
376 setAlpha(calculatedLevel);
377 }
378
379 if (activeAbsorption_)
380 {
381 const scalarField activeLevel(this->waterLevel());
382
383 forAll(U_, facei)
384 {
385 const label paddlei = faceToPaddle_[facei];
386
387 if (zMin_[facei] - zMin0_ < activeLevel[paddlei])
388 {
389 scalar UCorr =
390 (calculatedLevel[paddlei] - activeLevel[paddlei])
391 *sqrt(mag(g_)/activeLevel[paddlei]);
392
393 U_[facei].x() += UCorr;
394 }
395 else
396 {
397 U_[facei].x() = 0;
398 }
399 }
400 }
401
402 // Transform velocity into global coordinate system
403 U_ = Rlg_ & U_;
404
405 currTimeIndex_ = mesh_.time().timeIndex();
406 }
407}
408
409
411{
412 return U_;
413}
414
415
417{
418 return alpha_;
419}
420
421
423{
424 os << "Wave model: patch " << patch_.name() << nl
425 << " Type : " << type() << nl
426 << " Velocity field name : " << UName_ << nl
427 << " Phase fraction field name : " << alphaName_ << nl
428 << " Transformation from local to global system : " << Rlg_ << nl
429 << " Number of paddles: " << nPaddle_ << nl
430 << " Reference water depth : " << waterDepthRef_ << nl
431 << " Active absorption: " << activeAbsorption_ << nl;
432}
433
434
435// ************************************************************************* //
scalar y
Foam::label startTime
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const word & modelName() const
Return const access to the name of the sub-model.
Definition: subModelBase.C:107
A class for managing temporary objects.
Definition: tmp.H:65
Base class for waveModels.
Definition: waveModel.H:61
virtual void initialiseGeometry()
Initialise.
Definition: waveModel.C:57
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
Definition: waveModel.C:416
virtual const vectorField & U() const
Return the latest wave velocity prediction.
Definition: waveModel.C:410
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:191
static const word dictName
Dictionary name.
Definition: waveModel.H:230
virtual tmp< scalarField > waterLevel() const
Water level.
Definition: waveModel.C:128
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
Definition: waveModel.C:163
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: waveModel.C:295
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & Cp
Definition: EEqn.H:7
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
OBJstream os(runTime.globalPath()/outputName)
word timeName
Definition: getTimeIndex.H:3
Different types of constants.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
messageStream Info
Information stream (stdout output on master, null elsewhere)
static const Identity< scalar > I
Definition: Identity.H:94
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Tensor< scalar > tensor
Definition: symmTensor.H:61
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
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:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
propsDict readIfPresent("fields", acceptFields)
const scalar xMin
Definition: createFields.H:34
const scalar xMax
Definition: createFields.H:35