streamFunctionWaveModel.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::streamFunction
29
30Description
31 streamFunction wave model
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef waveModels_streamFunction_H
36#define waveModels_streamFunction_H
37
38#include "regularWaveModel.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace waveModels
45{
46
47/*---------------------------------------------------------------------------*\
48 Class streamFunction Declaration
49\*---------------------------------------------------------------------------*/
52:
53 public regularWaveModel
54{
55private:
56
57 // Private Member Functions
58
59 //- Wave height
60 virtual scalar eta
61 (
62 const scalar h,
63 const scalar kx,
64 const scalar ky,
65 const scalar T,
66 const scalar x,
67 const scalar y,
68 const scalar omega,
69 const scalar t,
70 const scalar phase
71 ) const;
72
73 //- Wave velocity
74 virtual vector Uf
75 (
76 const scalar d,
77 const scalar kx,
78 const scalar ky,
79 const scalar T,
80 const scalar x,
81 const scalar y,
82 const scalar omega,
83 const scalar t,
84 const scalar phase,
85 const scalar z
86 ) const;
87
88
89protected:
90
91 // Protected data
92
93 //- Mean fluid speed in frame of reference (stream function)
94 scalar uMean_;
95
96 //- Stream Function Bj coefficients
98
99 //- Stream Function Ej coefficients
101
102
103 // Protected Member Functions
104
105
106 //- Set the water level
107 virtual void setLevel
108 (
109 const scalar t,
110 const scalar tCoeff,
111 scalarField& level
112 ) const;
113
114 //- Calculate the wave model velocity
115 virtual void setVelocity
116 (
117 const scalar t,
118 const scalar tCoeff,
119 const scalarField& level
120 );
121
122
123public:
124
125 //- Runtime type information
126 TypeName("streamFunction");
127
128 //- Constructor
130 (
131 const dictionary& dict,
132 const fvMesh& mesh,
133 const polyPatch& patch,
134 const bool readFields = true
135 );
136
137 //- Destructor
138 virtual ~streamFunction() = default;
139
140
141 // Public Member Functions
142
143 //- Read from dictionary
144 virtual bool readDict(const dictionary& overrideDict);
145
146 //- Info
147 virtual void info(Ostream& os) const;
148};
149
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153} // End namespace waveModels
154} // End namespace Foam
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158#endif
159
160// ************************************************************************* //
scalar y
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
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
virtual ~streamFunction()=default
Destructor.
scalarList Ejs_
Stream Function Ej coefficients.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
scalar uMean_
Mean fluid speed in frame of reference (stream function)
scalarList Bjs_
Stream Function Bj coefficients.
TypeName("streamFunction")
Runtime type information.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
const volScalarField & T
dynamicFvMesh & mesh
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
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73