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-2019 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 {
35  defineTypeNameAndDebug(Stokes5, 0);
36  addToRunTimeSelectionTable(waveModel, Stokes5, objectRegistry);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const objectRegistry& db,
46  const dictionary& dict
47 )
48 :
49  Stokes2(db, dict)
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54 
56 {}
57 
58 
59 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
60 
62 (
63  const scalar t,
64  const scalarField& x
65 ) const
66 {
67  static const scalar kdGreat = log(great);
68  const scalar kd = min(max(k()*depth(), - kdGreat), kdGreat);
69  const scalar ka = k()*amplitude(t);
70 
71  const scalar S = deep() ? 0 : 1/cosh(2*kd), T = deep() ? 1 : tanh(kd);
72 
73  const scalar B31 =
74  - 3.0/8/pow3(1 - S)
75  *(
76  1 + 3*S + 3*sqr(S) + 2*pow3(S)
77  );
78 
79  const scalar B42 =
80  1.0/6/T/(3 + 2*S)/pow4(1 - S)
81  *(
82  6 - 26*S - 182*sqr(S) - 204*pow3(S) - 25*pow4(S) + 26*pow5(S)
83  );
84 
85  const scalar B44 =
86  1.0/24/T/(3 + 2*S)/pow4(1 - S)
87  *(
88  24 + 92*S + 122*sqr(S) + 66*pow3(S) + 67*pow4(S) + 34*pow5(S)
89  );
90 
91  const scalar B53 =
92  9.0/128/(3 + 2*S)/(4 + S)/pow6(1 - S)
93  *(
94  132 + 17*S - 2216*sqr(S) - 5897*pow3(S) - 6292*pow4(S)
95  - 2687*pow5(S) + 194*pow6(S) + 467*S*pow6(S) + 82*sqr(pow4(S))
96  );
97 
98  const scalar B55 =
99  5.0/384/(3 + 2*S)/(4 + S)/pow6(1 - S)
100  *(
101  300 + 1579*S + 3176*sqr(S) + 2949*pow3(S) + 1188*pow4(S)
102  + 675*pow5(S) + 1326*pow6(S) + 827*S*pow6(S) + 130*sqr(pow4(S))
103  );
104 
105  if (debug)
106  {
107  Info<< "B31 = " << B31 << endl
108  << "B42 = " << B42 << endl
109  << "B44 = " << B44 << endl
110  << "B53 = " << B53 << endl
111  << "B55 = " << B55 << endl;
112  }
113 
114  const scalarField phi(angle(t, x));
115 
116  return
117  Stokes2::elevation(t, x)
118  + (1/k())
119  *(
120  pow3(ka)*B31*(cos(phi) - cos(3*phi))
121  + pow4(ka)*(B42*cos(2*phi) + B44*cos(4*phi))
122  + pow5(ka)*(- (B53 + B55)*cos(phi) + B53*cos(3*phi) + B55*cos(5*phi))
123  );
124 }
125 
126 
128 (
129  const scalar t,
130  const vector2DField& xz
131 ) const
132 {
133  static const scalar kdGreat = log(great);
134  const scalar kd = min(max(k()*depth(), - kdGreat), kdGreat);
135  const scalar ka = k()*amplitude(t);
136 
137  const scalar S = deep() ? 0 : 1/cosh(2*kd);
138  const scalar SByA11 = deep() ? 0 : S*sinh(kd);
139 
140  const scalar A31ByA11 =
141  1.0/8/pow3(1 - S)
142  *(
143  - 4 - 20*S + 10*sqr(S) - 13*pow3(S)
144  );
145 
146  const scalar A33ByA11 =
147  1.0/8/pow3(1 - S)
148  *(
149  - 2*sqr(S) + 11*pow3(S)
150  );
151 
152  const scalar A42ByA11 =
153  SByA11/24/pow5(1 - S)
154  *(
155  12 - 14*S - 264*sqr(S) - 45*pow3(S) - 13*pow4(S)
156  );
157 
158  const scalar A44ByA11 =
159  SByA11/48/(3 + 2*S)/pow5(1 - S)
160  *(
161  10*sqr(S) - 174*pow3(S) + 291*pow4(S) + 278*pow5(S)
162  );
163 
164  const scalar A51ByA11 =
165  1.0/64/(3 + 2*S)/(4 + S)/pow6(1 - S)
166  *(
167  - 1184 + 32*S + 13232*sqr(S) + 21712*pow3(S) + 20940*pow4(S)
168  + 12554*pow5(S) - 500*pow6(S) - 3341*S*pow6(S) - 670*sqr(pow4(S))
169  );
170 
171  const scalar A53ByA11 =
172  1.0/32/(3 + 2*S)/pow6(1 - S)
173  *(
174  4*S + 105*sqr(S) + 198*pow3(S) - 1376*pow4(S) - 1302*pow5(S)
175  - 117*pow6(S) + 58*S*pow6(S)
176  );
177 
178  const scalar A55ByA11 =
179  1.0/64/(3 + 2*S)/(4 + S)/pow6(1 - S)
180  *(
181  - 6*pow3(S) + 272*pow4(S) - 1552*pow5(S) + 852*pow6(S)
182  + 2029*S*pow6(S) + 430*sqr(pow4(S))
183  );
184 
185  if (debug)
186  {
187  const scalar A11 = 1/sinh(kd);
188  Info<< "A31 = " << A31ByA11*A11 << endl
189  << "A33 = " << A33ByA11*A11 << endl
190  << "A42 = " << A42ByA11*A11 << endl
191  << "A44 = " << A44ByA11*A11 << endl
192  << "A51 = " << A51ByA11*A11 << endl
193  << "A53 = " << A53ByA11*A11 << endl
194  << "A55 = " << A55ByA11*A11 << endl;
195  }
196 
197  const vector2DField v1(vi(1, t, xz)), v3(vi(3, t, xz));
198 
199  return
200  Stokes2::velocity(t, xz)
201  + celerity()
202  *(
203  pow3(ka)*(A31ByA11*v1 + A33ByA11*v3)
204  + pow4(ka)*(A42ByA11*vi(2, t, xz) + A44ByA11*vi(4, t, xz))
205  + pow5(ka)*(A51ByA11*v1 + A53ByA11*v3 + A55ByA11*vi(5, t, xz))
206  );
207 }
208 
209 
210 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
defineTypeNameAndDebug(Airy, 0)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label k
Boltzmann constant.
virtual ~Stokes5()
Destructor.
Definition: Stokes5.C:55
dimensionedScalar pow5(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
phi
Definition: pEqn.H:104
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:62
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:87
dimensionedScalar cos(const dimensionedScalar &ds)
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:128
Stokes5(const objectRegistry &db, const dictionary &dict)
Construct from a database and a dictionary.
Definition: Stokes5.C:44
addToRunTimeSelectionTable(waveModel, Airy, objectRegistry)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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:62
dimensionedScalar pow3(const dimensionedScalar &ds)
Second-order wave model.
Definition: Stokes2.H:59
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
messageStream Info
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
Namespace for OpenFOAM.