34 bool Foam::lumpedMassTemperatureFvPatchScalarField::closed()
const
36 return mag(
gSum(patch().Sf()))/
gSum(patch().magSf()) < rootSmall;
39 Foam::scalar Foam::lumpedMassTemperatureFvPatchScalarField::V()
const
41 return -
gSum(patch().Sf() & patch().Cf())
42 /patch().boundaryMesh().mesh().nSolutionD();
55 fixedValueFvPatchScalarField(
p, iF),
56 rho_(
dict.lookup<scalar>(
"rho")),
57 Cv_(
dict.lookup<scalar>(
"Cv")),
62 "T_" + patch().
name(),
74 db().time().userUnits(),
82 if (!
dict.readIfPresent(
"volume", V_))
87 Info<<
"Volume for the thermal mass, enclosed by patch '"
88 << patch().name() <<
"', = " << V_;
93 <<
"Patch '" << patch().name()
94 <<
"' corresponding to a thermal mass is not closed." <<
nl
95 <<
"Please specify the volume with the optional 'volume' entry."
113 fixedValueFvPatchScalarField(ptf,
p, iF, mapper),
129 fixedValueFvPatchScalarField(ptf, iF),
146 fixedValueFvPatchScalarField::map(ptf, mapper);
155 fixedValueFvPatchScalarField::reset(ptf);
172 internalField().group()
180 const scalar Hs = rho_*Cv_*V_/db().time().deltaTValue();
182 const scalar t = db().time().userTimeValue();
186 Q_->value(t) +
gSum(Hf*patchInternalField())
187 + Hs*T_.oldTime().value()
192 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 keyword definitions, which are a keyword followed by any number of values (e....
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)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimPower
const dimensionSet dimTemperature
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)