34 template<
class BasicSol
idThermo,
class MixtureType>
41 scalarField& rhoCells = this->rho_.primitiveFieldRef();
42 scalarField& alphaCells = this->alpha_.primitiveFieldRef();
46 const typename MixtureType::thermoType& mixture_ =
47 this->cellMixture(celli);
49 const typename MixtureType::thermoType& volMixture_ =
50 this->cellVolMixture(pCells[celli], TCells[celli], celli);
54 TCells[celli] = mixture_.THE
62 rhoCells[celli] = volMixture_.rho(pCells[celli], TCells[celli]);
65 volMixture_.kappa(pCells[celli], TCells[celli])
67 mixture_.Cpv(pCells[celli], TCells[celli]);
70 volScalarField::Boundary& pBf = this->p_.boundaryFieldRef();
71 volScalarField::Boundary& TBf = this->T_.boundaryFieldRef();
72 volScalarField::Boundary& rhoBf = this->rho_.boundaryFieldRef();
73 volScalarField::Boundary& heBf = this->
he().boundaryFieldRef();
74 volScalarField::Boundary& alphaBf = this->alpha_.boundaryFieldRef();
76 forAll(this->T_.boundaryField(), patchi)
88 const typename MixtureType::thermoType& mixture_ =
89 this->patchFaceMixture(patchi, facei);
91 const typename MixtureType::thermoType& volMixture_ =
92 this->patchFaceVolMixture
100 phe[facei] = mixture_.HE(
pp[facei], pT[facei]);
101 prho[facei] = volMixture_.rho(
pp[facei], pT[facei]);
104 volMixture_.kappa(
pp[facei], pT[facei])
105 / mixture_.Cpv(
pp[facei], pT[facei]);
112 const typename MixtureType::thermoType& mixture_ =
113 this->patchFaceMixture(patchi, facei);
115 const typename MixtureType::thermoType& volMixture_ =
116 this->patchFaceVolMixture
126 pT[facei] = mixture_.THE(phe[facei],
pp[facei] ,pT[facei]);
129 prho[facei] = volMixture_.rho(
pp[facei], pT[facei]);
132 volMixture_.kappa(
pp[facei], pT[facei])
133 / mixture_.Cpv(
pp[facei], pT[facei]);
138 this->alpha_.correctBoundaryConditions();
143 template<
class BasicSol
idThermo,
class MixtureType>
148 const word& phaseName
151 heThermo<BasicSolidThermo, MixtureType>(
mesh, phaseName)
159 template<
class BasicSol
idThermo,
class MixtureType>
165 const word& phaseName
176 template<
class BasicSol
idThermo,
class MixtureType>
181 const word& phaseName,
198 template<
class BasicSol
idThermo,
class MixtureType>
205 template<
class BasicSol
idThermo,
class MixtureType>
216 template<
class BasicSol
idThermo,
class MixtureType>
229 auto& Kappa = tKappa.ref();
231 vectorField& KappaCells = Kappa.primitiveFieldRef();
243 ).Kappa(pCells[celli], TCells[celli]);
251 const scalarField& pT = this->T_.boundaryField()[patchi];
257 this->patchFaceVolMixture
263 ).Kappa(
pp[facei], pT[facei]);
271 template<
class BasicSol
idThermo,
class MixtureType>
279 const scalarField& Tp = this->T_.boundaryField()[patchi];
282 auto& Kappap = tKappa.ref();
287 this->patchFaceVolMixture
293 ).Kappa(
pp[facei], Tp[facei]);
A list of keyword definitions, which are a keyword followed by a number of values (eg...
constexpr char nl
The newline '\n' character (0x0a)
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
Type of boundary fields.
#define forAll(list, i)
Loop across all elements in list.
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
#define DebugInFunction
Report an information message using Foam::Info.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
#define DebugInfo
Report an information message using Foam::Info.
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const dimensionSet dimEnergy
Energy for a solid mixture.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Enthalpy/Internal energy for a mixture.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< vector > vectorField
Specialisation of Field<T> for vector.
A class for managing temporary objects.
virtual ~heSolidThermo()
Destructor.
Do not request registration (bool: false)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static constexpr const zero Zero
Global zero (0)