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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::waveModels::McCowan
29 
30 Description
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 
42 namespace Foam
43 {
44 namespace waveModels
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class McCowan Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class McCowan
52 :
53  public solitaryWaveModel
54 {
55 protected:
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 
124 public:
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 // ************************************************************************* //
Foam::waveModels::McCowan::Uf
virtual vector Uf(const scalar H, const scalar h, const scalar x, const scalar y, const scalar theta, const scalar t, const scalar X0, const scalar z) const
Wave velocity.
Definition: McCowanWaveModel.C:208
Foam::waveModels::solitaryWaveModel
Definition: solitaryWaveModel.H:49
Foam::waveModels::McCowan::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: McCowanWaveModel.C:290
Foam::waveModels::McCowan::setLevel
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
Definition: McCowanWaveModel.C:244
X0
scalarList X0(nSpecie, Zero)
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::waveModels::McCowan
McCowan wave model.
Definition: McCowanWaveModel.H:50
solitaryWaveModel.H
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
Foam::waveModels::McCowan::newtonRapsonF1
virtual scalar newtonRapsonF1(const scalar x0, const scalar H, const scalar h) const
Definition: McCowanWaveModel.C:95
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< scalar >
Foam::waveModels::McCowan::mn
virtual vector mn(const scalar H, const scalar h) const
Definition: McCowanWaveModel.C:77
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::waveModels::McCowan::eta
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.
Definition: McCowanWaveModel.C:52
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::waveModels::McCowan::newtonRapsonF2
virtual scalar newtonRapsonF2(const scalar x0, const scalar H, const scalar h, const scalar xa, const scalar m, const scalar n) const
Definition: McCowanWaveModel.C:151
Foam::waveModels::McCowan::~McCowan
virtual ~McCowan()=default
Destructor.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::waveModels::McCowan::McCowan
McCowan(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: McCowanWaveModel.C:272
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::waveModels::McCowan::setVelocity
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
Definition: McCowanWaveModel.C:302
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::waveModels::McCowan::TypeName
TypeName("McCowan")
Runtime type information.
y
scalar y
Definition: LISASMDCalcMethod1.H:14