31#include "phaseCompressibleTurbulenceModel.H"
42#include "surfaceInterpolate.H"
47template<
class BasePhaseModel>
53 IOobject::groupName(
"phi", this->
name()),
54 U.mesh().time().timeName(),
61 if (
io.typeHeaderOk<surfaceScalarField>(
true))
63 Info<<
"Reading face flux field " <<
io.name() <<
endl;
65 return tmp<surfaceScalarField>::New(
io,
U.mesh());
69 Info<<
"Calculating face flux field " <<
io.name() <<
endl;
71 io.readOpt(IOobject::NO_READ);
75 U.boundaryField().size(),
76 fvsPatchFieldBase::calculatedType()
83 isA<fixedValueFvPatchVectorField>(
U.boundaryField()[i])
84 || isA<slipFvPatchVectorField>(
U.boundaryField()[i])
85 || isA<partialSlipFvPatchVectorField>(
U.boundaryField()[i])
88 patchTypes[i] = fixedValueFvPatchScalarField::typeName;
92 return tmp<surfaceScalarField>::New
104template<
class BasePhaseModel>
107 const phaseSystem&
fluid,
108 const word& phaseName,
112 BasePhaseModel(
fluid, phaseName, index),
117 IOobject::groupName(
"U", this->
name()),
121 IOobject::AUTO_WRITE,
131 IOobject::groupName(
"alphaPhi", this->
name()),
136 dimensionedScalar(dimensionSet(0, 3, -1, 0, 0))
142 IOobject::groupName(
"alphaRhoPhi", this->
name()),
147 dimensionedScalar(dimensionSet(1, 0, -1, 0, 0))
154 phaseCompressibleTurbulenceModel::New
168 IOobject::groupName(
"continuityErrorFlow", this->
name()),
173 dimensionedScalar(dimDensity/dimTime)
175 continuityErrorSources_
179 IOobject::groupName(
"continuityErrorSources", this->
name()),
184 dimensionedScalar(dimDensity/dimTime)
196template<
class BasePhaseModel>
199 BasePhaseModel::correct();
211template<
class BasePhaseModel>
214 BasePhaseModel::correctKinematics();
235template<
class BasePhaseModel>
238 BasePhaseModel::correctTurbulence();
244template<
class BasePhaseModel>
247 BasePhaseModel::correctEnergyTransport();
253template<
class BasePhaseModel>
260template<
class BasePhaseModel>
272 +
fvm::SuSp(- this->continuityError(), U_)
279template<
class BasePhaseModel>
291 +
fvm::SuSp(- this->continuityErrorSources(), U_)
293 + turbulence_->divDevRhoReff(U_)
298template<
class BasePhaseModel>
306template<
class BasePhaseModel>
314template<
class BasePhaseModel>
322template<
class BasePhaseModel>
330template<
class BasePhaseModel>
338template<
class BasePhaseModel>
346template<
class BasePhaseModel>
354template<
class BasePhaseModel>
362template<
class BasePhaseModel>
375template<
class BasePhaseModel>
381 DUDtf_ =
byDt(phi_ - phi_.oldTime());
388template<
class BasePhaseModel>
396template<
class BasePhaseModel>
404template<
class BasePhaseModel>
412template<
class BasePhaseModel>
429template<
class BasePhaseModel>
437template<
class BasePhaseModel>
444template<
class BasePhaseModel>
452template<
class BasePhaseModel>
460template<
class BasePhaseModel>
468template<
class BasePhaseModel>
476template<
class BasePhaseModel>
484template<
class BasePhaseModel>
492template<
class BasePhaseModel>
500template<
class BasePhaseModel>
508template<
class BasePhaseModel>
516template<
class BasePhaseModel>
520 return turbulence_->pPrime();
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
writeOption writeOpt() const noexcept
Get the write option.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual void correctEnergyTransport()
Correct the energy transport e.g. alphat.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
tmp< volScalarField > divU_
Dilatation rate.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volVectorField > U() const
Access const reference to U.
tmp< surfaceScalarField > DUDtf_
Lagrangian acceleration field on the faces (needed for virtual-mass).
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
volScalarField continuityErrorFlow_
Continuity error due to the flow.
autoPtr< phaseCompressibleTurbulenceModel > turbulence_
Turbulence model.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
surfaceScalarField alphaRhoPhi_
Mass flux.
virtual void correctTurbulence()
Correct the turbulence.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
surfaceScalarField phi_
Flux.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi)).
tmp< volVectorField > DUDt_
Lagrangian acceleration field (needed for virtual-mass).
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual bool stationary() const
Return whether the phase is stationary.
tmp< volScalarField > K_
Kinetic Energy.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
volScalarField continuityErrorSources_
Continuity error due to any sources.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Class to represent a system of phases and model interfacial transfers between them.
const IOMRFZoneList & MRF() const
Return MRF zones.
fv::options & fvOptions() const
Access the fvOptions.
virtual tmp< volScalarField > rho() const
Density [kg/m^3] - uses current value of pressure.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the finiteVolume matrix for implicit and explicit sources.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
List< word > wordList
List of word.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
tmp< volScalarField > byDt(const volScalarField &vf)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
wordList patchTypes(nPatches)
tmp< volScalarField > trho
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.