pcEqn.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());
7 volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
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)
25  )
26  );
27 
28  MRF.makeRelative(fvc::interpolate(psi), phid);
29 
31  (
32  "phic",
33  fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf()
34  );
35 
36  HbyA -= (rAU - rAtU)*fvc::grad(p);
37 
38  volScalarField rhorAtU("rhorAtU", rho*rAtU);
39 
40  while (pimple.correctNonOrthogonal())
41  {
42  fvScalarMatrix pEqn
43  (
44  fvm::ddt(psi, p)
45  + fvm::div(phid, p)
46  + fvc::div(phic)
47  - fvm::laplacian(rhorAtU, p)
48  ==
49  fvOptions(psi, p, rho.name())
50  );
51 
52  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
53 
54  if (pimple.finalNonOrthogonalIter())
55  {
56  phi == phic + pEqn.flux();
57  }
58  }
59 }
60 else
61 {
63  (
64  "phiHbyA",
65  (
66  fvc::flux(rho*HbyA)
68  )
69  );
70 
71  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
72 
73  phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
74  HbyA -= (rAU - rAtU)*fvc::grad(p);
75 
76  volScalarField rhorAtU("rhorAtU", rho*rAtU);
77 
78  // Update the pressure BCs to ensure flux consistency
79  constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
80 
81  while (pimple.correctNonOrthogonal())
82  {
83  fvScalarMatrix pEqn
84  (
85  fvm::ddt(psi, p)
86  + fvc::div(phiHbyA)
88  ==
89  fvOptions(psi, p, rho.name())
90  );
91 
92  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
93 
94  if (pimple.finalNonOrthogonalIter())
95  {
96  phi = phiHbyA + pEqn.flux();
97  }
98  }
99 }
100 
101 #include "rhoEqn.H"
102 #include "compressibleContinuityErrs.H"
103 
104 // Explicitly relax pressure for momentum corrector
105 p.relax();
106 
107 U = HbyA - rAtU*fvc::grad(p);
108 U.correctBoundaryConditions();
109 fvOptions.correct(U);
110 K = 0.5*magSqr(U);
111 
112 if (thermo.dpdt())
113 {
114  dpdt = fvc::ddt(p);
115 }
116 
117 // Recalculate density from the relaxed pressure
118 rho = thermo.rho();
119 rho = max(rho, rhoMin);
120 rho = min(rho, rhoMax);
121 
122 if (!pimple.transonic())
123 {
124  rho.relax();
125 }
126 
127 Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;
surfaceScalarField & phi
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
tmp< surfaceScalarField > interpolate(const RhoType &rho)
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:74
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:73
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
surfaceScalarField phid("phid", fvc::interpolate(psi)*(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
dimensioned< scalar > magSqr(const dimensioned< Type > &)
U
Definition: pcEqn.H:107
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField rhorAtU("rhorAtU", rho *rAtU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
volScalarField rAtU(1.0/(1.0/rAU-UEqn.H1()))
messageStream Info
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