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()));
8 volVectorField HbyA("HbyA", U);
9 HbyA = rAU*UEqn().H();
10 
11 if (pimple.nCorrPISO() <= 1)
12 {
13  UEqn.clear();
14 }
15 
16 if (pimple.transonic())
17 {
19  (
20  "phid",
22  *(
23  (fvc::interpolate(HbyA) & mesh.Sf())
26  )
27  );
28 
29  MRF.makeRelative(fvc::interpolate(psi), phid);
30 
32  (
33  "phic",
35  );
36 
37  HbyA -= (rAU - rAtU)*fvc::grad(p);
38 
39  volScalarField rhorAtU("rhorAtU", rho*rAtU);
40 
41  while (pimple.correctNonOrthogonal())
42  {
43  fvScalarMatrix pEqn
44  (
45  fvm::ddt(psi, p)
46  + fvm::div(phid, p)
47  + fvc::div(phic)
49  ==
50  fvOptions(psi, p, rho.name())
51  );
52 
53  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
54 
55  if (pimple.finalNonOrthogonalIter())
56  {
57  phi == phic + pEqn.flux();
58  }
59  }
60 }
61 else
62 {
64  (
65  "phiHbyA",
66  (
67  (fvc::interpolate(rho*HbyA) & mesh.Sf())
69  )
70  );
71 
72  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
73 
75  HbyA -= (rAU - rAtU)*fvc::grad(p);
76 
77  volScalarField rhorAtU("rhorAtU", rho*rAtU);
78 
79  while (pimple.correctNonOrthogonal())
80  {
81  fvScalarMatrix pEqn
82  (
83  fvm::ddt(psi, p)
84  + fvc::div(phiHbyA)
86  ==
87  fvOptions(psi, p, rho.name())
88  );
89 
90  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
91 
92  if (pimple.finalNonOrthogonalIter())
93  {
94  phi = phiHbyA + pEqn.flux();
95  }
96  }
97 }
98 
99 #include "rhoEqn.H"
100 #include "compressibleContinuityErrs.H"
101 
102 // Explicitly relax pressure for momentum corrector
103 p.relax();
104 
106 U.correctBoundaryConditions();
107 fvOptions.correct(U);
108 K = 0.5*magSqr(U);
109 
110 if (thermo.dpdt())
111 {
112  dpdt = fvc::ddt(p);
113 }
114 
115 // Recalculate density from the relaxed pressure
116 rho = thermo.rho();
117 rho = max(rho, rhoMin);
118 rho = min(rho, rhoMax);
119 
120 if (!pimple.transonic())
121 {
122  rho.relax();
123 }
124 
125 Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
phid
surfaceScalarField phid("phid", fvc::interpolate(psi) *((mesh.Sf() &fvc::interpolate(HbyA))+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
rAU
volScalarField rAU(1.0/UEqn().A())
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
HbyA
HbyA
Definition: pcEqn.H:9
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
dpdt
volScalarField & dpdt
Definition: setRegionFluidFields.H:12
A
simpleMatrix< scalar > A(Nc)
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
rhoMin
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
phiHbyA
phiHbyA
Definition: pcEqn.H:74
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
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
rho
rho
Definition: pcEqn.H:1
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
U
U
Definition: pcEqn.H:105
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:15
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
phic
phic
Definition: alphaEqnsSubCycle.H:5
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
rhorAtU
volScalarField rhorAtU("rhorAtU", rho *rAtU)
rhoMax
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
rAtU
volScalarField rAtU(1.0/(1.0/rAU - UEqn().H1()))
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16