75 const scalar
xMin = bb.min().x();
76 const scalar
xMax = bb.max().x();
77 const scalar yMin = bb.min().y();
78 const scalar yMax = bb.max().y();
79 zSpan_ = bb.max().z() - bb.min().z();
82 xPaddle_.setSize(nPaddle_, 0);
83 yPaddle_.setSize(nPaddle_, 0);
85 const scalar paddleDy = (yMax - yMin)/scalar(nPaddle_);
86 for (label paddlei = 0; paddlei < nPaddle_; ++paddlei)
88 xPaddle_[paddlei] = xMid;
89 yPaddle_[paddlei] = paddlei*paddleDy + yMin + 0.5*paddleDy;
93 const vectorField CfLocal(Rgl_ & patch_.faceCentres());
94 z_ = CfLocal.component(2);
97 zMin_.setSize(patch_.size());
98 zMax_.setSize(patch_.size());
99 const faceList& faces = patch_.localFaces();
102 const face&
f = faces[facei];
103 const label nPoint =
f.size();
104 zMin_[facei] = CpLocal[
f[0]].z();
105 zMax_[facei] = CpLocal[
f[0]].z();
107 for (label fpi = 1; fpi < nPoint; ++fpi)
109 const label pointi =
f[fpi];
110 zMin_[facei] =
min(zMin_[facei], CpLocal[pointi].z());
111 zMax_[facei] =
max(zMax_[facei], CpLocal[pointi].z());
116 zMin0_ =
gMin(zMin_);
119 faceToPaddle_.setSize(patch_.size(), -1);
122 faceToPaddle_[facei] = floor((CfLocal[facei].
y() - yMin)/paddleDy);
131 auto& level = tlevel.ref();
136 const scalarField alphac(alphap.patchInternalField());
144 label paddlei = faceToPaddle_[facei];
145 paddleMagSf[paddlei] += magSf[facei];
146 paddleWettedMagSf[paddlei] += magSf[facei]*alphac[facei];
149 forAll(paddleMagSf, paddlei)
154 paddleWettedMagSf[paddlei]*zSpan_
155 /(paddleMagSf[paddlei] + ROOTVSMALL);
166 const label paddlei = faceToPaddle_[facei];
167 const scalar paddleCalc = level[paddlei];
169 const scalar zMin0 = zMin_[facei] - zMin0_;
170 const scalar zMax0 = zMax_[facei] - zMin0_;
172 if (zMax0 < paddleCalc)
176 else if (zMin0 > paddleCalc)
182 scalar dz = paddleCalc - zMin0;
183 alpha_[facei] = dz/(zMax0 - zMin0);
197 const label paddlei = faceToPaddle_[facei];
198 const scalar paddleCalc = level[paddlei];
199 const scalar paddleHeight =
min(paddleCalc, waterDepthRef_);
200 const scalar zMin = zMin_[facei] - zMin0_;
201 const scalar zMax = zMax_[facei] - zMin0_;
206 if (zMax < paddleHeight)
208 z = z_[facei] - zMin0_;
210 else if (zMin > paddleCalc)
216 if (paddleCalc < waterDepthRef_)
218 if ((zMax > paddleCalc) && (zMin < paddleCalc))
220 scalar dz = paddleCalc - zMin;
221 fraction = dz/(zMax - zMin);
222 z = z_[facei] - zMin0_;
227 if (zMax < paddleCalc)
231 else if ((zMax > paddleCalc) && (zMin < paddleCalc))
233 scalar dz = paddleCalc - zMin;
234 fraction = dz/(zMax - zMin);
281 activeAbsorption_(false),
299 IOdictionary::regIOobject::read();
307 readEntry(
"nPaddle", nPaddle_);
311 <<
"Number of paddles must be greater than zero. Supplied"
312 <<
" value nPaddles = " << nPaddle_
319 initialiseGeometry();
324 scalar waterDepth = 0;
327 waterDepthRef_ = waterDepth;
334 waterDepthRef_ = level.
first();
339 waterDepthRef_ += SMALL;
351 if (mesh_.time().timeIndex() != currTimeIndex_)
353 Info<<
"Updating " <<
type() <<
" wave model for patch "
354 << patch_.name() <<
endl;
357 const scalar tCoeff = timeCoeff(t);
369 setLevel(t, tCoeff, calculatedLevel);
372 setVelocity(t, tCoeff, calculatedLevel);
375 setAlpha(calculatedLevel);
378 if (activeAbsorption_)
384 const label paddlei = faceToPaddle_[facei];
386 if (zMin_[facei] - zMin0_ < activeLevel[paddlei])
389 (calculatedLevel[paddlei] - activeLevel[paddlei])
390 *
sqrt(
mag(g_)/activeLevel[paddlei]);
392 U_[facei].x() += UCorr;
423 os <<
"Wave model: patch " << patch_.name() <<
nl
424 <<
" Type : " <<
type() <<
nl
425 <<
" Velocity field name : " << UName_ <<
nl
426 <<
" Phase fraction field name : " << alphaName_ <<
nl
427 <<
" Transformation from local to global system : " << Rlg_ <<
nl
428 <<
" Number of paddles: " << nPaddle_ <<
nl
429 <<
" Reference water depth : " << waterDepthRef_ <<
nl
430 <<
" Active absorption: " << activeAbsorption_ <<
nl;
propsDict readIfPresent("fields", acceptFields)
label size() const noexcept
The number of elements in list.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
readOption readOpt() const noexcept
Get the read option.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Time & time() const noexcept
Return Time associated with the objectRegistry.
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
T & first()
Access first element of the list, position [0].
void size(const label n)
Older name for setAddressableSize.
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
const word & name() const
Name function is needed to disambiguate those inherited from regIOobject and dictionary.
A bounding box defined in terms of min/max extrema points.
const point & max() const noexcept
Maximum describing the bounding box.
const point & min() const noexcept
Minimum describing the bounding box.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
bool merge(const dictionary &dict)
Merge entries from the given dictionary.
A face is a list of labels corresponding to mesh vertices.
Mesh data needed to do the Finite Volume discretisation.
A patch is a list of labels that address the faces in the global face list.
bool headerOk()
Read and check header info. Does not check the headerClassName.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Base class for waveModels.
virtual void initialiseGeometry()
Initialise.
scalarField zMin_
Minimum z (point) height per patch face / [m].
scalarField z_
Patch face centre z coordinates / [m].
virtual scalar timeCoeff(const scalar t) const =0
Return the time scaling coefficient.
const polyPatch & patch_
Reference to the polyPatch.
tensor Rlg_
Rotation tensor from local to global system.
const fvMesh & mesh_
Reference to the mesh database.
scalar zSpan_
Overall (point) span in z-direction / [m].
word UName_
Name of velocity field, default = "U".
const vector & g_
Gravity.
virtual const scalarField & alpha() const
Return the latest wave indicator field prediction.
virtual void correct(const scalar t)
Correct the model for time, t[s].
scalar zMin0_
Minimum z reference level / [m].
scalarField alpha_
Wave indicator field.
bool activeAbsorption_
Active wave absorption switch.
virtual const vectorField & U() const
Return the latest wave velocity prediction.
tensor Rgl_
Rotation tensor from global to local system.
vectorField U_
Velocity field.
scalar waterDepthRef_
Reference water depth / [m].
virtual void setPaddlePropeties(const scalarField &level, const label facei, scalar &fraction, scalar &z) const
Set the paddle coverage fraction and reference height.
scalarField yPaddle_
Paddle y coordinates / [m].
static const word dictName
Dictionary name.
scalarField xPaddle_
Paddle x coordinates / [m].
virtual void setVelocity(const scalar t, const scalar tCoeff, const scalarField &level)=0
Calculate the wave model velocity.
waveModel(const dictionary &dict, const fvMesh &mesh, const polyPatch &patch, const bool readFields=true)
Constructor.
static autoPtr< waveModel > New(const word &dictName, const fvMesh &mesh, const polyPatch &patch)
Return a reference to the selected wave model.
scalar initialDepth_
Initial depth / [m].
scalarField zMax_
Maximum z (point) height per patch face / [m].
virtual void setLevel(const scalar t, const scalar tCoeff, scalarField &level) const =0
Set the water level.
static word modelName(const word &patchName)
Utility function to construct the model name.
label currTimeIndex_
Time index used for updating.
labelList faceToPaddle_
Addressing from patch face index to paddle index.
virtual tmp< scalarField > waterLevel() const
Water level.
word alphaName_
Name of phase fraction field, default = "alpha".
virtual void setAlpha(const scalarField &level)
Set the alpha field based on the water level.
label nPaddle_
Number of paddles.
virtual bool readDict(const dictionary &overrideDict)
Read from dictionary.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & Cp
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
OBJstream os(runTime.globalPath()/outputName)
Different types of constants.
const std::string patch
OpenFOAM patch number as a std::string.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
static const Identity< scalar > I
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< scalar > scalarList
List of scalar.
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.