streamFunctionWaveModel.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 
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(streamFunction, 0);
43  (
44  waveModel,
45  streamFunction,
46  patch
47  );
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
54 Foam::scalar Foam::waveModels::streamFunction::eta
55 (
56  const scalar h,
57  const scalar kx,
58  const scalar ky,
59  const scalar T,
60  const scalar x,
61  const scalar y,
62  const scalar omega,
63  const scalar t,
64  const scalar phase
65 ) const
66 {
67 
68  const scalar k = sqrt(kx*kx + ky*ky);
69  scalar strfnAux = 0.0;
70  forAll(Ejs_, iterSF)
71  {
72  strfnAux +=
73  Ejs_[iterSF]*cos((iterSF + 1) *(kx*x + ky*y - omega*t + phase));
74  }
75 
76  return (1/k)*strfnAux;
77 }
78 
79 Foam::vector Foam::waveModels::streamFunction::Uf
80 (
81  const scalar h,
82  const scalar kx,
83  const scalar ky,
84  const scalar T,
85  const scalar x,
86  const scalar y,
87  const scalar omega,
88  const scalar t,
89  const scalar phase,
90  const scalar z
91 ) const
92 {
93  const scalar k = sqrt(kx*kx + ky*ky);
94  const scalar phaseTot = kx*x + ky*y - omega*t + phase;
95 
96  scalar u = 0.0;
97  scalar w = 0.0;
98 
99  forAll(Bjs_, iterSF2)
100  {
101  u +=
102  (iterSF2 + 1)*Bjs_[iterSF2]*cosh((iterSF2 + 1)*k*z)
103  /cosh((iterSF2 + 1)*k*h)*cos((iterSF2 + 1)*phaseTot);
104 
105  w +=
106  (iterSF2 + 1)*Bjs_[iterSF2]*sinh((iterSF2 + 1)*k*z)
107  /cosh((iterSF2 + 1)*k*h)*sin((iterSF2 + 1)*phaseTot);
108  }
109 
110  u = waveLength_/T - uMean_ + sqrt(mag(g_)/k)*u;
111  w = sqrt(mag(g_)/k)*w;
112 
113  scalar v = u*sin(waveAngle_);
114  u *= cos(waveAngle_);
115 
116  return vector(u, v, w);
117 }
118 
119 
120 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
121 
123 (
124  const scalar t,
125  const scalar tCoeff,
126  scalarField& level
127 ) const
128 {
129  const scalar waveOmega = mathematical::twoPi/wavePeriod_;
130  const scalar waveK = mathematical::twoPi/waveLength_;
131  const scalar waveKx = waveK*cos(waveAngle_);
132  const scalar waveKy = waveK*sin(waveAngle_);
133 
134  forAll(level, paddlei)
135  {
136  const scalar eta =
137  this->eta
138  (
139  waterDepthRef_,
140  waveKx,
141  waveKy,
142  wavePeriod_,
143  xPaddle_[paddlei],
144  yPaddle_[paddlei],
145  waveOmega,
146  t,
147  wavePhase_
148  );
150  level[paddlei] = waterDepthRef_ + tCoeff*eta;
151  }
152 }
153 
154 
156 (
157  const scalar t,
158  const scalar tCoeff,
159  const scalarField& level
160 )
161 {
162  const scalar waveOmega = mathematical::twoPi/wavePeriod_;
163  const scalar waveK = mathematical::twoPi/waveLength_;
164  const scalar waveKx = waveK*cos(waveAngle_);
165  const scalar waveKy = waveK*sin(waveAngle_);
166 
167  forAll(U_, facei)
168  {
169  // Fraction of geometry represented by paddle - to be set
170  scalar fraction = 1;
171 
172  // Height - to be set
173  scalar z = 0;
174 
175  setPaddlePropeties(level, facei, fraction, z);
176 
177  if (fraction > 0)
178  {
179  const label paddlei = faceToPaddle_[facei];
180 
181  const vector Uf = this->Uf
182  (
183  waterDepthRef_,
184  waveKx,
185  waveKy,
186  wavePeriod_,
187  xPaddle_[paddlei],
188  yPaddle_[paddlei],
189  waveOmega,
190  t,
191  wavePhase_,
192  z
193  );
194 
195  U_[facei] = fraction*Uf*tCoeff;
196  }
197  }
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202 
204 (
205  const dictionary& dict,
206  const fvMesh& mesh,
207  const polyPatch& patch,
208  const bool readFields
209 )
210 :
211  regularWaveModel(dict, mesh, patch, false),
212  uMean_(0),
213  Bjs_(),
214  Ejs_()
215 {
216  if (readFields)
217  {
219  }
220 }
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
226 {
227  if (regularWaveModel::readDict(overrideDict))
228  {
229  overrideDict.readEntry("uMean", uMean_);
230  overrideDict.readEntry("waveLength", waveLength_);
231  overrideDict.readEntry("Bjs", Bjs_);
232  overrideDict.readEntry("Ejs", Ejs_);
233 
234  return true;
235  }
236 
237  return false;
238 }
239 
240 
242 {
244 
245  os << " uMean : " << uMean_ << nl
246  << " Stream function wavelength : " << waveLength_ << nl
247  << " Bj coefficients : " << Bjs_ << nl
248  << " Ej coefficients : " << Ejs_ << nl;
249 }
250 
251 
252 // ************************************************************************* //
Different types of constants.
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
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
streamFunction(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
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
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
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)
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
Vector< scalar > vector
Definition: vector.H:57
autoPtr< surfaceVectorField > Uf
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)
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
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)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const std::string patch
OpenFOAM patch number as a std::string.
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
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)