Go to the documentation of this file.
28 PtrList<surfaceScalarField>
AFfs(
fluid.AFfs());
30 PtrList<surfaceScalarField>
rAUfs;
35 IOobject::groupName(
"rAUf",
phase1.name()),
48 IOobject::groupName(
"rAUf",
phase2.name()),
84 IOobject::groupName(
"alpharAUf",
phase1.name()),
90 IOobject::groupName(
"alpharAUf",
phase2.name()),
125 IOobject::groupName(
"phiHbyA",
phase1.name()),
139 IOobject::groupName(
"phiHbyA",
phase2.name()),
152 PtrList<surfaceScalarField> phiKdPhifs(
fluid.phiKdPhifs(
rAUfs));
177 p_rgh.boundaryFieldRef(),
184 )/(
mesh.magSf().boundaryField()*
rAUf.boundaryField())
188 tmp<fvScalarMatrix> pEqnComp1, pEqnComp2;
190 if (
phase1.compressible())
196 IOobject::groupName(
"phid",
phase1.name()),
214 pEqnComp1.ref().relax();
227 if (
phase2.compressible())
233 IOobject::groupName(
"phid",
phase2.name()),
251 pEqnComp2.ref().relax();
269 if (pEqnComp1.valid())
271 pEqnComp1.ref() -= (optEqn1 &
rho1)/
rho1;
281 if (pEqnComp2.valid())
283 pEqnComp2.ref() -= (optEqn2 &
rho2)/
rho2;
295 PtrList<volScalarField> dmdts(
fluid.dmdts());
298 if (pEqnComp1.valid())
300 pEqnComp1.ref() -= dmdts[0]/
rho1;
309 if (pEqnComp2.valid())
311 pEqnComp2.ref() -= dmdts[1]/
rho2;
323 while (
pimple.correctNonOrthogonal())
334 if (pEqnComp1.valid())
339 if (pEqnComp2.valid())
351 if (
pimple.finalNonOrthogonalIter())
377 U1.correctBoundaryConditions();
381 U2.correctBoundaryConditions();
385 if (pEqnComp1.valid())
389 if (pEqnComp2.valid())
409 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)
const surfaceScalarField Kdf("Kdf", fluid.Kdf())
const surfaceScalarField & ghf
const dimensionedScalar & pMin
const volScalarField & alpha2
const surfaceScalarField alphaf2("alphaf2", scalar(1) - alphaf1)
surfaceScalarField & phi1
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
const volScalarField & alpha1
const volScalarField & gh
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
const volScalarField & psi2
const surfaceScalarField & rAUf1
surfaceScalarField & phi2
surfaceScalarField phir(fvc::flux(UdmModel.Udm()))
fvMatrix< scalar > fvScalarMatrix
const surfaceScalarField alphaf1("alphaf1", fvc::interpolate(alpha1))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label max(const labelHashSet &set, label maxValue=labelMin)
PtrList< surfaceScalarField > AFfs(fluid.AFfs())
const uniformDimensionedVectorField & g
const surfaceScalarField & alphaPhi1
const volScalarField & psi1
const surfaceScalarField & phiFf2
surfaceScalarField alphaRhof10("alphaRhof10", fvc::interpolate(max(alpha1.oldTime(), phase1.residualAlpha()) *rho1.oldTime()))
volScalarField p_rgh_0(p_rgh)
PtrList< surfaceScalarField > phiFfs(fluid.phiFfs(rAUfs))
dimensionedSymmTensor sqr(const dimensionedVector &dv)
surfaceScalarField alphaRhof20("alphaRhof20", fvc::interpolate(max(alpha2.oldTime(), phase2.residualAlpha()) *rho2.oldTime()))
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryFieldRef(),(phiHbyA.boundaryField() - MRF.relative(phib))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
PtrList< surfaceScalarField > rAUfs
const surfaceScalarField alpharAUf2(fvc::interpolate(max(alpha2, phase2.residualAlpha()) *rAU2))
const surfaceScalarField & rAUf2
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
tmp< surfaceScalarField > flux(const volVectorField &vvf)
tmp< volScalarField > byDt(const volScalarField &vf)
const surfaceScalarField alpharAUf1(fvc::interpolate(max(alpha1, phase1.residualAlpha()) *rAU1))
const surfaceScalarField & alphaPhi2
const surfaceScalarField & phiFf1
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))