Go to the documentation of this file.
60 rhoOutlet_(
dict.getOrDefault<scalar>(
"rhoOutlet", -VGREAT))
62 if (
dict.found(
"volumetricFlowRate"))
69 else if (
dict.found(
"massFlowRate"))
72 flowRate_ = Function1<scalar>::New(
"massFlowRate",
dict, &db());
73 rhoName_ =
dict.getOrDefault<word>(
"rho",
"rho");
78 <<
"Please supply either 'volumetricFlowRate' or"
83 if (
dict.found(
"value"))
100 const flowRateOutletVelocityFvPatchVectorField& ptf,
102 const DimensionedField<vector, volMesh>& iF,
103 const fvPatchFieldMapper& mapper
106 fixedValueFvPatchField<
vector>(ptf,
p, iF, mapper),
107 flowRate_(ptf.flowRate_.clone()),
108 volumetric_(ptf.volumetric_),
109 rhoName_(ptf.rhoName_),
110 rhoOutlet_(ptf.rhoOutlet_)
121 flowRate_(ptf.flowRate_.clone()),
122 volumetric_(ptf.volumetric_),
123 rhoName_(ptf.rhoName_),
124 rhoOutlet_(ptf.rhoOutlet_)
136 flowRate_(ptf.flowRate_.clone()),
137 volumetric_(ptf.volumetric_),
138 rhoName_(ptf.rhoName_),
139 rhoOutlet_(ptf.rhoOutlet_)
145 template<
class RhoType>
146 void Foam::flowRateOutletVelocityFvPatchVectorField::updateValues
151 const scalar t = db().time().timeOutputValue();
165 nUp =
max(nUp, scalar(0));
167 const scalar flowRate = flowRate_->value(t);
168 const scalar estimatedFlowRate =
gSum(
rho*(this->
patch().magSf()*nUp));
170 if (estimatedFlowRate > 0.5*flowRate)
172 nUp *= (
mag(flowRate)/
mag(estimatedFlowRate));
176 nUp += ((flowRate - estimatedFlowRate)/
gSum(
rho*
patch().magSf()));
194 if (volumetric_ || rhoName_ ==
"none")
201 if (db().foundObject<volScalarField>(rhoName_))
203 const fvPatchField<scalar>& rhop =
214 <<
"Did not find registered density field " << rhoName_
215 <<
" and no constant density 'rhoOutlet' specified"
219 updateValues(rhoOutlet_);
223 fixedValueFvPatchVectorField::updateCoeffs();
230 flowRate_->writeData(
os);
236 writeEntry(
"value",
os);
247 flowRateOutletVelocityFvPatchVectorField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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.
A class for handling words, derived from Foam::string.
Velocity outlet boundary condition which corrects the extrapolated velocity to match the specified fl...
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Type gSum(const FieldField< Field, Type > &f)
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Field< vector > vectorField
Specialisation of Field<T> for vector.
Generic templated field type.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
flowRateOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
label max(const labelHashSet &set, label maxValue=labelMin)
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.
Vector< scalar > vector
A scalar version of the templated Vector.
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
fvPatchField< vector > fvPatchVectorField
#define FatalErrorInFunction
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
#define FatalIOErrorInFunction(ios)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual void updateCoeffs()
virtual void write(Ostream &) const
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)