waveSuperposition.H
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-2020 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 Class
25  Foam::waveSuperposition
26 
27 Description
28  A wrapper around a list of wave models. Superimposes the modelled values of
29  elevation and velocity. The New method looks up or or constructs an
30  instance of this class on demand and returns a reference. Properties are
31  read from a dictionary in constant.
32 
33 Usage
34  \table
35  Property | Description | Req'd? | Default
36  origin | origin of the wave coordinate system | yes |
37  direction | direction of the wave coordinate system | yes |
38  waves | list of wave models to superimpose | yes |
39  UMean | velocity of the mean flow | yes |
40  scale | scale factor in the direction | no | None
41  crossScale | scale factor perpendicular to the direction | no | None
42  heightAboveWave | use the height above the wave as the vertical \\
43  coordinate | no | false
44  \endtable
45 
46  Example specification:
47  \verbatim
48  origin (0 25 0);
49  direction (1 0 0);
50  waves
51  (
52  Airy
53  {
54  length 40;
55  amplitude 0.5;
56  phase 0;
57  angle 0;
58  }
59  Airy
60  {
61  length 20;
62  amplitude 0.25;
63  phase 1.5708;
64  angle 0;
65  }
66  );
67  UMean (2 0 0);
68  scale table ((100 1) (200 0));
69  crossScale constant 1;
70  heightAboveWave no;
71  \endverbatim
72 
73 SourceFiles
74  waveSuperposition.C
75  waveSuperpositionNew.C
76 
77 \*---------------------------------------------------------------------------*/
78 
79 #ifndef waveSuperposition_H
80 #define waveSuperposition_H
81 
82 #include "waveModel.H"
83 #include "IOdictionary.H"
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 namespace Foam
88 {
89 
90 /*---------------------------------------------------------------------------*\
91  Class waveSuperposition Declaration
92 \*---------------------------------------------------------------------------*/
93 
94 class waveSuperposition
95 :
96  public IOdictionary
97 {
98 protected:
99 
100  // Protected Data
101 
102  //- The origin of the wave coordinate system
103  const vector origin_;
104 
105  //- The direction of the wave coordinate system
106  const vector direction_;
107 
108  //- Wave models to superimpose
109  PtrList<waveModel> waveModels_;
110 
111  //- The angle relative to the direction at which the waves propagate
113 
114  //- Mean velocity
115  const autoPtr<Function1<vector>> UMean_;
116 
117  //- Scaling in the local x-direction
118  const autoPtr<Function1<scalar>> scale_;
119 
120  //- Scaling perpendicular to the local x-direction
121  const autoPtr<Function1<scalar>> crossScale_;
122 
123  //- Calculate wave properties using the height above the wave (true) or
124  // the height above the origin (false)?
125  const Switch heightAboveWave_;
126 
127 
128  // Protected Member Functions
129 
130  //- Get the transformation to actual coordinates
131  void transformation
132  (
133  const scalar t,
134  const vectorField& p,
135  tensor& axes,
136  vectorField& xyz
137  ) const;
138 
139  //- Get the wave elevation relative to the mean at a given time and
140  // local coordinates. Local x is aligned with the direction, and y is
141  // perpendicular to both x and gravity.
143  (
144  const scalar t,
145  const vector2DField& xy
146  ) const;
148  //- Get the wave velocity at a given time and local coordinates. Local
149  // x is aligned with the direction, z with negative gravity, and y is
150  // perpendicular to both.
151  tmp<vectorField> velocity(const scalar t, const vectorField& xyz) const;
152 
153  //- Get the wave pressure at a given time and local coordinates. Local
154  // x is aligned with the direction, z with negative gravity, and y is
155  // perpendicular to both.
156  tmp<scalarField> pressure(const scalar t, const vectorField& xyz) const;
157 
158  //- Get the scaling factor, calculated from the optional scaling
159  // functions. X and y are the same as for the elevation method.
160  tmp<scalarField> scale(const vector2DField& xy) const;
161 
162 
163 public:
164 
165  //- Runtime type information
166  TypeName("wave");
167 
168 
169  // Declare runtime construction
171  (
172  autoPtr,
175  (const objectRegistry& db),
176  (db)
177  );
178 
179 
180  // Static Data
181 
182  //- The name of the dictionary
183  static const word dictName;
184 
185 
186  // Static Member Functions
187 
188  //- Return a reference to the wave model on the given database,
189  // constructing if it doesn't exist
190  static const waveSuperposition& New(const objectRegistry& db);
191 
192 
193  // Constructors
194 
195  //- Construct from a database
197 
198 
199  //- Destructor
201 
202 
203  // Member Functions
204 
205  //- Get the height above the waves at a given time and global positions
206  virtual tmp<scalarField> height
207  (
208  const scalar t,
209  const vectorField& p
210  ) const;
211 
212  //- Get the liquid velocity at a given time and global positions
213  virtual tmp<vectorField> ULiquid
214  (
215  const scalar t,
216  const vectorField& p
217  ) const;
218 
219  //- Get the gas velocity at a given time and global positions
220  virtual tmp<vectorField> UGas
221  (
222  const scalar t,
223  const vectorField& p
224  ) const;
225 
226  //- Get the liquid pressure at a given time and global positions
227  virtual tmp<scalarField> pLiquid
228  (
229  const scalar t,
230  const vectorField& p
231  ) const;
232 
233  //- Get the gas pressure at a given time and global positions
234  virtual tmp<scalarField> pGas
235  (
236  const scalar t,
237  const vectorField& p
238  ) const;
239 
240  //- Inherit write from regIOobject
241  using regIOobject::write;
242 
243  //- Write
244  void write(Ostream&) const;
245 };
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace Foam
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 #endif
255 
256 // ************************************************************************* //
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
virtual tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
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.
const vector direction_
The direction of the wave coordinate system.
const autoPtr< Function1< vector > > UMean_
Mean velocity.
scalarList waveAngles_
The angle relative to the direction at which the waves propagate.
~waveSuperposition()
Destructor.
virtual tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
const Switch heightAboveWave_
Calculate wave properties using the height above the wave (true) or.
tmp< scalarField > elevation(const scalar t, const vector2DField &xy) const
Get the wave elevation relative to the mean at a given time and.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
const autoPtr< Function1< scalar > > crossScale_
Scaling perpendicular to the local x-direction.
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
virtual tmp< scalarField > pGas(const scalar t, const vectorField &p) const
Get the gas pressure at a given time and global positions.
const autoPtr< Function1< scalar > > scale_
Scaling in the local x-direction.
A class for handling words, derived from string.
Definition: word.H:59
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
declareRunTimeSelectionTable(autoPtr, waveSuperposition, objectRegistry,(const objectRegistry &db),(db))
PtrList< waveModel > waveModels_
Wave models to superimpose.
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
TypeName("wave")
Runtime type information.
tmp< vectorField > velocity(const scalar t, const vectorField &xyz) const
Get the wave velocity at a given time and local coordinates. Local.
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.
virtual bool write(const bool write=true) const
Write using setting from DB.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
Namespace for OpenFOAM.
waveSuperposition(const objectRegistry &db)
Construct from a database.