Go to the documentation of this file.
42 mixedFvPatchVectorField(
p, iF),
47 refGrad() = vector::zero;
48 valueFraction() = 0.0;
61 mixedFvPatchVectorField(ptf,
p, iF, mapper),
75 mixedFvPatchVectorField(
p, iF),
76 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
77 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho"))
81 refGrad() = vector::zero;
82 valueFraction() = 0.0;
92 mixedFvPatchVectorField(pivpvf),
105 mixedFvPatchVectorField(pivpvf, iF),
131 refValue() =
n*phip/magS;
138 refValue() =
n*phip/(rhop*magS);
143 <<
"dimensions of phi are not correct"
144 <<
"\n on patch " << this->patch().name()
150 valueFraction() = 1.0 -
pos(phip);
152 mixedFvPatchVectorField::updateCoeffs();
164 writeEntry(
"value", os);
170 void Foam::pressureNormalInletOutletVelocityFvPatchVectorField::operator=
177 valueFraction()*(patch().nf()*(patch().nf() & pvf))
178 + (1 - valueFraction())*pvf
virtual void write(Ostream &) const
Write.
This velocity inlet/outlet boundary condition is applied to patches where the pressure is specified....
A class for handling words, derived from string.
pressureNormalInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
rDeltaT dimensionedInternalField()
word rhoName_
Density field name.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const dimensionSet dimArea(sqr(dimLength))
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volVectorField vectorField(fieldObject, mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
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.
word phiName_
Flux field name.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void write(Ostream &) const
Write.
dimensionedScalar pos(const dimensionedScalar &ds)