McCowanWaveModel.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) 2017 IH-Cantabria
9 Copyright (C) 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::McCowan
29
30Description
31 McCowan wave model
32
33\*---------------------------------------------------------------------------*/
34
35#ifndef waveModels_McCowan_H
36#define waveModels_McCowan_H
37
38#include "solitaryWaveModel.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42namespace Foam
43{
44namespace waveModels
45{
46
47/*---------------------------------------------------------------------------*\
48 Class McCowan Declaration
49\*---------------------------------------------------------------------------*/
51class McCowan
52:
54{
55protected:
56
57 // Protected Member Functions
58
59 //- Wave height
60 virtual scalar eta
61 (
62 const scalar H,
63 const scalar h,
64 const scalar x,
65 const scalar y,
66 const scalar theta,
67 const scalar t,
68 const scalar X0
69 ) const;
70
71 virtual vector mn
72 (
73 const scalar H,
74 const scalar h
75 ) const;
76
77 virtual scalar newtonRapsonF1
78 (
79 const scalar x0,
80 const scalar H,
81 const scalar h
82 ) const;
83
84 virtual scalar newtonRapsonF2
85 (
86 const scalar x0,
87 const scalar H,
88 const scalar h,
89 const scalar xa,
90 const scalar m,
91 const scalar n
92 ) const;
93
94 //- Wave velocity
95 virtual vector Uf
96 (
97 const scalar H,
98 const scalar h,
99 const scalar x,
100 const scalar y,
101 const scalar theta,
102 const scalar t,
103 const scalar X0,
104 const scalar z
105 ) const;
106
107 //- Set the water level
108 virtual void setLevel
109 (
110 const scalar t,
111 const scalar tCoeff,
112 scalarField& level
113 ) const;
114
115 //- Calculate the wave model velocity
116 virtual void setVelocity
117 (
118 const scalar t,
119 const scalar tCoeff,
120 const scalarField& level
121 );
122
123
124public:
125
126 //- Runtime type information
127 TypeName("McCowan");
128
129 //- Constructor
130 McCowan
131 (
132 const dictionary& dict,
133 const fvMesh& mesh,
134 const polyPatch& patch,
135 const bool readFields = true
136 );
137
138 //- Destructor
139 virtual ~McCowan() = default;
140
141
142 // Public Member Functions
143
144 //- Read from dictionary
145 virtual bool readDict(const dictionary& overrideDict);
146
147 //- Info
148 virtual void info(Ostream& os) const;
149};
150
151// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152
153} // End namespace waveModels
154} // End namespace Foam
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158#endif
159
160// ************************************************************************* //
scalar y
label n
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
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
McCowan wave model.
virtual scalar newtonRapsonF2(const scalar x0, const scalar H, const scalar h, const scalar xa, const scalar m, const scalar n) const
virtual scalar newtonRapsonF1(const scalar x0, const scalar H, const scalar h) const
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
virtual ~McCowan()=default
Destructor.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
virtual vector mn(const scalar H, const scalar h) const
virtual scalar eta(const scalar H, const scalar h, const scalar x, const scalar y, const scalar theta, const scalar t, const scalar X0) const
Wave height.
TypeName("McCowan")
Runtime type information.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
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.
scalarList X0(nSpecie, Zero)
dictionary dict
volScalarField & h
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73