pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 rho = max(rho, rhoMin);
3 rho = min(rho, rhoMax);
4 rho.relax();
5 
6 volScalarField rAU(1.0/UEqn.A());
9 
10 if (pimple.nCorrPISO() <= 1)
11 {
12  tUEqn.clear();
13 }
14 
15 if (pimple.transonic())
16 {
18  (
19  "phid",
21  *(
22  fvc::flux(HbyA)
23  + MRF.zeroFilter
24  (
26  )
27  )
28  );
29 
30  MRF.makeRelative(fvc::interpolate(psi), phid);
31 
32  while (pimple.correctNonOrthogonal())
33  {
34  fvScalarMatrix pEqn
35  (
36  fvm::ddt(psi, p)
37  + fvm::div(phid, p)
39  ==
40  parcels.Srho()
41  + fvOptions(psi, p, rho.name())
42  );
43 
44  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
45 
46  if (pimple.finalNonOrthogonalIter())
47  {
48  phi == pEqn.flux();
49  }
50  }
51 }
52 else
53 {
55  (
56  "phiHbyA",
57  (
58  fvc::flux(rho*HbyA)
60  )
61  );
62 
63  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
64 
65  // Update the pressure BCs to ensure flux consistency
66  constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF);
67 
68  while (pimple.correctNonOrthogonal())
69  {
70  fvScalarMatrix pEqn
71  (
72  fvm::ddt(psi, p)
73  + fvc::div(phiHbyA)
75  ==
76  parcels.Srho()
77  + fvOptions(psi, p, rho.name())
78  );
79 
80  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
81 
82  if (pimple.finalNonOrthogonalIter())
83  {
84  phi = phiHbyA + pEqn.flux();
85  }
86  }
87 }
88 
89 #include "rhoEqn.H"
90 #include "compressibleContinuityErrs.H"
91 
92 // Explicitly relax pressure for momentum corrector
93 p.relax();
94 
95 // Recalculate density from the relaxed pressure
96 rho = thermo.rho();
97 rho = max(rho, rhoMin);
98 rho = min(rho, rhoMax);
99 rho.relax();
100 Info<< "rho max/min : " << max(rho).value()
101  << " " << min(rho).value() << endl;
102 
103 U = HbyA - rAU*fvc::grad(p);
104 U.correctBoundaryConditions();
105 fvOptions.correct(U);
106 K = 0.5*magSqr(U);
107 
108 if (thermo.dpdt())
109 {
110  dpdt = fvc::ddt(p);
111 }
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
rho
Definition: pEqn.H:1
pimpleNoLoopControl & pimple
p
Definition: pEqn.H:50
volScalarField rAU(1.0/UEqn.A())
IOMRFZoneList & MRF
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
phiHbyA
Definition: pEqn.H:20
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
rhoReactionThermo & thermo
Definition: createFields.H:28
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
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
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
volScalarField & dpdt
dynamicFvMesh & mesh
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
fvVectorMatrix & UEqn
Definition: UEqn.H:13
messageStream Info
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
phi
Definition: pEqn.H:18
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
Solve the continuity for density.
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))