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
32  (
33  pressureReference.refCell(),
34  pressureReference.refValue()
35  );
36 
37  pEqn.solve();
38 
39  if (simple.finalNonOrthogonalIter())
40  {
41  phi = phiHbyA - pEqn.flux();
42  }
43  }
44 
45  #include "continuityErrs.H"
46 
47  // Explicitly relax pressure for momentum corrector
48  p.relax();
49 
50  // Momentum corrector
51  U = HbyA - rAtU()*fvc::grad(p);
52  U.correctBoundaryConditions();
54 }
tmp< volScalarField > rAtU(pimple.consistent() ? volScalarField::New("rAtU", 1.0/(1.0/rAU - UEqn.H1())) :tmp< volScalarField >(nullptr))
rAU
Definition: pEqn.H:1
pressureReference & pressureReference
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
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
phi
Definition: pEqn.H:99
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
Foam::fvConstraints & fvConstraints
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
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.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
phiHbyA
Definition: pEqn.H:32
adjustPhi(phiHbyA, Urel, p)
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
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
p
Definition: pEqn.H:101
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
simpleControl simple(mesh)
fvVectorMatrix & UEqn
Definition: UEqn.H:13