pEqnRhoSimpleC.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 volScalarField rAU(1.0 / UEqn.A());
4 volScalarField rAtU(1.0 / (1.0 / rAU - UEqn.H1()));
5 //volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
6 //***************** NOTE *******************
7 // constrainHbyA has been used since OpenFOAM-v1606; however, We do NOT use the constrainHbyA
8 // function in DAFoam because we found it significantly degrades the accuracy of shape derivatives.
9 // Basically, we should not constrain any variable because it will create discontinuity.
10 // Instead, we use the old implementation used in OpenFOAM-3.0+ and before
11 volVectorField HbyA("HbyA", U);
12 HbyA = rAU * UEqn.H();
13 tUEqn.clear();
14 
15 bool closedVolume = false;
16 
17 surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho) * fvc::flux(HbyA));
18 
19 // NOTE: we don't support transonic = false
20 
21 volScalarField rhorAtU("rhorAtU", rho* rAtU);
22 
23 surfaceScalarField phid(
24  "phid",
25  (fvc::interpolate(psi) / fvc::interpolate(rho)) * phiHbyA);
26 
28  fvc::interpolate(rho * (rAtU - rAU)) * fvc::snGrad(p) * mesh.magSf()
29  - fvc::interpolate(psi * p) * phiHbyA / fvc::interpolate(rho);
30 
31 HbyA -= (rAU - rAtU) * fvc::grad(p);
32 
33 while (simple.correctNonOrthogonal())
34 {
35  fvScalarMatrix pEqn(
36  fvc::div(phiHbyA)
37  + fvm::div(phid, p)
38  - fvm::laplacian(rhorAtU, p));
39 
40  // Relax the pressure equation to maintain diagonal dominance
41  pEqn.relax();
42 
43  pEqn.setReference(
44  pressureControl.refCell(),
45  pressureControl.refValue());
46 
47  // get the solver performance info such as initial
48  // and final residuals
49  SolverPerformance<scalar> solverP = pEqn.solve();
50 
51  this->primalResidualControl<scalar>(solverP, printToScreen, printInterval, "p");
52 
53  if (simple.finalNonOrthogonalIter())
54  {
55  phi = phiHbyA + pEqn.flux();
56  }
57 }
58 
59 if (printToScreen)
60 {
61 #include "continuityErrsPython.H"
62 }
63 
64 // Explicitly relax pressure for momentum corrector
65 p.relax();
66 
67 // bound p
68 DAUtility::boundVar(allOptions, p, printToScreen);
69 
70 U = HbyA - rAtU * fvc::grad(p);
71 // bound U
72 DAUtility::boundVar(allOptions, U, printToScreen);
73 U.correctBoundaryConditions();
74 
75 bool pLimited = pressureControl.limit(p);
76 
77 // For closed-volume cases adjust the pressure and density levels
78 // to obey overall mass continuity
80 {
81  p += (initialMass_ - fvc::domainIntegrate(psi * p))
82  / fvc::domainIntegrate(psi);
83 }
84 
86 {
87  p.correctBoundaryConditions();
88 }
89 
90 rho = thermo.rho();
91 
92 // bound rho
93 DAUtility::boundVar(allOptions, rho, printToScreen);
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
pressureControl
pressureControl & pressureControl
Definition: createRefsRhoSimpleC.H:12
simple
simpleControl & simple
Definition: createRefsRhoSimpleC.H:5
initialMass_
initialMass_
Definition: createRefsRhoSimpleC.H:13
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
closedVolume
bool closedVolume
Definition: pEqnRhoSimpleC.H:15
rho
rho
Definition: pEqnRhoSimpleC.H:1
HbyA
HbyA
Definition: pEqnRhoSimpleC.H:12
p
volScalarField & p
Definition: createRefsPimple.H:6
pLimited
bool pLimited
Definition: pEqnRhoSimpleC.H:75
continuityErrsPython.H
U
U
Definition: pEqnRhoSimpleC.H:70
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
primalResidualControl< scalar >
this primalResidualControl< scalar >(solverE, printToScreen, printInterval, "he")
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
rAU
volScalarField rAU(1.0/UEqn.A())
phiHbyA
phiHbyA
Definition: pEqnRhoSimpleC.H:27
thermo
fluidThermo & thermo
Definition: createRefsRhoSimpleC.H:6
rhorAtU
volScalarField rhorAtU("rhorAtU", rho *rAtU)
psi
const volScalarField & psi
Definition: createRefsRhoSimpleC.H:14
rAtU
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))