52 const volScalarField::GeometricBoundaryField& tbf =
53 this->T_.boundaryField();
59 if (isA<fixedJumpFvPatchScalarField>(tbf[
patchi]))
61 const fixedJumpFvPatchScalarField& pf =
62 dynamic_cast<const fixedJumpFvPatchScalarField&
>(tbf[
patchi]);
64 hbt[
patchi] = pf.interfaceFieldType();
66 else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
68 const fixedJumpAMIFvPatchScalarField& pf =
69 dynamic_cast<const fixedJumpAMIFvPatchScalarField&
> 74 hbt[
patchi] = pf.interfaceFieldType();
84 const volScalarField::GeometricBoundaryField& tbf =
85 this->T_.boundaryField();
91 if (isA<fixedValueFvPatchScalarField>(tbf[
patchi]))
93 hbt[
patchi] = fixedEnergyFvPatchScalarField::typeName;
97 isA<zeroGradientFvPatchScalarField>(tbf[patchi])
98 || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
101 hbt[
patchi] = gradientEnergyFvPatchScalarField::typeName;
103 else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
105 hbt[
patchi] = mixedEnergyFvPatchScalarField::typeName;
107 else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
111 else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
115 else if (tbf[patchi].
type() ==
"energyRegionCoupledFvPatchScalarField")
117 hbt[
patchi] =
"energyRegionCoupledFvPatchScalarField";
165 const word& phaseName
172 phasePropertyName(
dictName, phaseName),
180 phaseName_(phaseName),
182 p_(lookupOrConstruct(mesh,
"p")),
188 phasePropertyName(
"T"),
201 phasePropertyName(
"thermo:alpha"),
211 dpdt_(lookupOrDefault<Switch>(
"dpdt",
true))
219 const word& phaseName
226 phasePropertyName(
dictName, phaseName),
235 phaseName_(phaseName),
237 p_(lookupOrConstruct(mesh,
"p")),
243 phasePropertyName(
"T"),
256 phasePropertyName(
"thermo:alpha"),
273 const word& phaseName
276 return New<basicThermo>(
mesh, phaseName);
305 iter != thermos.
end();
330 if (!(
he().
name() == phasePropertyName(a)))
333 <<
"Supported energy type is " << phasePropertyName(a)
334 <<
", thermodynamics package provides " <<
he().
name()
349 he().
name() == phasePropertyName(a)
350 ||
he().
name() == phasePropertyName(b)
355 <<
"Supported energy types are " << phasePropertyName(a)
356 <<
" and " << phasePropertyName(b)
357 <<
", thermodynamics package provides " <<
he().
name()
373 he().
name() == phasePropertyName(a)
374 ||
he().
name() == phasePropertyName(b)
375 ||
he().
name() == phasePropertyName(c)
380 <<
"Supported energy types are " << phasePropertyName(a)
381 <<
", " << phasePropertyName(b)
382 <<
" and " << phasePropertyName(c)
383 <<
", thermodynamics package provides " <<
he().
name()
400 he().
name() == phasePropertyName(a)
401 ||
he().
name() == phasePropertyName(b)
402 ||
he().
name() == phasePropertyName(c)
403 ||
he().
name() == phasePropertyName(d)
408 <<
"Supported energy types are " << phasePropertyName(a)
409 <<
", " << phasePropertyName(b)
410 <<
", " << phasePropertyName(c)
411 <<
" and " << phasePropertyName(d)
412 <<
", thermodynamics package provides " <<
he().
name()
420 const word& thermoName,
431 (endb = thermoName.find(
'<', beg)) != string::npos
432 || (endc = thermoName.find(
',', beg)) != string::npos
435 if (endb == string::npos)
439 else if ((endc = thermoName.find(
',', beg)) != string::npos)
441 end =
min(endb, endc);
450 cmpts[i] = thermoName.substr(beg, end-beg);
451 cmpts[i++].replaceAll(
">",
"");
456 if (beg < thermoName.size())
458 cmpts[i] = thermoName.substr(beg, string::npos);
459 cmpts[i++].replaceAll(
">",
"");
static autoPtr< Thermo > New(const fvMesh &, const word &phaseName=word::null)
Generic New for each of the related thermodynamics packages.
static const char *const typeName
Mesh data needed to do the Finite Volume discretisation.
virtual const volScalarField & T() const
Temperature [K].
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
word name(const complex &)
Return a string representation of a complex.
volScalarField & lookupOrConstruct(const fvMesh &mesh, const char *name) const
PtrList< solidThermo > thermos(solidRegions.size())
bool foundObject(const word &name) const
Is the named Type found?
An STL-conforming iterator.
An STL-conforming hash table.
basicThermo(const basicThermo &)
Construct as copy (not implemented)
Abstract base-class for fluid and solid thermodynamic properties.
static Table::iterator lookupThermo(const dictionary &thermoDict, Table *tablePtr)
Generic lookup for each of the related thermodynamics packages.
A class for handling words, derived from string.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
virtual volScalarField & p()
Pressure [Pa].
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
const word & constant() const
Return constant name.
virtual ~basicThermo()
Destructor.
void validate(const string &app, const word &) const
Check that the thermodynamics package is consistent.
const word dictName("particleTrackDict")
const Time & time() const
Return the top-level database.
graph_traits< Graph >::vertices_size_type size_type
virtual bool read()
Read thermophysical properties dictionary.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Dimension set for the base types.
void store()
Transfer ownership of this object to its registry.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
virtual bool read()
Read object.
static wordList splitThermoName(const word &thermoName, const int nCmpt)
Split name of thermo package into a list of the components names.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const word & name() const
Return name.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
rDeltaT dimensionedInternalField()
wordList heBoundaryTypes()
Return the enthalpy/internal energy field boundary types.
const objectRegistry & db() const
Return local objectRegistry.
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
iterator begin()
Iterator set to the beginning of the HashTable.
const word dictName() const
Return the local dictionary name (final part of scoped name)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
wordList heBoundaryBaseTypes()
Return the enthalpy/internal energy field boundary base types.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
static const word null
An empty word.
defineTypeNameAndDebug(combustionModel, 0)