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-2022 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 {
38  defineTypeNameAndDebug(waveSuperposition, 0);
39  defineRunTimeSelectionTable(waveSuperposition, objectRegistry);
41  (
42  waveSuperposition,
43  waveSuperposition,
44  objectRegistry
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 {
61  db().lookupObject<uniformDimensionedVectorField>("g");
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("origin")),
185  direction_(lookup("direction")),
186  waveModels_(),
187  waveAngles_(),
188  UMean_(Function1<vector>::New("UMean", *this)),
189  scale_
190  (
191  found("scale")
192  ? Function1<scalar>::New("scale", *this)
193  : autoPtr<Function1<scalar>>()
194  ),
195  crossScale_
196  (
197  found("crossScale")
198  ? Function1<scalar>::New("crossScale", *this)
199  : autoPtr<Function1<scalar>>()
200  ),
201  heightAboveWave_(lookupOrDefault<Switch>("heightAboveWave", false))
202 {
204  {
207  (
208  IOobject
209  (
210  "g",
211  db.time().constant(),
212  db,
215  ),
217  );
218 
219  gPtr->store();
220  }
221 
224 
225  const PtrList<entry> waveEntries(lookup("waves"));
226 
227  waveModels_.setSize(waveEntries.size());
228  waveAngles_.setSize(waveEntries.size());
229 
230  forAll(waveEntries, wavei)
231  {
232  const dictionary waveDict = waveEntries[wavei].dict();
233  waveModels_.set
234  (
235  wavei,
236  waveModel::New(waveDict.dictName(), waveDict, mag(g.value()))
237  );
238  waveAngles_[wavei] = waveDict.lookup<scalar>("angle");
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
244 
246 {}
247 
248 
249 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
250 
252 (
253  const scalar t,
254  const vectorField& p
255 ) const
256 {
257  tensor axes;
258  vector drift;
259  vectorField xyz(p.size());
260  transformation(t, p, axes, drift, xyz);
261 
262  return
263  xyz.component(2) + drift.z()
264  - elevation
265  (
266  t,
267  vector2D(drift.x(), drift.y()),
268  zip(xyz.component(0), xyz.component(1))
269  );
270 }
271 
272 
274 (
275  const scalar t,
276  const vectorField& p
277 ) const
278 {
279  tensor axes;
280  vector drift;
281  vectorField xyz(p.size());
282  transformation(t, p, axes, drift, xyz);
283 
284  if (heightAboveWave_)
285  {
286  xyz.replace(2, height(t, p));
287  }
288 
289  return UMean_->value(t) + (velocity(t, drift, xyz) & axes);
290 }
291 
292 
294 (
295  const scalar t,
296  const vectorField& p
297 ) const
298 {
299  tensor axes;
300  vector drift;
301  vectorField xyz(p.size());
302  transformation(t, p, axes, drift, xyz);
303 
304  axes = tensor(- axes.x(), - axes.y(), axes.z());
305 
306  if (heightAboveWave_)
307  {
308  xyz.replace(2, height(t, p));
309  }
310 
311  xyz.replace(2, - xyz.component(2));
312 
313  return UMean_->value(t) + (velocity(t, drift, xyz) & axes);
314 }
315 
316 
318 {
319  writeEntry(os, "origin", origin_);
320  writeEntry(os, "direction", direction_);
321  writeKeyword(os, "waves") << nl << token::BEGIN_LIST << nl << incrIndent;
322  forAll(waveModels_, wavei)
323  {
324  writeKeyword(os, waveModels_[wavei].type()) << nl << indent
326  waveModels_[wavei].write(os);
327  writeKeyword(os, "angle") << waveAngles_[wavei] << token::END_STATEMENT
328  << nl << decrIndent << indent << token::END_BLOCK << nl;
329  }
331  writeEntry(os, UMean_());
332  if (scale_.valid())
333  {
334  writeEntry(os, scale_());
335  }
336  if (crossScale_.valid())
337  {
338  writeEntry(os, crossScale_());
339  }
340  if (heightAboveWave_)
341  {
342  writeEntry(os, "heightAboveWave", heightAboveWave_);
343  }
344 }
345 
346 
347 // ************************************************************************* //
tmp< vector2DField > zip(const tmp< scalarField > &x, const tmp< scalarField > &y)
Definition: vector2DField.C:31
Run-time selectable general function of one variable.
Definition: Function1.H:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
const vector origin_
The origin of the wave coordinate system.
void write(Ostream &) const
Write.
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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.
const vector direction_
The direction of the wave coordinate system.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const autoPtr< Function1< vector > > UMean_
Mean velocity.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
scalarList waveAngles_
The angle relative to the direction at which the waves propagate.
Vector< Cmpt > x() const
Definition: TensorI.H:289
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Definition: vector2D.H:49
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
~waveSuperposition()
Destructor.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
bool foundObject(const word &name) const
Is the named Type found?
virtual tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
Vector< Cmpt > y() const
Definition: TensorI.H:296
const Cmpt & z() const
Definition: VectorI.H:87
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:467
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
const Switch heightAboveWave_
Calculate wave properties using the height above the wave (true) or.
const Cmpt & x() const
Definition: Vector2DI.H:68
Macros for easy insertion into run-time selection tables.
const Cmpt & y() const
Definition: VectorI.H:81
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:121
Vector< Cmpt > z() const
Definition: TensorI.H:303
scalar y
const autoPtr< Function1< scalar > > crossScale_
Scaling perpendicular to the local x-direction.
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void transformation(const scalar t, const vectorField &p, tensor &axes, vector &drift, vectorField &xyz) const
Get the transformation to actual coordinates.
stressControl lookup("compactNormalStress") >> compactNormalStress
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
dimensionedScalar cos(const dimensionedScalar &ds)
const dimensionSet dimAcceleration
const Cmpt & y() const
Definition: Vector2DI.H:74
const autoPtr< Function1< scalar > > scale_
Scaling in the local x-direction.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
const Type & value() const
Return const reference to value.
static const zero Zero
Definition: zero.H:97
const Cmpt & x() const
Definition: VectorI.H:75
PtrList< waveModel > waveModels_
Wave models to superimpose.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:455
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:260
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void setSize(const label)
Reset size of List.
Definition: List.C:281
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dimensioned< scalar > mag(const dimensioned< Type > &)
static autoPtr< waveModel > New(const dictionary &dict, const scalar g)
Select.
Definition: waveModelNew.C:31
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
tmp< scalarField > scale(const vector2DField &xy) const
Get the scaling factor, calculated from the optional scaling.
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
const dimensionedVector & g
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
bool found
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.
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
const word dictName("noiseDict")
Namespace for OpenFOAM.
waveSuperposition(const objectRegistry &db)
Construct from a database.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864