pcEqn.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()));
6 tUEqn.clear();
7 
8 bool closedVolume = false;
9 
11 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
12 
13 volScalarField rhorAtU("rhorAtU", rho*rAtU);
14 
15 // Update the pressure BCs to ensure flux consistency
17 
18 if (simple.transonic())
19 {
21  (
22  "phid",
24  );
25 
26  phiHbyA +=
29 
30  HbyA -= (rAU - rAtU)*fvc::grad(p);
31 
32  while (simple.correctNonOrthogonal())
33  {
34  fvScalarMatrix pEqn
35  (
37  + fvm::div(phid, p)
39  ==
40  fvOptions(psi, p, rho.name())
41  );
42 
43  // Relax the pressure equation to maintain diagonal dominance
44  pEqn.relax();
45 
46  pEqn.setReference
47  (
48  pressureControl.refCell(),
49  pressureControl.refValue()
50  );
51 
52  pEqn.solve();
53 
54  if (simple.finalNonOrthogonalIter())
55  {
56  phi = phiHbyA + pEqn.flux();
57  }
58  }
59 }
60 else
61 {
63 
65  HbyA -= (rAU - rAtU)*fvc::grad(p);
66 
67  while (simple.correctNonOrthogonal())
68  {
69  fvScalarMatrix pEqn
70  (
73  ==
74  fvOptions(psi, p, rho.name())
75  );
76 
77  pEqn.setReference
78  (
79  pressureControl.refCell(),
80  pressureControl.refValue()
81  );
82 
83  pEqn.solve();
84 
85  if (simple.finalNonOrthogonalIter())
86  {
87  phi = phiHbyA + pEqn.flux();
88  }
89  }
90 }
91 
92 // The incompressible form of the continuity error check is appropriate for
93 // steady-state compressible also.
95 
96 // Explicitly relax pressure for momentum corrector
97 p.relax();
98 
100 U.correctBoundaryConditions();
101 fvOptions.correct(U);
102 
103 bool pLimited = pressureControl.limit(p);
104 
105 // For closed-volume cases adjust the pressure and density levels
106 // to obey overall mass continuity
108 {
111 }
112 
114 {
115  p.correctBoundaryConditions();
116 }
117 
118 rho = thermo.rho();
119 
120 if (!simple.transonic())
121 {
122  rho.relax();
123 }
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:35
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:47
rAtU
volScalarField rAtU(1.0/(1.0/rAU - UEqn.H1()))
phid
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+MRF.zeroFilter(rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho))))
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
rAU
volScalarField rAU(1.0/UEqn.A())
constrainPressure
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF)
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
HbyA
HbyA
Definition: pcEqn.H:74
continuityErrs.H
Calculates and prints the continuity errors.
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
phiHbyA
phiHbyA
Definition: pcEqn.H:73
pLimited
bool pLimited
Definition: pcEqn.H:103
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
rho
rho
Definition: pcEqn.H:1
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
U
U
Definition: pcEqn.H:107
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
pressureControl
const pressureControl & pressureControl
Definition: setRegionFluidFields.H:69
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
rhorAtU
volScalarField rhorAtU("rhorAtU", rho *rAtU)
initialMass
dimensionedScalar initialMass
Definition: createFields.H:57
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
closedVolume
bool closedVolume
Definition: pcEqn.H:8
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.