waveMakerPointPatchVectorField.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) 2018-2019 IH-Cantabria
9 Copyright (C) 2018-2019 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
31#include "pointPatchFields.H"
33#include "Time.H"
34#include "gravityMeshObject.H"
35
36#include "polyMesh.H"
37#include "surfaceFields.H"
38#include "volFields.H"
39
40// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
41
44({
45 { motionTypes::piston, "piston" },
46 { motionTypes::flap, "flap" },
47 { motionTypes::solitary, "solitary" }
48});
49
50
51// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
52
54{
56
57 if (mag(gf.value()) < SMALL)
58 {
60 << "Gravity vector is not set. Please update "
61 << gf.uniformDimensionedVectorField::path()
62 << exit(FatalError);
63 }
64
65 return gf.value();
66}
67
68
70(
71 const scalar h,
72 const scalar T
73)
74{
75 const scalar L0 = mag(g())*T*T/(constant::mathematical::twoPi);
76 scalar L = L0;
77
78 for (label i=1; i<=100; ++i)
79 {
81 }
82
83 return L;
84}
85
86
88(
89 const scalar t
90) const
91{
92 return max(0, min(t/rampTime_, 1));
93}
94
95
97{
98 // Global patch extents
99 const vectorField& Cp = this->patch().localPoints();
100 const vectorField CpLocal(Cp);
101 boundBox bb(CpLocal, true);
102
103 const scalar xMin = bb.min().x();
104 const scalar xMax = bb.max().x();
105 const scalar yMin = bb.min().y();
106 const scalar yMax = bb.max().y();
107 zSpan_ = bb.max().z() - bb.min().z();
108
109 zMinGb_ = bb.min().z();
110 reduce(zMinGb_, minOp<scalar>());
111
112 // Global x, y positions of the paddle centres
113 xPaddle_.setSize(nPaddle_, 0);
114 yPaddle_.setSize(nPaddle_, 0);
115 const scalar xMid = xMin + 0.5*(xMax - xMin);
116 const scalar paddleDy = (yMax - yMin)/scalar(nPaddle_);
117
118 for (label paddlei = 0; paddlei < nPaddle_; ++paddlei)
119 {
120 xPaddle_[paddlei] = xMid;
121 yPaddle_[paddlei] = paddlei*paddleDy + yMin + 0.5*paddleDy;
122 }
123
124 // Local face centres
125 x_ = this->patch().localPoints().component(0);
126 y_ = this->patch().localPoints().component(1);
127 z_ = this->patch().localPoints().component(2);
128
129 // Local point-to-paddle addressing
130 pointToPaddle_.setSize(this->patch().size(), -1);
131
132 forAll(pointToPaddle_, ppi)
133 {
134 pointToPaddle_[ppi] = floor((y_[ppi] - yMin)/(paddleDy+0.01*paddleDy));
135 }
136}
137
138// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139
141(
142 const pointPatch& p,
144)
145:
147 motionType_(motionTypes::piston),
148 n_(Zero),
149 gHat_(Zero),
150 initialDepth_(0),
151 wavePeriod_(0),
152 waveHeight_(0),
153 wavePhase_(0),
154 waveAngle_(0),
155 startTime_(0),
156 rampTime_(1),
157 secondOrder_(false),
158 nPaddle_(0)
159{}
160
161
163(
164 const pointPatch& p,
166 const dictionary& dict
167)
168:
170 motionType_(motionTypeNames.get("motionType", dict)),
171 n_(dict.get<vector>("n")),
172 gHat_(Zero),
173 initialDepth_(dict.get<scalar>("initialDepth")),
174 wavePeriod_(dict.get<scalar>("wavePeriod")),
175 waveHeight_(dict.get<scalar>("waveHeight")),
176 wavePhase_(dict.get<scalar>("wavePhase")),
177 waveAngle_(dict.getOrDefault<scalar>("waveAngle", 0)),
178 startTime_
179 (
180 dict.getOrDefault<scalar>
181 (
182 "startTime",
183 db().time().startTime().value()
184 )
185 ),
186 rampTime_(dict.get<scalar>("rampTime")),
187 secondOrder_(dict.getOrDefault<bool>("secondOrder", false)),
188 nPaddle_(dict.getOrDefault<label>("nPaddle", 1))
189{
190 // Create the co-ordinate system
191 if (mag(n_) < SMALL)
192 {
194 << "Patch normal direction vector is not set. 'n' = " << n_
195 << exit(FatalIOError);
196 }
197 n_.normalise();
198
199 gHat_ = (g() - n_*(n_&g()));
200 if (mag(gHat_) < SMALL)
201 {
203 << "Patch normal and gravity directions must not be aligned. "
204 << "'n' = " << n_ << " 'g' = " << g()
205 << exit(FatalIOError);
206 }
207 gHat_.normalise();
208
210
212
214
215 if (!dict.found("value"))
216 {
217 updateCoeffs();
218 }
219}
220
221
223(
225 const pointPatch& p,
227 const pointPatchFieldMapper& mapper
228)
229:
230 fixedValuePointPatchField<vector>(ptf, p, iF, mapper),
231 motionType_(ptf.motionType_),
232 n_(ptf.n_),
233 gHat_(ptf.gHat_),
234 initialDepth_(ptf.initialDepth_),
235 wavePeriod_(ptf.wavePeriod_),
236 waveHeight_(ptf.waveHeight_),
237 wavePhase_(ptf.wavePhase_),
238 waveAngle_(ptf.waveAngle_),
239 startTime_(ptf.startTime_),
240 rampTime_(ptf.rampTime_),
241 secondOrder_(ptf.secondOrder_),
242 nPaddle_(ptf.nPaddle_)
243{}
244
245
247(
250)
251:
253 motionType_(ptf.motionType_),
254 n_(ptf.n_),
255 gHat_(ptf.gHat_),
256 initialDepth_(ptf.initialDepth_),
257 wavePeriod_(ptf.wavePeriod_),
258 waveHeight_(ptf.waveHeight_),
259 wavePhase_(ptf.wavePhase_),
260 waveAngle_(ptf.waveAngle_),
261 startTime_(ptf.startTime_),
262 rampTime_(ptf.rampTime_),
263 secondOrder_(ptf.secondOrder_),
264 nPaddle_(ptf.nPaddle_)
265{}
266
267
268// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
269
271{
272 if (this->updated())
273 {
274 return;
275 }
276
277 if (firstTime == 0)
278 {
279 // Set the reference water depth
280 if (initialDepth_ != 0 )
281 {
282 forAll(waterDepthRef_, paddlei)
283 {
284 waterDepthRef_[paddlei] = initialDepth_;
285 }
286 }
287 else
288 {
290 << "initialDepth is not set. Please update "
291 << abort(FatalError);
292 }
293
294
295 Info<< " WaterDepth at the wavepaddles = " << waterDepthRef_ << endl;
296 firstTime = 1;
297 }
298
299 const scalar t = db().time().value() - startTime_;
300
301 scalarField waveLength_(nPaddle_, -1);
302
303 scalarField waveK(nPaddle_, -1);
304 scalarField waveKx(nPaddle_, -1);
305 scalarField waveKy(nPaddle_, -1);
306
307 forAll(waveK, padddlei)
308 {
309 waveLength_[padddlei] =
310 waveLength(waterDepthRef_[padddlei], wavePeriod_);
311
312 waveK[padddlei] = constant::mathematical::twoPi/waveLength_[padddlei];
313 waveKx[padddlei] = waveK[padddlei]*cos(waveAngle_);
314 waveKy[padddlei] = waveK[padddlei]*sin(waveAngle_);
315 }
316 const scalar sigma = 2*constant::mathematical::pi/wavePeriod_;
317
318 switch (motionType_)
319 {
320 case motionTypes::flap:
321 {
322 const pointField& points = patch().localPoints();
323 scalarField motionX(patch().localPoints().size(), -1);
324
325 forAll(points, pointi)
326 {
327 const label paddlei = pointToPaddle_[pointi];
328
329 const scalar phaseTot =
330 waveKx[paddlei]*xPaddle_[paddlei]
331 + waveKy[paddlei]*yPaddle_[paddlei];
332
333 const scalar depthRef = waterDepthRef_[paddlei];
334 const scalar kh = waveK[paddlei]*depthRef;
335 const scalar pz = points[pointi].component(2);
336
337 const scalar m1 =
338 (4*sinh(kh)/(sinh(2*kh) + 2*kh))
339 * (sinh(kh) + 1/kh*(1 - cosh(kh)));
340
341 const scalar boardStroke = waveHeight_/m1;
342
343 motionX[pointi] = 0.5*boardStroke*sin(phaseTot - sigma*t);
344
345 if (secondOrder_)
346 {
347 motionX[pointi] +=
348 sqr(waveHeight_)/(16*depthRef)
349 * (3*cosh(kh)/pow3(sinh(kh)) - 2/m1)
350 * sin(phaseTot - 2*sigma*t);
351
352 }
353
354 motionX[pointi] *= 1.0 + (pz - zMinGb_ - depthRef)/depthRef;
355
356 }
357
358 Field<vector>::operator=(timeCoeff(t)*n_*motionX);
359
360 break;
361 }
362 case motionTypes::piston:
363 {
364 const pointField& points = patch().localPoints();
365 scalarField motionX(patch().localPoints().size(), -1);
366
367 forAll(points, pointi)
368 {
369 const label paddlei = pointToPaddle_[pointi];
370
371 const scalar phaseTot =
372 waveKx[paddlei]*xPaddle_[paddlei]
373 + waveKy[paddlei]*yPaddle_[paddlei];
374
375 const scalar depthRef = waterDepthRef_[paddlei];
376 const scalar kh = waveK[paddlei]*depthRef;
377 const scalar m1 = 2*(cosh(2*kh) - 1.0)/(sinh(2*kh) + 2*kh);
378 const scalar boardStroke = waveHeight_/m1;
379
380 motionX[pointi] = 0.5*boardStroke*sin(phaseTot - sigma*t);
381
382 if (secondOrder_)
383 {
384 motionX[pointi] +=
385 + sqr(waveHeight_)
386 / (32*depthRef)*(3*cosh(kh)/pow3(sinh(kh)) - 2.0/m1)
387 * sin(phaseTot - 2*sigma*t);
388 }
389 }
390
391 Field<vector>::operator=(timeCoeff(t)*n_*motionX);
392
393 break;
394 }
395 case motionTypes::solitary:
396 {
397 const pointField& points = patch().localPoints();
398 scalarField motionX(patch().localPoints().size(), -1);
399 const scalar magG = mag(g());
400
401 forAll(points, pointi)
402 {
403 const label paddlei = pointToPaddle_[pointi];
404 const scalar depthRef = waterDepthRef_[paddlei];
405
406 const scalar kappa = sqrt(0.75*waveHeight_/pow3(depthRef));
407 const scalar celerity = sqrt(magG*(depthRef + waveHeight_));
408 const scalar stroke = sqrt(16*waveHeight_*depthRef/3.0);
409 const scalar hr = waveHeight_/depthRef;
410 wavePeriod_ = 2.0/(kappa*celerity)*(3.8 + hr);
411 const scalar tSolitary = -0.5*wavePeriod_ + t;
412
413 // Newton-Raphson
414 scalar theta1 = 0;
415 scalar theta2 = 0;
416 scalar er = 10000;
417 const scalar error = 0.001;
418 while (er > error)
419 {
420 theta2 =
421 theta1
422 - (theta1 - kappa*celerity*tSolitary + hr*tanh(theta1))
423 /(1.0 + hr*(1.0/cosh(theta1))*(1.0/cosh(theta1)));
424
425 er = mag(theta1 - theta2);
426 theta1 = theta2;
427 }
428
429 motionX[pointi] =
430 waveHeight_/(kappa*depthRef)*tanh(theta1) + 0.5*stroke;
431 }
432
433 Field<vector>::operator=(n_*motionX);
434
435 break;
436 }
437 default:
438 {
440 << "Unhandled enumeration " << motionTypeNames[motionType_]
441 << abort(FatalError);
442 }
443 }
444
446}
447
448
450{
452 os.writeEntry("motionType", motionTypeNames[motionType_]);
453 os.writeEntry("n", n_);
454 os.writeEntry("initialDepth", initialDepth_);
455 os.writeEntry("wavePeriod", wavePeriod_);
456 os.writeEntry("waveHeight", waveHeight_);
457 os.writeEntry("wavePhase", wavePhase_);
458 os.writeEntry("waveAngle", waveAngle_);
459 os.writeEntry("startTime", startTime_);
460 os.writeEntry("rampTime", rampTime_);
461 os.writeEntry("secondOrder", secondOrder_);
462 os.writeEntry("nPaddle", nPaddle_);
463 writeEntry("value", os);
464}
465
466
467// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
468
469namespace Foam
470{
472 (
475 );
476}
477
478// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Foam::label startTime
const uniformDimensionedVectorField & g
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:641
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
const Type & value() const
Return const reference to value.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:77
A FixedValue boundary condition for pointField.
virtual bool write()
Write the output fields.
Gravitational acceleration vector Although termed a MeshObject it is registered on Time only and thus...
Foam::pointPatchFieldMapper.
const objectRegistry & db() const
Return local objectRegistry.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
Point motion boundary condition to generate waves based on either piston or flap motions.
static const Enum< motionTypes > motionTypeNames
Names for motion types.
scalarField waterDepthRef_
Calculated water depth at the patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual scalar timeCoeff(const scalar t) const
Return the time scaling coefficient.
virtual scalar waveLength(const scalar h, const scalar T)
Dispersion equation.
const vector & g()
Return the gravitational acceleration.
volScalarField & p
const volScalarField & T
const volScalarField & Cp
Definition: EEqn.H:7
bool
Definition: EEqn.H:20
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
constexpr scalar pi(M_PI)
constexpr scalar twoPi(2 *M_PI)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar sinh(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensionedScalar cos(const dimensionedScalar &ds)
#define makePointPatchTypeField(PatchTypeField, typePatchTypeField)
dictionary dict
volScalarField & h
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.
const scalar xMin
Definition: createFields.H:34
const scalar xMax
Definition: createFields.H:35
const vector L(dict.get< vector >("L"))