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