34 bool Foam::lumpedMassTemperatureFvPatchScalarField::closed()
const
36 return mag(
gSum(patch().Sf()))/
gSum(patch().magSf()) < rootSmall;
40 Foam::scalar Foam::lumpedMassTemperatureFvPatchScalarField::V()
const
43 -
gSum(patch().Sf() & patch().Cf())
44 /patch().boundaryMesh().mesh().nSolutionD();
58 fixedValueFvPatchScalarField(
p, iF),
59 rho_(
dict.lookup<scalar>(
"rho")),
60 Cv_(
dict.lookup<scalar>(
"Cv")),
65 "T_" + patch().
name(),
77 db().time().userUnits(),
85 if (!
dict.readIfPresent(
"volume", V_))
90 Info<<
"Volume for the thermal mass, enclosed by patch '"
91 << patch().name() <<
"', = " << V_;
96 <<
"Patch '" << patch().name()
97 <<
"' corresponding to a thermal mass is not closed." <<
nl
98 <<
"Please specify the volume with the optional 'volume' entry."
116 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
132 fixedValueFvPatchScalarField(ptf, iF),
149 fixedValueFvPatchScalarField::map(ptf, mapper);
158 fixedValueFvPatchScalarField::reset(ptf);
175 internalField().group()
180 ttm.
kappaEff(patch().index())*patch().magSf()*patch().deltaCoeffs()
183 const scalar Hs = rho_*Cv_*V_/db().time().deltaTValue();
187 Q_->value(db().time().value())
188 +
gSum(Hf*patchInternalField())
189 + Hs*T_.oldTime().value()
194 fixedValueFvPatchScalarField::updateCoeffs();
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...
Run-time selectable general function of one variable.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Abstract base class for field mapping.
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 finiteVolume patch using a polyPatch and a fvBoundaryMesh.
This boundary condition is applied to a patch which bounds a solid body, wholly or partially....
lumpedMassTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void write(Ostream &) const
Write.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Abstract base class for all fluid and solid thermophysical transport models.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimPower
const dimensionSet dimTemperature
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)