Go to the documentation of this file.
41 directionMixedFvPatchVectorField(
p, iF),
44 refValue() = vector::zero;
45 refGrad() = vector::zero;
46 valueFraction() = symmTensor::zero;
59 directionMixedFvPatchVectorField(ptf,
p, iF, mapper),
77 directionMixedFvPatchVectorField(
p, iF),
78 phiName_(
dict.lookupOrDefault<
word>(
"phi",
"phi"))
82 if (
dict.found(
"tangentialVelocity"))
91 refValue() = vector::zero;
94 refGrad() = vector::zero;
95 valueFraction() = symmTensor::zero;
105 directionMixedFvPatchVectorField(pivpvf),
118 directionMixedFvPatchVectorField(pivpvf, iF),
140 directionMixedFvPatchVectorField::autoMap(m);
141 if (tangentialVelocity_.size())
143 tangentialVelocity_.autoMap(m);
154 directionMixedFvPatchVectorField::rmap(ptf, addr);
156 if (tangentialVelocity_.size())
159 refCast<const pressureInletOutletVelocityFvPatchVectorField>(ptf);
176 valueFraction() =
neg(phip)*(
I -
sqr(patch().nf()));
178 directionMixedFvPatchVectorField::updateCoeffs();
179 directionMixedFvPatchVectorField::evaluate();
190 writeEntryIfDifferent<word>(os,
"phi",
"phi", phiName_);
191 if (tangentialVelocity_.size())
193 tangentialVelocity_.writeEntry(
"tangentialVelocity", os);
195 writeEntry(
"value", os);
201 void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
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 operator=(const UList< Type > &)
A class for handling words, derived from string.
A class for managing temporary objects.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual void write(Ostream &) const
Write.
const vectorField & tangentialVelocity() const
Return the tangential velocity.
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionSet transform(const dimensionSet &)
Pre-declare SubField and related Field type.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
volVectorField vectorField(fieldObject, mesh)
static const sphericalTensor I(1)
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 rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
word phiName_
Flux field name.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
void setTangentialVelocity(const vectorField &tangentialVelocity)
Reset the tangential velocity.
vectorField tangentialVelocity_
Optional tangential velocity component.
dimensionedScalar neg(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...