waveSuperposition.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "waveSuperposition.H"
27 #include "dimensionedTypes.H"
28 #include "Time.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 const Foam::word Foam::waveSuperposition::dictName("waveProperties");
35 
36 namespace Foam
37 {
41  (
45  );
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 (
53  const scalar t,
54  const vectorField& p,
55  tensor& axes,
56  vector& drift,
57  vectorField& xyz
58 ) const
59 {
62  const scalar magG = mag(g.value());
63  const vector gHat = g.value()/magG;
64 
65  const vector dSurf = direction_ - gHat*(gHat & direction_);
66  const scalar magDSurf = mag(dSurf);
67  const vector dSurfHat = direction_/magDSurf;
68 
69  axes = tensor(dSurfHat, - gHat ^ dSurfHat, - gHat);
70 
71  drift = - axes & UMean_->integral(0, t);
72 
73  xyz = axes & (p - origin_);
74 }
75 
76 
78 (
79  const scalar t,
80  const vector2D& drift,
81  const vector2DField& xy
82 ) const
83 {
84  scalarField result(xy.size(), 0);
85 
86  forAll(waveModels_, wavei)
87  {
88  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
89  result += waveModels_[wavei].elevation(t, d & (drift + xy));
90  }
91 
92  return scale(xy)*result;
93 }
94 
95 
97 (
98  const scalar t,
99  const vector& drift,
100  const vectorField& xyz
101 ) const
102 {
103  vectorField result(xyz.size(), vector::zero);
104 
105  forAll(waveModels_, wavei)
106  {
107  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
108  const vector2DField xz
109  (
110  zip
111  (
112  d
113  & zip
114  (
115  drift.x() + xyz.component(0),
116  drift.y() + xyz.component(1)
117  ),
118  tmp<scalarField>(drift.z() + xyz.component(2))
119  )
120  );
121  const vector2DField uw
122  (
123  waveModels_[wavei].velocity(t, xz)
124  );
125  result += zip
126  (
127  d.x()*uw.component(0),
128  d.y()*uw.component(0),
129  uw.component(1)
130  );
131  }
132 
133  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
134 
135  return s*result;
136 }
137 
138 
140 (
141  const vector2DField& xy
142 ) const
143 {
144  tmp<scalarField> tResult(new scalarField(xy.size(), 1));
145  scalarField& result = tResult.ref();
146 
147  if (scale_.valid())
148  {
149  const scalarField x(xy.component(0));
150  forAll(result, i)
151  {
152  result[i] *= scale_->value(x[i]);
153  }
154  }
155 
156  if (crossScale_.valid())
157  {
158  const scalarField y(xy.component(1));
159  forAll(result, i)
160  {
161  result[i] *= crossScale_->value(y[i]);
162  }
163  }
164 
165  return tResult;
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
172 :
174  (
175  IOobject
176  (
177  dictName,
178  db.time().constant(),
179  db,
180  IOobject::MUST_READ,
181  IOobject::NO_WRITE
182  )
183  ),
184  origin_(lookup<vector>("origin", dimLength)),
185  direction_(lookup<vector>("direction", dimless)),
186  waveModels_(),
187  waveAngles_(),
188  UMean_
189  (
191  (
192  "UMean",
193  db.time().userUnits(),
194  dimVelocity,
195  *this
196  )
197  ),
198  scale_
199  (
200  found("scale")
201  ? Function1<scalar>::New
202  (
203  "scale",
204  db.time().userUnits(),
205  dimless,
206  *this
207  )
208  : autoPtr<Function1<scalar>>()
209  ),
210  crossScale_
211  (
212  found("crossScale")
213  ? Function1<scalar>::New
214  (
215  "crossScale",
216  db.time().userUnits(),
217  dimless,
218  *this
219  )
220  : autoPtr<Function1<scalar>>()
221  ),
222  heightAboveWave_(lookupOrDefault<Switch>("heightAboveWave", false))
223 {
225  {
228  (
229  IOobject
230  (
231  "g",
232  db.time().constant(),
233  db,
236  ),
238  );
239 
240  gPtr->store();
241  }
242 
245 
246  const PtrList<entry> waveEntries(lookup("waves"));
247 
248  waveModels_.setSize(waveEntries.size());
249  waveAngles_.setSize(waveEntries.size());
250 
251  forAll(waveEntries, wavei)
252  {
253  const dictionary waveDict = waveEntries[wavei].dict();
254  waveModels_.set
255  (
256  wavei,
257  waveModel::New(waveDict.dictName(), waveDict, mag(g.value()))
258  );
259  waveAngles_[wavei] = waveDict.lookup<scalar>("angle");
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
265 
267 {}
268 
269 
270 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
271 
272 Foam::scalar Foam::waveSuperposition::maxWaveSpeed(const scalar t) const
273 {
274  scalar maxCelerity = 0;
275  forAll(waveModels_, wavei)
276  {
277  maxCelerity = max(waveModels_[wavei].celerity(), maxCelerity);
278  }
279 
280  return mag(maxCelerity + (direction_ & UMean_->value(t)));
281 }
282 
283 
285 (
286  const scalar t,
287  const vectorField& p
288 ) const
289 {
290  tensor axes;
291  vector drift;
292  vectorField xyz(p.size());
293  transformation(t, p, axes, drift, xyz);
294 
295  return
296  xyz.component(2) + drift.z()
297  - elevation
298  (
299  t,
300  vector2D(drift.x(), drift.y()),
301  zip(xyz.component(0), xyz.component(1))
302  );
303 }
304 
305 
307 (
308  const scalar t,
309  const vectorField& p
310 ) const
311 {
312  tensor axes;
313  vector drift;
314  vectorField xyz(p.size());
315  transformation(t, p, axes, drift, xyz);
316 
317  if (heightAboveWave_)
318  {
319  xyz.replace(2, height(t, p));
320  }
321 
322  return UMean_->value(t) + (velocity(t, drift, xyz) & axes);
323 }
324 
325 
327 (
328  const scalar t,
329  const vectorField& p
330 ) const
331 {
332  tensor axes;
333  vector drift;
334  vectorField xyz(p.size());
335  transformation(t, p, axes, drift, xyz);
336 
337  axes = tensor(- axes.x(), - axes.y(), axes.z());
338 
339  if (heightAboveWave_)
340  {
341  xyz.replace(2, height(t, p));
342  }
343 
344  xyz.replace(2, - xyz.component(2));
345 
346  return UMean_->value(t) + (velocity(t, drift, xyz) & axes);
347 }
348 
349 
351 {
352  writeEntry(os, "origin", origin_);
353  writeEntry(os, "direction", direction_);
354  writeKeyword(os, "waves") << nl << token::BEGIN_LIST << nl << incrIndent;
355  forAll(waveModels_, wavei)
356  {
357  writeKeyword(os, waveModels_[wavei].type()) << nl << indent
359  waveModels_[wavei].write(os);
360  writeKeyword(os, "angle") << waveAngles_[wavei] << token::END_STATEMENT
361  << nl << decrIndent << indent << token::END_BLOCK << nl;
362  }
364  writeEntry(os, db().time().userUnits(), dimVelocity, UMean_());
365  if (scale_.valid())
366  {
367  writeEntry(os, db().time().userUnits(), dimless, scale_());
368  }
369  if (crossScale_.valid())
370  {
371  writeEntry(os, db().time().userUnits(), dimless, crossScale_());
372  }
373  if (heightAboveWave_)
374  {
375  writeEntry(os, "heightAboveWave", heightAboveWave_);
376  }
377 }
378 
379 
380 // ************************************************************************* //
scalar y
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:491
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:479
Run-time selectable general function of one variable.
Definition: Function1.H:125
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual Ostream & write(const char)=0
Write character.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Vector< Cmpt > z() const
Definition: TensorI.H:303
Vector< Cmpt > y() const
Definition: TensorI.H:296
Vector< Cmpt > x() const
Definition: TensorI.H:289
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
Type & value()
Return a reference to the value.
const Cmpt & y() const
Definition: Vector2DI.H:74
const Cmpt & x() const
Definition: Vector2DI.H:68
static const Form zero
Definition: VectorSpace.H:118
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Registry of regIOobjects.
const Time & time() const
Return time.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
bool foundObject(const word &name) const
Is the named Type in registry.
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
@ END_STATEMENT
Definition: token.H:108
@ BEGIN_LIST
Definition: token.H:109
@ END_LIST
Definition: token.H:110
static autoPtr< waveModel > New(const dictionary &dict, const scalar g)
Select.
Definition: waveModelNew.C:31
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
tmp< scalarField > elevation(const scalar t, const vector2D &drift, const vector2DField &xy) const
Get the wave elevation relative to the mean at a given time and.
const autoPtr< Function1< vector > > UMean_
Mean velocity.
const vector origin_
The origin of the wave coordinate system.
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
virtual tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
virtual tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
waveSuperposition(const objectRegistry &db)
Construct from a database.
static const word dictName
The name of the dictionary.
tmp< scalarField > scale(const vector2DField &xy) const
Get the scaling factor, calculated from the optional scaling.
void transformation(const scalar t, const vectorField &p, tensor &axes, vector &drift, vectorField &xyz) const
Get the transformation to actual coordinates.
scalarList waveAngles_
The angle relative to the direction at which the waves propagate.
const vector direction_
The direction of the wave coordinate system.
virtual bool write(const bool write=true) const
Inherit write from regIOobject.
PtrList< waveModel > waveModels_
Wave models to superimpose.
scalar maxWaveSpeed(const scalar t) const
Return the maximum wave speed for the given time t.
tmp< vectorField > velocity(const scalar t, const vector &drift, const vectorField &xyz) const
Get the wave velocity at a given time and local coordinates. Local.
A class for handling words, derived from string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Definition: vector2D.H:49
static const zero Zero
Definition: zero.H:97
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimLength
const dimensionSet dimAcceleration
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
UniformDimensionedField< vector > uniformDimensionedVectorField
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVelocity
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
tmp< vector2DField > zip(const tmp< scalarField > &x, const tmp< scalarField > &y)
Definition: vector2DField.C:31
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField & p