pEqnSimpleT.H
Go to the documentation of this file.
1 {
2  volScalarField rAU(1.0 / UEqn.A());
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  surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA));
12 
13  MRF.makeRelative(phiHbyA);
14 
15  adjustPhi(phiHbyA, U, p);
16 
17  tmp<volScalarField> rAtU(rAU);
18 
19  if (simple.consistent())
20  {
21  rAtU = 1.0 / (1.0 / rAU - UEqn.H1());
22  phiHbyA +=
23  fvc::interpolate(rAtU() - rAU) * fvc::snGrad(p) * mesh.magSf();
24  HbyA -= (rAU - rAtU()) * fvc::grad(p);
25  }
26 
27  tUEqn.clear();
28 
29  // Update the pressure BCs to ensure flux consistency
31 
32  // Non-orthogonal pressure corrector loop
33  while (simple.correctNonOrthogonal())
34  {
35  fvScalarMatrix pEqn(
36  fvm::laplacian(rAtU(), p) == fvc::div(phiHbyA));
37 
38  pEqn.setReference(pRefCell, pRefValue);
39 
40  // get the solver performance info such as initial
41  // and final residuals
42  SolverPerformance<scalar> solverP = pEqn.solve();
43 
44  this->primalResidualControl<scalar>(solverP, printToScreen, printInterval, "p");
45 
46  if (simple.finalNonOrthogonalIter())
47  {
48  phi = phiHbyA - pEqn.flux();
49  }
50  }
51 
52  if (printToScreen)
53  {
54 #include "continuityErrsPython.H"
55  }
56 
57  // Explicitly relax pressure for momentum corrector
58  p.relax();
59 
60  // Momentum corrector
61  U = HbyA - rAtU() * fvc::grad(p);
62  U.correctBoundaryConditions();
63  fvOptions.correct(U);
64 }
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(U)+fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
simple
simpleControl & simple
Definition: createRefsRhoSimpleC.H:5
pRefValue
scalar & pRefValue
Definition: createRefsPimple.H:12
UEqn
fvVectorMatrix & UEqn
Definition: UEqnPimpleDyM.H:13
pRefCell
label & pRefCell
Definition: createRefsPimple.H:11
MRF
IOMRFZoneListDF & MRF
Definition: createRefsRhoSimple.H:18
p
volScalarField & p
Definition: createRefsPimple.H:6
continuityErrsPython.H
mesh
fvMesh & mesh
Definition: createRefsHeatTransfer.H:4
primalResidualControl< scalar >
this primalResidualControl< scalar >(solverE, printToScreen, printInterval, "he")
phi
surfaceScalarField & phi
Definition: createRefsPimple.H:8
U
U
Definition: pEqnSimpleT.H:61
constrainPressure
constrainPressure(p, U, phiHbyA, rAtU(), MRF)
rAtU
tmp< volScalarField > rAtU(rAU)
adjustPhi
adjustPhi(phiHbyA, U, p)
rAU
volScalarField rAU(1.0/UEqn.A())
HbyA
HbyA
Definition: pEqnSimpleT.H:10
phiHbyA
surfaceScalarField phiHbyA("phiHbyA", fvc::flux(HbyA))