pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
5 
7 
9  (
10  "phiHbyA",
12  + MRF.zeroFilter(rAUf*fvc::ddtCorr(U, phi))
13  + phig
14  );
15 
16  MRF.makeRelative(phiHbyA);
17 
18  // Update the pressure BCs to ensure flux consistency
20 
21  while (pimple.correctNonOrthogonal())
22  {
23  fvScalarMatrix p_rghEqn
24  (
26  );
27 
28  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
29 
30  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
31 
32  if (pimple.finalNonOrthogonalIter())
33  {
34  // Calculate the conservative fluxes
35  phi = phiHbyA - p_rghEqn.flux();
36 
37  // Explicitly relax pressure for momentum corrector
38  p_rgh.relax();
39 
40  // Correct the momentum source with the pressure gradient flux
41  // calculated from the relaxed pressure
42  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rAUf);
43  U.correctBoundaryConditions();
44  fvOptions.correct(U);
45  }
46  }
47 
48  #include "continuityErrs.H"
49 
50  p = p_rgh + rhok*gh;
51 
52  if (p_rgh.needReference())
53  {
55  (
56  "p",
57  p.dimensions(),
59  );
60  p_rgh = p - rhok*gh;
61  }
62 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
pimpleNoLoopControl & pimple
p
Definition: pEqn.H:50
volScalarField rAU(1.0/UEqn.A())
IOMRFZoneList & MRF
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
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:136
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
p_rgh
Definition: pEqn.H:152
const surfaceScalarField & ghf
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volVectorField & HbyA
Definition: pEqn.H:13
label pRefCell
Definition: createFields.H:106
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
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
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
rhok
Definition: TEqn.H:27
scalar pRefValue
Definition: createFields.H:107