Go to the documentation of this file.
18 IOobject::groupName(
"rAU", phase.name()),
22 +
max(phase.residualAlpha() -
alpha, scalar(0))
23 *phase.rho()/runTime.deltaT()
40 autoPtr<PtrList<volVectorField> > Fs =
fluid.Fs();
53 IOobject::groupName(
"phiF", phase.name()),
61 autoPtr<PtrList<surfaceScalarField> > phiDs =
fluid.phiDs(
rAUs);
80 IOobject::groupName(
"phiF", phase.name()),
103 MRF.correctBoundaryFlux(phase.U(), phase.phi());
110 IOobject::groupName(
"HbyA", phase.name()),
119 +
max(phase.residualAlpha() -
alpha, scalar(0))
120 *phase.rho()*phase.U().oldTime()/runTime.deltaT()
133 PtrList<surfaceScalarField> phigFs(
phases.size());
146 -
fluid.surfaceTension(phase)*
mesh.magSf()
191 || isA<cyclicAMIFvPatch>(
mesh.boundary()[
patchi])
194 phiCorrCoeff.boundaryField()[
patchi] = 0;
203 IOobject::groupName(
"phiHbyA", phase.name()),
211 MRF.absolute(phase.phi().oldTime())
220 phaseSystem::KdTable,
227 const phasePair& pair(
fluid.phasePairs()[KdIter.key()]);
229 const phaseModel*
phase1 = &pair.phase1();
230 const phaseModel*
phase2 = &pair.phase2();
274 surfaceScalarField::GeometricBoundaryField
phib(
phi.boundaryField());
284 p_rgh.boundaryField(),
287 )/(
mesh.magSf().boundaryField()*
rAUf.boundaryField())
291 PtrList<fvScalarMatrix> pEqnComps(
phases.size());
296 if (phase.compressible())
305 IOobject::groupName(
"phid", phase.name()),
314 phase.continuityError()
332 pEqnComps[
phasei].faceFluxCorrectionPtr()
334 pEqnComps[
phasei].relax();
343 phase.continuityError()
356 if (
fluid.transfersMass(phase))
358 if (pEqnComps.set(
phasei))
378 while (
pimple.correctNonOrthogonal())
392 if (pEqnComps.set(
phasei))
394 pEqn += pEqnComps[
phasei];
406 if (
pimple.finalNonOrthogonalIter())
419 if (phase.compressible())
428 mSfGradp = pEqnIncomp.flux()/
rAUf;
441 phase.U().correctBoundaryConditions();
457 phase.rho()() += phase.thermo().psi()*(
p_rgh -
p_rgh_0);
463 p_rgh.correctBoundaryConditions();
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
surfaceScalarField phid("phid", fvc::interpolate(psi) *((mesh.Sf() &fvc::interpolate(HbyA))+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
const surfaceScalarField & ghf
const dimensionSet dimVelocity
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
dimensionedScalar pMin("pMin", dimPressure, fluid)
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionedVector & g
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryField(),(phiHbyA.boundaryField() - MRF.relative(mesh.Sf().boundaryField() &U.boundaryField()))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
CGAL::Exact_predicates_exact_constructions_kernel K
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void deleteDemandDrivenData(DataPtr &dataPtr)
const dimensionSet dimArea(sqr(dimLength))
fvMatrix< scalar > fvScalarMatrix
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
PtrList< volVectorField > HbyAs(fluid.phases().size())
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
const volScalarField & gh
volScalarField p_rgh_0(p_rgh)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
PtrList< surfaceScalarField > phiFs(phases.size())
multiphaseSystem::phaseModelList & phases
PtrList< surfaceScalarField > alpharAUfs(phases.size())
PtrList< surfaceScalarField > alphafs(phases.size())
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
const dictionary & pimple
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
PtrList< volScalarField > rAUs(fluid.phases().size())
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
word name(const complex &)
Return a string representation of a complex.
dimensionedScalar pos(const dimensionedScalar &ds)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))