pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 const volScalarField rAU("rAU", 1.0/UEqn.A());
6 
7 tUEqn.clear();
8 
10 (
11  "phiHbyA",
13 );
14 
15 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
16 
17 bool closedVolume = simple.transonic() ? false : adjustPhi(phiHbyA, U, p_rgh);
18 
20 
22 
23 // Update the pressure BCs to ensure flux consistency
25 
27 
28 if (simple.transonic())
29 {
31  (
32  "phid",
34  );
35 
36  phiHbyA -= fvc::interpolate(psi*p_rgh)*phiHbyA/fvc::interpolate(rho);
37 
38  while (simple.correctNonOrthogonal())
39  {
40  p_rghEqn =
41  fvc::div(phiHbyA) + fvm::div(phid, p_rgh)
42  - fvm::laplacian(rhorAUf, p_rgh)
43  ==
44  fvOptions(psi, p_rgh, rho.name());
45 
46  // Relax the pressure equation to ensure diagonal-dominance
47  p_rghEqn.relax();
48 
49  p_rghEqn.setReference
50  (
51  pressureControl.refCell(),
52  pressureControl.refValue()
53  );
54 
55  p_rghEqn.solve();
56  }
57 }
58 else
59 {
60  while (simple.correctNonOrthogonal())
61  {
62  p_rghEqn =
65  ==
66  fvOptions(psi, p_rgh, rho.name());
67 
68  p_rghEqn.setReference
69  (
70  pressureControl.refCell(),
71  pressureControl.refValue()
72  );
73 
74  p_rghEqn.solve();
75  }
76 }
77 
78 phi = phiHbyA + p_rghEqn.flux();
79 
80 p = p_rgh + rho*gh + pRef;
81 
83 
84 // Explicitly relax pressure for momentum corrector
85 p_rgh.relax();
86 
87 // Correct the momentum source with the pressure gradient flux
88 // calculated from the relaxed pressure
90 U.correctBoundaryConditions();
91 fvOptions.correct(U);
92 
93 pressureControl.limit(p);
94 
95 // For closed-volume compressible cases adjust the pressure level
96 // to obey overall mass continuity
97 if (closedVolume && !thermo.incompressible())
98 {
101  p_rgh = p - rho*gh - pRef;
102  p.correctBoundaryConditions();
103 }
104 
105 rho = thermo.rho();
106 
107 if (!simple.transonic())
108 {
109  rho.relax();
110 }
rAU
Definition: pEqn.H:1
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
rho
Definition: pEqn.H:1
dimensionedScalar initialMass
Definition: createFields.H:69
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())
p_rgh
Definition: pEqn.H:140
Calculates and prints the continuity errors.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
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
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
pressureControl & pressureControl
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
const surfaceScalarField & ghf
scalar pRef
Definition: createFields.H:19
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvOptions(rho, U))
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
adjustPhi(phiHbyA, Urel, p)
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
const volScalarField & gh
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
const volScalarField & psi
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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
simpleControl simple(mesh)