StokesVWaveModel.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::waveModels::StokesV
29
30Description
31 Stokes V wave model
32
33 Reference
34 \verbatim
35 Skjelbreia, L., Hendrickson, J. (1960).
36 Fifth order gravity wave theory.
37 Proc. Coastal Engineering Conf., pp. 184-196.
38 \endverbatim
39
40\*---------------------------------------------------------------------------*/
41
42#ifndef waveModels_StokesV_H
43#define waveModels_StokesV_H
44
45#include "StokesIWaveModel.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51namespace waveModels
52{
53
54/*---------------------------------------------------------------------------*\
55 Class StokesV Declaration
56\*---------------------------------------------------------------------------*/
58class StokesV
59:
60 public StokesI
61{
62private:
63
64 // Private Member Functions
65
66 // Model coefficients
67
68 scalar A11(const scalar h, const scalar k) const;
69
70 scalar A13(const scalar h, const scalar k) const;
71
72 scalar A15(const scalar h, const scalar k) const;
73
74 scalar A22(const scalar h, const scalar k) const;
75
76 scalar A24(const scalar h, const scalar k) const;
77
78 scalar A33(const scalar h, const scalar k) const;
79
80 scalar A35(const scalar h, const scalar k) const;
81
82 scalar A44(const scalar h, const scalar k) const;
83
84 scalar A55(const scalar h, const scalar k) const;
85
86 scalar B22(const scalar h, const scalar k) const;
87
88 scalar B24(const scalar h, const scalar k) const;
89
90 scalar B33(const scalar h, const scalar k) const;
91
92 scalar B33k(const scalar h, const scalar k) const;
93
94 scalar B35(const scalar h, const scalar k) const;
95
96 scalar B35k(const scalar h, const scalar k) const;
97
98 scalar B44(const scalar h, const scalar k) const;
99
100 scalar B55(const scalar h, const scalar k) const;
101
102 scalar B55k(const scalar h, const scalar k) const;
103
104 scalar C1(const scalar h, const scalar k) const;
105
106 scalar C1k(const scalar h, const scalar k) const;
107
108 scalar C2(const scalar h, const scalar k) const;
109
110 scalar C2k(const scalar h, const scalar k) const;
111
112 scalar C3(const scalar h, const scalar k) const;
113
114 scalar C4(const scalar h, const scalar k) const;
115
116
117 //- Model initialisation
118 void initialise
119 (
120 const scalar H,
121 const scalar d,
122 const scalar T,
123 scalar& kOut,
124 scalar& LambdaOut,
125 scalar& f1Out,
126 scalar& f2Out
127 ) const;
128
129 //- Wave height
130 scalar eta
131 (
132 const scalar h,
133 const scalar kx,
134 const scalar ky,
135 const scalar lambda,
136 const scalar T,
137 const scalar x,
138 const scalar y,
139 const scalar t,
140 const scalar phase
141 ) const;
142
143 //- Wave velocity
144 vector Uf
145 (
146 const scalar d,
147 const scalar kx,
148 const scalar ky,
149 const scalar lambda,
150 const scalar T,
151 const scalar x,
152 const scalar y,
153 const scalar t,
154 const scalar phase,
155 const scalar z
156 ) const;
157
158
159protected:
160
161 // Protected Data
162
163 //-
164 scalar lambda_;
165
166
167 // Protected Member Functions
168
169 //- Set the water level
170 virtual void setLevel
171 (
172 const scalar t,
173 const scalar tCoeff,
174 scalarField& level
175 ) const;
176
177 //- Calculate the wave model velocity
178 virtual void setVelocity
179 (
180 const scalar t,
181 const scalar tCoeff,
182 const scalarField& level
183 );
184
185
186public:
187
188 //- Runtime type information
189 TypeName("StokesV");
190
191 //- Constructor
192 StokesV
193 (
194 const dictionary& dict,
195 const fvMesh& mesh,
196 const polyPatch& patch,
197 const bool readFields = true
198 );
199
200 //- Destructor
201 virtual ~StokesV() = default;
202
203
204 // Public Member Functions
205
206 //- Read from dictionary
207 virtual bool readDict(const dictionary& overrideDict);
208
209 //- Info
210 virtual void info(Ostream& os) const;
211};
212
213
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215
216} // End namespace waveModels
217} // End namespace Foam
218
219// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220
221#endif
222
223// ************************************************************************* //
scalar y
label k
InfoProxy< IOobject > info() const
Return info proxy, for printing information to a stream.
Definition: IOobject.H:690
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Stokes I wave model.
Stokes V wave model.
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
virtual ~StokesV()=default
Destructor.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
TypeName("StokesV")
Runtime type information.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
const volScalarField & T
dynamicFvMesh & mesh
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
autoPtr< surfaceVectorField > Uf
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
dictionary dict
volScalarField & h
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73