Go to the documentation of this file.
70 fixedValueFvPatchScalarField(
p, iF),
85 fixedValueFvPatchScalarField(
p, iF),
86 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
87 zetaName_(
dict.lookupOrDefault<
word>(
"zeta",
"zeta")),
88 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho"))
106 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
119 fixedValueFvPatchScalarField(wspsf),
133 fixedValueFvPatchScalarField(wspsf, iF),
149 const label patchI = patch().index();
151 const scalar dt = db().time().deltaTValue();
162 const word ddtSchemeName(zeta.mesh().ddtScheme(zeta.name()));
173 vectorField dZetap(dt*nf()*
phi.boundaryField()[patchI]/patch().magSf());
188 zetap = zeta.
oldTime().boundaryField()[patchI] + dZetap;
194 scalar dt0 = db().time().deltaT0Value();
196 scalar
c = 1.0 + dt/(dt + dt0);
197 scalar c00 = dt*dt/(dt0*(dt + dt0));
202 c0*zeta.
oldTime().boundaryField()[patchI]
203 - c00*zeta.
oldTime().oldTime().boundaryField()[patchI]
212 << ddtSchemeName <<
nl
213 <<
" on patch " << this->patch().name()
221 Info<<
"min/max zetap = " <<
gMin(zetap & nf()) <<
", "
230 fixedValueFvPatchScalarField::updateCoeffs();
237 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
238 writeEntryIfDifferent<word>(os,
"zeta",
"zeta", zetaName_);
239 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
240 writeEntry(
"value", os);
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
A class for handling words, derived from string.
A class for managing temporary objects.
static const NamedEnum< ddtSchemeType, 3 > ddtSchemeTypeNames_
Time scheme type names.
const dimensionSet dimVelocity
const dimensionSet dimDensity
rDeltaT dimensionedInternalField()
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionedVector & g
Ostream & endl(Ostream &os)
Add newline and flush stream.
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
word zetaName_
Wave height field name.
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word phiName_
Flux field name.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
virtual void write(Ostream &) const
Write.
errorManip< error > abort(error &err)
ddtSchemeType
Enumeration defining the available ddt schemes.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Second-order backward-differencing ddt using the current and two previous time-step values.
word rhoName_
Density field for mass-based flux evaluations.
volScalarField scalarField(fieldObject, mesh)
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Type gMin(const FieldField< Field, Type > &f)
This is a pressure boundary condition, whose value is calculated as the hydrostatic pressure based on...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Generic GeometricField class.
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Type gMax(const FieldField< Field, Type > &f)
Initialise the NamedEnum HashTable from the static list of names.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...