pEqn.H
Go to the documentation of this file.
1 if (!mesh.steady() && !pimple.simpleRho())
2 {
3  rho = thermo.rho();
4 }
5 
6 // Thermodynamic density needs to be updated by psi*d(p) after the
7 // pressure solution
8 const volScalarField psip0(psi*p);
9 
10 const volScalarField rAU("rAU", 1.0/UEqn.A());
11 const surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
13 
14 if (pimple.nCorrPiso() <= 1)
15 {
16  tUEqn.clear();
17 }
18 
20 (
21  "phiHbyA",
23  + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi, rhoUf))
24 );
25 
26 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
27 
28 bool adjustMass = mesh.steady() && adjustPhi(phiHbyA, U, p_rgh);
29 
30 surfaceScalarField phig(-rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf());
31 
33 
34 // Update the pressure BCs to ensure flux consistency
35 constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF);
36 
38 
40 
41 if (pimple.transonic())
42 {
44  (
45  "phid",
47  );
48 
49  phiHbyA -= fvc::interpolate(psi*p_rgh)*phiHbyA/fvc::interpolate(rho);
50 
51  fvScalarMatrix p_rghDDtEqn
52  (
54  + fvc::div(phiHbyA) + fvm::div(phid, p_rgh)
55  ==
56  fvModels.source(psi, p_rgh, rho.name())
57  );
58 
59  while (pimple.correctNonOrthogonal())
60  {
61  p_rghEqn = p_rghDDtEqn - fvm::laplacian(rhorAUf, p_rgh);
62 
63  // Relax the pressure equation to ensure diagonal-dominance
64  p_rghEqn.relax();
65 
66  p_rghEqn.setReference
67  (
68  pressureReference.refCell(),
69  pressureReference.refValue()
70  );
71 
72  p_rghEqn.solve();
73  }
74 }
75 else
76 {
77  fvScalarMatrix p_rghDDtEqn
78  (
80  + fvc::div(phiHbyA)
81  ==
82  fvModels.source(psi, p_rgh, rho.name())
83  );
84 
85  while (pimple.correctNonOrthogonal())
86  {
87  p_rghEqn = p_rghDDtEqn - fvm::laplacian(rhorAUf, p_rgh);
88 
89  p_rghEqn.setReference
90  (
91  pressureReference.refCell(),
92  pressureReference.refValue()
93  );
94 
95  p_rghEqn.solve();
96  }
97 }
98 
99 phi = phiHbyA + p_rghEqn.flux();
100 
101 p = p_rgh + rho*gh + pRef;
102 
103 if (mesh.steady())
104 {
106 }
107 else
108 {
109  const bool constrained = fvConstraints.constrain(p);
110 
111  // Thermodynamic density update
112  thermo.correctRho(psi*p - psip0);
113 
114  if (constrained)
115  {
116  rho = thermo.rho();
117  }
118 
119  #include "rhoEqn.H"
121 }
122 
123 // Explicitly relax pressure for momentum corrector
124 p_rgh.relax();
125 
126 // Correct the momentum source with the pressure gradient flux
127 // calculated from the relaxed pressure
128 U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf);
129 U.correctBoundaryConditions();
131 K = 0.5*magSqr(U);
132 
133 if (mesh.steady())
134 {
136 }
137 
138 // For steady closed-volume compressible cases adjust the pressure level
139 // to obey overall mass continuity
140 if (adjustMass && !thermo.incompressible())
141 {
142  p += (initialMass - fvc::domainIntegrate(thermo.rho()))
144  p_rgh = p - rho*gh - pRef;
145  p.correctBoundaryConditions();
146 }
147 
148 if (mesh.steady() || pimple.simpleRho() || adjustMass)
149 {
150  rho = thermo.rho();
151 }
152 
153 if ((mesh.steady() || pimple.simpleRho()) && !pimple.transonic())
154 {
155  rho.relax();
156 }
157 
158 // Correct rhoUf if the mesh is moving
160 
161 if (thermo.dpdt())
162 {
163  dpdt = fvc::ddt(p);
164 
165  if (mesh.moving())
166  {
167  dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
168  }
169 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
rAU
Definition: pEqn.H:1
pressureReference & pressureReference
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
dimensionedScalar initialMass
Definition: createFields.H:68
p_rgh
Definition: pEqn.H:159
IOMRFZoneList & MRF
rhoUf
Definition: pEqn.H:94
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
fvScalarMatrix p_rghEqn(p_rgh, dimMass/dimTime)
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
Calculates and prints the continuity errors.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
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
CGAL::Exact_predicates_exact_constructions_kernel K
phi
Definition: pEqn.H:99
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
volScalarField & dpdt
const dimensionSet dimTime
dynamicFvMesh & mesh
rho
Definition: pEqn.H:1
Foam::fvConstraints & fvConstraints
bool adjustMass
Definition: pEqn.H:47
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
const volScalarField psip0(psi *p)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
MRF makeRelative(phiHbyA)
const surfaceScalarField & ghf
scalar pRef
Definition: createFields.H:6
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
const dimensionSet dimMass
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 phig(-rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
phiHbyA
Definition: pEqn.H:32
Foam::fvModels & fvModels
adjustPhi(phiHbyA, Urel, p)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:279
Calculates and prints the continuity errors.
const volScalarField & gh
const volScalarField & psi
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
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
fvVectorMatrix & UEqn
Definition: UEqn.H:13