pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
4 
5  volVectorField HbyA("HbyA", U);
6  HbyA = rAU*UEqn.H();
7 
9  (
10  "phiHbyA",
11  (fvc::interpolate(HbyA) & mesh.Sf())
13  );
15 
17  (
18  (
19  interface.surfaceTensionForce()
21  )*rAUf*mesh.magSf()
22  );
23 
25 
26  // Update the fixedFluxPressure BCs to ensure flux consistency
28  (
29  p_rgh.boundaryField(),
30  (
31  phiHbyA.boundaryField()
32  - (mesh.Sf().boundaryField() & U.boundaryField())
33  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
34  );
35 
36  Pair<tmp<volScalarField> > vDotP = mixture->vDotP();
37  const volScalarField& vDotcP = vDotP[0]();
38  const volScalarField& vDotvP = vDotP[1]();
39 
40  while (pimple.correctNonOrthogonal())
41  {
42  fvScalarMatrix p_rghEqn
43  (
45  - (vDotvP - vDotcP)*(pSat - rho*gh) + fvm::Sp(vDotvP - vDotcP, p_rgh)
46  );
47 
48  p_rghEqn.setReference(pRefCell, pRefValue);
49 
50  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
51 
52  if (pimple.finalNonOrthogonalIter())
53  {
54  phi = phiHbyA + p_rghEqn.flux();
55 
56  U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf);
57  U.correctBoundaryConditions();
58  fvOptions.correct(U);
59  }
60  }
61 
62  p == p_rgh + rho*gh;
63 
64  if (p_rgh.needReference())
65  {
67  (
68  "p",
69  p.dimensions(),
71  );
72  p_rgh = p - rho*gh;
73  }
74 }
phi
Definition: pEqn.H:20
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryField(),( phiHbyA.boundaryField() -MRF.relative(mesh.Sf().boundaryField()&U.boundaryField()) *rho.boundaryField() )/(mesh.magSf().boundaryField()*rhorAUf.boundaryField()))
const dimensionedScalar & pSat
Definition: createFields.H:41
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
volScalarField & p_rgh
tmp< surfaceScalarField > interpolate(const RhoType &rho)
dynamicFvMesh & mesh
surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
const dictionary & pimple
p
Definition: pEqn.H:59
fv::IOoptionList & fvOptions
const scalar pRefValue
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fvVectorMatrix UEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
phiHbyA
Definition: pEqn.H:21
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
adjustPhi(phiHbyA, U, p_rgh)
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:154
volScalarField rAU(1.0/UEqn.A())
const label pRefCell
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
rho
Definition: pEqn.H:1
interfaceProperties interface(alpha1, U, mixture())
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
U
Definition: pEqn.H:82
const surfaceScalarField & ghf
HbyA
Definition: pEqn.H:7
const volScalarField & gh