waveModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015 IH-Cantabria
9  Copyright (C) 2016-2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::waveModel
29 
30 Description
31  Base class for waveModels
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef waveModel_H
36 #define waveModel_H
37 
38 #include "autoPtr.H"
39 #include "runTimeSelectionTables.H"
40 
41 #include "IOdictionary.H"
42 #include "dictionary.H"
43 #include "scalarField.H"
44 #include "vectorField.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 class fvMesh;
52 class polyPatch;
53 
54 /*---------------------------------------------------------------------------*\
55  Class waveModel Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class waveModel
59 :
60  public refCount,
61  public IOdictionary
62 {
63 protected:
64 
65  // Protected data
66 
67  //- Reference to the mesh database
68  const fvMesh& mesh_;
69 
70  //- Reference to the polyPatch
71  const polyPatch& patch_;
72 
73  //- Gravity
74  const vector& g_;
75 
76  //- Name of velocity field, default = "U"
77  word UName_;
78 
79  //- Name of phase fraction field, default = "alpha"
81 
82  //- Rotation tensor from global to local system
83  tensor Rgl_;
84 
85  //- Rotation tensor from local to global system
86  tensor Rlg_;
87 
88  //- Number of paddles
89  label nPaddle_;
90 
91  //- Paddle x coordinates / [m]
93 
94  //- Paddle y coordinates / [m]
96 
97  //- Addressing from patch face index to paddle index
99 
100  //- Patch face centre z coordinates / [m]
101  scalarField z_;
102 
103  //- Overall (point) span in z-direction / [m]
104  scalar zSpan_;
105 
106  //- Minimum z (point) height per patch face / [m]
108 
109  //- Maximum z (point) height per patch face / [m]
111 
112  //- Minimum z reference level / [m]
113  scalar zMin0_;
114 
115  //- Reference water depth / [m]
116  scalar waterDepthRef_;
117 
118  //- Initial depth / [m]
119  scalar initialDepth_;
120 
121  //- Time index used for updating
122  label currTimeIndex_;
123 
124  //- Active wave absorption switch
125  bool activeAbsorption_;
126 
127 
128  // Current values
129 
130  //- Velocity field
131  vectorField U_;
132 
133  //- Wave indicator field
135 
136 
137  // Protected Member Functions
138 
139  //- Initialise
140  virtual void initialiseGeometry();
141 
142  //- Water level
143  virtual tmp<scalarField> waterLevel() const;
144 
145  //- Return the time scaling coefficient
146  virtual scalar timeCoeff(const scalar t) const = 0;
147 
148  //- Set the water level
149  virtual void setLevel
150  (
151  const scalar t,
152  const scalar tCoeff,
153  scalarField& level
154  ) const = 0;
155 
156  //- Calculate the wave model velocity
157  virtual void setVelocity
158  (
159  const scalar t,
160  const scalar tCoeff,
161  const scalarField& level
162  ) = 0;
163 
164  //- Set the alpha field based on the water level
165  virtual void setAlpha(const scalarField& level);
166 
167  //- Set the paddle coverage fraction and reference height
168  virtual void setPaddlePropeties
169  (
170  const scalarField& level,
171  const label facei,
172  scalar& fraction,
173  scalar& z
174  ) const;
175 
176 
177 public:
178 
179  //- Runtime type information
180  TypeName("waveModel");
181 
182  // Declare runtime constructor selection table
183 
185  (
186  autoPtr,
187  waveModel,
188  patch,
189  (
190  const dictionary& dict,
191  const fvMesh& mesh,
192  const polyPatch& patch
193  ),
194  (dict, mesh, patch)
195  );
196 
197 
198  //- Constructor
199  waveModel
200  (
201  const dictionary& dict,
202  const fvMesh& mesh,
203  const polyPatch& patch,
204  const bool readFields = true
205  );
206 
207 
208  // Selectors
209 
210  //- Return a reference to the selected wave model
211  static autoPtr<waveModel> New
212  (
213  const word& dictName,
214  const fvMesh& mesh,
215  const polyPatch& patch
216  );
217 
218  //- Lookup waveModel from database, or create new
220  (
221  const polyPatch& patch,
222  const fvMesh& mesh,
223  const word& waveDictName
224  );
225 
226 
227  //- Destructor
228  virtual ~waveModel() = default;
229 
230  //- Dictionary name
231  static const word dictName;
232 
233 
234  // Public Member Functions
235 
236  //- Utility function to construct the model name
237  static word modelName(const word& patchName);
238 
239  //- Read from dictionary
240  virtual bool readDict(const dictionary& overrideDict);
241 
242  //- Return the latest wave velocity prediction
243  virtual const vectorField& U() const;
244 
245  //- Return the latest wave indicator field prediction
246  virtual const scalarField& alpha() const;
247 
248  //- Correct the model for time, t[s]
249  virtual void correct(const scalar t);
250 
251  //- Info
252  virtual void info(Ostream& os) const;
253 };
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 } // End namespace Foam
259 
260 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
261 
262 #endif
263 
264 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::Tensor< scalar >
Foam::waveModel::U_
vectorField U_
Velocity field.
Definition: waveModel.H:130
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::waveModel::initialiseGeometry
virtual void initialiseGeometry()
Initialise.
Definition: waveModel.C:57
Foam::waveModel
Base class for waveModels.
Definition: waveModel.H:57
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
scalarField.H
Foam::refCount
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
Foam::waveModel::zMax_
scalarField zMax_
Maximum z (point) height per patch face / [m].
Definition: waveModel.H:109
Foam::waveModel::z_
scalarField z_
Patch face centre z coordinates / [m].
Definition: waveModel.H:100
Foam::waveModel::currTimeIndex_
label currTimeIndex_
Time index used for updating.
Definition: waveModel.H:121
Foam::waveModel::activeAbsorption_
bool activeAbsorption_
Active wave absorption switch.
Definition: waveModel.H:124
Foam::waveModel::alpha
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
Definition: waveModel.C:416
Foam::waveModel::setVelocity
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)=0
Calculate the wave model velocity.
Foam::waveModel::Rgl_
tensor Rgl_
Rotation tensor from global to local system.
Definition: waveModel.H:82
Foam::waveModel::~waveModel
virtual ~waveModel()=default
Destructor.
Foam::waveModel::zMin0_
scalar zMin0_
Minimum z reference level / [m].
Definition: waveModel.H:112
Foam::waveModel::waveModel
waveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: waveModel.C:246
Foam::waveModel::alphaName_
word alphaName_
Name of phase fraction field, default = "alpha".
Definition: waveModel.H:79
Foam::waveModel::New
static autoPtr< waveModel > New(const word &dictName, const fvMesh &mesh, const polyPatch &patch)
Return a reference to the selected wave model.
Definition: waveModelNew.C:33
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
Foam::waveModel::Rlg_
tensor Rlg_
Rotation tensor from local to global system.
Definition: waveModel.H:85
Foam::waveModel::waterLevel
virtual tmp< scalarField > waterLevel() const
Water level.
Definition: waveModel.C:128
Foam::waveModel::U
virtual const vectorField & U() const
Return the latest wave velocity prediction.
Definition: waveModel.C:410
Foam::waveModel::xPaddle_
scalarField xPaddle_
Paddle x coordinates / [m].
Definition: waveModel.H:91
Foam::waveModel::yPaddle_
scalarField yPaddle_
Paddle y coordinates / [m].
Definition: waveModel.H:94
Foam::waveModel::setAlpha
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
Definition: waveModel.C:163
Foam::Field< scalar >
Foam::waveModel::nPaddle_
label nPaddle_
Number of paddles.
Definition: waveModel.H:88
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::waveModel::TypeName
TypeName("waveModel")
Runtime type information.
Foam::waveModel::modelName
static word modelName(const word &patchName)
Utility function to construct the model name.
Definition: waveModel.C:49
Foam::waveModel::setPaddlePropeties
virtual void setPaddlePropeties(const scalarField &level, const label facei, scalar &fraction, scalar &z) const
Set the paddle coverage fraction and reference height.
Definition: waveModel.C:191
Foam::waveModel::UName_
word UName_
Name of velocity field, default = "U".
Definition: waveModel.H:76
Foam::waveModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, waveModel, patch,(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch),(dict, mesh, patch))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::waveModel::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: waveModel.H:67
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::waveModel::patch_
const polyPatch & patch_
Reference to the polyPatch.
Definition: waveModel.H:70
Foam::waveModel::zMin_
scalarField zMin_
Minimum z (point) height per patch face / [m].
Definition: waveModel.H:106
Foam::waveModel::faceToPaddle_
labelList faceToPaddle_
Addressing from patch face index to paddle index.
Definition: waveModel.H:97
Foam::waveModel::initialDepth_
scalar initialDepth_
Initial depth / [m].
Definition: waveModel.H:118
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::waveModel::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: waveModel.C:295
IOdictionary.H
Foam::waveModel::setLevel
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const =0
Set the water level.
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::waveModel::alpha_
scalarField alpha_
Wave indicator field.
Definition: waveModel.H:133
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::Vector< scalar >
Foam::List< label >
vectorField.H
dictionary.H
Foam::waveModel::dictName
static const word dictName
Dictionary name.
Definition: waveModel.H:230
Foam::waveModel::timeCoeff
virtual scalar timeCoeff(const scalar t) const =0
Return the time scaling coefficient.
Foam::waveModel::waterDepthRef_
scalar waterDepthRef_
Reference water depth / [m].
Definition: waveModel.H:115
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::waveModel::correct
virtual void correct(const scalar t)
Correct the model for time, t[s].
Definition: waveModel.C:350
Foam::waveModel::zSpan_
scalar zSpan_
Overall (point) span in z-direction / [m].
Definition: waveModel.H:103
Foam::waveModel::g_
const vector & g_
Gravity.
Definition: waveModel.H:73
autoPtr.H
Foam::waveModel::lookupOrCreate
static tmp< waveModel > lookupOrCreate(const polyPatch &patch, const fvMesh &mesh, const word &waveDictName)
Lookup waveModel from database, or create new.
Definition: waveModelNew.C:86
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60