BoussinesqWaveModel.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 "BoussinesqWaveModel.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace waveModels
37{
40 (
43 patch
44 );
45}
46}
47
48
49// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50
52(
53 const scalar H,
54 const scalar h,
55 const scalar x,
56 const scalar y,
57 const scalar theta,
58 const scalar t,
59 const scalar X0
60) const
61{
62 scalar C = sqrt(mag(g_)*(H + h));
63 scalar ts = 3.5*h/sqrt(H/h);
64 scalar aux = sqrt(3.0*H/(4.0*h))/h;
65 scalar Xa = -C*t + ts - X0 + x*cos(theta) + y*sin(theta);
66
67 return H*1.0/sqr(cosh(aux*Xa));
68}
69
70
71Foam::vector Foam::waveModels::Boussinesq::Deta
72(
73 const scalar H,
74 const scalar h,
75 const scalar x,
76 const scalar y,
77 const scalar theta,
78 const scalar t,
79 const scalar X0
80) const
81{
82 vector deta(Zero);
83
84 scalar C = sqrt(mag(g_)*(H + h));
85 scalar ts = 3.5*h/sqrt(H/h);
86 scalar a = sqrt(3*H/(4*h))/h;
87 scalar Xa = -C*t + ts - X0 + x*cos(theta) + y*sin(theta);
88 scalar expTerm = exp(2*a*Xa);
89 scalar b = 8*a*h*expTerm;
90
91 deta[0] =
92 b*(1 - expTerm)
93 /pow3(1 + expTerm);
94
95 deta[1] =
96 2*a*b*(exp(4*a*Xa) - 4*expTerm + 1)
97 /pow4(1 + expTerm);
98
99 deta[2] =
100 -4*sqr(a)*b*(exp(6*a*Xa) - 11*exp(4*a*Xa) + 11*expTerm - 1)
101 /pow5(1 + expTerm);
102
103 return deta;
104}
105
106
108(
109 const scalar H,
110 const scalar h,
111 const scalar x,
112 const scalar y,
113 const scalar theta,
114 const scalar t,
115 const scalar X0,
116 const scalar z
117) const
118{
119 scalar C = sqrt(mag(g_)*(H + h));
120 scalar eta = this->eta(H, h, x, y, theta, t, X0);
121 vector Deta = this->Deta(H, h, x, y, theta, t, X0);
122
123 scalar u =
124 C*eta/h
125 *(
126 1.0
127 - eta/(4.0*h)
128 + sqr(h)/(3.0*eta)*(1.0 - 3.0/2.0*sqr(z/h))*Deta[1]
129 );
130
131 scalar w =
132 -C*z/h
133 *(
134 (1.0 - eta/(2.0*h))*Deta[0]
135 + sqr(h)/3.0*(1.0 - 1.0/2.0*sqr(z/h))*Deta[2]
136 );
137
138 scalar v = u*sin(waveAngle_);
139 u *= cos(waveAngle_);
140
141 return vector(u, v, w);
142}
143
144
145// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
146
148(
149 const scalar t,
150 const scalar tCoeff,
151 scalarField& level
152) const
153{
154 forAll(level, paddlei)
155 {
156 const scalar eta =
157 this->eta
158 (
159 waveHeight_,
160 waterDepthRef_,
161 xPaddle_[paddlei],
162 yPaddle_[paddlei],
163 waveAngle_,
164 t,
165 x0_
166 );
167
168 level[paddlei] = waterDepthRef_ + tCoeff*eta;
169 }
170}
171
172
174(
175 const scalar t,
176 const scalar tCoeff,
177 const scalarField& level
178)
179{
180 forAll(U_, facei)
181 {
182 // Fraction of geometry represented by paddle - to be set
183 scalar fraction = 1;
184
185 // Height - to be set
186 scalar z = 0;
187
188 setPaddlePropeties(level, facei, fraction, z);
189
190 if (fraction > 0)
191 {
192 const label paddlei = faceToPaddle_[facei];
193
194 const vector Uf = this->Uf
195 (
196 waveHeight_,
197 waterDepthRef_,
198 xPaddle_[paddlei],
199 yPaddle_[paddlei],
200 waveAngle_,
201 t,
202 x0_,
203 z
204 );
205
206 U_[facei] = fraction*Uf*tCoeff;
207 }
208 }
209}
210
211
212// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
213
215(
216 const dictionary& dict,
217 const fvMesh& mesh,
218 const polyPatch& patch,
219 const bool readFields
220)
221:
222 solitaryWaveModel(dict, mesh, patch, false)
223{
224 if (readFields)
225 {
226 readDict(dict);
227 }
228}
229
230
231// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
232
234{
235 if (solitaryWaveModel::readDict(overrideDict))
236 {
237 return true;
238 }
239
240 return false;
241}
242
243
245{
247}
248
249
250// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Incompressible gas equation of state using the Boussinesq approximation for the density as a function...
Definition: Boussinesq.H:96
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
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
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
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)
Namespace for OpenFOAM.
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar cosh(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)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Vector< scalar > vector
Definition: vector.H:61
dimensionedScalar cos(const dimensionedScalar &ds)
scalarList X0(nSpecie, Zero)
volScalarField & C
dictionary dict
volScalarField & h
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333