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