31 template<
class BasicSol
idThermo,
class MixtureType>
38 scalarField& rhoCells = this->rho_.primitiveFieldRef();
39 scalarField& alphaCells = this->alpha_.primitiveFieldRef();
43 const typename MixtureType::thermoType& mixture_ =
44 this->cellMixture(celli);
46 rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
48 TCells[celli] = mixture_.THE
55 const typename MixtureType::thermoType& volMixture_ =
56 this->cellVolMixture(pCells[celli], TCells[celli], celli);
59 volMixture_.kappa(pCells[celli], TCells[celli])
60 /mixture_.Cpv(pCells[celli], TCells[celli]);
63 volScalarField::Boundary& pBf =
64 this->p_.boundaryFieldRef();
66 volScalarField::Boundary& TBf =
67 this->T_.boundaryFieldRef();
69 volScalarField::Boundary& rhoBf =
70 this->rho_.boundaryFieldRef();
72 volScalarField::Boundary& heBf =
73 this->
he().boundaryFieldRef();
75 volScalarField::Boundary& alphaBf =
76 this->alpha_.boundaryFieldRef();
90 const typename MixtureType::thermoType& mixture_ =
91 this->patchFaceMixture(patchi, facei);
93 const typename MixtureType::thermoType& volMixture_ =
94 this->patchFaceVolMixture
103 phe[facei] = mixture_.HE(pp[facei], pT[facei]);
104 prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
107 volMixture_.kappa(pp[facei], pT[facei])
108 / mixture_.Cpv(pp[facei], pT[facei]);
115 const typename MixtureType::thermoType& mixture_ =
116 this->patchFaceMixture(patchi, facei);
118 const typename MixtureType::thermoType& volMixture_ =
119 this->patchFaceVolMixture
127 pT[facei] = mixture_.THE(phe[facei], pp[facei] ,pT[facei]);
128 prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
131 volMixture_.kappa(pp[facei], pT[facei])
132 / mixture_.Cpv(pp[facei], pT[facei]);
137 this->alpha_.correctBoundaryConditions();
143 template<
class BasicSol
idThermo,
class MixtureType>
148 const word& phaseName
157 template<
class BasicSol
idThermo,
class MixtureType>
163 const word& phaseName
174 template<
class BasicSol
idThermo,
class MixtureType>
181 template<
class BasicSol
idThermo,
class MixtureType>
198 template<
class BasicSol
idThermo,
class MixtureType>
202 const fvMesh& mesh = this->T_.mesh();
227 ).Kappa(pCells[celli], TCells[celli]);
241 this->patchFaceVolMixture
247 ).Kappa(pp[facei], pT[facei]);
255 template<
class BasicSol
idThermo,
class MixtureType>
271 this->patchFaceVolMixture
277 ).Kappa(pp[facei], Tp[facei]);
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A list of keyword definitions, which are a keyword followed by any number of values (e...
T & ref() const
Return non-const reference or generate a fatal error.
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
volVectorField vectorField(fieldObject, mesh)
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
fvPatchField< scalar > fvPatchScalarField
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
volScalarField scalarField(fieldObject, mesh)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const dimensionSet dimEnergy
Energy for a solid mixture.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
A class for managing temporary objects.
virtual ~heSolidThermo()
Destructor.
#define InfoInFunction
Report an information message using Foam::Info.