Go to the documentation of this file.
55 void Foam::waveModels::cnoidal::initialise
64 const scalar mTolerance = 0.0001;
65 scalar mElliptic = 0.5;
67 scalar phaseSpeed = 0;
70 scalar mMinError = GREAT;
72 while (mElliptic < 1.0)
74 scalar KElliptic, EElliptic;
77 LElliptic = KElliptic*
sqrt(16.0*
pow3(d)*mElliptic/(3.0*
H));
81 *(1.0 -
H/d/2.0 +
H/d/mElliptic*(1.0 - 3.0/2.0*EElliptic/KElliptic));
83 mError =
mag(
T - LElliptic/phaseSpeed);
85 if (mError <= mMinError)
92 mElliptic += mTolerance;
97 Foam::scalar Foam::waveModels::cnoidal::eta
112 const scalar uCnoidal =
118 return H*((1.0 - E/
K)/m - 1.0 +
sqr(cn));
122 Foam::scalar Foam::waveModels::cnoidal::eta1D
133 const scalar uCnoidal = -2.0*
K*(t/
T);
138 return H*((1.0 - E/
K)/m - 1.0 +
sqr(cn));
142 Foam::scalar Foam::waveModels::cnoidal::etaMeanSq
152 for (
int i=0; i<1000; i++)
154 eta = eta1D(
H, m, i*
T/(1000.0),
T);
167 const scalar uCnoidal,
173 const scalar dudx = 2.0*
K/
L;
174 const scalar dudxx = 2.0*
K/
L*dudx;
175 const scalar dudxxx = 2.0*
K/
L*dudxx;
180 scalar d1 = -2.0*
H*cn*dn*sn*dudx;
181 scalar d2 = 2.0*
H*(dn*dn*sn*sn - cn*cn*dn*dn + m*cn*cn*sn*sn)*dudxx;
185 cn*sn*dn*dn*dn*(-4.0 - 2.0*m)
186 + 4.0*m*cn*sn*sn*sn*dn
187 - 2.0*m*cn*cn*cn*sn*dn
191 return vector(d1, d2, d3);
212 const scalar uCnoidal =
214 const scalar
k =
sqrt(kx*kx + ky*ky);
216 const scalar
c =
L/
T;
218 const scalar etaCN = eta(
H, m, kx, ky,
T,
x,
y, t);
219 const vector etaX = this->dEtaDx(
H, m, uCnoidal,
L,
K, E);
220 const scalar etaMS = etaMeanSq(
H, m,
T);
224 -
c*(etaCN*etaCN/
h/
h + etaMS*etaMS/
h/
h)
225 + 1.0/2.0*
c*
h*(1.0/3.0 - z*z/
h/
h)*etaX[1];
228 -
c*z*(etaX[0]/
h*(1.0 - 2.0*etaCN/
h) + 1.0/6.0*
h*(1.0 - z*z/
h/
h)*etaX[2]);
230 scalar v = u*
sin(waveAngle_);
231 u *=
cos(waveAngle_);
247 const scalar waveKx = waveK*
cos(waveAngle_);
248 const scalar waveKy = waveK*
sin(waveAngle_);
265 level[paddlei] = waterDepthRef_ + tCoeff*eta;
278 const scalar waveKx = waveK*
cos(waveAngle_);
279 const scalar waveKy = waveK*
sin(waveAngle_);
289 setPaddlePropeties(level, facei, fraction, z);
293 const label paddlei = faceToPaddle_[facei];
309 U_[facei] = fraction*
Uf*tCoeff;
362 os <<
" Cnoidal m parameter : " << m_ <<
nl;
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
const vector L(dict.get< vector >("L"))
virtual bool readDict(const dictionary &overrideDict)
Different types of constants.
dimensionedScalar sin(const dimensionedScalar &ds)
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
autoPtr< surfaceVectorField > Uf
InfoProxy< IOobject > info() const
CGAL::Exact_predicates_exact_constructions_kernel K
void JacobiSnCnDn(const scalar u, const scalar m, scalar &Sn, scalar &Cn, scalar &Dn)
Generic templated field type.
virtual bool readDict(const dictionary &overrideDict)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionedScalar h
A patch is a list of labels that address the faces in the global face list.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr scalar twoPi(2 *M_PI)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
void ellipticIntegralsKE(const scalar m, scalar &K, scalar &E)
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
Vector< scalar > vector
A scalar version of the templated Vector.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr scalar pi(M_PI)
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
cnoidal(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
const dimensionedScalar c
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar cos(const dimensionedScalar &ds)