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