Go to the documentation of this file.
28 #include "twoPhaseSystem.H"
29 #include "dragModel.H"
30 #include "virtualMassModel.H"
44 template<
class BasicTurbulenceModel>
53 const word& propertiesName,
69 liquidTurbulencePtr_(NULL),
140 this->runTime_.timeName(),
152 this->runTime_.timeName(),
160 bound(k_, this->kMin_);
161 bound(epsilon_, this->epsilonMin_);
163 if (
type == typeName)
165 this->printCoeffs(
type);
170 template<
class BasicTurbulenceModel>
182 if (isA<fixedValueFvPatchScalarField>(ebf[
patchi]))
184 ebt[
patchi] = fixedValueFvPatchScalarField::typeName;
192 template<
class BasicTurbulenceModel>
207 isA<inletOutletFvPatchScalarField>(bf[
patchi])
208 && isA<inletOutletFvPatchScalarField>(refBf[
patchi])
211 refCast<inletOutletFvPatchScalarField>
213 refCast<const inletOutletFvPatchScalarField>
214 (refBf[
patchi]).refValue();
220 template<
class BasicTurbulenceModel>
223 if (rhom_.valid())
return;
236 this->runTime_.timeName(this->runTime_.startTime().value())
287 correctInletOutlet(km_(), kl);
301 mix(epsilonl, epsilong),
302 epsilonBoundaryTypes(epsilonl)
305 correctInletOutlet(epsilonm_(), epsilonl);
311 template<
class BasicTurbulenceModel>
316 Cmu_.readIfPresent(this->coeffDict());
317 C1_.readIfPresent(this->coeffDict());
318 C2_.readIfPresent(this->coeffDict());
319 C3_.readIfPresent(this->coeffDict());
320 Cp_.readIfPresent(this->coeffDict());
321 sigmak_.readIfPresent(this->coeffDict());
322 sigmaEps_.readIfPresent(this->coeffDict());
333 template<
class BasicTurbulenceModel>
336 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
337 this->nut_.correctBoundaryConditions();
341 template<
class BasicTurbulenceModel>
345 if (!liquidTurbulencePtr_)
351 refCast<const twoPhaseSystem>(gas.fluid());
354 liquidTurbulencePtr_ =
368 return *liquidTurbulencePtr_;
372 template<
class BasicTurbulenceModel>
376 this->liquidTurbulence();
388 (6*this->Cmu_/(4*
sqrt(3.0/2.0)))
389 *
fluid.drag(gas).K()/liquid.rho()
390 *(liquidTurbulence.
k_/liquidTurbulence.
epsilon_)
393 volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
395 return sqr(1 + (Ct0 - 1)*
exp(-fAlphad));
399 template<
class BasicTurbulenceModel>
404 return fluid.otherPhase(gas).rho();
408 template<
class BasicTurbulenceModel>
415 +
fluid.virtualMass(gas).Cvm()*
fluid.otherPhase(gas).rho();
419 template<
class BasicTurbulenceModel>
425 return alphal*rholEff() + alphag*rhogEff();
429 template<
class BasicTurbulenceModel>
439 return (
alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
443 template<
class BasicTurbulenceModel>
454 (
alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
455 /(
alphal*rholEff() + alphag*rhogEff()*Ct2_());
459 template<
class BasicTurbulenceModel>
481 template<
class BasicTurbulenceModel>
485 this->liquidTurbulence();
500 +
pow(
fluid.drag(gas).CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
517 template<
class BasicTurbulenceModel>
520 return fvm::Su(bubbleG()/rhom_(), km_());
524 template<
class BasicTurbulenceModel>
527 return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
531 template<
class BasicTurbulenceModel>
538 if (&gas != &
fluid.phase1())
542 this->liquidTurbulence();
547 if (!this->turbulence_)
565 this->liquidTurbulence();
644 bound(km, this->kMin_);
645 epsilonm ==
mix(epsilonl, epsilong);
646 bound(epsilonm, this->epsilonMin_);
657 -
fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
658 -
fvm::Sp(C2_*epsilonm/km, epsilonm)
667 bound(epsilonm, this->epsilonMin_);
686 bound(km, this->kMin_);
692 epsilonl = Cc2*epsilonm;
699 epsilong = Ct2_()*epsilonl;
701 nutg = Ct2_()*(liquidTurbulence.nu()/this->
nu())*nutl;
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static word groupName(Name name, const word &group)
BasicTurbulenceModel::rhoField rhoField
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
void updateCoeffs()
Update the boundary condition coefficients.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from string.
Class which solves the volume fraction equations for two phases.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Templated abstract base class for RAS turbulence models.
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
surfaceScalarField phig(-rAUf *ghf *fvc::snGrad(rhok) *mesh.magSf())
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
tmp< volScalarField > bubbleG() const
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.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
static const word propertiesName
Default name of the turbulence properties dictionary.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
virtual void correctNut()
dimensionedScalar exp(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
BasicTurbulenceModel::transportModel transportModel
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > rhom() const
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedScalar pow3(const dimensionedScalar &ds)
mixtureKEpsilon< BasicTurbulenceModel > & liquidTurbulence() const
Return the turbulence model for the other phase.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Bound the given scalar field if it has gone unbounded.
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< volScalarField > rholEff() const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Generic dimensioned Type class.
tmp< volScalarField > Ct2() const
virtual bool read()
Re-read model coefficients if they have changed.
Calculate the matrix for implicit and explicit sources.
void correctBoundaryConditions()
Correct boundary field.
wordList types() const
Return a list of the patch types.
tmp< volScalarField > rhogEff() const
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< fvScalarMatrix > kSource() const
Eddy viscosity turbulence model base class.
BasicTurbulenceModel::alphaField alphaField
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void clear() const
If object pointer points to valid object:
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
mixtureKEpsilon(const mixtureKEpsilon &)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)