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-2018 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"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::waveSuperposition::transformation
32 (
33  const vectorField& p,
34  tensor& axes,
35  scalar& u,
36  vectorField& xyz
37 ) const
38 {
41  const scalar magG = mag(g.value());
42  const vector gHat = g.value()/magG;
43 
44  const vector dSurf = direction_ - gHat*(gHat & direction_);
45  const scalar magDSurf = mag(dSurf);
46  const vector dSurfHat = direction_/magDSurf;
47 
48  axes = tensor(dSurfHat, - gHat ^ dSurfHat, - gHat);
49 
50  u = speed_*magDSurf;
51 
52  xyz = axes & (p - origin_);
53 }
54 
55 
56 Foam::tmp<Foam::scalarField> Foam::waveSuperposition::elevation
57 (
58  const scalar t,
59  const vector2DField& xy
60 ) const
61 {
62  scalarField result(xy.size(), 0);
63 
64  forAll(waveModels_, wavei)
65  {
66  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
67  result += waveModels_[wavei].elevation(t, d.x()*speed_, d & xy);
68  }
69 
70  return scale(xy)*result;
71 }
72 
73 
74 Foam::tmp<Foam::vectorField> Foam::waveSuperposition::velocity
75 (
76  const scalar t,
77  const vectorField& xyz
78 ) const
79 {
80  vectorField result(xyz.size(), vector::zero);
81 
82  forAll(waveModels_, wavei)
83  {
84  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
85  const vector2DField xz
86  (
87  zip
88  (
89  d & zip(xyz.component(0), xyz.component(1)),
90  tmp<scalarField>(xyz.component(2))
91  )
92  );
93  const vector2DField uw
94  (
95  waveModels_[wavei].velocity(t, d.x()*speed_, xz)
96  );
97  result += zip
98  (
99  d.x()*uw.component(0),
100  d.y()*uw.component(0),
101  uw.component(1)
102  );
103  }
104 
105  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
106 
107  return s*result;
108 }
109 
110 
111 Foam::tmp<Foam::scalarField> Foam::waveSuperposition::pressure
112 (
113  const scalar t,
114  const vectorField& xyz
115 ) const
116 {
117  scalarField result(xyz.size(), 0);
118 
119  forAll(waveModels_, wavei)
120  {
121  const vector2D d(cos(waveAngles_[wavei]), sin(waveAngles_[wavei]));
122  const vector2DField xz
123  (
124  zip
125  (
126  d & zip(xyz.component(0), xyz.component(1)),
127  tmp<scalarField>(xyz.component(2))
128  )
129  );
130  const vector2DField uw
131  (
132  waveModels_[wavei].velocity(t, d.x()*speed_, xz)
133  );
134  result += waveModels_[wavei].pressure(t, d.x()*speed_, xz);
135  }
136 
137  tmp<scalarField> s = scale(zip(xyz.component(0), xyz.component(1)));
138 
139  return s*result;
140 }
141 
142 
143 Foam::tmp<Foam::scalarField> Foam::waveSuperposition::scale
144 (
145  const vector2DField& xy
146 ) const
147 {
148  tmp<scalarField> tResult(new scalarField(xy.size(), 1));
149  scalarField& result = tResult.ref();
150 
151  if (scale_.valid())
152  {
153  const scalarField x(xy.component(0));
154  forAll(result, i)
155  {
156  result[i] *= scale_->value(x[i]);
157  }
158  }
159 
160  if (crossScale_.valid())
161  {
162  const scalarField y(xy.component(1));
163  forAll(result, i)
164  {
165  result[i] *= crossScale_->value(y[i]);
166  }
167  }
168 
169  return tResult;
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
174 
176 :
177  db_(db),
178  origin_(vector::zero),
179  direction_(vector(1, 0, 0)),
180  speed_(0),
181  waveModels_(),
182  waveAngles_(),
183  ramp_(),
184  scale_(),
185  crossScale_(),
186  heightAboveWave_(false)
187 {}
188 
189 
191 :
192  db_(waves.db_),
193  origin_(waves.origin_),
194  direction_(waves.direction_),
195  speed_(waves.speed_),
196  waveModels_(waves.waveModels_),
197  waveAngles_(waves.waveAngles_),
198  ramp_(waves.ramp_, false),
199  scale_(waves.scale_, false),
200  crossScale_(waves.crossScale_, false),
201  heightAboveWave_(waves.heightAboveWave_)
202 {}
203 
204 
206 (
207  const objectRegistry& db,
208  const dictionary& dict
209 )
210 :
211  db_(db),
212  origin_(dict.lookup("origin")),
213  direction_(dict.lookup("direction")),
214  speed_(readScalar(dict.lookup("speed"))),
215  waveModels_(),
216  waveAngles_(),
217  ramp_
218  (
219  dict.found("ramp")
220  ? Function1<scalar>::New("ramp", dict)
221  : autoPtr<Function1<scalar>>()
222  ),
223  scale_
224  (
225  dict.found("scale")
226  ? Function1<scalar>::New("scale", dict)
227  : autoPtr<Function1<scalar>>()
228  ),
229  crossScale_
230  (
231  dict.found("crossScale")
232  ? Function1<scalar>::New("crossScale", dict)
233  : autoPtr<Function1<scalar>>()
234  ),
235  heightAboveWave_(dict.lookupOrDefault<Switch>("heightAboveWave", false))
236 {
237  const PtrList<entry> waveEntries(dict.lookup("waves"));
238 
239  waveModels_.setSize(waveEntries.size());
240  waveAngles_.setSize(waveEntries.size());
241 
242  forAll(waveEntries, wavei)
243  {
244  const dictionary waveDict = waveEntries[wavei].dict();
245  waveModels_.set
246  (
247  wavei,
248  waveModel::New(waveDict.dictName(), db, waveDict)
249  );
250  waveAngles_[wavei] = readScalar(waveDict.lookup("angle"));
251  }
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
256 
258 {}
259 
260 
261 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
262 
264 (
265  const scalar t,
266  const vectorField& p
267 ) const
268 {
269  tensor axes;
270  scalar u;
271  vectorField xyz(p.size());
272  transformation(p, axes, u, xyz);
273 
274  return
275  xyz.component(2)
276  - elevation(t, zip(xyz.component(0), xyz.component(1)));
277 }
278 
279 
281 (
282  const scalar t,
283  const vectorField& p
284 ) const
285 {
286  tensor axes;
287  scalar u;
288  vectorField xyz(p.size());
289  transformation(p, axes, u, xyz);
290 
291  if (heightAboveWave_)
292  {
293  xyz.replace(2, height(t, p));
294  }
295 
296  return UMean(t) + (velocity(t, xyz) & axes);
297 }
298 
299 
301 (
302  const scalar t,
303  const vectorField& p
304 ) const
305 {
306  tensor axes;
307  scalar u;
308  vectorField xyz(p.size());
309  transformation(p, axes, u, xyz);
310 
311  axes = tensor(- axes.x(), - axes.y(), axes.z());
312 
313  if (heightAboveWave_)
314  {
315  xyz.replace(2, height(t, p));
316  }
317 
318  return UMean(t) + (velocity(t, xyz) & axes);
319 }
320 
321 
323 (
324  const scalar t,
325  const vectorField& p
326 ) const
327 {
328  tensor axes;
329  scalar u;
330  vectorField xyz(p.size());
331  transformation(p, axes, u, xyz);
332 
333  if (heightAboveWave_)
334  {
335  xyz.replace(2, height(t, p));
336  }
337 
338  return pressure(t, xyz);
339 }
340 
341 
343 (
344  const scalar t,
345  const vectorField& p
346 ) const
347 {
348  return - pLiquid(t, p);
349 }
350 
351 
353 {
354  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
355  os.writeKeyword("direction") << direction_ << token::END_STATEMENT << nl;
356  os.writeKeyword("speed") << speed_ << token::END_STATEMENT << nl;
357  os.writeKeyword("waves") << nl << token::BEGIN_LIST << nl << incrIndent;
358  forAll(waveModels_, wavei)
359  {
360  os.writeKeyword(waveModels_[wavei].type()) << nl << indent
361  << token::BEGIN_BLOCK << nl << incrIndent;
362  waveModels_[wavei].write(os);
363  os.writeKeyword("angle") << waveAngles_[wavei] << token::END_STATEMENT
364  << nl << decrIndent << indent << token::END_BLOCK << nl;
365  }
367  if (ramp_.valid())
368  {
369  ramp_->writeData(os);
370  }
371  if (scale_.valid())
372  {
373  scale_->writeData(os);
374  }
375  if (crossScale_.valid())
376  {
377  crossScale_->writeData(os);
378  }
379  if (heightAboveWave_)
380  {
381  os.writeKeyword("heightAboveWave") << heightAboveWave_
382  << token::END_STATEMENT << nl;
383  }
384 }
385 
386 
387 // ************************************************************************* //
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...
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
void write(Ostream &) const
Write.
Field< vector2D > vector2DField
Forward declarations of the specialisation of Field<T> for vector2D.
tmp< scalarField > pLiquid(const scalar t, const vectorField &p) const
Get the liquid pressure at a given time and global positions.
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Vector< Cmpt > x() const
Definition: TensorI.H:279
Vector2D< scalar > vector2D
vector2D obtained from generic Vector2D
Definition: vector2D.H:49
~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.
Definition: Switch.H:60
tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
volVectorField vectorField(fieldObject, mesh)
Vector< Cmpt > y() const
Definition: TensorI.H:286
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void replace(const direction, const Cmpt &)
Definition: VectorSpaceI.H:149
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
Vector< Cmpt > z() const
Definition: TensorI.H:293
scalar y
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)
tmp< scalarField > pGas(const scalar t, const vectorField &p) const
Get the gas pressure at a given time and global positions.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar sin(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:265
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
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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
vector UMean(const scalar t) const
Get the mean flow velocity.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:233
bool found
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
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:576