Go to the documentation of this file.
41 mixedFvPatchVectorField(
p, iF),
47 refGrad() = vector::zero;
48 valueFraction() = 0.0;
61 mixedFvPatchVectorField(ptf,
p, iF, mapper),
76 mixedFvPatchVectorField(
p, iF),
77 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi")),
78 rhoName_(
dict.lookupOrDefault<
word>(
"rho",
"rho")),
79 inletDir_(
"inletDirection",
dict,
p.size())
83 refGrad() = vector::zero;
84 valueFraction() = 0.0;
94 mixedFvPatchVectorField(pivpvf),
108 mixedFvPatchVectorField(pivpvf, iF),
122 mixedFvPatchVectorField::autoMap(m);
123 inletDir_.autoMap(m);
133 mixedFvPatchVectorField::rmap(ptf, addr);
136 refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
168 refValue() =
inletDir_*phip/(rhop*ndmagS);
173 <<
"dimensions of phi are not correct"
174 <<
"\n on patch " << this->patch().name()
180 valueFraction() = 1.0 -
pos(phip);
182 mixedFvPatchVectorField::updateCoeffs();
192 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
193 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
194 inletDir_.writeEntry(
"inletDirection", os);
195 writeEntry(
"value", os);
201 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
208 valueFraction()*(inletDir_*(inletDir_ & pvf))
209 + (1 - valueFraction())*pvf
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A class for handling words, derived from string.
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionSet dimDensity
rDeltaT dimensionedInternalField()
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const dimensionSet dimArea(sqr(dimLength))
word rhoName_
Density field name.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volVectorField vectorField(fieldObject, mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
vectorField inletDir_
Inlet direction.
word phiName_
Flux field name.
pressureDirectedInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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...
dimensionedScalar pos(const dimensionedScalar &ds)