pEqn.H
Go to the documentation of this file.
1 bool closedVolume = p_rgh.needReference();
3 bool compressible = (compressibility.value() > SMALL);
4 
5 rho = thermo.rho();
6 
7 // Thermodynamic density needs to be updated by psi*d(p) after the
8 // pressure solution
10 
11 volScalarField rAU("rAU", 1.0/UEqn.A());
14 
16 
18 (
19  "phiHbyA",
20  (
22  + MRF.zeroFilter(rhorAUf*fvc::ddtCorr(rho, U, phi))
23  )
24  + phig
25 );
26 
27 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
28 
29 // Update the pressure BCs to ensure flux consistency
31 
32 {
34  (
36  + fvc::div(phiHbyA)
37  );
38 
39  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
40  {
41  fvScalarMatrix p_rghEqn
42  (
45  );
46 
47  p_rghEqn.setReference
48  (
49  pRefCell,
51  );
52 
53  p_rghEqn.solve
54  (
55  p_rgh.select
56  (
57  (
58  oCorr == nOuterCorr-1
59  && corr == nCorr-1
60  && nonOrth == nNonOrthCorr
61  )
62  )
63  );
64 
65  if (nonOrth == nNonOrthCorr)
66  {
67  phi = phiHbyA + p_rghEqn.flux();
68 
69  p_rgh.relax();
70 
71  U = HbyA
72  + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf);
73  U.correctBoundaryConditions();
74  fvOptions.correct(U);
75  K = 0.5*magSqr(U);
76  }
77  }
78 
79  p = p_rgh + rho*gh;
80 
81 }
82 
83 pressureControl.limit(p);
84 
85 // For closed-volume cases adjust the pressure and density levels
86 // to obey overall mass continuity
88 {
89  if (!compressible)
90  {
92  (
93  "p",
94  p.dimensions(),
96  );
97  }
98  else
99  {
102  thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
103  rho = thermo.rho();
104  p_rgh = p - rho*gh;
105  p_rgh.correctBoundaryConditions();
106  }
107 }
108 else
109 {
110  thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
111 }
112 
113 #include "rhoEqn.H"
115 
116 rho = thermo.rho();
117 
118 // Update pressure time derivative if needed
119 if (thermo.dpdt())
120 {
121  dpdt = fvc::ddt(p);
122 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
rho
Definition: pEqn.H:1
dimensionedScalar initialMass
Definition: createFields.H:57
p
Definition: pEqn.H:50
volScalarField rAU(1.0/UEqn.A())
phiHbyA
Definition: pEqn.H:20
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
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:165
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:127
CGAL::Exact_predicates_exact_constructions_kernel K
const dimensionedScalar rhoMin
const int nNonOrthCorr
const pressureControl & pressureControl
fv::options & fvOptions
dimensionedScalar compressibility
Definition: pEqn.H:1
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
dynamicFvMesh & mesh
IOMRFZoneList & MRF
const int nOuterCorr
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
volScalarField & dpdt
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:28
HbyA
Definition: pEqn.H:4
p_rgh
Definition: pEqn.H:139
bool compressible
Definition: pEqn.H:2
const scalar pRefValue
const label pRefCell
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
U
Definition: pEqn.H:72
fvScalarMatrix p_rghDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p_rgh))+fvc::div(phiHbyA)==fvOptions(psi, p_rgh, rho.name()))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
mesh interpolate(rAU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
phi
Definition: pEqn.H:18
const volScalarField & psi
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
const volScalarField psip0(psi *p)
const dimensionedScalar rhoMax
bool closedVolume
Definition: pEqn.H:10