pEqn.H
Go to the documentation of this file.
1 rAU = 1.0/UEqn.A();
4 
5 if (pimple.nCorrPiso() <= 1)
6 {
7  tUEqn.clear();
8 }
9 
11 (
12  "phiHbyA",
14  + MRF.zeroFilter(rAUf*fvc::ddtCorr(U, phi, Uf))
15 );
16 
17 MRF.makeRelative(phiHbyA);
18 
19 if (p_gh.needReference())
20 {
22  adjustPhi(phiHbyA, U, p_gh);
24 }
25 
26 // Update the pressure BCs to ensure flux consistency
28 
29 // Non-orthogonal pressure corrector loop
30 while (pimple.correctNonOrthogonal())
31 {
32  fvScalarMatrix p_ghEqn
33  (
35  );
36 
37  p_ghEqn.setReference(p_ghRefCell, p_ghRefValue);
38 
39  p_ghEqn.solve();
40 
41  if (pimple.finalNonOrthogonalIter())
42  {
43  phi = phiHbyA - p_ghEqn.flux();
44 
45  // Explicitly relax pressure for momentum corrector
46  p_gh.relax();
47 
48  U = HbyA - rAU*fvc::grad(p_gh);
49  U.correctBoundaryConditions();
50  fvOptions.correct(U);
51  }
52 }
53 
54 #include "continuityErrs.H"
55 
56 // Correct Uf if the mesh is moving
58 
59 // Make the fluxes relative to the mesh motion
61 
62 p = p_gh + (g & mesh.C());
rAU
Definition: pEqn.H:1
phiHbyA
Definition: pEqn.H:22
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
pimpleNoLoopControl & pimple
p
Definition: pEqn.H:50
IOMRFZoneList & MRF
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
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:170
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
dynamicFvMesh & mesh
void correctUf(autoPtr< surfaceVectorField > &Uf, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:223
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
autoPtr< surfaceVectorField > Uf
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.
adjustPhi(phiHbyA, Urel, p)
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:114
phi
Definition: pEqn.H:18
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
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
const dimensionedVector & g
MRF makeRelative(fvc::interpolate(rho), phiHbyA)