pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU(1.0/UEqn.A());
5  MRF.makeRelative(phiHbyA);
6  adjustPhi(phiHbyA, U, p);
7 
8  tmp<volScalarField> rAtU(rAU);
9 
10  if (simple.consistent())
11  {
12  rAtU = 1.0/(1.0/rAU - UEqn.H1());
13  phiHbyA +=
14  fvc::interpolate(rAtU() - rAU)*fvc::snGrad(p)*mesh.magSf();
15  HbyA -= (rAU - rAtU())*fvc::grad(p);
16  }
17 
18  tUEqn.clear();
19 
20  // Update the pressure BCs to ensure flux consistency
22 
23  // Non-orthogonal pressure corrector loop
24  while (simple.correctNonOrthogonal())
25  {
26  fvScalarMatrix pEqn
27  (
29  );
30 
31  pEqn.setReference(pRefCell, pRefValue);
32 
33  pEqn.solve();
34 
35  if (simple.finalNonOrthogonalIter())
36  {
37  phi = phiHbyA - pEqn.flux();
38  }
39  }
40 
41  #include "continuityErrs.H"
42 
43  // Explicitly relax pressure for momentum corrector
44  p.relax();
45 
46  // Momentum corrector
47  U = HbyA - rAtU()*fvc::grad(p);
48  U.correctBoundaryConditions();
49  fvOptions.correct(U);
50 }
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
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
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
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
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
adjustPhi(phiHbyA, U, p_rgh)
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.
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
label pRefCell
Definition: createFields.H:106
tmp< volScalarField > rAtU(rAU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
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
scalar pRefValue
Definition: createFields.H:107
simpleControl simple(mesh)