Airy.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 "Airy.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace waveModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const scalar depth,
48  const scalar amplitude,
49  const scalar g,
50  scalar (*modelCelerity)(scalar, scalar, scalar, scalar)
51 )
52 {
53  const bool haveLength = dict.found("length");
54  const bool havePeriod = dict.found("period");
55 
56  if (haveLength == havePeriod)
57  {
59  << "Exactly one of either length or period must be specified"
60  << exit(FatalIOError);
61  }
62 
63  if (haveLength)
64  {
65  return dict.lookup<scalar>("length");
66  }
67  else
68  {
69  const scalar period = dict.lookup<scalar>("period");
70 
71  scalar length0 = 0;
72  scalar length1 = 1.5*g*sqr(period)/(2*constant::mathematical::pi);
73 
74  // Bisect down to round off error to find the wave length
75  for (label i = 0; i < ceil(std::log2(1/small)); ++ i)
76  {
77  const scalar length = (length0 + length1)/2;
78 
79  const scalar t = length/modelCelerity(depth, amplitude, length, g);
80 
81  (t < period ? length0 : length1) = length;
82  }
83 
84  return (length0 + length1)/2;
85  }
86 }
87 
88 
89 // * * * * * * * * * * Static Protected Member Functions * * * * * * ** * * //
90 
91 Foam::scalar Foam::waveModels::Airy::k(const scalar length)
92 {
93  return 2*Foam::constant::mathematical::pi/length;
94 }
95 
96 
97 bool Foam::waveModels::Airy::deep(const scalar depth, const scalar length)
98 {
99  return depth*k(length) > log(great);
100 }
101 
102 
104 (
105  const scalar depth,
106  const scalar amplitude,
107  const scalar length,
108  const scalar g
109 )
110 {
111  return sqrt(g/k(length)*tanh(k(length)*depth));
112 }
113 
114 
115 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
116 
118 (
119  const scalar t,
120  const scalarField& x
121 ) const
122 {
123  return phase_ + k()*(x - celerity()*t);
124 }
125 
126 
128 (
129  const label i,
130  const scalar t,
131  const vector2DField& xz
132 ) const
133 {
134  const scalarField x(xz.component(0));
135  const scalarField z(xz.component(1));
136 
137  const scalarField phi(angle(t, x));
138 
139  const scalar kzGreat = log(i*great);
140  const scalarField kz(min(max(k()*z, - kzGreat), kzGreat));
141 
142  if (deep())
143  {
144  return i*exp(kz)*zip(cos(i*phi), sin(i*phi));
145  }
146  else
147  {
148  const scalar kd = k()*depth();
149  const scalarField kdz(max(scalar(0), kd + kz));
150 
151  return i*zip(cosh(i*kdz)*cos(i*phi), sinh(i*kdz)*sin(i*phi))/sinh(kd);
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
159 :
160  waveModel(wave),
161  depth_(wave.depth_),
162  amplitude_(wave.amplitude_, false),
163  length_(wave.length_),
164  phase_(wave.phase_)
165 {}
166 
167 
169 (
170  const dictionary& dict,
171  const scalar g,
172  const word& modelName,
173  scalar (*modelCelerity)(scalar, scalar, scalar, scalar)
174 )
175 :
176  waveModel(dict, g),
177  depth_(dict.lookupOrDefault<scalar>("depth", great)),
178  amplitude_(Function1<scalar>::New("amplitude", dict)),
179  length_(length(dict, depth_, amplitude(), g, modelCelerity)),
180  phase_(dict.lookup<scalar>("phase"))
181 {
182  const scalar c = modelCelerity(depth(), amplitude(), length(), g);
183 
184  Info<< waveModel::typeName << ": " << modelName
185  << ": period = " << length_/c
186  << ", length = " << length_ << endl;;
187 }
188 
189 
190 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 
193 {}
194 
195 
196 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
197 
199 (
200  const scalar t,
201  const scalarField& x
202 ) const
203 {
204  return amplitude(t)*cos(angle(t, x));
205 }
206 
207 
209 (
210  const scalar t,
211  const vector2DField& xz
212 ) const
213 {
214  const scalar ka = k()*amplitude(t);
215 
216  return Airy::celerity()*ka*vi(1, t, xz);
217 }
218 
219 
221 {
222  waveModel::write(os);
223 
224  if (!deep())
225  {
226  writeEntry(os, "depth", depth_);
227  }
228  writeEntry(os, amplitude_());
229  writeEntry(os, "length", length_);
230  writeEntry(os, "phase", phase_);
231 }
232 
233 
234 // ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:455
Run-time selectable general function of one variable.
Definition: Function1.H:64
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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
scalar g() const
Get the value of gravity.
Definition: waveModel.H:119
virtual void write(Ostream &os) const
Write.
Definition: waveModel.C:63
First-order wave model.
Definition: Airy.H:62
bool deep() const
Return whether shallow and intermediate effects are to be omitted.
Definition: Airy.H:124
virtual ~Airy()
Destructor.
Definition: Airy.C:192
scalar length() const
Get the length.
Definition: Airy.H:206
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
scalar k() const
The angular wavenumber [rad/m].
Definition: Airy.H:118
virtual void write(Ostream &os) const
Write.
Definition: Airy.C:220
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
scalar amplitude() const
Get the amplitude at steady state.
Definition: Airy.H:200
virtual scalar celerity() const
The wave celerity [m/s].
Definition: Airy.H:130
scalar depth() const
Get the depth.
Definition: Airy.H:188
tmp< scalarField > angle(const scalar t, const scalarField &x) const
Angle of the oscillation [rad].
Definition: Airy.C:118
tmp< vector2DField > vi(const label i, const scalar t, const vector2DField &xz) const
Return the non-dimensionalised i-th harmonic of the velocity.
Definition: Airy.C:128
Airy(const Airy &wave)
Construct a copy.
Definition: Airy.C:158
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
const dimensionedScalar c
Speed of light in a vacuum.
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
defineTypeNameAndDebug(Airy, 0)
addToRunTimeSelectionTable(waveModel, Airy, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
label log2(label i)
Return the log base 2 by successive bit-shifting of the given label.
Definition: label.H:79
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:251
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(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)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
tmp< vector2DField > zip(const tmp< scalarField > &x, const tmp< scalarField > &y)
Definition: vector2DField.C:31
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict