Go to the documentation of this file.
45 Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeTypeNames_
48 ddtSchemeType::tsEuler,
49 fv::EulerDdtScheme<scalar>::typeName_()
52 ddtSchemeType::tsCrankNicolson,
53 fv::CrankNicolsonDdtScheme<scalar>::typeName_()
56 ddtSchemeType::tsBackward,
57 fv::backwardDdtScheme<scalar>::typeName_()
71 fixedValueFvPatchScalarField(
p, iF),
86 fixedValueFvPatchScalarField(
p, iF,
dict),
87 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
88 zetaName_(
dict.getOrDefault<
word>(
"zeta",
"zeta")),
89 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho"))
102 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
103 phiName_(ptf.phiName_),
104 zetaName_(ptf.zetaName_),
105 rhoName_(ptf.rhoName_)
115 fixedValueFvPatchScalarField(wspsf),
116 phiName_(wspsf.phiName_),
117 zetaName_(wspsf.zetaName_),
118 rhoName_(wspsf.rhoName_)
129 fixedValueFvPatchScalarField(wspsf, iF),
130 phiName_(wspsf.phiName_),
131 zetaName_(wspsf.zetaName_),
132 rhoName_(wspsf.rhoName_)
145 const label patchi =
patch().index();
147 const scalar dt = db().time().deltaTValue();
151 vectorField& zetap = zeta.boundaryFieldRef()[patchi];
154 const word ddtSchemeName(zeta.mesh().ddtScheme(zeta.name()));
162 tmp<vectorField> nf(
patch().nf());
182 zetap = zeta0.boundaryField()[patchi] + dZetap;
188 scalar dt0 = db().time().deltaT0Value();
190 scalar
c = 1.0 + dt/(dt + dt0);
191 scalar c00 = dt*dt/(dt0*(dt + dt0));
196 c0*zeta0.boundaryField()[patchi]
197 - c00*zeta0.oldTime().boundaryField()[patchi]
206 << ddtSchemeName <<
nl
207 <<
" on patch " << this->
patch().name()
208 <<
" of field " << this->internalField().name()
209 <<
" in file " << this->internalField().objectPath()
215 Info<<
"min/max zetap = " <<
gMin(zetap & nf()) <<
", "
224 fixedValueFvPatchScalarField::updateCoeffs();
234 writeEntry(
"value",
os);
245 waveSurfacePressureFvPatchScalarField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
fvPatchField< scalar > fvPatchScalarField
virtual void write(Ostream &) const
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Ostream & endl(Ostream &os)
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimArea(sqr(dimLength))
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
UniformDimensionedField< vector > uniformDimensionedVectorField
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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 void write(Ostream &) const
const uniformDimensionedVectorField & g
errorManip< error > abort(error &err)
GeometricField< vector, fvPatchField, volMesh > volVectorField
#define FatalErrorInFunction
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
const dimensionedScalar c
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
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()
Generic GeometricField class.
Type gMax(const FieldField< Field, Type > &f)
static const gravity & New(const Time &runTime)
const Boundary & boundaryField() const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...