StokesIWaveModel.C
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 \*---------------------------------------------------------------------------*/
28 
29 #include "StokesIWaveModel.H"
30 #include "mathematicalConstants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace waveModels
40 {
41  defineTypeNameAndDebug(StokesI, 0);
43  (
44  waveModel,
45  StokesI,
46  patch
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
54 Foam::scalar Foam::waveModels::StokesI::eta
55 (
56  const scalar H,
57  const scalar Kx,
58  const scalar x,
59  const scalar Ky,
60  const scalar y,
61  const scalar omega,
62  const scalar t,
63  const scalar phase
64 ) const
65 {
66  scalar phaseTot = Kx*x + Ky*y - omega*t + phase;
67 
68  return H*0.5*cos(phaseTot);
69 }
70 
71 
72 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
73 
75 (
76  const scalar h,
77  const scalar T
78 ) const
79 {
80  scalar L0 = mag(g_)*T*T/(2.0*mathematical::pi);
81  scalar L = L0;
82 
83  for (int i=1; i<=100; i++)
84  {
85  L = L0*tanh(2.0*mathematical::pi*h/L);
86  }
87 
88  return L;
89 }
90 
91 
93 (
94  const scalar H,
95  const scalar h,
96  const scalar Kx,
97  const scalar x,
98  const scalar Ky,
99  const scalar y,
100  const scalar omega,
101  const scalar t,
102  const scalar phase,
103  const scalar z
104 ) const
105 {
106  scalar k = sqrt(Kx*Kx + Ky*Ky);
107  scalar phaseTot = Kx*x + Ky*y - omega*t + phase;
108 
109  scalar u = H*0.5*omega*cos(phaseTot)*cosh(k*z)/sinh(k*h);
110  scalar w = H*0.5*omega*sin(phaseTot)*sinh(k*z)/sinh(k*h);
111  scalar v = u*sin(waveAngle_);
112  u *= cos(waveAngle_);
113 
114  return vector(u, v, w);
115 }
116 
117 
119 (
120  const scalar t,
121  const scalar tCoeff,
122  scalarField& level
123 ) const
124 {
125  const scalar waveOmega = mathematical::twoPi/wavePeriod_;
126  const scalar waveK = mathematical::twoPi/waveLength_;
127  const scalar waveKx = waveK*cos(waveAngle_);
128  const scalar waveKy = waveK*sin(waveAngle_);
129 
130  forAll(level, paddlei)
131  {
132  const scalar eta =
133  this->eta
134  (
135  waveHeight_,
136  waveKx,
137  xPaddle_[paddlei],
138  waveKy,
139  yPaddle_[paddlei],
140  waveOmega,
141  t,
142  wavePhase_
143  );
145  level[paddlei] = waterDepthRef_ + tCoeff*eta;
146  }
147 }
148 
149 
151 (
152  const scalar t,
153  const scalar tCoeff,
154  const scalarField& level
155 )
156 {
157  const scalar waveOmega = mathematical::twoPi/wavePeriod_;
158  const scalar waveK = mathematical::twoPi/waveLength_;
159  const scalar waveKx = waveK*cos(waveAngle_);
160  const scalar waveKy = waveK*sin(waveAngle_);
161 
162  forAll(U_, facei)
163  {
164  // Fraction of geometry represented by paddle - to be set
165  scalar fraction = 1;
166 
167  // Height - to be set
168  scalar z = 0;
169 
170  setPaddlePropeties(level, facei, fraction, z);
171 
172  if (fraction > 0)
173  {
174  const label paddlei = faceToPaddle_[facei];
175 
176  const vector Uf = UfBase
177  (
178  waveHeight_,
179  waterDepthRef_,
180  waveKx,
181  xPaddle_[paddlei],
182  waveKy,
183  yPaddle_[paddlei],
184  waveOmega,
185  t,
186  wavePhase_,
187  z
188  );
189 
190  U_[facei] = fraction*Uf*tCoeff;
191  }
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
197 
199 (
200  const dictionary& dict,
201  const fvMesh& mesh,
202  const polyPatch& patch,
203  const bool readFields
204 )
205 :
206  regularWaveModel(dict, mesh, patch, false)
207 {
208  if (readFields)
209  {
211  }
212 }
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 bool Foam::waveModels::StokesI::readDict(const dictionary& overrideDict)
218 {
219  if (regularWaveModel::readDict(overrideDict))
220  {
221  waveLength_ = waveLength(waterDepthRef_, wavePeriod_);
222 
223  return true;
224  }
225 
226  return false;
227 }
228 
229 
231 {
233 
234  os << " Wave type: " << waveType() << nl;
235 }
236 
237 
238 // ************************************************************************* //
Different types of constants.
dimensionedScalar tanh(const dimensionedScalar &ds)
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
scalar y
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
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.
constexpr scalar twoPi(2 *M_PI)
constexpr scalar pi(M_PI)
virtual scalar waveLength(const scalar h, const scalar T) const
Return the wavelength.
Vector< scalar > vector
Definition: vector.H:57
autoPtr< surfaceVectorField > Uf
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.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition: IOobject.H:955
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar sin(const dimensionedScalar &ds)
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionedScalar h
Planck constant.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
dimensionedScalar sinh(const dimensionedScalar &ds)
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
dimensionedScalar cosh(const dimensionedScalar &ds)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
StokesI(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)