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)
24  )
25  );
26 
27  MRF.makeRelative(fvc::interpolate(psi), phid);
28 
29  while (pimple.correctNonOrthogonal())
30  {
31  fvScalarMatrix pEqn
32  (
33  fvm::ddt(psi, p)
34  + fvm::div(phid, p)
36  ==
37  fvOptions(psi, p, rho.name())
38  );
39 
40  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
41 
42  if (pimple.finalNonOrthogonalIter())
43  {
44  phi == pEqn.flux();
45  }
46  }
47 }
48 else
49 {
51  (
52  "phiHbyA",
53  (
54  fvc::flux(rho*HbyA)
56  )
57  );
58 
59  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
60 
61  // Update the pressure BCs to ensure flux consistency
62  constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF);
63 
64  while (pimple.correctNonOrthogonal())
65  {
66  fvScalarMatrix pEqn
67  (
68  fvm::ddt(psi, p)
69  + fvc::div(phiHbyA)
71  ==
72  fvOptions(psi, p, rho.name())
73  );
74 
75  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
76 
77  if (pimple.finalNonOrthogonalIter())
78  {
79  phi = phiHbyA + pEqn.flux();
80  }
81  }
82 }
83 
84 #include "rhoEqn.H"
85 #include "compressibleContinuityErrs.H"
86 
87 // Explicitly relax pressure for momentum corrector
88 p.relax();
89 
90 // Recalculate density from the relaxed pressure
91 rho = thermo.rho();
92 rho = max(rho, rhoMin);
93 rho = min(rho, rhoMax);
94 rho.relax();
95 Info<< "rho max/min : " << max(rho).value()
96  << " " << min(rho).value() << endl;
97 
98 U = HbyA - rAU*fvc::grad(p);
99 U.correctBoundaryConditions();
100 fvOptions.correct(U);
101 K = 0.5*magSqr(U);
102 
103 if (thermo.dpdt())
104 {
105  dpdt = fvc::ddt(p);
106 }
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
U
Definition: pEqn.H:83
p
Definition: pEqn.H:50
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
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
fv::options & fvOptions
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
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
surfaceScalarField phid("phid", fvc::interpolate(psi)*(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
messageStream Info
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
rho
Definition: pEqn.H:1
volScalarField rAU(1.0/UEqn.A())