Go to the documentation of this file.
26 #include "phasePressureModel.H"
27 #include "twoPhaseSystem.H"
38 const transportModel& phase,
39 const word& propertiesName,
66 dimensionSet(1, -1, -2, 0, 0),
93 RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel> >
97 coeffDict().lookup(
"alphaMax") >> alphaMax_;
98 coeffDict().lookup(
"preAlphaExp") >> preAlphaExp_;
99 coeffDict().lookup(
"expMax") >> expMax_;
100 g0_.readIfPresent(coeffDict());
130 return tmp<volSymmTensorField>
143 dimensioned<symmTensor>
146 dimensionSet(0, 2, -2, 0, 0),
157 tmp<volScalarField> tpPrime
162 exp(preAlphaExp_*(alpha_ - alphaMax_)),
167 volScalarField::GeometricBoundaryField& bpPrime = tpPrime().boundaryField();
171 if (!bpPrime[
patchi].coupled())
184 tmp<surfaceScalarField> tpPrime
194 surfaceScalarField::GeometricBoundaryField& bpPrime =
195 tpPrime().boundaryField();
199 if (!bpPrime[
patchi].coupled())
212 return tmp<volSymmTensorField>
225 dimensioned<symmTensor>
228 rho_.dimensions()*dimensionSet(0, 2, -2, 0, 0),
242 return tmp<fvVectorMatrix>
247 rho_.dimensions()*dimensionSet(0, 4, -2, 0, 0)
static word groupName(Name name, const word &group)
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
dimensionedScalar exp(const dimensionedScalar &ds)
phasePressureModel(const phasePressureModel &)
Disallow default bitwise copy construct.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual ~phasePressureModel()
Destructor.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
fvMatrix< vector > fvVectorMatrix
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const SymmTensor zero
virtual bool read()
Re-read model coefficients if they have changed.
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);singlePhaseTransportModel laminarTransport(U, phi);autoPtr< incompressible::RASModel > RASModel(incompressible::New< incompressible::RASModel >(U, phi, laminarTransport))
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
stressControl lookup("compactNormalStress") >> compactNormalStress