pEqn.H
Go to the documentation of this file.
2 bool compressible = (compressibility.value() > SMALL);
3 
4 rho = thermo.rho();
5 
6 // Thermodynamic density needs to be updated by psi*d(p) after the
7 // pressure solution
9 
10 volScalarField rAU(1.0/UEqn.A());
13 
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
30 
32 
34 (
36  + fvc::div(phiHbyA)
37  ==
38  fvOptions(psi, p_rgh, rho.name())
39 );
40 
41 while (pimple.correctNonOrthogonal())
42 {
43  fvScalarMatrix p_rghEqn
44  (
47  );
48 
49  p_rghEqn.setReference
50  (
51  pRefCell,
53  );
54 
55  p_rghEqn.solve(p_rgh.select(pimple.finalInnerIter()));
56 
57  if (pimple.finalNonOrthogonalIter())
58  {
59  // Calculate the conservative fluxes
60  phi = phiHbyA + p_rghEqn.flux();
61 
62  // Explicitly relax pressure for momentum corrector
63  p_rgh.relax();
64 
65  // Correct the momentum source with the pressure gradient flux
66  // calculated from the relaxed pressure
67  U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf);
68  U.correctBoundaryConditions();
69  fvOptions.correct(U);
70  K = 0.5*magSqr(U);
71  }
72 }
73 
74 p = p_rgh + rho*gh;
75 
76 pressureControl.limit(p);
77 
78 if (p_rgh.needReference())
79 {
80  if (!compressible)
81  {
83  (
84  "p",
85  p.dimensions(),
87  );
88  }
89  else
90  {
93  thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
94  rho = thermo.rho();
95  p_rgh = p - rho*gh;
96  p_rgh.correctBoundaryConditions();
97  }
98 }
99 else
100 {
101  thermo.correctRho(psi*p - psip0, rhoMin, rhoMax);
102 }
103 
104 #include "rhoEqn.H"
105 #include "compressibleContinuityErrs.H"
106 
107 rho = thermo.rho();
108 
109 // Correct rhoUf if the mesh is moving
111 
112 if (thermo.dpdt())
113 {
114  dpdt = fvc::ddt(p);
115 
116  if (mesh.moving())
117  {
118  dpdt -= fvc::div(fvc::meshPhi(rho, U), p);
119  }
120 }
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 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
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
pimpleControl & pimple
const scalar pRefValue
const label pRefCell
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
rhoUf
Definition: pEqn.H:89
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
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:29
const volScalarField & gh
void correctRhoUf(autoPtr< surfaceVectorField > &rhoUf, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi)
Definition: fvcMeshPhi.C:236
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))
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
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