Stokes2.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-2021 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 "Stokes2.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace waveModels
34 {
35  defineTypeNameAndDebug(Stokes2, 0);
36  addToRunTimeSelectionTable(waveModel, Stokes2, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * Static Protected Member Functions * * * * * * ** * * //
42 
44 (
45  const scalar depth,
46  const scalar amplitude,
47  const scalar length,
48  const scalar g
49 )
50 {
51  static const scalar kdGreat = log(great);
52  const scalar kd = min(max(k(length)*depth, - kdGreat), kdGreat);
53  const scalar ka = k(length)*amplitude;
54 
55  const scalar S = deep(depth, length) ? 0 : 1/cosh(2*kd);
56 
57  const scalar C0 = Airy::celerity(depth, amplitude, length, g);
58  const scalar C2ByC0 = (2 + 7*sqr(S))/4/sqr(1 - S);
59 
60  if (debug)
61  {
62  Info<< "C2 = " << C2ByC0*C0 << endl;
63  }
64 
65  return Airy::celerity(depth, amplitude, length, g) + sqr(ka)*C2ByC0*C0;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const dictionary& dict,
74  const scalar g,
75  const word& modelName,
76  scalar (*modelCelerity)(scalar, scalar, scalar, scalar)
77 )
78 :
79  Airy(dict, g, modelName, modelCelerity)
80 {}
81 
82 
83 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
84 
86 {}
87 
88 
89 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
90 
92 (
93  const scalar t,
94  const scalarField& x
95 ) const
96 {
97  static const scalar kdGreat = log(great);
98  const scalar kd = min(max(k()*depth(), - kdGreat), kdGreat);
99  const scalar ka = k()*amplitude(t);
100 
101  const scalar T = deep() ? 1 : tanh(kd);
102 
103  const scalar B22 = (3/sqr(T) - 1)/T/4;
104 
105  if (debug)
106  {
107  Info<< "B22 = " << B22 << endl;
108  }
109 
110  return
111  Airy::elevation(t, x)
112  + (1/k())*sqr(ka)*B22*cos(2*angle(t, x));
113 }
114 
115 
117 (
118  const scalar t,
119  const vector2DField& xz
120 ) const
121 {
122  static const scalar kdGreat = log(great);
123  const scalar kd = min(max(k()*depth(), - kdGreat), kdGreat);
124  const scalar ka = k()*amplitude(t);
125 
126  const scalar A22ByA11 = deep() ? 0 : 0.375/pow3(sinh(kd));
127 
128  if (debug)
129  {
130  const scalar A11 = 1/sinh(kd);
131  Info<< "A22 = " << A22ByA11*A11 << endl;
132  }
133 
134  return
135  Airy::velocity(t, xz)
136  + Airy::celerity()*sqr(ka)*A22ByA11*vi(2, t, xz);
137 }
138 
139 
140 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
defineTypeNameAndDebug(Airy, 0)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< scalarField > elevation(const scalar t, const scalarField &x) const
Get the wave elevation at a given time and local coordinates. Local.
Definition: Airy.C:199
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
First-order wave model.
Definition: Airy.H:59
virtual tmp< vector2DField > velocity(const scalar t, const vector2DField &xz) const
Get the wave velocity at a given time and local coordinates. Local.
Definition: Airy.C:209
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Airy.H:130
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 ~Stokes2()
Destructor.
Definition: Stokes2.C:85
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:117
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensionedScalar pow3(const dimensionedScalar &ds)
addToRunTimeSelectionTable(waveModel, Airy, dictionary)
dimensionedScalar sinh(const dimensionedScalar &ds)
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Stokes2.H:77
messageStream Info
dimensionedScalar cosh(const dimensionedScalar &ds)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
Stokes2(const dictionary &dict, const scalar g, const word &modelName=Stokes2::typeName, scalar(*modelCelerity)(scalar, scalar, scalar, scalar)=&Stokes2::celerity)
Construct from a dictionary and gravity.
Definition: Stokes2.C:72