StokesIWaveModel.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::StokesI
29
30Description
31 Stokes I wave model
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef waveModels_StokesI_H
36#define waveModels_StokesI_H
37
38#include "regularWaveModel.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace waveModels
45{
46
47/*---------------------------------------------------------------------------*\
48 Class StokesI Declaration
49\*---------------------------------------------------------------------------*/
51class StokesI
52:
53 public regularWaveModel
54{
55private:
56
57 // Private Member Functions
58
59 //- Wave height
60 scalar eta
61 (
62 const scalar H,
63 const scalar Kx,
64 const scalar x,
65 const scalar Ky,
66 const scalar y,
67 const scalar omega,
68 const scalar t,
69 const scalar phase
70 ) const;
71
72
73protected:
74
75 // Protected Member Functions
76
77 //- Return the wavelength
78 virtual scalar waveLength(const scalar h, const scalar T) const;
79
80 //- Wave velocity
81 virtual vector UfBase
82 (
83 const scalar H,
84 const scalar h,
85 const scalar Kx,
86 const scalar x,
87 const scalar Ky,
88 const scalar y,
89 const scalar omega,
90 const scalar t,
91 const scalar phase,
92 const scalar z
93 ) const;
94
95 //- Set the water level
96 virtual void setLevel
97 (
98 const scalar t,
99 const scalar tCoeff,
100 scalarField& level
101 ) const;
102
103 //- Calculate the wave model velocity
104 virtual void setVelocity
105 (
106 const scalar t,
107 const scalar tCoeff,
108 const scalarField& level
109 );
110
111
112public:
113
114 //- Runtime type information
115 TypeName("StokesI");
116
117 //- Constructor
118 StokesI
119 (
120 const dictionary& dict,
121 const fvMesh& mesh,
122 const polyPatch& patch,
123 const bool readFields = true
124 );
125
126 //- Destructor
127 virtual ~StokesI() = default;
128
129
130 // Public Member Functions
131
132 //- Read from dictionary
133 virtual bool readDict(const dictionary& overrideDict);
134
135 //- Info
136 virtual void info(Ostream& os) const;
137};
138
139
140// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141
142} // End namespace waveModels
143} // End namespace Foam
144
145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146
147#endif
148
149// ************************************************************************* //
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
Stokes I wave model.
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
TypeName("StokesI")
Runtime type information.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
virtual ~StokesI()=default
Destructor.
virtual scalar waveLength(const scalar h, const scalar T) const
Return the wavelength.
virtual vector UfBase(const scalar H, const scalar h, const scalar Kx, const scalar x, const scalar Ky, const scalar y, const scalar omega, const scalar t, const scalar phase, const scalar z) const
Wave velocity.
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))
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