35 void Foam::flowRateInletVelocityFvPatchVectorField::setWallDist()
41 patch().patch().boundaryMesh().findIndices<wallPolyPatch>()
44 const patchPatchDist pwd(patch().patch(), otherPatchIDs);
49 area_ =
gSum(patch().magSf());
54 Foam::flowRateInletVelocityFvPatchVectorField::profile()
58 return profile_->value(y_);
62 return tmp<scalarField>(
new scalarField(size(), scalar(1)));
67 template<
class ScaleType,
class AlphaType,
class RhoType>
68 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
70 const ScaleType& scale,
71 const AlphaType&
alpha,
79 *flowRate_->value(db().time().value())
86 template<
class AlphaType>
87 void Foam::flowRateInletVelocityFvPatchVectorField::updateValues
89 const AlphaType&
alpha
94 updateValues(area_,
alpha, one());
96 else if (volumetric_ || rhoName_ ==
"none")
98 updateValues(one(),
alpha, one());
103 if (db().foundObject<volScalarField>(rhoName_))
105 const fvPatchField<scalar>& rhop =
108 updateValues(one(),
alpha, rhop);
116 <<
"Did not find registered density field " << rhoName_
117 <<
" and no constant density 'rhoInlet' specified"
121 updateValues(one(),
alpha, rhoInlet_);
127 bool Foam::flowRateInletVelocityFvPatchVectorField::canEvaluate()
131 || !patch().boundaryMesh().mesh().time().processorCase();
151 rhoInlet_(
dict.lookupOrDefault<scalar>(
"rhoInlet",
dimDensity, -vGreat)),
152 alphaName_(
dict.lookupOrDefault<
word>(
"alpha",
word::null)),
156 if (
dict.found(
"meanVelocity"))
158 meanVelocity_ =
true;
164 db().time().userUnits(),
169 else if (
dict.found(
"volumetricFlowRate"))
171 meanVelocity_ =
false;
176 "volumetricFlowRate",
177 db().time().userUnits(),
182 else if (
dict.found(
"massFlowRate"))
184 meanVelocity_ =
false;
190 db().time().userUnits(),
194 rhoName_ =
word(
dict.lookupOrDefault<
word>(
"rho",
"rho"));
199 <<
"Please supply 'meanVelocity', 'volumetricFlowRate' or"
203 if (
dict.found(
"profile"))
213 if (!canEvaluate() ||
dict.found(
"value"))
237 flowRate_(ptf.flowRate_, false),
238 meanVelocity_(ptf.meanVelocity_),
239 volumetric_(ptf.volumetric_),
240 profile_(ptf.profile_, false),
241 rhoName_(ptf.rhoName_),
242 rhoInlet_(ptf.rhoInlet_),
243 alphaName_(ptf.alphaName_),
246 profile_.
valid() && canEvaluate()
262 flowRate_(ptf.flowRate_, false),
263 meanVelocity_(ptf.meanVelocity_),
264 volumetric_(ptf.volumetric_),
265 profile_(ptf.profile_, false),
266 rhoName_(ptf.rhoName_),
267 rhoInlet_(ptf.rhoInlet_),
268 alphaName_(ptf.alphaName_),
282 fixedValueFvPatchVectorField::map(ptf, mapper);
285 refCast<const flowRateInletVelocityFvPatchVectorField>(ptf);
287 if (profile_.valid() && canEvaluate())
289 mapper(y_, tiptf.y_);
299 fixedValueFvPatchVectorField::reset(ptf);
302 refCast<const flowRateInletVelocityFvPatchVectorField>(ptf);
304 if (profile_.valid() && canEvaluate())
321 <<
"Cannot evaluate flow rate on a non-parallel processor case"
330 updateValues(alphap);
337 fixedValueFvPatchVectorField::updateCoeffs();
345 if (profile_.valid())
351 writeEntryIfDifferent<word>(os,
"rho",
"rho", rhoName_);
352 writeEntryIfDifferent<scalar>(os,
"rhoInlet", -vGreat, rhoInlet_);
354 writeEntryIfDifferent<word>(os,
"alpha",
word::null, alphaName_);
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static bool & parRun()
Is this a parallel run?
A list of keyword definitions, which are a keyword followed by any number of values (e....
Abstract base class for field mapping.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Velocity inlet boundary condition creating a velocity field with optionally specified profile normal ...
flowRateInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchVectorField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
A class for managing temporary objects.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
bool valid(const PtrList< ModelType > &l)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gSum(const FieldField< Field, Type > &f)
const unitConversion unitAny
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
const dimensionSet dimDensity
VolField< scalar > volScalarField
Field< vector > vectorField
Specialisation of Field<T> for vector.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
const dimensionSet dimVelocity
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Type gMax(const FieldField< Field, Type > &f)