Loading...
Searching...
No Matches
McCowanWaveModel.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) 2017 IH-Cantabria
9 Copyright (C) 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 "McCowanWaveModel.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace waveModels
37{
40 (
42 McCowan,
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 vector vec = this->mn(H, h);
63 scalar mm = vec[0];
64 scalar nn = vec[1];
65
66 scalar C = sqrt(((mag(g_)*h)/mm)*tan(mm));
67 scalar ts = 3.5*h/sqrt(H/h);
68 scalar Xa = -C*t + ts - X0 + x*cos(theta) + y*sin(theta);
70 scalar xin = 0.5*H;
71 scalar etas = newtonRapsonF2(xin, H, h, Xa, mm, nn);
72 return etas;
73}
74
75
77(
78 const scalar H,
79 const scalar h
80) const
81{
82 // m
83 scalar xin = 1;
84 scalar m = newtonRapsonF1(xin, H, h);
85
86 // n
87 scalar c1 = sin(m + (1.0 + (2.0*H/(3.0*h))));
88 scalar n = (2.0/3.0)*sqr(c1);
89
90 return vector(m, n, n);
91}
92
93
95(
96 const scalar x0,
97 const scalar H,
98 const scalar h
99) const
100{
101 label N = 10000;
102 scalar eps = 1.e-5;
103 scalar maxval = 10000.0;
104
105 label iter = 1;
106 scalar x = x0;
107 scalar residual = 0;
108 while (iter <= N)
109 {
110 // f
111 scalar a = x + 1.0 + 2.0*H/(3.0*h);
112 scalar b = 0.5*x*(1.0 + H/h);
113 scalar c = 0.5*x*(1.0 + h/H);
114 scalar c1 = sin(a);
115 scalar fx = (2.0/3.0)*sqr(c1) - x*H/(h*tan(b));
116
117 residual = mag(fx);
118
119 if (residual < eps)
120 {
121 return x;
122 }
123 else if ((iter > 1) && (residual > maxval))
124 {
126 << "Newton-Raphson iterations diverging: "
127 << "iterations = " << iter
128 << ", residual = " << residual
129 << exit(FatalError);
130 }
131
132 // f-prime
133 scalar c2 = 1.0/tan(c);
134 scalar c3 = 1.0/sin(b);
135
136 scalar fprime = (4.0/3.0)*c1*cos(a) - c2*h/H - b*sqr(c3);
137
138 x -= fx/fprime;
139 iter++;
140 }
141
143 << "Failed to converge in " << iter << " iterations. Residual = "
144 << residual << nl << endl;
145
146 return x;
147}
148
149
151(
152 const scalar x0,
153 const scalar H,
154 const scalar h,
155 const scalar xa,
156 const scalar m,
157 const scalar n
158) const
159{
160 label N = 10000;
161 scalar eps = 1.e-5;
162 scalar maxval = 10000;
163
164 label iter = 1;
165 scalar x = x0;
166 scalar residual = 0;
167 while (iter <= N)
168 {
169 // f
170 scalar a = m*(1.0 + x/h);
171 scalar c1 = cos(a);
172 scalar c2 = sin(a);
173
174 scalar fx = x - (h*n/m*(c2/(c1 + cosh(m*xa/h))));
175
176 residual = mag(fx);
177
178 if (residual < eps)
179 {
180 return x;
181 }
182 else if ((iter > 1) && (residual > maxval))
183 {
185 << "Newton-Raphson iterations diverging: "
186 << "iterations = " << iter
187 << ", residual = " << residual
188 << exit(FatalError);
189 }
190
191 // f-prime
192 scalar c3 = cosh(xa*m/h) + c1;
193 scalar fprime = 1 - n/c3*(c1 - sqr(c2)/c3);
194
195 x -= fx/fprime;
196 iter++;
197 }
198
200 << "Failed to converge in " << iter << " iterations. Residual = "
201 << residual << nl << endl;
202
203 return x;
204}
205
206
208(
209 const scalar H,
210 const scalar h,
211 const scalar x,
212 const scalar y,
213 const scalar theta,
214 const scalar t,
215 const scalar X0,
216 const scalar z
217) const
218{
219 const vector vec = this->mn(H, h);
220 const scalar mm = vec[0];
221 const scalar nn = vec[1];
222
223 const scalar C = sqrt((mag(g_)*h)/mm*tan(mm));
224 const scalar ts = 3.5*h/sqrt(H/h);
225 const scalar Xa = -C*t + ts - X0 + x*cos(theta) + y*sin(theta);
226
227 scalar outa = C*nn*(1.0 + cos(mm*z/h)*cosh(mm*Xa/h));
228 scalar outb = sqr(cos(mm*z/h) + cosh(mm*Xa/h));
229
230 scalar u = outa/outb;
231
232 outa = C*nn*sin(mm*z/h)*sinh(mm*Xa/h);
233
234 const scalar w = outa/outb;
235
236 const scalar v = u*sin(waveAngle_);
237 u *= cos(waveAngle_);
238
239 return vector(u, v, w);
240}
241
242
244(
245 const scalar t,
246 const scalar tCoeff,
247 scalarField& level
248) const
249{
250 forAll(level, paddlei)
251 {
252 const scalar eta =
253 this->eta
254 (
255 waveHeight_,
256 waterDepthRef_,
257 xPaddle_[paddlei],
258 yPaddle_[paddlei],
259 waveAngle_,
260 t,
261 x0_
262 );
263
264 level[paddlei] = waterDepthRef_ + tCoeff*eta;
265 }
266}
267
268
269// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
270
272(
273 const dictionary& dict,
274 const fvMesh& mesh,
275 const polyPatch& patch,
276 const bool readFields
277)
278:
279 solitaryWaveModel(dict, mesh, patch, false)
280{
281 if (readFields)
282 {
284 }
285}
286
287
288// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289
290bool Foam::waveModels::McCowan::readDict(const dictionary& overrideDict)
291{
292 if (solitaryWaveModel::readDict(overrideDict))
293 {
294 return true;
295 }
296
297 return false;
298}
299
300
302(
303 const scalar t,
304 const scalar tCoeff,
305 const scalarField& level
306)
307{
308 forAll(U_, facei)
309 {
310 // Fraction of geometry represented by paddle - to be set
311 scalar fraction = 1;
312
313 // Height - to be set
314 scalar z = 0;
315
316 setPaddlePropeties(level, facei, fraction, z);
317
318 if (fraction > 0)
319 {
320 const label paddlei = faceToPaddle_[facei];
321
322 const vector Uf = this->Uf
323 (
324 waveHeight_,
325 waterDepthRef_,
326 xPaddle_[paddlei],
327 yPaddle_[paddlei],
328 waveAngle_,
329 t,
330 x0_,
331 z
332 );
334 U_[facei] = fraction*Uf*tCoeff;
335 }
336 }
337}
338
339
340void Foam::waveModels::McCowan::info(Ostream& os) const
341{
343}
344
345
346// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Graphite solid properties.
Definition C.H:49
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Base class for waveModels.
Definition waveModel.H:55
const vector & g_
Gravity.
Definition waveModel.H:73
vectorField U_
Velocity field.
Definition waveModel.H:166
scalar waterDepthRef_
Reference water depth / [m].
Definition waveModel.H:143
virtual void setPaddlePropeties(const scalarField &level, const label facei, scalar &fraction, scalar &z) const
Set the paddle coverage fraction and reference height.
Definition waveModel.C:183
scalarField yPaddle_
Paddle y coordinates / [m].
Definition waveModel.H:108
scalarField xPaddle_
Paddle x coordinates / [m].
Definition waveModel.H:103
labelList faceToPaddle_
Addressing from patch face index to paddle index.
Definition waveModel.H:113
virtual scalar newtonRapsonF2(const scalar x0, const scalar H, const scalar h, const scalar xa, const scalar m, const scalar n) const
virtual scalar newtonRapsonF1(const scalar x0, const scalar H, const scalar h) const
McCowan(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
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 vector mn(const scalar H, const scalar h) const
virtual scalar eta(const scalar H, const scalar h, const scalar x, const scalar y, const scalar theta, const scalar t, const scalar X0) const
Wave height.
virtual vector Uf(const scalar H, const scalar h, const scalar x, const scalar y, const scalar theta, const scalar t, const scalar X0, const scalar z) const
Wave velocity.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
solitaryWaveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
scalar waveHeight_
Wave height / [m].
scalar waveAngle_
Wave angle / [rad] (read in degrees).
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
scalarList X0(nSpecie, Zero)
dictionary dict
volScalarField & h
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const Vector< label > N(dict.get< Vector< label > >("N"))