45template<
class BasePhaseModel,
class phaseThermo>
53 BasePhaseModel(
fluid, phaseName),
71 <<
" The selected thermo is pure. Use a multicomponent thermo."
80 thermoPtr_().template getOrDefault<bool>(
"addDiffusion",
false);
97 fluid.mesh().time().timeName(),
115template<
class BasePhaseModel,
class phaseThermo>
135 X_[i].correctBoundaryConditions();
145 const scalarField& ybf =
Y()[i].boundaryField()[patchi];
149 xbf[facei] = Wbf[facei]*ybf[facei]/Wi.
value();
155 X_[inertIndex_] = 1.0 - Xtotal;
160template<
class BasePhaseModel,
class phaseThermo>
165 for(label i=1; i< species_.size(); i++)
174 Info<<
Y()[i].name() <<
" mass fraction = "
175 <<
" Min(Y) = " <<
min(
Y()[i]).value()
176 <<
" Max(Y) = " <<
max(
Y()[i]).value()
182template<
class BasePhaseModel,
class phaseThermo>
190template<
class BasePhaseModel,
class phaseThermo>
198template<
class BasePhaseModel,
class phaseThermo>
205template<
class BasePhaseModel,
class phaseThermo>
218 scalar cAlpha(MULEScontrols.
get<scalar>(
"cYi"));
250 if (!phirp.coupled())
260 if (inertIndex_ != i)
269 "phi" + Yi.name() +
"Corr",
274 "div(phi," + Yi.name() +
')'
305 forAll(phiYiCorr.boundaryField(), patchi)
308 phiYiCorr.boundaryFieldRef()[patchi];
310 if (!phiYiCorrp.coupled())
313 const scalarField& Yip = Yi.boundaryField()[patchi];
317 if (phi1p[facei] < 0)
319 phiYiCorrp[facei] = Yip[facei]*phi1p[facei];
350 if (inertIndex_ != i)
372 phiYiCorr += YiEqn.flux();
374 if (nYiSubCycles > 1)
379 !(++YiSubCycle).
end();
423 YiDiffEqn.solve(
"diffusion" + Yi.name());
430 X_[inertIndex_] = scalar(1) -
Yt;
437template<
class BasePhaseModel,
class phaseThermo>
445template<
class BasePhaseModel,
class phaseThermo>
453template<
class BasePhaseModel,
class phaseThermo>
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
volScalarField Yt(0.0 *Y[0])
const volScalarField & alpha1
const volScalarField & alpha2
void clamp_min(const Type &lower)
Impose lower (floor) clamp on the field values (in-place).
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.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which represents a phase with multiple species. Returns the species' mass fractions,...
void calculateMassFractions()
Transfor volume fraction into mass fractions.
virtual void solveYi(PtrList< volScalarField::Internal > &, PtrList< volScalarField::Internal > &)
Solve species fraction equation.
virtual void correct()
Correct phase thermo.
label inertIndex() const
Return inert species index.
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
autoPtr< phaseThermo > thermoPtr_
Thermophysical model.
PtrList< volScalarField > X_
Ptr list of volumetric fractions for species.
scalar Sct_
Schmidt number.
virtual const phaseThermo & thermo() const
Access to thermo.
void calculateVolumeFractions()
Transfor mass fraction into volume fractions.
MultiComponentPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
label inertIndex_
Inert species index.
bool addDiffusion_
Add diffusion transport on Yi's Eq.
hashedWordList species_
Species table.
virtual const PtrList< volScalarField > & Y() const
Constant access the species mass fractions.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
scalar deltaTValue() const noexcept
Return time step value.
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
static word phasePropertyName(const word &name, const word &phaseName)
The phase property name as property.phase.
static const word dictName
The dictionary name ("thermophysicalProperties").
virtual void correct()=0
Update properties.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const Type & value() const noexcept
Return const reference to value.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
virtual bool fixesValue() const
True if the patch field fixes a value.
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Basic second-order convection using face-gradients and Gauss' theorem.
virtual bool coupled() const
True if the patch field is coupled.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
HashTable< autoPtr< multiphaseInter::phaseModel > > phaseModelTable
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
bool end() const
Return true if the number of sub-cycles has been reached.
Perform a subCycleTime on a field.
const surfaceScalarField & phi() const
Return the mixture flux.
Upwind differencing scheme class.
A class for handling words, derived from Foam::string.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
basicSpecieMixture & composition
PtrList< volScalarField > & Y
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
Calculate the substantive (total) derivative.
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 matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
fvsPatchField< scalar > fvsPatchScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
fvPatchField< scalar > fvPatchScalarField
psiReactionThermo & thermo
multiphaseSystem::phaseModelList & phases
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object.