Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::waveModel
29
30Description
31 Base class for waveModels
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef waveModel_H
36#define waveModel_H
37
38#include "autoPtr.H"
40
41#include "IOdictionary.H"
42#include "dictionary.H"
43#include "scalarField.H"
44#include "vectorField.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51class fvMesh;
52class polyPatch;
53
54/*---------------------------------------------------------------------------*\
55 Class waveModel Declaration
56\*---------------------------------------------------------------------------*/
57
58class waveModel
59:
60 public refCount,
61 public IOdictionary
62{
63protected:
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"
79 //- Name of phase fraction field, default = "alpha"
81
82 //- Rotation tensor from global to local system
84
85 //- Rotation tensor from local to global system
87
88 //- Number of paddles
89 label nPaddle_;
90
91 //- Paddle x coordinates / [m]
94 //- Paddle y coordinates / [m]
96
97 //- Addressing from patch face index to paddle index
99
100 //- Patch face centre z coordinates / [m]
102
103 //- Overall (point) span in z-direction / [m]
104 scalar zSpan_;
105
106 //- Minimum z (point) height per patch face / [m]
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_;
124 //- Active wave absorption switch
126
127
128 // Current values
129
130 //- Velocity field
132
133 //- Wave indicator field
135
136
137 // Protected Member Functions
139 //- Initialise
140 virtual void initialiseGeometry();
141
142 //- Water level
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,
154 ) const = 0;
155
156 //- Calculate the wave model velocity
157 virtual void setVelocity
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);
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
177public:
178
179 //- Runtime type information
180 TypeName("waveModel");
181
182 // Declare runtime constructor selection table
183
185 (
186 autoPtr,
187 waveModel,
188 patch,
190 const dictionary& dict,
191 const fvMesh& mesh,
192 const polyPatch& patch
193 ),
194 (dict, mesh, patch)
195 );
196
197
198 //- Constructor
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
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
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// ************************************************************************* //
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
constexpr refCount() noexcept
Default construct, initializing count to 0.
Definition refCount.H:63
A class for managing temporary objects.
Definition tmp.H:75
virtual void initialiseGeometry()
Initialise.
Definition waveModel.C:50
scalarField zMin_
Minimum z (point) height per patch face / [m].
Definition waveModel.H:128
scalarField z_
Patch face centre z coordinates / [m].
Definition waveModel.H:118
virtual scalar timeCoeff(const scalar t) const =0
Return the time scaling coefficient.
const polyPatch & patch_
Reference to the polyPatch.
Definition waveModel.H:68
tensor Rlg_
Rotation tensor from local to global system.
Definition waveModel.H:93
const fvMesh & mesh_
Reference to the mesh database.
Definition waveModel.H:63
scalar zSpan_
Overall (point) span in z-direction / [m].
Definition waveModel.H:123
word UName_
Name of velocity field, default = "U".
Definition waveModel.H:78
const vector & g_
Gravity.
Definition waveModel.H:73
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
Definition waveModel.C:408
scalar zMin0_
Minimum z reference level / [m].
Definition waveModel.H:138
scalarField alpha_
Wave indicator field.
Definition waveModel.H:171
bool activeAbsorption_
Active wave absorption switch.
Definition waveModel.H:158
virtual const vectorField & U() const
Return the latest wave velocity prediction.
Definition waveModel.C:402
tensor Rgl_
Rotation tensor from global to local system.
Definition waveModel.H:88
vectorField U_
Velocity field.
Definition waveModel.H:166
scalar waterDepthRef_
Reference water depth / [m].
Definition waveModel.H:143
declareRunTimeSelectionTable(autoPtr, waveModel, patch,(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch),(dict, mesh, patch))
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:183
TypeName("waveModel")
Runtime type information.
scalarField yPaddle_
Paddle y coordinates / [m].
Definition waveModel.H:108
static const word dictName
Dictionary name.
Definition waveModel.H:294
scalarField xPaddle_
Paddle x coordinates / [m].
Definition waveModel.H:103
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)=0
Calculate the wave model velocity.
waveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition waveModel.C:238
static autoPtr< waveModel > New(const word &dictName, const fvMesh &mesh, const polyPatch &patch)
Return a reference to the selected wave model.
virtual ~waveModel()=default
Destructor.
scalar initialDepth_
Initial depth / [m].
Definition waveModel.H:148
scalarField zMax_
Maximum z (point) height per patch face / [m].
Definition waveModel.H:133
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const =0
Set the water level.
static word modelName(const word &patchName)
Utility function to construct the model name.
Definition waveModel.C:42
label currTimeIndex_
Time index used for updating.
Definition waveModel.H:153
labelList faceToPaddle_
Addressing from patch face index to paddle index.
Definition waveModel.H:113
virtual tmp< scalarField > waterLevel() const
Water level.
Definition waveModel.C:120
word alphaName_
Name of phase fraction field, default = "alpha".
Definition waveModel.H:83
static tmp< waveModel > lookupOrCreate(const polyPatch &patch, const fvMesh &mesh, const word &waveDictName)
Lookup waveModel from database, or create new.
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
Definition waveModel.C:155
label nPaddle_
Number of paddles.
Definition waveModel.H:98
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition waveModel.C:287
A class for handling words, derived from Foam::string.
Definition word.H:66
thermo correct()
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
List< label > labelList
A List of labels.
Definition List.H:62
Tensor< scalar > tensor
Definition symmTensor.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Vector< scalar > vector
Definition vector.H:57
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68