pcEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 volScalarField rAU(1.0/UEqn.A());
4 volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
6 
7 if (pimple.nCorrPISO() <= 1)
8 {
9  tUEqn.clear();
10 }
11 
12 if (pimple.transonic())
13 {
15  (
16  "phid",
18  *(
19  fvc::flux(HbyA)
22  )
23  );
24 
25  MRF.makeRelative(fvc::interpolate(psi), phid);
26 
28  (
29  "phic",
30  fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
31  );
32 
33  HbyA -= (rAU - rAtU)*fvc::grad(p);
34 
35  volScalarField rhorAtU("rhorAtU", rho*rAtU);
36 
37  while (pimple.correctNonOrthogonal())
38  {
39  fvScalarMatrix pEqn
40  (
41  fvm::ddt(psi, p)
42  + fvm::div(phid, p)
43  + fvc::div(phic)
44  - fvm::laplacian(rhorAtU, p)
45  ==
46  fvOptions(psi, p, rho.name())
47  );
48 
49  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
50 
51  if (pimple.finalNonOrthogonalIter())
52  {
53  phi == phic + pEqn.flux();
54  }
55  }
56 }
57 else
58 {
60  (
61  "phiHbyA",
62  (
63  fvc::flux(rho*HbyA)
65  )
66  );
67 
68  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
69 
70  phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
71  HbyA -= (rAU - rAtU)*fvc::grad(p);
72 
73  volScalarField rhorAtU("rhorAtU", rho*rAtU);
74 
75  // Update the pressure BCs to ensure flux consistency
76  constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
77 
78  while (pimple.correctNonOrthogonal())
79  {
80  fvScalarMatrix pEqn
81  (
82  fvm::ddt(psi, p)
83  + fvc::div(phiHbyA)
85  ==
86  fvOptions(psi, p, rho.name())
87  );
88 
89  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
90 
91  if (pimple.finalNonOrthogonalIter())
92  {
93  phi = phiHbyA + pEqn.flux();
94  }
95  }
96 }
97 
98 #include "rhoEqn.H"
99 #include "compressibleContinuityErrs.H"
100 
101 // Explicitly relax pressure for momentum corrector
102 p.relax();
103 
104 U = HbyA - rAtU*fvc::grad(p);
105 U.correctBoundaryConditions();
106 fvOptions.correct(U);
107 K = 0.5*magSqr(U);
108 
109 if (pressureControl.limit(p))
110 {
111  p.correctBoundaryConditions();
112  rho = thermo.rho();
113 }
114 
115 if (thermo.dpdt())
116 {
117  dpdt = fvc::ddt(p);
118 }
surfaceScalarField & phi
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
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF)
volScalarField rAU(1.0/UEqn.A())
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
const dictionary & pimple
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
CGAL::Exact_predicates_exact_constructions_kernel K
HbyA
Definition: pcEqn.H:71
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:70
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:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
IOMRFZoneList & MRF
volScalarField & dpdt
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
rho
Definition: pcEqn.H:1
dimensioned< scalar > magSqr(const dimensioned< Type > &)
U
Definition: pcEqn.H:104
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 phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
volScalarField rhorAtU("rhorAtU", rho *rAtU)
surfaceScalarField phic(mixture.cAlpha() *mag(phi/mesh.magSf()))
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & psi
volScalarField & p
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
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))