cnoidalWaveModel.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-------------------------------------------------------------------------------
11License
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 "cnoidalWaveModel.H"
32#include "Elliptic.H"
33
34using namespace Foam::constant;
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace waveModels
41{
44 (
46 cnoidal,
47 patch
48 );
49}
50}
51
52
53// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54
56(
57 const scalar H,
58 const scalar d,
59 const scalar T,
60 scalar& mOut,
61 scalar& LOut
62) const
63{
64 const scalar mTolerance = 0.0001;
65 scalar mElliptic = 0.5;
66 scalar LElliptic = 0;
67 scalar phaseSpeed = 0;
68
69 scalar mError = 0.0;
70 scalar mMinError = GREAT;
71
72 while (mElliptic < 1.0)
73 {
74 scalar KElliptic, EElliptic;
75 Elliptic::ellipticIntegralsKE(mElliptic, KElliptic, EElliptic);
76
77 LElliptic = KElliptic*sqrt(16.0*pow3(d)*mElliptic/(3.0*H));
78
79 phaseSpeed =
80 sqrt(mag(g_)*d)
81 *(1.0 - H/d/2.0 + H/d/mElliptic*(1.0 - 3.0/2.0*EElliptic/KElliptic));
82
83 mError = mag(T - LElliptic/phaseSpeed);
84
85 if (mError <= mMinError)
86 {
87 mOut = mElliptic;
88 LOut = LElliptic;
89 mMinError = mError;
90 }
91
92 mElliptic += mTolerance;
93 }
94}
95
96
98(
99 const scalar H,
100 const scalar m,
101 const scalar kx,
102 const scalar ky,
103 const scalar T,
104 const scalar x,
105 const scalar y,
106 const scalar t
107) const
108{
109 scalar K, E;
111
112 const scalar uCnoidal =
113 K/mathematical::pi*(kx*x + ky*y - 2.0*mathematical::pi*t/T);
114
115 scalar sn, cn, dn;
116 Elliptic::JacobiSnCnDn(uCnoidal, m, sn, cn, dn);
117
118 return H*((1.0 - E/K)/m - 1.0 + sqr(cn));
119}
120
121
122Foam::scalar Foam::waveModels::cnoidal::eta1D
123(
124 const scalar H,
125 const scalar m,
126 const scalar t,
127 const scalar T
128) const
129{
130 scalar K, E;
132
133 const scalar uCnoidal = -2.0*K*(t/T);
134
135 scalar sn, cn, dn;
136 Elliptic::JacobiSnCnDn(uCnoidal, m, sn, cn, dn);
137
138 return H*((1.0 - E/K)/m - 1.0 + sqr(cn));
139}
140
141
142Foam::scalar Foam::waveModels::cnoidal::etaMeanSq
143(
144 const scalar H,
145 const scalar m,
146 const scalar T
147) const
148{
149 scalar eta = 0;
150 scalar etaSumSq = 0;
151
152 for (int i=0; i<1000; i++)
153 {
154 eta = eta1D(H, m, i*T/(1000.0), T);
155 etaSumSq += eta*eta;
156 }
157
158 etaSumSq /= 1000.0;
159 return etaSumSq;
160}
161
162
163Foam::vector Foam::waveModels::cnoidal::dEtaDx
164(
165 const scalar H,
166 const scalar m,
167 const scalar uCnoidal,
168 const scalar L,
169 const scalar K,
170 const scalar E
171) const
172{
173 const scalar dudx = 2.0*K/L;
174 const scalar dudxx = 2.0*K/L*dudx;
175 const scalar dudxxx = 2.0*K/L*dudxx;
176
177 scalar sn, cn, dn;
178 Elliptic::JacobiSnCnDn(uCnoidal, m, sn, cn, dn);
179
180 scalar d1 = -2.0*H*cn*dn*sn*dudx;
181 scalar d2 = 2.0*H*(dn*dn*sn*sn - cn*cn*dn*dn + m*cn*cn*sn*sn)*dudxx;
182 scalar d3 =
183 8.0*H
184 *(
185 cn*sn*dn*dn*dn*(-4.0 - 2.0*m)
186 + 4.0*m*cn*sn*sn*sn*dn
187 - 2.0*m*cn*cn*cn*sn*dn
188 )
189 *dudxxx;
190
191 return vector(d1, d2, d3);
192}
193
194
196(
197 const scalar H,
198 const scalar h,
199 const scalar m,
200 const scalar kx,
201 const scalar ky,
202 const scalar T,
203 const scalar x,
204 const scalar y,
205 const scalar t,
206 const scalar z
207) const
208{
209 scalar K, E;
211
212 const scalar uCnoidal =
213 K/mathematical::pi*(kx*x + ky*y - 2.0*mathematical::pi*t/T);
214 const scalar k = sqrt(kx*kx + ky*ky);
215 const scalar L = 2.0*mathematical::pi/k;
216 const scalar c = L/T;
217
218 const scalar etaCN = eta(H, m, kx, ky, T, x, y, t);
219 const vector etaX = this->dEtaDx(H, m, uCnoidal, L, K, E);
220 const scalar etaMS = etaMeanSq(H, m, T);
221
222 scalar u =
223 c*etaCN/h
224 - c*(etaCN*etaCN/h/h + etaMS*etaMS/h/h)
225 + 1.0/2.0*c*h*(1.0/3.0 - z*z/h/h)*etaX[1];
226
227 scalar w =
228 -c*z*(etaX[0]/h*(1.0 - 2.0*etaCN/h) + 1.0/6.0*h*(1.0 - z*z/h/h)*etaX[2]);
229
230 scalar v = u*sin(waveAngle_);
231 u *= cos(waveAngle_);
232
233 return vector(u, v, w);
234}
235
236
237// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
238
240(
241 const scalar t,
242 const scalar tCoeff,
243 scalarField& level
244) const
245{
246 const scalar waveK = mathematical::twoPi/waveLength_;
247 const scalar waveKx = waveK*cos(waveAngle_);
248 const scalar waveKy = waveK*sin(waveAngle_);
249
250 forAll(level, paddlei)
251 {
252 const scalar eta =
253 this->eta
254 (
255 waveHeight_,
256 m_,
257 waveKx,
258 waveKy,
259 wavePeriod_,
260 xPaddle_[paddlei],
261 yPaddle_[paddlei],
262 t
263 );
264
265 level[paddlei] = waterDepthRef_ + tCoeff*eta;
266 }
267}
268
269
271(
272 const scalar t,
273 const scalar tCoeff,
274 const scalarField& level
275)
276{
277 const scalar waveK = mathematical::twoPi/waveLength_;
278 const scalar waveKx = waveK*cos(waveAngle_);
279 const scalar waveKy = waveK*sin(waveAngle_);
280
281 forAll(U_, facei)
282 {
283 // Fraction of geometry represented by paddle - to be set
284 scalar fraction = 1;
285
286 // Height - to be set
287 scalar z = 0;
288
289 setPaddlePropeties(level, facei, fraction, z);
290
291 if (fraction > 0)
292 {
293 const label paddlei = faceToPaddle_[facei];
294
295 const vector Uf = this->Uf
296 (
297 waveHeight_,
298 waterDepthRef_,
299 m_,
300 waveKx,
301 waveKy,
302 wavePeriod_,
303 xPaddle_[paddlei],
304 yPaddle_[paddlei],
305 t,
306 z
307 );
308
309 U_[facei] = fraction*Uf*tCoeff;
310 }
311 }
312}
313
314
315// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
316
318(
319 const dictionary& dict,
320 const fvMesh& mesh,
321 const polyPatch& patch,
322 const bool readFields
323)
324:
325 regularWaveModel(dict, mesh, patch, false),
326 m_(0)
327{
328 if (readFields)
329 {
330 readDict(dict);
331 }
332}
333
334
335// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
336
338{
339 if (regularWaveModel::readDict(overrideDict))
340 {
341 // Initialise m parameter and wavelength
342 initialise
343 (
344 waveHeight_,
345 waterDepthRef_,
346 wavePeriod_,
347 m_,
348 waveLength_
349 );
350
351 return true;
352 }
353
354 return false;
355}
356
357
359{
361
362 os << " Cnoidal m parameter : " << m_ << nl;
363}
364
365
366// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
scalar y
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
InfoProxy< ensightCells > info() const
Return info proxy.
Definition: ensightCells.H:254
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
void initialise()
Initialise integral-scale box properties.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Base class for waveModels.
Definition: waveModel.H:61
const vector & g_
Gravity.
Definition: waveModel.H:73
Cnoidal wave model.
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
autoPtr< surfaceVectorField > Uf
OBJstream os(runTime.globalPath()/outputName)
void JacobiSnCnDn(const scalar u, const scalar m, scalar &Sn, scalar &Cn, scalar &Dn)
Definition: Elliptic.H:134
void ellipticIntegralsKE(const scalar m, scalar &K, scalar &E)
Definition: Elliptic.H:47
constexpr scalar pi(M_PI)
constexpr scalar twoPi(2 *M_PI)
const dimensionedScalar c
Speed of light in a vacuum.
Different types of constants.
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Vector< scalar > vector
Definition: vector.H:61
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
volScalarField & h
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const vector L(dict.get< vector >("L"))