32#include "twoPhaseSystem.H"
34#include "virtualMassModel.H"
48template<
class BasicTurbulenceModel>
49mixtureKEpsilon<BasicTurbulenceModel>::mixtureKEpsilon
57 const word& propertiesName,
73 liquidTurbulencePtr_(
nullptr),
144 this->runTime_.timeName(),
156 this->runTime_.timeName(),
169 this->printCoeffs(
type);
174template<
class BasicTurbulenceModel>
188 ebt[patchi] = fixedValueFvPatchScalarField::typeName;
196template<
class BasicTurbulenceModel>
216 (bf[patchi]).refValue() =
218 (refBf[patchi]).refValue();
224template<
class BasicTurbulenceModel>
234 mixtureKEpsilon<BasicTurbulenceModel>& turbc = this->liquidTurbulence();
240 this->runTime_.timeName(this->runTime_.startTime().value())
288 kl.boundaryField().types()
291 correctInletOutlet(km_(), kl);
305 mix(epsilonl, epsilong),
306 epsilonBoundaryTypes(epsilonl)
315template<
class BasicTurbulenceModel>
320 Cmu_.readIfPresent(this->coeffDict());
321 C1_.readIfPresent(this->coeffDict());
322 C2_.readIfPresent(this->coeffDict());
323 C3_.readIfPresent(this->coeffDict());
324 Cp_.readIfPresent(this->coeffDict());
325 sigmak_.readIfPresent(this->coeffDict());
326 sigmaEps_.readIfPresent(this->coeffDict());
335template<
class BasicTurbulenceModel>
338 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
339 this->nut_.correctBoundaryConditions();
342 BasicTurbulenceModel::correctNut();
346template<
class BasicTurbulenceModel>
347mixtureKEpsilon<BasicTurbulenceModel>&
348mixtureKEpsilon<BasicTurbulenceModel>::liquidTurbulence()
const
350 if (!liquidTurbulencePtr_)
354 const transportModel& gas = this->transport();
355 const twoPhaseSystem&
fluid =
357 const transportModel& liquid =
fluid.otherPhase(gas);
359 liquidTurbulencePtr_ =
360 &
const_cast<mixtureKEpsilon<BasicTurbulenceModel>&
>
362 U.db().lookupObject<mixtureKEpsilon<BasicTurbulenceModel>>
373 return *liquidTurbulencePtr_;
377template<
class BasicTurbulenceModel>
380 const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
381 this->liquidTurbulence();
393 (6*this->Cmu_/(4*
sqrt(3.0/2.0)))
395 *(liquidTurbulence.k_/liquidTurbulence.epsilon_)
400 return sqr(1 + (Ct0 - 1)*
exp(-fAlphad));
404template<
class BasicTurbulenceModel>
409 return fluid.otherPhase(gas).rho();
413template<
class BasicTurbulenceModel>
422 + virtualMass.
Cvm()*
fluid.otherPhase(gas).rho();
426template<
class BasicTurbulenceModel>
436template<
class BasicTurbulenceModel>
450template<
class BasicTurbulenceModel>
466template<
class BasicTurbulenceModel>
488template<
class BasicTurbulenceModel>
492 this->liquidTurbulence();
526template<
class BasicTurbulenceModel>
533template<
class BasicTurbulenceModel>
540template<
class BasicTurbulenceModel>
547 if (&gas != &
fluid.phase1())
551 this->liquidTurbulence();
556 if (!this->turbulence_)
565 tmp<surfaceScalarField>
phig = this->
phi();
573 mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
574 this->liquidTurbulence();
575 tmp<surfaceScalarField> phil = liquidTurbulence.phi();
589 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
607 tmp<volScalarField> Gc;
609 tmp<volTensorField> tgradUl =
fvc::grad(Ul);
610 Gc = tmp<volScalarField>
621 kl.boundaryFieldRef().updateCoeffs();
623 kl.boundaryFieldRef().evaluateCoupled<coupledFvPatch>();
625 epsilonl.boundaryFieldRef().updateCoeffs();
626 epsilonl.boundaryFieldRef().evaluateCoupled<coupledFvPatch>();
631 tmp<volScalarField> Gd;
633 tmp<volTensorField> tgradUg =
fvc::grad(Ug);
634 Gd = tmp<volScalarField>
645 kg.boundaryFieldRef().updateCoeffs();
647 kg.boundaryFieldRef().evaluateCoupled<coupledFvPatch>();
648 epsilong.boundaryFieldRef().updateCoeffs();
650 epsilong.boundaryFieldRef().evaluateCoupled<coupledFvPatch>();
663 bound(km, this->kMin_);
664 epsilonm == mix(epsilonl, epsilong);
665 bound(epsilonm, this->epsilonMin_);
668 tmp<fvScalarMatrix> epsEqn
676 -
fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
677 -
fvm::Sp(C2_*epsilonm/km, epsilonm)
682 epsEqn.ref().relax();
684 epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
687 bound(epsilonm, this->epsilonMin_);
691 tmp<fvScalarMatrix> kmEqn
709 bound(km, this->kMin_);
710 km.correctBoundaryConditions();
714 kl.correctBoundaryConditions();
715 epsilonl = Cc2*epsilonm;
716 epsilonl.correctBoundaryConditions();
717 liquidTurbulence.correctNut();
721 kg.correctBoundaryConditions();
722 epsilong = Ct2_()*epsilonl;
723 epsilong.correctBoundaryConditions();
724 nutg = Ct2_()*(liquidTurbulence.nu()/this->
nu())*nutl;
Bound the given scalar field if it has gone unbounded.
wordList types() const
Return a list of the patch types.
void evaluateCoupled(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions on coupled patches of the given type, using specified or default comms.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word group(const word &name)
Return group (extension part of name).
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Templated abstract base class for RAS turbulence models.
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > DepsilonEff(const volScalarField &nutm) const
Return the effective diffusivity for epsilon.
dimensionedScalar sigmak_
BasicTurbulenceModel::rhoField rhoField
tmp< volScalarField > Ct2() const
tmp< volScalarField > rhogEff() const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
dimensionedScalar sigmaEps_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > rhom() const
tmp< volScalarField > bubbleG() const
autoPtr< volScalarField > rhom_
tmp< volScalarField > DkEff(const volScalarField &nutm) const
Return the effective diffusivity for k.
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
autoPtr< volScalarField > km_
autoPtr< volScalarField > Ct2_
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
tmp< volScalarField > rholEff() const
autoPtr< volScalarField > epsilonm_
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual bool read()
Re-read model coefficients if they have changed.
An abstract base class for patches that couple regions of the computational domain e....
Generic dimensioned Type class.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
static const word propertiesName
Default name of the turbulence properties dictionary.
Class which solves the volume fraction equations for two phases.
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
A class for handling words, derived from Foam::string.
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
Calculate the finiteVolume matrix for implicit and explicit sources.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
List< word > wordList
List of word.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.