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-2019 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 "Time.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 const Foam::word Foam::waveSuperposition::dictName("waveProperties");
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(waveSuperposition, 0);
38  defineRunTimeSelectionTable(waveSuperposition, objectRegistry);
40  (
41  waveSuperposition,
42  waveSuperposition,
43  objectRegistry
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 (
52  const scalar t,
53  const vectorField& p,
54  tensor& axes,
55  vectorField& xyz
56 ) const
57 {
59  db().lookupObject<uniformDimensionedVectorField>("g");
60  const scalar magG = mag(g.value());
61  const vector gHat = g.value()/magG;
62 
63  const vector dSurf = direction_ - gHat*(gHat & direction_);
64  const scalar magDSurf = mag(dSurf);
65  const vector dSurfHat = direction_/magDSurf;
66 
67  axes = tensor(dSurfHat, - gHat ^ dSurfHat, - gHat);
68 
69  xyz = axes & (p - origin_ - UMean_->integrate(0, t));
70 }
71 
72 
74 (
75  const scalar t,
76  const vector2DField& xy
77 ) const
78 {
79  scalarField result(xy.size(), 0);
80 
81  forAll(waveModels_, wavei)
82  {
83  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
84  result += waveModels_[wavei].elevation(t, d & xy);
85  }
86 
87  return scale(xy)*result;
88 }
89 
90 
92 (
93  const scalar t,
94  const vectorField& xyz
95 ) const
96 {
97  vectorField result(xyz.size(), vector::zero);
98 
99  forAll(waveModels_, wavei)
100  {
101  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
102  const vector2DField xz
103  (
104  zip
105  (
106  d & zip(xyz.component(0), xyz.component(1)),
108  )
109  );
110  const vector2DField uw
111  (
112  waveModels_[wavei].velocity(t, xz)
113  );
114  result += zip
115  (
116  d.x()*uw.component(0),
117  d.y()*uw.component(0),
118  uw.component(1)
119  );
120  }
121 
122  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
123 
124  return s*result;
125 }
126 
127 
129 (
130  const scalar t,
131  const vectorField& xyz
132 ) const
133 {
134  scalarField result(xyz.size(), 0);
135 
136  forAll(waveModels_, wavei)
137  {
138  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
139  const vector2DField xz
140  (
141  zip
142  (
143  d & zip(xyz.component(0), xyz.component(1)),
145  )
146  );
147  const vector2DField uw
148  (
149  waveModels_[wavei].velocity(t, xz)
150  );
151  result += waveModels_[wavei].pressure(t, xz);
152  }
153 
154  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
155 
156  return s*result;
157 }
158 
159 
161 (
162  const vector2DField& xy
163 ) const
164 {
165  tmp<scalarField> tResult(new scalarField(xy.size(), 1));
166  scalarField& result = tResult.ref();
167 
168  if (scale_.valid())
169  {
170  const scalarField x(xy.component(0));
171  forAll(result, i)
172  {
173  result[i] *= scale_->value(x[i]);
174  }
175  }
176 
177  if (crossScale_.valid())
178  {
179  const scalarField y(xy.component(1));
180  forAll(result, i)
181  {
182  result[i] *= crossScale_->value(y[i]);
183  }
184  }
185 
186  return tResult;
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191 
193 :
195  (
196  IOobject
197  (
198  dictName,
199  db.time().constant(),
200  db,
201  IOobject::MUST_READ,
202  IOobject::NO_WRITE
203  )
204  ),
205  origin_(lookup("origin")),
206  direction_(lookup("direction")),
207  waveModels_(),
208  waveAngles_(),
209  UMean_(Function1<vector>::New("UMean", *this)),
210  scale_
211  (
212  found("scale")
213  ? Function1<scalar>::New("scale", *this)
214  : autoPtr<Function1<scalar>>()
215  ),
216  crossScale_
217  (
218  found("crossScale")
219  ? Function1<scalar>::New("crossScale", *this)
220  : autoPtr<Function1<scalar>>()
221  ),
222  heightAboveWave_(lookupOrDefault<Switch>("heightAboveWave", false))
223 {
224  const PtrList<entry> waveEntries(lookup("waves"));
225 
226  waveModels_.setSize(waveEntries.size());
227  waveAngles_.setSize(waveEntries.size());
228 
229  forAll(waveEntries, wavei)
230  {
231  const dictionary waveDict = waveEntries[wavei].dict();
232  waveModels_.set
233  (
234  wavei,
235  waveModel::New(waveDict.dictName(), db, waveDict)
236  );
237  waveAngles_[wavei] = waveDict.lookup<scalar>("angle");
238  }
239 }
240 
241 
242 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
243 
245 {}
246 
247 
248 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
249 
251 (
252  const scalar t,
253  const vectorField& p
254 ) const
255 {
256  tensor axes;
257  vectorField xyz(p.size());
258  transformation(t, p, axes, xyz);
259 
260  return
261  xyz.component(2)
262  - elevation(t, zip(xyz.component(0), xyz.component(1)));
263 }
264 
265 
267 (
268  const scalar t,
269  const vectorField& p
270 ) const
271 {
272  tensor axes;
273  vectorField xyz(p.size());
274  transformation(t, p, axes, xyz);
275 
276  if (heightAboveWave_)
277  {
278  xyz.replace(2, height(t, p));
279  }
280 
281  return UMean_->value(t) + (velocity(t, xyz) & axes);
282 }
283 
284 
286 (
287  const scalar t,
288  const vectorField& p
289 ) const
290 {
291  tensor axes;
292  vectorField xyz(p.size());
293  transformation(t, p, axes, xyz);
294 
295  axes = tensor(- axes.x(), - axes.y(), axes.z());
296 
297  if (heightAboveWave_)
298  {
299  xyz.replace(2, height(t, p));
300  }
301 
302  xyz.replace(2, - xyz.component(2));
303 
304  return UMean_->value(t) + (velocity(t, xyz) & axes);
305 }
306 
307 
309 (
310  const scalar t,
311  const vectorField& p
312 ) const
313 {
314  tensor axes;
315  vectorField xyz(p.size());
316  transformation(t, p, axes, xyz);
317 
318  if (heightAboveWave_)
319  {
320  xyz.replace(2, height(t, p));
321  }
322 
323  return pressure(t, xyz);
324 }
325 
326 
328 (
329  const scalar t,
330  const vectorField& p
331 ) const
332 {
333  tensor axes;
334  vectorField xyz(p.size());
335  transformation(t, p, axes, 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 pressure(t, xyz);
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, UMean_());
365  if (scale_.valid())
366  {
367  writeEntry(os, scale_());
368  }
369  if (crossScale_.valid())
370  {
371  writeEntry(os, crossScale_());
372  }
373  if (heightAboveWave_)
374  {
375  writeEntry(os, "heightAboveWave", heightAboveWave_);
376  }
377 }
378 
379 
380 // ************************************************************************* //
tmp< vector2DField > zip(const tmp< scalarField > &x, const tmp< scalarField > &y)
Definition: vector2DField.C:31
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
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
virtual tmp< scalarField > pLiquid(const scalar t, const vectorField &p) const
Get the liquid pressure at a given time and global positions.
tmp< scalarField > pressure(const scalar t, const vectorField &xyz) const
Get the wave pressure at a given time and local coordinates. Local.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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:174
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
~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
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
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:472
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.
tmp< scalarField > elevation(const scalar t, const vector2DField &xy) const
Get the wave elevation relative to the mean at a given time and.
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:52
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
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))
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)
virtual tmp< scalarField > pGas(const scalar t, const vectorField &p) const
Get the gas pressure at a given time and global positions.
const Cmpt & y() const
Definition: Vector2DI.H:74
virtual Type value(const scalar x) const =0
Return value as a function of (scalar) independent variable.
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 Type & value() const
Return const reference to value.
const word dictName("particleTrackDict")
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
PtrList< waveModel > waveModels_
Wave models to superimpose.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:460
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)
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
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
static autoPtr< waveModel > New(const objectRegistry &db, const dictionary &dict)
Select.
Definition: waveModelNew.C:31
tmp< vectorField > velocity(const scalar t, const vectorField &xyz) const
Get the wave velocity at a given time and local coordinates. Local.
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 > &)
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.
void transformation(const scalar t, const vectorField &p, tensor &axes, vectorField &xyz) const
Get the transformation to actual coordinates.
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:322
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:92
bool found
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
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:812