waveSuperposition.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2017 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::scale
112 (
113  const vector2DField& xy
114 ) const
115 {
116  tmp<scalarField> tResult(new scalarField(xy.size(), 1));
117  scalarField& result = tResult.ref();
118 
119  if (scale_.valid())
120  {
121  const scalarField x(xy.component(0));
122  forAll(result, i)
123  {
124  result[i] *= scale_->value(x[i]);
125  }
126  }
127 
128  if (crossScale_.valid())
129  {
130  const scalarField y(xy.component(1));
131  forAll(result, i)
132  {
133  result[i] *= crossScale_->value(y[i]);
134  }
135  }
136 
137  return tResult;
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
144 :
145  db_(db),
146  origin_(vector::zero),
147  direction_(vector(1, 0, 0)),
148  speed_(0),
149  waveModels_(),
150  waveAngles_(),
151  scale_(),
152  crossScale_()
153 {}
154 
155 
157 :
158  db_(waves.db_),
159  origin_(waves.origin_),
160  direction_(waves.direction_),
161  speed_(waves.speed_),
162  waveModels_(waves.waveModels_),
163  waveAngles_(waves.waveAngles_),
164  scale_(waves.scale_, false),
165  crossScale_(waves.crossScale_, false)
166 {}
167 
168 
170 (
171  const objectRegistry& db,
172  const dictionary& dict
173 )
174 :
175  db_(db),
176  origin_(dict.lookup("origin")),
177  direction_(dict.lookup("direction")),
178  speed_(readScalar(dict.lookup("speed"))),
179  waveModels_(),
180  waveAngles_(),
181  scale_
182  (
183  dict.found("scale")
184  ? Function1<scalar>::New("scale", dict)
185  : autoPtr<Function1<scalar>>()
186  ),
187  crossScale_
188  (
189  dict.found("crossScale")
190  ? Function1<scalar>::New("crossScale", dict)
191  : autoPtr<Function1<scalar>>()
192  )
193 {
194  const PtrList<entry> waveEntries(dict.lookup("waves"));
195 
196  waveModels_.setSize(waveEntries.size());
197  waveAngles_.setSize(waveEntries.size());
198 
199  forAll(waveEntries, wavei)
200  {
201  const dictionary waveDict = waveEntries[wavei].dict();
202  waveModels_.set
203  (
204  wavei,
205  waveModel::New(waveDict.dictName(), db, waveDict)
206  );
207  waveAngles_[wavei] = readScalar(waveDict.lookup("angle"));
208  }
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213 
215 {}
216 
217 
218 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
219 
221 (
222  const scalar t,
223  const vectorField& p
224 ) const
225 {
226  tensor axes;
227  scalar u;
228  vectorField xyz(p.size());
229  transformation(p, axes, u, xyz);
230 
231  return
232  xyz.component(2)
233  - elevation(t, zip(xyz.component(0), xyz.component(1)));
234 }
235 
236 
238 (
239  const scalar t,
240  const vectorField& p
241 ) const
242 {
243  tensor axes;
244  scalar u;
245  vectorField xyz(p.size());
246  transformation(p, axes, u, xyz);
247 
248  return UMean() + (velocity(t, xyz) & axes);
249 }
250 
251 
253 (
254  const scalar t,
255  const vectorField& p
256 ) const
257 {
258  tensor axes;
259  scalar u;
260  vectorField xyz(p.size());
261  transformation(p, axes, u, xyz);
262 
263  axes = tensor(- axes.x(), - axes.y(), axes.z());
264 
265  return UMean() + (velocity(t, xyz) & axes);
266 }
267 
268 
270 {
271  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
272  os.writeKeyword("direction") << direction_ << token::END_STATEMENT << nl;
273  os.writeKeyword("speed") << speed_ << token::END_STATEMENT << nl;
274  os.writeKeyword("waves") << nl << token::BEGIN_LIST << nl << incrIndent;
275  forAll(waveModels_, wavei)
276  {
277  os.writeKeyword(waveModels_[wavei].type()) << nl << indent
278  << token::BEGIN_BLOCK << nl << incrIndent;
279  waveModels_[wavei].write(os);
280  os.writeKeyword("angle") << waveAngles_[wavei] << token::END_STATEMENT
281  << nl << decrIndent << indent << token::END_BLOCK << nl;
282  }
284  if (scale_.valid())
285  {
286  scale_->writeData(os);
287  }
288  if (crossScale_.valid())
289  {
290  crossScale_->writeData(os);
291  }
292 }
293 
294 
295 // ************************************************************************* //
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
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:223
void write(Ostream &) const
Write.
Field< vector2D > vector2DField
Forward declarations of the specialisation of Field<T> for vector2D.
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.
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
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)
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 succesful.
Definition: doubleScalar.H:63
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:262
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
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:485
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
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.
vector UMean() const
Get the mean flow velocity.
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:230
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