pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
4 
6 
8  (
9  "phiHbyA",
12  );
13  MRF.makeRelative(phiHbyA);
15 
17  (
18  (
19  mixture.surfaceTensionForce()
21  )*rAUf*mesh.magSf()
22  );
23 
25 
26  // Update the pressure BCs to ensure flux consistency
28 
29  while (pimple.correctNonOrthogonal())
30  {
31  fvScalarMatrix p_rghEqn
32  (
34  );
35 
36  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
37 
38  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
39 
40  if (pimple.finalNonOrthogonalIter())
41  {
42  phi = phiHbyA - p_rghEqn.flux();
43 
44  p_rgh.relax();
45 
46  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rAUf);
47  U.correctBoundaryConditions();
48  fvOptions.correct(U);
49  }
50  }
51 
52  #include "continuityErrs.H"
53 
54  p == p_rgh + rho*gh;
55 
56  if (p_rgh.needReference())
57  {
59  (
60  "p",
61  p.dimensions(),
63  );
64  p_rgh = p - rho*gh;
65  }
66 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
p
Definition: pEqn.H:50
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
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
const dictionary & pimple
tmp< surfaceScalarField > interpolate(const RhoType &rho)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:130
fv::options & fvOptions
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
IOMRFZoneList & MRF
surfaceScalarField phig("phig",-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
p_rgh
Definition: pEqn.H:120
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
adjustPhi(phiHbyA, U, p_rgh)
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
const scalar pRefValue
const label pRefCell
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
phi
Definition: pEqn.H:18
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
rho
Definition: pEqn.H:1
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
volScalarField rAU(1.0/UEqn.A())