Loading...
Searching...
No Matches
StokesIIWaveModel.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 "StokesIIWaveModel.H"
32
33using namespace Foam::constant;
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace waveModels
40{
43 (
46 patch
47 );
48}
49}
50
51
52// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53
54Foam::scalar Foam::waveModels::StokesII::eta
55(
56 const scalar H,
57 const scalar h,
58 const scalar Kx,
59 const scalar x,
60 const scalar Ky,
61 const scalar y,
62 const scalar omega,
63 const scalar t,
64 const scalar phase
65) const
66{
67 const scalar k = sqrt(Kx*Kx + Ky*Ky);
68 const scalar sigma = tanh(k*h);
69 const scalar phaseTot = Kx*x + Ky*y - omega*t + phase;
70
71 return
72 H*0.5*cos(phaseTot)
73 + k*H*H/4.0*(3.0 - sigma*sigma)/(4.0*pow3(sigma))*cos(2.0*phaseTot);
74}
75
76
77// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
78
80(
81 const scalar H,
82 const scalar h,
83 const scalar Kx,
84 const scalar x,
85 const scalar Ky,
86 const scalar y,
87 const scalar omega,
88 const scalar t,
89 const scalar phase,
90 const scalar z
91) const
92{
93 const scalar k = sqrt(Kx*Kx + Ky*Ky);
94 const scalar phaseTot = Kx*x + Ky*y - omega*t + phase;
95
96 scalar u =
97 H*0.5*omega*cos(phaseTot)*cosh(k*z)/sinh(k*h)
98 + 3.0/4.0*H*H/4.0*omega*k*cosh(2.0*k*z)/pow4(sinh(k*h))*cos(2.0*phaseTot);
99 scalar w =
100 H*0.5*omega*sin(phaseTot)*sinh(k*z)/sinh(k*h)
101 + 3.0/4.0*H*H/4.0*omega*k*sinh(2.0*k*z)/pow4(sinh(k*h))*sin(2.0*phaseTot);
102 scalar v = u*sin(waveAngle_);
103 u *= cos(waveAngle_);
104
105 return vector(u, v, w);
106}
107
108
110(
111 const scalar t,
112 const scalar tCoeff,
113 scalarField& level
114) const
115{
116 const scalar waveOmega = mathematical::twoPi/wavePeriod_;
117 const scalar waveK = mathematical::twoPi/waveLength_;
118 const scalar waveKx = waveK*cos(waveAngle_);
119 const scalar waveKy = waveK*sin(waveAngle_);
120
121 forAll(level, paddlei)
122 {
123 const scalar eta =
124 this->eta
125 (
126 waveHeight_,
127 waterDepthRef_,
128 waveKx,
129 xPaddle_[paddlei],
130 waveKy,
131 yPaddle_[paddlei],
132 waveOmega,
133 t,
134 wavePhase_
135 );
136
137 level[paddlei] = waterDepthRef_ + tCoeff*eta;
138 }
139}
140
141
142// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143
145(
146 const dictionary& dict,
147 const fvMesh& mesh,
148 const polyPatch& patch,
149 const bool readFields
150)
151:
152 StokesI(dict, mesh, patch, false)
153{
154 if (readFields)
155 {
157 }
158}
159
160
161// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
162
163bool Foam::waveModels::StokesII::readDict(const dictionary& overrideDict)
164{
165 if (StokesI::readDict(overrideDict))
166 {
167 return true;
168 }
169
170 return false;
171}
172
173
175{
177}
178
179
180// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
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
scalar waterDepthRef_
Reference water depth / [m].
Definition waveModel.H:143
scalarField yPaddle_
Paddle y coordinates / [m].
Definition waveModel.H:108
scalarField xPaddle_
Paddle x coordinates / [m].
Definition waveModel.H:103
Stokes II wave model.
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
Set the water level.
virtual vector UfBase(const scalar H, const scalar h, const scalar Kx, const scalar x, const scalar Ky, const scalar y, const scalar omega, const scalar t, const scalar phase, const scalar z) const
Wave velocity.
StokesII(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
StokesI(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
scalar waveHeight_
Wave height / [m].
scalar waveAngle_
Wave angle / [rad] (read in degrees).
#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))
OBJstream os(runTime.globalPath()/outputName)
constexpr scalar twoPi(2 *M_PI)
Different types of constants.
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 pow3(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
dictionary dict
volScalarField & h
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299