pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
5  //tUEqn.clear();
6 
8 
10  (
11  "phiHbyA",
13  );
14 
15  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
16 
18 
20 
21  // Update the pressure BCs to ensure flux consistency
23 
25  bool compressible = (compressibility.value() > SMALL);
26 
27  // Solve pressure
28  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
29  {
30  fvScalarMatrix p_rghEqn
31  (
33  );
34 
35  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
36 
37  p_rghEqn.solve();
38 
39  if (nonOrth == nNonOrthCorr)
40  {
41  // Calculate the conservative fluxes
42  phi = phiHbyA - p_rghEqn.flux();
43 
44  // Explicitly relax pressure for momentum corrector
45  p_rgh.relax();
46 
47  // Correct the momentum source with the pressure gradient flux
48  // calculated from the relaxed pressure
49  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rhorAUf);
50  U.correctBoundaryConditions();
51  fvOptions.correct(U);
52  }
53  }
54 
55  p = p_rgh + rho*gh;
56 
57  #include "continuityErrs.H"
58 
59  // For closed-volume cases adjust the pressure level
60  // to obey overall mass continuity
62  {
63  if (!compressible)
64  {
66  (
67  "p",
68  p.dimensions(),
70  );
71  }
72  else
73  {
76  }
77  p_rgh = p - rho*gh;
78  }
79 
80  rho = thermo.rho();
81  rho = max(rho, rhoMin[i]);
82  rho = min(rho, rhoMax[i]);
83  rho.relax();
84 
85  Info<< "Min/max rho:" << min(rho).value() << ' '
86  << max(rho).value() << endl;
87 }
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:56
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
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:18
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
phiHbyA
phiHbyA
Definition: pEqn.H:20
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:17
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
HbyA
HbyA
Definition: pEqn.H:4
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
rhoMin
const dimensionedScalar rhoMin
Definition: setRegionFluidFields.H:67
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
compressible
bool compressible
Definition: pEqn.H:2
p
p
Definition: pEqn.H:50
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
rAU
volScalarField rAU(1.0/UEqn.A())
nNonOrthCorr
const int nNonOrthCorr
Definition: readFluidMultiRegionSIMPLEControls.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
phig
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
Foam::getRefCellValue
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:134
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
rhoMax
const dimensionedScalar rhoMax
Definition: setRegionFluidFields.H:66
phi
phi
Definition: pEqn.H:18
initialMass
dimensionedScalar initialMass
Definition: createFields.H:57
closedVolume
bool closedVolume
Definition: pEqn.H:10
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
compressibility
dimensionedScalar compressibility
Definition: pEqn.H:1
p_rgh
p_rgh
Definition: pEqn.H:139
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
adjustPhi
adjustPhi(phiHbyA, U, p_rgh)
rho
rho
Definition: pEqn.H:1
interpolate
mesh interpolate(rAU)
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:36
constrainPressure
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:37