Go to the documentation of this file.
44 mixedFvPatchVectorField(
p, iF),
51 valueFraction() = 0.0;
64 mixedFvPatchVectorField(ptf,
p, iF, mapper),
65 phiName_(ptf.phiName_),
66 rhoName_(ptf.rhoName_),
67 inletDir_(ptf.inletDir_, mapper)
79 mixedFvPatchVectorField(
p, iF),
80 phiName_(
dict.getOrDefault<
word>(
"phi",
"phi")),
81 rhoName_(
dict.getOrDefault<
word>(
"rho",
"rho")),
82 inletDir_(
"inletDirection",
dict,
p.size())
88 valueFraction() = 0.0;
98 mixedFvPatchVectorField(pivpvf),
99 phiName_(pivpvf.phiName_),
100 rhoName_(pivpvf.rhoName_),
101 inletDir_(pivpvf.inletDir_)
112 mixedFvPatchVectorField(pivpvf, iF),
113 phiName_(pivpvf.phiName_),
114 rhoName_(pivpvf.rhoName_),
115 inletDir_(pivpvf.inletDir_)
126 mixedFvPatchVectorField::autoMap(m);
127 inletDir_.autoMap(m);
137 mixedFvPatchVectorField::rmap(ptf, addr);
140 refCast<const pressureDirectedInletOutletVelocityFvPatchVectorField>
143 inletDir_.rmap(tiptf.inletDir_, addr);
157 const fvsPatchField<scalar>& phip =
160 tmp<vectorField>
n =
patch().nf();
161 tmp<scalarField> ndmagS = (
n & inletDir_)*
patch().magSf();
165 refValue() = inletDir_*phip/ndmagS;
169 const fvPatchField<scalar>& rhop =
172 refValue() = inletDir_*phip/(rhop*ndmagS);
177 <<
"dimensions of phi are not correct"
178 <<
"\n on patch " << this->
patch().name()
179 <<
" of field " << this->internalField().name()
180 <<
" in file " << this->internalField().objectPath()
184 valueFraction() = 1.0 -
pos0(phip);
186 mixedFvPatchVectorField::updateCoeffs();
198 inletDir_.writeEntry(
"inletDirection",
os);
199 writeEntry(
"value",
os);
205 void Foam::pressureDirectedInletOutletVelocityFvPatchVectorField::operator=
212 valueFraction()*(inletDir_*(inletDir_ & pvf))
213 + (1 - valueFraction())*pvf
225 pressureDirectedInletOutletVelocityFvPatchVectorField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
virtual void autoMap(const fvPatchFieldMapper &)
virtual void operator=(const UList< Type > &)
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
static constexpr const zero Zero
const dimensionSet dimVelocity
const dimensionSet dimDensity
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 &)
dimensionedScalar pos0(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimArea(sqr(dimLength))
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual void updateCoeffs()
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
This velocity inlet/outlet boundary condition is applied to velocity boundaries where the pressure is...
errorManipArg< error, int > exit(error &err, const int errNo=1)
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
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...
pressureDirectedInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Generic GeometricField class.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...