Stokes5.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Stokes5.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace waveModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * Static Protected Member Functions * * * * * * ** * * //
42 
44 {
45  static const scalar kdGreat = log(great);
46  const scalar kd = min(max(coeffs.k()*coeffs.depth, - kdGreat), kdGreat);
47  const scalar ka = coeffs.k()*coeffs.amplitude;
48 
49  const scalar S = coeffs.deep() ? 0 : 1/cosh(2*kd);
50 
51  const scalar C0 = coeffs.celerity();
52  const scalar C4ByC0 =
53  1.0/32/pow5(1 - S)
54  *(4 + 32*S - 116*sqr(S) - 400*pow3(S) - 71*pow4(S) + 146*pow5(S));
55 
56  if (debug)
57  {
58  Info<< "C4 = " << C4ByC0*C0 << endl;
59  }
60 
61  return Stokes2::celerity(coeffs) + pow4(ka)*C4ByC0*C0;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
66 
68 (
69  const dictionary& dict,
70  const scalar g,
71  const word& modelName,
72  scalar (*celerityPtr)(const AiryCoeffs&)
73 )
74 :
75  Stokes2(dict, g, modelName, celerityPtr)
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
80 
82 {}
83 
84 
85 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
86 
88 {
89  return celerity(coeffs());
90 }
91 
92 
94 (
95  const scalar t,
96  const scalarField& x
97 ) const
98 {
99  const AiryCoeffs coeffs = this->coeffs();
100 
101  static const scalar kdGreat = log(great);
102  const scalar kd = min(max(coeffs.k()*depth(), - kdGreat), kdGreat);
103  const scalar ka = coeffs.k()*amplitude(t);
104 
105  const scalar S = coeffs.deep() ? 0 : 1/cosh(2*kd);
106  const scalar T = coeffs.deep() ? 1 : tanh(kd);
107 
108  const scalar B31 =
109  - 3.0/8/pow3(1 - S)
110  *(
111  1 + 3*S + 3*sqr(S) + 2*pow3(S)
112  );
113 
114  const scalar B42 =
115  1.0/6/T/(3 + 2*S)/pow4(1 - S)
116  *(
117  6 - 26*S - 182*sqr(S) - 204*pow3(S) - 25*pow4(S) + 26*pow5(S)
118  );
119 
120  const scalar B44 =
121  1.0/24/T/(3 + 2*S)/pow4(1 - S)
122  *(
123  24 + 92*S + 122*sqr(S) + 66*pow3(S) + 67*pow4(S) + 34*pow5(S)
124  );
125 
126  const scalar B53 =
127  9.0/128/(3 + 2*S)/(4 + S)/pow6(1 - S)
128  *(
129  132 + 17*S - 2216*sqr(S) - 5897*pow3(S) - 6292*pow4(S)
130  - 2687*pow5(S) + 194*pow6(S) + 467*S*pow6(S) + 82*sqr(pow4(S))
131  );
132 
133  const scalar B55 =
134  5.0/384/(3 + 2*S)/(4 + S)/pow6(1 - S)
135  *(
136  300 + 1579*S + 3176*sqr(S) + 2949*pow3(S) + 1188*pow4(S)
137  + 675*pow5(S) + 1326*pow6(S) + 827*S*pow6(S) + 130*sqr(pow4(S))
138  );
139 
140  if (debug)
141  {
142  Info<< "B31 = " << B31 << endl
143  << "B42 = " << B42 << endl
144  << "B44 = " << B44 << endl
145  << "B53 = " << B53 << endl
146  << "B55 = " << B55 << endl;
147  }
148 
149  const scalarField phi(coeffs.angle(phase(), t, x));
150 
151  return
153  + (1/coeffs.k())
154  *(
155  pow3(ka)*B31*(cos(phi) - cos(3*phi))
156  + pow4(ka)*(B42*cos(2*phi) + B44*cos(4*phi))
157  + pow5(ka)*(- (B53 + B55)*cos(phi) + B53*cos(3*phi) + B55*cos(5*phi))
158  );
159 }
160 
161 
163 (
164  const scalar t,
165  const vector2DField& xz
166 ) const
167 {
168  const AiryCoeffs coeffs = this->coeffs();
169 
170  static const scalar kdGreat = log(great);
171  const scalar kd = min(max(coeffs.k()*depth(), - kdGreat), kdGreat);
172  const scalar ka = coeffs.k()*amplitude(t);
173 
174  const scalar S = coeffs.deep() ? 0 : 1/cosh(2*kd);
175  const scalar SByA11 = coeffs.deep() ? 0 : S*sinh(kd);
176 
177  const scalar A31ByA11 =
178  1.0/8/pow3(1 - S)
179  *(
180  - 4 - 20*S + 10*sqr(S) - 13*pow3(S)
181  );
182 
183  const scalar A33ByA11 =
184  1.0/8/pow3(1 - S)
185  *(
186  - 2*sqr(S) + 11*pow3(S)
187  );
188 
189  const scalar A42ByA11 =
190  SByA11/24/pow5(1 - S)
191  *(
192  12 - 14*S - 264*sqr(S) - 45*pow3(S) - 13*pow4(S)
193  );
194 
195  const scalar A44ByA11 =
196  SByA11/48/(3 + 2*S)/pow5(1 - S)
197  *(
198  10*sqr(S) - 174*pow3(S) + 291*pow4(S) + 278*pow5(S)
199  );
200 
201  const scalar A51ByA11 =
202  1.0/64/(3 + 2*S)/(4 + S)/pow6(1 - S)
203  *(
204  - 1184 + 32*S + 13232*sqr(S) + 21712*pow3(S) + 20940*pow4(S)
205  + 12554*pow5(S) - 500*pow6(S) - 3341*S*pow6(S) - 670*sqr(pow4(S))
206  );
207 
208  const scalar A53ByA11 =
209  1.0/32/(3 + 2*S)/pow6(1 - S)
210  *(
211  4*S + 105*sqr(S) + 198*pow3(S) - 1376*pow4(S) - 1302*pow5(S)
212  - 117*pow6(S) + 58*S*pow6(S)
213  );
214 
215  const scalar A55ByA11 =
216  1.0/64/(3 + 2*S)/(4 + S)/pow6(1 - S)
217  *(
218  - 6*pow3(S) + 272*pow4(S) - 1552*pow5(S) + 852*pow6(S)
219  + 2029*S*pow6(S) + 430*sqr(pow4(S))
220  );
221 
222  if (debug)
223  {
224  const scalar A11 = 1/sinh(kd);
225  Info<< "A31 = " << A31ByA11*A11 << endl
226  << "A33 = " << A33ByA11*A11 << endl
227  << "A42 = " << A42ByA11*A11 << endl
228  << "A44 = " << A44ByA11*A11 << endl
229  << "A51 = " << A51ByA11*A11 << endl
230  << "A53 = " << A53ByA11*A11 << endl
231  << "A55 = " << A55ByA11*A11 << endl;
232  }
233 
234  auto v = [&](const label i) { return coeffs.vi(i, phase(), t, xz); };
235 
236  const vector2DField v1(v(1)), v3(v(3));
237 
238  return
239  Stokes2::velocity(t, xz)
240  + Airy::celerity()
241  *(
242  pow3(ka)*(A31ByA11*v1 + A33ByA11*v3)
243  + pow4(ka)*(A42ByA11*v(2) + A44ByA11*v(4))
244  + pow5(ka)*(A51ByA11*v1 + A53ByA11*v3 + A55ByA11*v(5))
245  );
246 }
247 
248 
249 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
Generic base class for waves. Derived classes must implement field functions which return the elevati...
Definition: waveModel.H:56
Calculation engine for the Airy wave model and other models that are a correction on top of the Airy ...
Definition: AiryCoeffs.H:52
const scalar depth
Depth [m].
Definition: AiryCoeffs.H:76
bool deep() const
Return whether shallow and intermediate effects are to be omitted.
Definition: AiryCoeffs.C:99
tmp< vector2DField > vi(const label i, const scalar phase, const scalar t, const vector2DField &xz) const
Return the non-dimensionalised i-th harmonic of the velocity.
Definition: AiryCoeffs.C:140
const scalar amplitude
Amplitude [m].
Definition: AiryCoeffs.H:79
tmp< scalarField > angle(const scalar phase, const scalar t, const scalarField &x) const
Angle of the oscillation [rad].
Definition: AiryCoeffs.C:118
static scalar celerity(const AiryCoeffs &coeffs)
The wave celerity [m/s].
Definition: AiryCoeffs.C:105
scalar k() const
The angular wavenumber [rad/m].
Definition: AiryCoeffs.C:93
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Airy.C:131
AiryCoeffs coeffs() const
Return the wave coefficients at steady state.
Definition: AiryI.H:39
Second-order wave model.
Definition: Stokes2.H:68
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Stokes2.C:85
virtual tmp< scalarField > elevation(const scalar t, const scalarField &x) const
Get the wave elevation at a given time and local coordinates. Local.
Definition: Stokes2.C:92
virtual tmp< vector2DField > velocity(const scalar t, const vector2DField &xz) const
Get the wave velocity at a given time and local coordinates. Local.
Definition: Stokes2.C:119
Fifth-order wave model.
Definition: Stokes5.H:67
virtual ~Stokes5()
Destructor.
Definition: Stokes5.C:81
Stokes5(const dictionary &dict, const scalar g, const word &modelName=Stokes5::typeName, scalar(*celerityPtr)(const AiryCoeffs &)=&Stokes5::celerity)
Construct from a dictionary and gravity.
Definition: Stokes5.C:68
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Stokes5.C:87
virtual tmp< scalarField > elevation(const scalar t, const scalarField &x) const
Get the wave elevation at a given time and local coordinates. Local.
Definition: Stokes5.C:94
virtual tmp< vector2DField > velocity(const scalar t, const vector2DField &xz) const
Get the wave velocity at a given time and local coordinates. Local.
Definition: Stokes5.C:163
A class for handling words, derived from string.
Definition: word.H:62
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
defineTypeNameAndDebug(Airy, 0)
addToRunTimeSelectionTable(waveModel, Airy, dictionary)
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict