pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 rho.relax();
3 
4 const volScalarField rAU("rAU", 1.0/UEqn.A());
6 
7 tmp<volScalarField> rAtU
8 (
9  simple.consistent()
10  ? volScalarField::New("rAtU", 1.0/(1.0/rAU - UEqn.H1()))
11  : tmp<volScalarField>(nullptr)
12 );
13 tmp<surfaceScalarField> rhorAtUf
14 (
15  simple.consistent()
17  : tmp<surfaceScalarField>(nullptr)
18 );
19 
20 const volScalarField& rAAtU = simple.consistent() ? rAtU() : rAU;
21 const surfaceScalarField& rhorAAtUf =
22  simple.consistent() ? rhorAtUf() : rhorAUf;
23 
25 
26 tUEqn.clear();
27 
29 (
30  "phiHbyA",
32 );
33 
34 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
35 
36 bool closedVolume = false;
37 
38 // Update the pressure BCs to ensure flux consistency
39 constrainPressure(p, rho, U, phiHbyA, rhorAAtUf, MRF);
40 
41 if (simple.transonic())
42 {
43  const surfaceScalarField phid
44  (
45  "phid",
47  );
48 
49  const fvScalarMatrix divPhidp(fvm::div(phid, p));
50  phiHbyA -= divPhidp.flux();
51 
52  if (simple.consistent())
53  {
54  phiHbyA += (rhorAAtUf - rhorAUf)*fvc::snGrad(p)*mesh.magSf();
55  HbyA += (rAAtU - rAU)*fvc::grad(p);
56  }
57 
58  while (simple.correctNonOrthogonal())
59  {
60  fvScalarMatrix pEqn
61  (
62  fvc::div(phiHbyA) + divPhidp
63  - fvm::laplacian(rhorAAtUf, p)
64  ==
65  fvModels.source(psi, p, rho.name())
66  );
67 
68  // Relax the pressure equation to ensure diagonal-dominance
69  pEqn.relax();
70 
71  pEqn.setReference
72  (
73  pressureReference.refCell(),
74  pressureReference.refValue()
75  );
76 
77  pEqn.solve();
78 
79  if (simple.finalNonOrthogonalIter())
80  {
81  phi = phiHbyA + pEqn.flux();
82  }
83  }
84 }
85 else
86 {
87  closedVolume = adjustPhi(phiHbyA, U, p);
88 
89  if (simple.consistent())
90  {
91  phiHbyA += (rhorAAtUf - rhorAUf)*fvc::snGrad(p)*mesh.magSf();
92  HbyA += (rAAtU - rAU)*fvc::grad(p);
93  }
94 
95  while (simple.correctNonOrthogonal())
96  {
97  fvScalarMatrix pEqn
98  (
100  - fvm::laplacian(rhorAAtUf, p)
101  ==
102  fvModels.source(psi, p, rho.name())
103  );
104 
105  pEqn.setReference
106  (
107  pressureReference.refCell(),
108  pressureReference.refValue()
109  );
110 
111  pEqn.solve();
112 
113  if (simple.finalNonOrthogonalIter())
114  {
115  phi = phiHbyA + pEqn.flux();
116  }
117  }
118 }
119 
121 
122 // Explicitly relax pressure for momentum corrector
123 p.relax();
124 
125 U = HbyA - rAAtU*fvc::grad(p);
126 U.correctBoundaryConditions();
128 
130 
131 // For closed-volume cases adjust the pressure and density levels
132 // to obey overall mass continuity
133 if (closedVolume && !thermo.incompressible())
134 {
137  p.correctBoundaryConditions();
138 }
139 
140 rho = thermo.rho();
141 rho.relax();
tmp< volScalarField > rAtU(pimple.consistent() ? volScalarField::New("rAtU", 1.0/(1.0/rAU - UEqn.H1())) :tmp< volScalarField >(nullptr))
rAU
Definition: pEqn.H:1
pressureReference & pressureReference
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fluidReactionThermo & thermo
Definition: createFields.H:28
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
dimensionedScalar initialMass
Definition: createFields.H:68
U
Definition: pEqn.H:72
mixture MRF().makeRelative(phiHbyA)
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > rhorAtUf(pimple.consistent() ? surfaceScalarField::New("rhoRAtUf", fvc::interpolate(rho *rAtU())) :tmp< surfaceScalarField >(nullptr))
Calculates and prints the continuity errors.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
p
Definition: pEqn.H:125
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
phiHbyA
Definition: pEqn.H:30
Foam::fvConstraints & fvConstraints
const volScalarField & psi
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(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.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
phi
Definition: pEqn.H:98
adjustPhi(phiHbyA, Urel, p)
volVectorField & HbyA
Definition: pEqn.H:13
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
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
rho
Definition: pEqn.H:1
simpleControl simple(mesh)
fvVectorMatrix & UEqn
Definition: UEqn.H:13