Go to the documentation of this file.
27 PtrList<surfaceScalarField>
alphafs(
fluid.phases().size());
28 PtrList<volVectorField>
HbyAs(
fluid.phases().size());
30 PtrList<volScalarField>
rAUs(
fluid.phases().size());
36 phaseModel& phase = iter();
38 MRF.makeAbsolute(phase.phi().oldTime());
39 MRF.makeAbsolute(phase.phi());
67 phaseModel& phase = iter();
82 zeroGradientFvPatchScalarField::typeName
84 dragCoeffi.correctBoundaryConditions();
104 MRF.makeRelative(phase.phi().oldTime());
105 MRF.makeRelative(phase.phi());
115 multiphaseSystem::dragModelTable::const_iterator dmIter =
116 fluid.dragModels().begin();
117 multiphaseSystem::dragCoeffFields::const_iterator dcIter =
126 const phaseModel *phase2Ptr = NULL;
128 if (&phase == &dmIter()->
phase1())
130 phase2Ptr = &dmIter()->phase2();
132 else if (&phase == &dmIter()->
phase2())
134 phase2Ptr = &dmIter()->phase1();
182 phaseModel& phase = iter();
190 surfaceScalarField::GeometricBoundaryField
phib(
phi.boundaryField());
195 phaseModel& phase = iter();
199 *(
mesh.Sf().boundaryField() & phase.U().boundaryField());
206 p_rgh.boundaryField(),
209 )/(
mesh.magSf().boundaryField()*
rAUf.boundaryField())
233 if (
pimple.finalNonOrthogonalIter())
241 phaseModel& phase = iter();
249 *mSfGradp/phase.rho();
265 mSfGradp = pEqnIncomp.flux()/
rAUf;
272 phaseModel& phase = iter();
289 phase.U().correctBoundaryConditions();
300 #include "continuityErrs.H"
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 & ghf
const dimensionSet dimVelocity
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
const word & name() const
Return const reference to name.
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()))
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
dimensioned< scalar > mag(const dimensioned< Type > &)
simpleMatrix< scalar > A(Nc)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
PtrList< surfaceScalarField > rAlphaAUfs(fluid.phases().size())
const dimensionSet dimArea(sqr(dimLength))
fvMatrix< scalar > fvScalarMatrix
PtrList< volVectorField > HbyAs(fluid.phases().size())
PtrList< fvVectorMatrix > UEqns(fluid.phases().size())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
PtrList< surfaceScalarField > phiHbyAs(fluid.phases().size())
const volScalarField & gh
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
PtrList< surfaceScalarField > alphafs(phases.size())
surfaceScalarField ghSnGradRho(ghf *fvc::snGrad(rho) *mesh.magSf())
const dictionary & pimple
PtrList< volScalarField > rAUs(fluid.phases().size())
forAllIter(PtrDictionary< phaseModel >, fluid.phases(), iter)
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))