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());
8 
9 volVectorField HbyA("HbyA", U);
10 HbyA = rAU*UEqn().H();
11 
12 if (pimple.nCorrPISO() <= 1)
13 {
14  UEqn.clear();
15 }
16 
17 if (pimple.transonic())
18 {
20  (
21  "phid",
23  *(
24  (fvc::interpolate(HbyA) & mesh.Sf())
26  )
27  );
28 
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  fvOptions(psi, p, rho.name())
41  );
42 
43  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
44 
45  if (pimple.finalNonOrthogonalIter())
46  {
47  phi == pEqn.flux();
48  }
49  }
50 }
51 else
52 {
54  (
55  "phiHbyA",
56  (fvc::interpolate(rho*HbyA) & mesh.Sf())
58  );
59 
61  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
62 
63  while (pimple.correctNonOrthogonal())
64  {
65  // Pressure corrector
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 {
105  surfaceVectorField n(mesh.Sf()/mesh.magSf());
106  rhoUf += n*(fvc::absolute(phi, rho, U)/mesh.magSf() - (n & rhoUf));
107 }
108 
109 if (thermo.dpdt())
110 {
111  dpdt = fvc::ddt(p);
112 
113  if (mesh.moving())
114  {
115  dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
116  }
117 }
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)))
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...
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::fvc::meshPhi
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:33
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
phi
phi
Definition: pEqn.H:20
A
simpleMatrix< scalar > A(Nc)
U
U
Definition: pEqn.H:46
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())
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
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
n
surfaceVectorField n(mesh.Sf()/mesh.magSf())
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:15
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
HbyA
HbyA
Definition: pEqn.H:4
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
makeRelative
MRF makeRelative(phiHbyA)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
rhoMax
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
phiHbyA
phiHbyA
Definition: pEqn.H:21
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
rhoUf
rhoUf
Definition: pEqn.H:90
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16