Go to the documentation of this file.
68 scalar Xa = -
C*t + ts -
X0 +
x*
cos(theta) +
y*
sin(theta);
84 scalar m = newtonRapsonF1(xin,
H,
h);
87 scalar
c1 =
sin(m + (1.0 + (2.0*
H/(3.0*
h))));
103 scalar maxval = 10000.0;
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);
123 else if ((iter > 1) && (residual > maxval))
126 <<
"Newton-Raphson iterations diverging: "
127 <<
"iterations = " << iter
128 <<
", residual = " << residual
134 scalar c3 = 1.0/
sin(
b);
143 <<
"Failed to converge in " << iter <<
" iterations. Residual = "
162 scalar maxval = 10000;
170 scalar a = m*(1.0 +
x/
h);
182 else if ((iter > 1) && (residual > maxval))
185 <<
"Newton-Raphson iterations diverging: "
186 <<
"iterations = " << iter
187 <<
", residual = " << residual
193 scalar fprime = 1 -
n/c3*(
c1 -
sqr(
c2)/c3);
200 <<
"Failed to converge in " << iter <<
" iterations. Residual = "
220 const scalar mm = vec[0];
221 const scalar nn = vec[1];
224 const scalar ts = 3.5*
h/
sqrt(
H/
h);
225 const scalar Xa = -
C*t + ts -
X0 +
x*
cos(theta) +
y*
sin(theta);
227 scalar outa =
C*nn*(1.0 +
cos(mm*z/
h)*
cosh(mm*Xa/
h));
230 scalar u = outa/outb;
234 const scalar w = outa/outb;
236 const scalar v = u*
sin(waveAngle_);
264 level[paddlei] = waterDepthRef_ + tCoeff*eta;
316 setPaddlePropeties(level, facei, fraction, z);
320 const label paddlei = faceToPaddle_[facei];
334 U_[facei] = fraction*
Uf*tCoeff;
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
dimensionedScalar tan(const dimensionedScalar &ds)
virtual bool readDict(const dictionary &overrideDict)
dimensionedScalar cosh(const dimensionedScalar &ds)
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const
scalarList X0(nSpecie, Zero)
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
addToRunTimeSelectionTable(waveModel, shallowWaterAbsorption, patch)
Ostream & endl(Ostream &os)
virtual bool readDict(const dictionary &overrideDict)
InfoProxy< IOobject > info() const
virtual scalar newtonRapsonF1(const scalar x0, const scalar H, const scalar h) const
const dimensionedScalar c1
const dimensionedScalar b
Generic templated field type.
virtual vector mn(const scalar H, const scalar h) const
const dimensionedScalar h
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
A patch is a list of labels that address the faces in the global face list.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
virtual scalar newtonRapsonF2(const scalar x0, const scalar H, const scalar h, const scalar xa, const scalar m, const scalar n) const
Mesh data needed to do the Finite Volume discretisation.
const dimensionedScalar c2
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
McCowan(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
#define FatalErrorInFunction
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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,...
const Vector< label > N(dict.get< Vector< label >>("N"))
Graphite solid properties.
defineTypeNameAndDebug(waveAbsorptionModel, 0)
#define WarningInFunction
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)