All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
pEqn.H
Go to the documentation of this file.
1 if (!pimple.simpleRho())
2 {
3  rho = thermo.rho();
4 }
5 
6 // Thermodynamic density needs to be updated by psi*d(p) after the
7 // pressure solution
8 const volScalarField psip0(psi*p);
9 
10 volScalarField rAU(1.0/UEqn.A());
13 
14 surfaceScalarField phig(-rhorAUf*ghf*fvc::snGrad(rho)*mesh.magSf());
15 
17 (
18  "phiHbyA",
19  (
21  + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi))
22  )
23  + phig
24 );
25 
26 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
27 
28 // Update the pressure BCs to ensure flux consistency
29 constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF);
30 
32 (
34  + fvc::div(phiHbyA)
35  ==
36  parcels.Srho()
37  + surfaceFilm.Srho()
38  + fvOptions(psi, p_rgh, rho.name())
39 );
40 
41 while (pimple.correctNonOrthogonal())
42 {
44  (
45  p_rghDDtEqn
46  - fvm::laplacian(rhorAUf, p_rgh)
47  );
48 
49  p_rghEqn.solve();
50 
51  if (pimple.finalNonOrthogonalIter())
52  {
53  phi = phiHbyA + p_rghEqn.flux();
54 
55  // Explicitly relax pressure for momentum corrector
56  p_rgh.relax();
57 
58  U = HbyA + rAU*fvc::reconstruct((p_rghEqn.flux() + phig)/rhorAUf);
59  U.correctBoundaryConditions();
60  fvOptions.correct(U);
61  K = 0.5*magSqr(U);
62  }
63 }
64 
65 p = p_rgh + rho*gh;
66 
67 bool limitedp = pressureControl.limit(p);
68 
69 if (limitedp)
70 {
71  p_rgh = p - rho*gh;
72 }
73 
74 // Thermodynamic density update
75 thermo.correctRho(psi*p - psip0);
76 
77 if (limitedp)
78 {
79  rho = thermo.rho();
80 }
81 
82 #include "rhoEqn.H"
84 
85 if (pimple.simpleRho())
86 {
87  rho = thermo.rho();
88 }
89 
90 if (thermo.dpdt())
91 {
92  dpdt = fvc::ddt(p);
93 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
rAU
Definition: pEqn.H:1
fvScalarMatrix p_rghDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p_rgh))+fvc::div(phiHbyA)==parcels.Srho()+surfaceFilm.Srho()+fvOptions(psi, p_rgh, rho.name()))
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
rho
Definition: pEqn.H:1
pimpleNoLoopControl & pimple
Solve the continuity for density.
p
Definition: pEqn.H:50
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
fvScalarMatrix p_rghEqn(p_rgh, dimMass/dimTime)
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
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
p_rgh
Definition: pEqn.H:140
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
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
phi
Definition: pEqn.H:104
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:57
volScalarField & dpdt
dynamicFvMesh & mesh
pressureControl & pressureControl
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
const volScalarField psip0(psi *p)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
regionModels::surfaceFilmModel & surfaceFilm
const surfaceScalarField & ghf
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.
phiHbyA
Definition: pEqn.H:32
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculates and prints the continuity errors.
const volScalarField & gh
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
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
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45