pcEqn.H
Go to the documentation of this file.
1 volScalarField rAU(1.0/UEqn.A());
2 volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()));
4 tUEqn.clear();
5 
6 bool closedVolume = false;
7 
8 if (simple.transonic())
9 {
11  (
12  "phid",
14  *fvc::flux(HbyA)
15  );
16 
17  MRF.makeRelative(fvc::interpolate(psi), phid);
18 
20  (
21  "phic",
23  );
24 
25  HbyA -= (rAU - rAtU)*fvc::grad(p);
26 
27  volScalarField rhorAtU("rhorAtU", rho*rAtU);
28 
29  while (simple.correctNonOrthogonal())
30  {
31  fvScalarMatrix pEqn
32  (
33  fvm::div(phid, p)
34  + fvc::div(phic)
35  - fvm::laplacian(rhorAtU, p)
36  ==
37  fvOptions(psi, p, rho.name())
38  );
39 
40  // Relax the pressure equation to maintain diagonal dominance
41  pEqn.relax();
42 
43  pEqn.setReference(pRefCell, pRefValue);
44 
45  pEqn.solve();
46 
47  if (simple.finalNonOrthogonalIter())
48  {
49  phi == phic + pEqn.flux();
50  }
51  }
52 }
53 else
54 {
55  surfaceScalarField phiHbyA("phiHbyA", fvc::flux(rho*HbyA));
56  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
57 
58  closedVolume = adjustPhi(phiHbyA, U, p);
59 
60  phiHbyA += fvc::interpolate(rho*(rAtU - rAU))*fvc::snGrad(p)*mesh.magSf();
61  HbyA -= (rAU - rAtU)*fvc::grad(p);
62 
63  volScalarField rhorAtU("rhorAtU", rho*rAtU);
64 
65  // Update the pressure BCs to ensure flux consistency
66  constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF);
67 
68  while (simple.correctNonOrthogonal())
69  {
70  fvScalarMatrix pEqn
71  (
72  fvc::div(phiHbyA)
74  ==
75  fvOptions(psi, p, rho.name())
76  );
77 
78  pEqn.setReference(pRefCell, pRefValue);
79 
80  pEqn.solve();
81 
82  if (simple.finalNonOrthogonalIter())
83  {
84  phi = phiHbyA + pEqn.flux();
85  }
86  }
87 }
88 
89 // The incompressibe form of the continuity error check is appropriate for
90 // steady-state compressible also.
92 
93 // Explicitly relax pressure for momentum corrector
94 p.relax();
95 
96 U = HbyA - rAtU*fvc::grad(p);
97 U.correctBoundaryConditions();
98 fvOptions.correct(U);
99 
100 // For closed-volume cases adjust the pressure and density levels
101 // to obey overall mass continuity
102 if (closedVolume)
103 {
106 }
107 
108 // Recalculate density from the relaxed pressure
109 rho = thermo.rho();
110 rho = max(rho, rhoMin);
111 rho = min(rho, rhoMax);
112 
113 if (!simple.transonic())
114 {
115  rho.relax();
116 }
117 
118 Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value() << endl;
surfaceScalarField & phi
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
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
dimensionedScalar initialMass
Definition: createFields.H:82
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF)
volScalarField rAU(1.0/UEqn.A())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
Calculates and prints the continuity errors.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
tmp< surfaceScalarField > interpolate(const RhoType &rho)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
HbyA
Definition: pcEqn.H:74
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
fv::options & fvOptions
phiHbyA
Definition: pcEqn.H:73
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
IOMRFZoneList & MRF
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
rho
Definition: pcEqn.H:1
surfaceScalarField phid("phid", fvc::interpolate(psi)*(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
U
Definition: pcEqn.H:107
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const scalar pRefValue
const label pRefCell
const dictionary & simple
volScalarField rhorAtU("rhorAtU", rho *rAtU)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
volScalarField rAtU(1.0/(1.0/rAU-UEqn.H1()))
messageStream Info
const volScalarField & psi
volScalarField & p
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