pEqnTurbo.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 // bound rho
3 DAUtility::boundVar(allOptions, rho, printToScreen);
4 
5 if (!simple.transonic())
6 {
7  rho.relax();
8 }
9 
10 volScalarField p0(p);
11 
12 volScalarField AU(UEqn.A());
13 volScalarField AtU(AU - UEqn.H1());
14 volVectorField HbyA("HbyA", U);
15 HbyA = UEqn.H() / AU;
16 
17 volScalarField rAU(1.0 / UEqn.A());
18 tUEqn.clear();
19 
20 bool closedVolume = false;
21 
22 if (simple.transonic())
23 {
24 
25  surfaceScalarField phid(
26  "phid",
27  fvc::interpolate(psi) * (fvc::interpolate(HbyA) & mesh.Sf()));
28 
29  MRF.makeRelative(fvc::interpolate(psi), phid);
30 
31  while (simple.correctNonOrthogonal())
32  {
33  fvScalarMatrix pEqn(
34  fvm::div(phid, p)
35  - fvm::laplacian(rho * rAU, p));
36 
37  // Relax the pressure equation to ensure diagonal-dominance
38  pEqn.relax();
39 
40  pEqn.setReference(
41  pressureControl.refCell(),
42  pressureControl.refValue());
43 
44  // get the solver performance info such as initial
45  // and final residuals
46  SolverPerformance<scalar> solverP = pEqn.solve();
47 
48  this->primalResidualControl<scalar>(solverP, printToScreen, printInterval, "p");
49 
50  if (simple.finalNonOrthogonalIter())
51  {
52  phi == pEqn.flux();
53  }
54  }
55 }
56 else
57 {
58  surfaceScalarField phiHbyA(
59  "phiHbyA",
60  fvc::interpolate(rho * HbyA) & mesh.Sf());
61 
62  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
63 
65  phiHbyA += fvc::interpolate(rho / AtU - rho / AU) * fvc::snGrad(p) * mesh.magSf();
66 
67  while (simple.correctNonOrthogonal())
68  {
69  fvScalarMatrix pEqn(
70  fvc::div(phiHbyA)
71  - fvm::laplacian(rho / AtU, p));
72 
73  pEqn.setReference(
74  pressureControl.refCell(),
75  pressureControl.refValue());
76 
77  // get the solver performance info such as initial
78  // and final residuals
79  SolverPerformance<scalar> solverP = pEqn.solve();
80 
81  this->primalResidualControl<scalar>(solverP, printToScreen, printInterval, "p");
82 
83  if (simple.finalNonOrthogonalIter())
84  {
85  phi = phiHbyA + pEqn.flux();
86  }
87  }
88 }
89 
90 if (printToScreen)
91 {
92 #include "continuityErrsPython.H"
93 }
94 
95 // Explicitly relax pressure for momentum corrector
96 p.relax();
97 
98 // bound p
99 DAUtility::boundVar(allOptions, p, printToScreen);
100 
101 U = HbyA - (fvc::grad(p0) * (1.0 / AU - 1.0 / AtU) + fvc::grad(p) / AtU);
102 
103 // bound U
104 DAUtility::boundVar(allOptions, U, printToScreen);
105 U.correctBoundaryConditions();
106 
107 // For closed-volume cases adjust the pressure and density levels
108 // to obey overall mass continuity
110 {
111  p += (initialMass_ - fvc::domainIntegrate(psi * p))
112  / fvc::domainIntegrate(psi);
113 }
114 
115 rho = thermo.rho();
116 DAUtility::boundVar(allOptions, rho, printToScreen);
117 
118 if (!simple.transonic())
119 {
120  rho.relax();
121 }
phid
surfaceScalarField phid("phid",(fvc::interpolate(psi)/fvc::interpolate(rho)) *phiHbyA)
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(U)+fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
allOptions
const dictionary & allOptions
Definition: createRefsRhoSimpleC.H:15
adjustPhi
adjustPhi(phiHbyA, U, p)
pressureControl
pressureControl & pressureControl
Definition: createRefsRhoSimpleC.H:12
simple
simpleControl & simple
Definition: createRefsRhoSimpleC.H:5
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
HbyA
HbyA
Definition: pEqnTurbo.H:15
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
MRF
IOMRFZoneListDF & MRF
Definition: createRefsRhoSimple.H:18
p0
volScalarField p0(p)
p
volScalarField & p
Definition: createRefsPimple.H:6
continuityErrsPython.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
rAU
volScalarField rAU(1.0/UEqn.A())
primalResidualControl< scalar >
this primalResidualControl< scalar >(solverE, printToScreen, printInterval, "he")
AU
volScalarField AU(UEqn.A())
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
rho
rho
Definition: pEqnTurbo.H:1
phiHbyA
phiHbyA
Definition: pEqnTurbo.H:65
closedVolume
bool closedVolume
Definition: pEqnTurbo.H:20
AtU
volScalarField AtU(AU - UEqn.H1())
U
U
Definition: pEqnTurbo.H:101
thermo
fluidThermo & thermo
Definition: createRefsRhoSimpleC.H:6
psi
const volScalarField & psi
Definition: createRefsRhoSimpleC.H:14