GrimshawWaveModel.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::Grimshaw
29 
30 Description
31  Grimshaw wave model
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef waveModels_Grimshaw_H
36 #define waveModels_Grimshaw_H
37 
38 #include "solitaryWaveModel.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace waveModels
45 {
46 
47 /*---------------------------------------------------------------------------*\
48  Class Grimshaw Declaration
49 \*---------------------------------------------------------------------------*/
50 
51 class Grimshaw
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  //-
72  virtual scalar alfa
73  (
74  const scalar H,
75  const scalar h
76  ) const;
77 
78  //- Wave velocity
79  virtual vector Uf
80  (
81  const scalar H,
82  const scalar h,
83  const scalar x,
84  const scalar y,
85  const scalar theta,
86  const scalar t,
87  const scalar X0,
88  const scalar z
89  ) const;
90 
91  //- Set the water level
92  virtual void setLevel
93  (
94  const scalar t,
95  const scalar tCoeff,
96  scalarField& level
97  ) const;
98 
99  //- Calculate the wave model velocity
100  virtual void setVelocity
101  (
102  const scalar t,
103  const scalar tCoeff,
104  const scalarField& level
105  );
106 
107 
108 public:
109 
110  //- Runtime type information
111  TypeName("Grimshaw");
112 
113  //- Constructor
114  Grimshaw
115  (
116  const dictionary& dict,
117  const fvMesh& mesh,
118  const polyPatch& patch,
119  const bool readFields = true
120  );
121 
122  //- Destructor
123  virtual ~Grimshaw() = default;
124 
125 
126  // Public Member Functions
127 
128  //- Read from dictionary
129  virtual bool readDict(const dictionary& overrideDict);
130 
131  //- Info
132  virtual void info(Ostream& os) const;
133 };
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace waveModels
138 } // End namespace Foam
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 #endif
143 
144 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
virtual ~Grimshaw()=default
Destructor.
scalar y
dynamicFvMesh & mesh
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject *> &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
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.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition: IOobject.H:963
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
const dimensionedScalar h
Planck constant.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
scalarList X0(nSpecie, Zero)
Grimshaw(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
TypeName("Grimshaw")
Runtime type information.
virtual scalar alfa(const scalar H, const scalar h) const
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
Grimshaw wave model.
Namespace for OpenFOAM.