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 -------------------------------------------------------------------------------
11 License
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 
36 using namespace Foam::constant;
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
44 }
45 
46 const 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 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);
120  forAll(faceToPaddle_, facei)
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 
162 void 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  }
237  }
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 :
252  IOdictionary
253  (
254  IOobject
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 
294 bool Foam::waveModel::readDict(const dictionary& overrideDict)
295 {
296  readOpt(IOobject::READ_IF_PRESENT);
297  if (headerOk())
298  {
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 
349 void 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 
421 void Foam::waveModel::info(Ostream& os) const
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 // ************************************************************************* //
Different types of constants.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
Definition: waveModel.C:408
dictionary dict
virtual void initialiseGeometry()
Initialise.
Definition: waveModel.C:50
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
Definition: waveModel.C:155
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
virtual void correct(const scalar t)
Correct the model for time, t[s].
Definition: waveModel.C:342
Type gMin(const FieldField< Field, Type > &f)
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition: OSstream.H:134
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual tmp< scalarField > waterLevel() const
Water level.
Definition: waveModel.C:120
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
dimensionedScalar sqrt(const dimensionedScalar &ds)
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const scalar xMin
Definition: createFields.H:34
static word modelName(const word &patchName)
Utility function to construct the model name.
Definition: waveModel.C:42
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
scalar y
virtual const vectorField & U() const
Return the latest wave velocity prediction.
Definition: waveModel.C:402
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Base class for waveModels.
Definition: waveModel.H:51
static const Identity< scalar > I
Definition: Identity.H:100
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.
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Reading is optional [identical to LAZY_READ].
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition: IOobject.H:963
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:608
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const volScalarField & Cp
Definition: EEqn.H:7
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
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
labelList f(nPoints)
const scalar xMax
Definition: createFields.H:35
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
Type gAverage(const FieldField< Field, Type > &f)
meshDefDict readIfPresent("polyMeshPatches", polyPatchNames)
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: waveModel.C:287
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Foam::label startTime
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
waveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: waveModel.C:238