Go to the documentation of this file.
43 template<
class BasePhaseModel>
47 word phiName(IOobject::groupName(
"phi", this->
name()));
52 U.mesh().time().timeName(),
57 if (phiHeader.headerOk())
59 Info<<
"Reading face flux field " << phiName <<
endl;
61 return tmp<surfaceScalarField>
68 U.mesh().time().timeName(),
79 Info<<
"Calculating face flux field " << phiName <<
endl;
83 U.boundaryField().size(),
84 calculatedFvPatchScalarField::typeName
91 isA<fixedValueFvPatchVectorField>(
U.boundaryField()[i])
92 || isA<slipFvPatchVectorField>(
U.boundaryField()[i])
93 || isA<partialSlipFvPatchVectorField>(
U.boundaryField()[i])
96 phiTypes[i] = fixedValueFvPatchScalarField::typeName;
100 return tmp<surfaceScalarField>
107 U.mesh().time().timeName(),
122 template<
class BasePhaseModel>
125 const phaseSystem&
fluid,
126 const word& phaseName,
130 BasePhaseModel(
fluid, phaseName, index),
135 IOobject::groupName(
"U", this->
name()),
148 IOobject::groupName(
"alphaPhi", this->
name()),
159 IOobject::groupName(
"alphaRhoPhi", this->
name()),
170 IOobject::groupName(
"DUDt", this->
name()),
194 IOobject::groupName(
"continuityError", this->
name()),
209 template<
class BasePhaseModel>
216 template<
class BasePhaseModel>
221 this->
fluid().MRF().correctBoundaryVelocity(U_);
227 - (this->
fluid().fvOptions()(*
this,
rho)&rho);
231 template<
class BasePhaseModel>
234 BasePhaseModel::correctKinematics();
240 template<
class BasePhaseModel>
243 BasePhaseModel::correctTurbulence();
245 turbulence_->correct();
249 template<
class BasePhaseModel>
252 BasePhaseModel::correctEnergyTransport();
253 turbulence_->correctEnergyTransport();
257 template<
class BasePhaseModel>
265 -
fvm::Sp(continuityError_, U_)
267 + turbulence_->divDevRhoReff(U_)
272 template<
class BasePhaseModel>
287 template<
class BasePhaseModel>
290 const tmp<surfaceScalarField>& DbyA
297 template<
class BasePhaseModel>
305 template<
class BasePhaseModel>
313 template<
class BasePhaseModel>
321 template<
class BasePhaseModel>
329 template<
class BasePhaseModel>
333 const tmp<volScalarField>&
divU
340 template<
class BasePhaseModel>
344 return continuityError_;
348 template<
class BasePhaseModel>
356 template<
class BasePhaseModel>
364 template<
class BasePhaseModel>
372 template<
class BasePhaseModel>
380 template<
class BasePhaseModel>
388 template<
class BasePhaseModel>
396 template<
class BasePhaseModel>
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
virtual tmp< volVectorField > U() const
Constant access the velocity.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionSet dimDensity
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
Calculate the divergence of the given field.
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
Calculate the matrix for the divergence of the given field and flux.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
virtual const surfaceScalarField & DbyA() const
Return the phase diffusivity divided by the momentum coefficient.
List< word > wordList
A List of words.
virtual void correctKinematics()
Correct the kinematics.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual void correctTurbulence()
Correct the turbulence.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual const phaseCompressibleTurbulenceModel & turbulence() const
Return the turbulence model.
Vector< scalar > vector
A scalar version of the templated Vector.
Calculate the matrix for implicit and explicit sources.
const dimensionSet dimAcceleration
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Constant access the mass flux of the phase.
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual const tmp< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual tmp< volScalarField > rho() const
Density [kg/m^3].
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
Calculate the first temporal derivative.
virtual tmp< volScalarField > continuityError() const
Constant access the continuity error.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
word name(const complex &)
Return a string representation of a complex.
Calulate the matrix for the first temporal derivative.
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)