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 -------------------------------------------------------------------------------
11 License
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"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace waveModels
37 {
38  defineTypeNameAndDebug(McCowan, 0);
40  (
41  waveModel,
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);
69 
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  {
283  readDict(dict);
284  }
285 }
286 
287 
288 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289 
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  );
333 
334  U_[facei] = fraction*Uf*tCoeff;
335  }
336  }
337 }
338 
339 
341 {
343 }
344 
345 
346 // ************************************************************************* //
Foam::waveModels::McCowan::Uf
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.
Definition: McCowanWaveModel.C:208
Foam::waveModels::solitaryWaveModel
Definition: solitaryWaveModel.H:49
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::waveModels::McCowan::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: McCowanWaveModel.C:290
Foam::cosh
dimensionedScalar cosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:271
Foam::waveModels::McCowan::setLevel
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
Definition: McCowanWaveModel.C:244
X0
scalarList X0(nSpecie, Zero)
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::waveModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(waveModel, shallowWaterAbsorption, patch)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::waveModels::solitaryWaveModel::readDict
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
Definition: solitaryWaveModel.C:87
Foam::IOobject::info
InfoProxy< IOobject > info() const
Return info proxy.
Definition: IOobject.H:689
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::waveModels::McCowan::newtonRapsonF1
virtual scalar newtonRapsonF1(const scalar x0, const scalar H, const scalar h) const
Definition: McCowanWaveModel.C:95
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< scalar >
Foam::waveModels::McCowan::mn
virtual vector mn(const scalar H, const scalar h) const
Definition: McCowanWaveModel.C:77
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::waveModels::McCowan::eta
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.
Definition: McCowanWaveModel.C:52
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::waveModels::McCowan::newtonRapsonF2
virtual scalar newtonRapsonF2(const scalar x0, const scalar H, const scalar h, const scalar xa, const scalar m, const scalar n) const
Definition: McCowanWaveModel.C:151
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::constant::physicoChemical::c2
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::waveModels::McCowan::McCowan
McCowan(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
Definition: McCowanWaveModel.C:272
McCowanWaveModel.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::waveModels::McCowan::setVelocity
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
Calculate the wave model velocity.
Definition: McCowanWaveModel.C:302
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
N
const Vector< label > N(dict.get< Vector< label >>("N"))
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::waveModels::defineTypeNameAndDebug
defineTypeNameAndDebug(waveAbsorptionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::sinh
dimensionedScalar sinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:270