33 template<
class BasicSol
idThermo,
class MixtureType>
37 const auto& pCells = this->p_;
40 scalarField& CpCells = this->Cp_.primitiveFieldRef();
41 scalarField& CvCells = this->Cv_.primitiveFieldRef();
42 scalarField& rhoCells = this->rho_.primitiveFieldRef();
43 scalarField& alphaCells = this->alpha_.primitiveFieldRef();
47 const typename MixtureType::thermoMixtureType& thermoMixture =
48 this->cellThermoMixture(celli);
50 const typename MixtureType::transportMixtureType& transportMixture =
51 this->cellTransportMixture(celli, thermoMixture);
53 TCells[celli] = thermoMixture.THE
60 CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
61 CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
62 rhoCells[celli] = thermoMixture.rho(pCells[celli], TCells[celli]);
65 transportMixture.kappa(pCells[celli], TCells[celli])
66 /thermoMixture.Cv(pCells[celli], TCells[celli]);
70 volScalarField::Boundary& heBf =
71 this->
he().boundaryFieldRef();
73 const auto& pBf = this->p_.boundaryField();
75 volScalarField::Boundary& TBf =
76 this->T_.boundaryFieldRef();
78 volScalarField::Boundary& CpBf =
79 this->Cp_.boundaryFieldRef();
81 volScalarField::Boundary& CvBf =
82 this->Cv_.boundaryFieldRef();
84 volScalarField::Boundary& rhoBf =
85 this->rho_.boundaryFieldRef();
87 volScalarField::Boundary& alphaBf =
88 this->alpha_.boundaryFieldRef();
93 const auto& pp = pBf[
patchi];
104 const typename MixtureType::thermoMixtureType&
105 thermoMixture = this->patchFaceThermoMixture(patchi, facei);
107 const typename MixtureType::transportMixtureType&
109 this->patchFaceTransportMixture
110 (patchi, facei, thermoMixture);
112 phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
114 prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
115 pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
116 pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
119 transportMixture.kappa(pp[facei], pT[facei])
120 /thermoMixture.Cv(pp[facei], pT[facei]);
127 const typename MixtureType::thermoMixtureType&
128 thermoMixture = this->patchFaceThermoMixture(patchi, facei);
130 const typename MixtureType::transportMixtureType&
132 this->patchFaceTransportMixture
133 (patchi, facei, thermoMixture);
135 pT[facei] = thermoMixture.THE(phe[facei], pp[facei] ,pT[facei]);
137 prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
138 pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
139 pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
142 transportMixture.kappa(pp[facei], pT[facei])
143 /thermoMixture.Cv(pp[facei], pT[facei]);
148 this->alpha_.correctBoundaryConditions();
154 template<
class BasicSol
idThermo,
class MixtureType>
159 const word& phaseName
168 template<
class BasicSol
idThermo,
class MixtureType>
174 const word& phaseName
185 template<
class BasicSol
idThermo,
class MixtureType>
192 template<
class BasicSol
idThermo,
class MixtureType>
209 template<
class BasicSol
idThermo,
class MixtureType>
213 const fvMesh& mesh = this->T_.mesh();
225 const auto& pCells = this->p_;
234 this->cellTransportMixture
237 ).Kappa(pCells[celli], TCells[celli]);
244 const auto& pp = this->p_.boundaryField()[
patchi];
252 this->patchFaceTransportMixture
256 ).Kappa(pp[facei], pT[facei]);
264 template<
class BasicSol
idThermo,
class MixtureType>
271 const auto& pp = this->p_.boundaryField()[
patchi];
280 this->patchFaceTransportMixture
284 ).Kappa(pp[patchi], Tp[facei]);
291 template<
class BasicSol
idThermo,
class MixtureType>
295 const fvMesh& mesh = this->T_.mesh();
316 KappaLocal.primitiveFieldRef() =
317 coordinates.
R(mesh.
C()).transformVector(Kappa);
321 KappaLocal.boundaryFieldRef()[
patchi] =
330 template<
class BasicSol
idThermo,
class MixtureType>
337 const fvMesh& mesh = this->T_.mesh();
346 .transformVector(Kappa(patchi));
350 template<
class BasicSol
idThermo,
class MixtureType>
354 const fvMesh& mesh = this->T_.mesh();
366 template<
class BasicSol
idThermo,
class MixtureType>
383 KappaLocal()/this->
Cv(),
385 "laplacian(" + this->
alpha().
name() +
",e)" Base class for other coordinate system specifications.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
scalar Cv(const scalar p, const scalar T) const
#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...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
void size(const label)
Override size to be inconsistent with allocated storage.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Calculate the matrix for the laplacian of the field.
void setFluxRequired(const word &name) const
volVectorField vectorField(fieldObject, mesh)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
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.
Calculate the laplacian of the given field.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
volScalarField scalarField(fieldObject, mesh)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
heSolidThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
const dimensionSet dimEnergy
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
word name(const complex &)
Return a string representation of a complex.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Internal & ref()
Return a reference to the dimensioned internal field.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Energy for a solid mixture.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
const volVectorField & C() const
Return cell centres as volVectorField.
A class for managing temporary objects.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
virtual ~heSolidThermo()
Destructor.
const dimensionSet dimTemperature
const coordinateRotation & R() const
Return const reference to co-ordinate rotation.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
#define InfoInFunction
Report an information message using Foam::Info.