irregularMultiDirectionalWaveModel.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 -------------------------------------------------------------------------------
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::irregularMultiDirectional
29 
30 Description
31  Multi-directional irregular wave model
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef waveModels_irregularMultiDirectional_H
36 #define waveModels_irregularMultiDirectional_H
37 
38 #include "irregularWaveModel.H"
39 #include "scalarList.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace waveModels
46 {
47 
48 /*---------------------------------------------------------------------------*\
49  Class irregularMultiDirectional Declaration
50 \*---------------------------------------------------------------------------*/
51 
53 :
54  public irregularWaveModel
55 {
56 private:
57 
58  // Proteced Data
59 
60  //- Wave periods for irregularMultiDirectional case [seconds]
61  List<scalarList> irregWavePeriods_;
62 
63  //- Wave heights for irregularMultiDirectional case [metres]
64  List<scalarList> irregWaveHeights_;
65 
66  //- Wave lengths for irregularMultiDirectional case [metres]
67  List<scalarList> irregWaveLengths_;
68 
69  //- Wave phases for irregularMultiDirectional case [radians]
70  List<scalarList> irregWavePhases_;
71 
72  //- Direction of propagation for irregularMultiDirectional case
73  //- [degrees] from X axis
74  List<scalarList> irregWaveDirs_;
75 
76 
77  // Private Member Functions
78 
79  //- First Order Wave height (same as StokesI)
80  virtual scalar eta
81  (
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;
91 
92 
93 protected:
94 
95  // Protected Member Functions
96 
97  //- Return the wavelength
98  virtual scalar waveLength
99  (
100  const scalar h,
101  const scalar T
102  ) const;
103 
104  //- Wave velocity
105  virtual vector Uf
106  (
107  const scalar d,
108  const scalar x,
109  const scalar y,
110  const scalar t,
111  const scalar z
112  ) const;
113 
114  //-
115  virtual vector uMultiDirec
116  (
117  const scalar irregH,
118  const scalar irregWaveOmega,
119  const scalar phaseTot,
120  const scalar irregWaveKs,
121  const scalar z,
122  const scalar h,
123  const scalar irregDir
124  ) const;
125 
126  //- Set the water level
127  virtual void setLevel
128  (
129  const scalar t,
130  const scalar tCoeff,
131  scalarField& level
132  ) const;
133 
134  //- Calculate the wave model velocity
135  virtual void setVelocity
136  (
137  const scalar t,
138  const scalar tCoeff,
139  const scalarField& level
140  );
141 
142 
143 public:
144 
145  //- Runtime type information
146  TypeName("irregularMultiDirectional");
147 
148  //- Constructor
150  (
151  const dictionary& dict,
152  const fvMesh& mesh,
153  const polyPatch& patch,
154  const bool readFields = true
155  );
156 
157  //- Destructor
158  virtual ~irregularMultiDirectional() = default;
159 
160 
161  // Public Member Functions
162 
163  //- Read from dictionary
164  virtual bool readDict(const dictionary& overrideDict);
165 
166  //- Info
167  virtual void info(Ostream& os) const;
168 };
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace waveModels
174 } // End namespace Foam
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 #endif
179 
180 // ************************************************************************* //
Foam::waveModels::irregularWaveModel
Definition: irregularWaveModel.H:49
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::waveModels::irregularMultiDirectional::setLevel
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
Definition: irregularMultiDirectionalWaveModel.C:139
Foam::waveModels::irregularMultiDirectional::waveLength
virtual scalar waveLength(const scalar h, const scalar T) const
Return the wavelength.
Definition: irregularMultiDirectionalWaveModel.C:74
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
Foam::waveModels::irregularMultiDirectional::uMultiDirec
virtual vector uMultiDirec(const scalar irregH, const scalar irregWaveOmega, const scalar phaseTot, const scalar irregWaveKs, const scalar z, const scalar h, const scalar irregDir) const
Definition: irregularMultiDirectionalWaveModel.C:178
scalarList.H
Foam::waveModels::irregularMultiDirectional::setVelocity
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
Definition: irregularMultiDirectionalWaveModel.C:205
Foam::Field< scalar >
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
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::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::waveModels::irregularMultiDirectional::Uf
virtual vector Uf(const scalar d, const scalar x, const scalar y, const scalar t, const scalar z) const
Wave velocity.
Definition: irregularMultiDirectionalWaveModel.C:92
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::irregularMultiDirectional::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: irregularMultiDirectionalWaveModel.C:262
Foam::waveModels::irregularMultiDirectional::~irregularMultiDirectional
virtual ~irregularMultiDirectional()=default
Destructor.
irregularWaveModel.H
Foam::waveModels::irregularMultiDirectional::TypeName
TypeName("irregularMultiDirectional")
Runtime type information.
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::List< scalarList >
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::waveModels::irregularMultiDirectional::irregularMultiDirectional
irregularMultiDirectional(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: irregularMultiDirectionalWaveModel.C:243
Foam::waveModels::irregularMultiDirectional
Multi-directional irregular wave model.
Definition: irregularMultiDirectionalWaveModel.H:51
y
scalar y
Definition: LISASMDCalcMethod1.H:14