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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 const Foam::word Foam::waveSuperposition::dictName("waveProperties");
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(waveSuperposition, 0);
37  defineRunTimeSelectionTable(waveSuperposition, objectRegistry);
39  (
40  waveSuperposition,
41  waveSuperposition,
42  objectRegistry
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 (
51  const scalar t,
52  const vectorField& p,
53  tensor& axes,
54  vectorField& xyz
55 ) const
56 {
58  db().lookupObject<uniformDimensionedVectorField>("g");
59  const scalar magG = mag(g.value());
60  const vector gHat = g.value()/magG;
61 
62  const vector dSurf = direction_ - gHat*(gHat & direction_);
63  const scalar magDSurf = mag(dSurf);
64  const vector dSurfHat = direction_/magDSurf;
65 
66  axes = tensor(dSurfHat, - gHat ^ dSurfHat, - gHat);
67 
68  xyz = axes & (p - origin_ - UMean_->integrate(0, t));
69 }
70 
71 
73 (
74  const scalar t,
75  const vector2DField& xy
76 ) const
77 {
78  scalarField result(xy.size(), 0);
79 
80  forAll(waveModels_, wavei)
81  {
82  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
83  result += waveModels_[wavei].elevation(t, d & xy);
84  }
85 
86  return scale(xy)*result;
87 }
88 
89 
91 (
92  const scalar t,
93  const vectorField& xyz
94 ) const
95 {
96  vectorField result(xyz.size(), vector::zero);
97 
98  forAll(waveModels_, wavei)
99  {
100  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
101  const vector2DField xz
102  (
103  zip
104  (
105  d & zip(xyz.component(0), xyz.component(1)),
107  )
108  );
109  const vector2DField uw
110  (
111  waveModels_[wavei].velocity(t, xz)
112  );
113  result += zip
114  (
115  d.x()*uw.component(0),
116  d.y()*uw.component(0),
117  uw.component(1)
118  );
119  }
120 
121  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
122 
123  return s*result;
124 }
125 
126 
128 (
129  const scalar t,
130  const vectorField& xyz
131 ) const
132 {
133  scalarField result(xyz.size(), 0);
134 
135  forAll(waveModels_, wavei)
136  {
137  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
138  const vector2DField xz
139  (
140  zip
141  (
142  d & zip(xyz.component(0), xyz.component(1)),
144  )
145  );
146  const vector2DField uw
147  (
148  waveModels_[wavei].velocity(t, xz)
149  );
150  result += waveModels_[wavei].pressure(t, xz);
151  }
152 
153  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
154 
155  return s*result;
156 }
157 
158 
160 (
161  const vector2DField& xy
162 ) const
163 {
164  tmp<scalarField> tResult(new scalarField(xy.size(), 1));
165  scalarField& result = tResult.ref();
166 
167  if (scale_.valid())
168  {
169  const scalarField x(xy.component(0));
170  forAll(result, i)
171  {
172  result[i] *= scale_->value(x[i]);
173  }
174  }
175 
176  if (crossScale_.valid())
177  {
178  const scalarField y(xy.component(1));
179  forAll(result, i)
180  {
181  result[i] *= crossScale_->value(y[i]);
182  }
183  }
184 
185  return tResult;
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
190 
192 :
194  (
195  IOobject
196  (
197  dictName,
198  db.time().constant(),
199  db,
200  IOobject::MUST_READ,
201  IOobject::NO_WRITE
202  )
203  ),
204  origin_(lookup("origin")),
205  direction_(lookup("direction")),
206  waveModels_(),
207  waveAngles_(),
208  UMean_(Function1<vector>::New("UMean", *this)),
209  scale_
210  (
211  found("scale")
212  ? Function1<scalar>::New("scale", *this)
213  : autoPtr<Function1<scalar>>()
214  ),
215  crossScale_
216  (
217  found("crossScale")
218  ? Function1<scalar>::New("crossScale", *this)
219  : autoPtr<Function1<scalar>>()
220  ),
221  heightAboveWave_(lookupOrDefault<Switch>("heightAboveWave", false))
222 {
223  const PtrList<entry> waveEntries(lookup("waves"));
224 
225  waveModels_.setSize(waveEntries.size());
226  waveAngles_.setSize(waveEntries.size());
227 
228  forAll(waveEntries, wavei)
229  {
230  const dictionary waveDict = waveEntries[wavei].dict();
231  waveModels_.set
232  (
233  wavei,
234  waveModel::New(waveDict.dictName(), db, waveDict)
235  );
236  waveAngles_[wavei] = readScalar(waveDict.lookup("angle"));
237  }
238 }
239 
240 
241 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
242 
244 {}
245 
246 
247 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
248 
250 (
251  const scalar t,
252  const vectorField& p
253 ) const
254 {
255  tensor axes;
256  vectorField xyz(p.size());
257  transformation(t, p, axes, xyz);
258 
259  return
260  xyz.component(2)
261  - elevation(t, zip(xyz.component(0), xyz.component(1)));
262 }
263 
264 
266 (
267  const scalar t,
268  const vectorField& p
269 ) const
270 {
271  tensor axes;
272  vectorField xyz(p.size());
273  transformation(t, p, axes, xyz);
274 
275  if (heightAboveWave_)
276  {
277  xyz.replace(2, height(t, p));
278  }
279 
280  return UMean_->value(t) + (velocity(t, xyz) & axes);
281 }
282 
283 
285 (
286  const scalar t,
287  const vectorField& p
288 ) const
289 {
290  tensor axes;
291  vectorField xyz(p.size());
292  transformation(t, p, axes, xyz);
293 
294  axes = tensor(- axes.x(), - axes.y(), axes.z());
295 
296  if (heightAboveWave_)
297  {
298  xyz.replace(2, height(t, p));
299  }
300 
301  xyz.replace(2, - xyz.component(2));
302 
303  return UMean_->value(t) + (velocity(t, xyz) & axes);
304 }
305 
306 
308 (
309  const scalar t,
310  const vectorField& p
311 ) const
312 {
313  tensor axes;
314  vectorField xyz(p.size());
315  transformation(t, p, axes, xyz);
316 
317  if (heightAboveWave_)
318  {
319  xyz.replace(2, height(t, p));
320  }
321 
322  return pressure(t, xyz);
323 }
324 
325 
327 (
328  const scalar t,
329  const vectorField& p
330 ) const
331 {
332  tensor axes;
333  vectorField xyz(p.size());
334  transformation(t, p, axes, xyz);
335 
336  axes = tensor(- axes.x(), - axes.y(), axes.z());
337 
338  if (heightAboveWave_)
339  {
340  xyz.replace(2, height(t, p));
341  }
342 
343  xyz.replace(2, - xyz.component(2));
344 
345  return pressure(t, xyz);
346 }
347 
348 
350 {
351  writeEntry(os, "origin", origin_);
352  writeEntry(os, "direction", direction_);
353  os.writeKeyword("waves") << nl << token::BEGIN_LIST << nl << incrIndent;
354  forAll(waveModels_, wavei)
355  {
356  os.writeKeyword(waveModels_[wavei].type()) << nl << indent
358  waveModels_[wavei].write(os);
359  os.writeKeyword("angle") << waveAngles_[wavei] << token::END_STATEMENT
360  << nl << decrIndent << indent << token::END_BLOCK << nl;
361  }
363  writeEntry(os, UMean_());
364  if (scale_.valid())
365  {
366  writeEntry(os, scale_());
367  }
368  if (crossScale_.valid())
369  {
370  writeEntry(os, crossScale_());
371  }
372  if (heightAboveWave_)
373  {
374  writeEntry(os, "heightAboveWave", heightAboveWave_);
375  }
376 }
377 
378 
379 // ************************************************************************* //
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:53
#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:226
const vector origin_
The origin of the wave coordinate system.
void write(Ostream &) const
Write.
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:163
scalarList waveAngles_
The angle relative to the direction at which the waves propagate.
Vector< Cmpt > x() const
Definition: TensorI.H:279
~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:286
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:293
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
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 readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
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:53
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:265
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:240
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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: newWaveModel.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:354
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
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:583