pEqn.H
Go to the documentation of this file.
1 {
2  rho = thermo.rho();
3  rho = max(rho, rhoMin[i]);
4  rho = min(rho, rhoMax[i]);
5  rho.relax();
6 
7  volScalarField rAU("rAU", 1.0/UEqn().A());
9 
10  volVectorField HbyA("HbyA", U);
11  HbyA = rAU*UEqn().H();
12  UEqn.clear();
13 
15 
17  (
18  "phiHbyA",
19  (fvc::interpolate(rho*HbyA) & mesh.Sf())
20  );
21 
22  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
23 
24  bool closedVolume = adjustPhi(phiHbyA, U, p_rgh);
25 
27 
28  // Update the fixedFluxPressure BCs to ensure flux consistency
30  (
31  p_rgh.boundaryField(),
32  (
33  phiHbyA.boundaryField()
34  - MRF.relative(mesh.Sf().boundaryField() & U.boundaryField())
35  *rho.boundaryField()
36  )/(mesh.magSf().boundaryField()*rhorAUf.boundaryField())
37  );
38 
40  bool compressible = (compressibility.value() > SMALL);
41 
42  // Solve pressure
43  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
44  {
45  fvScalarMatrix p_rghEqn
46  (
48  );
49 
50  p_rghEqn.setReference
51  (
52  pRefCell,
53  compressible ? getRefCellValue(p_rgh, pRefCell) : pRefValue
54  );
55 
56  p_rghEqn.solve();
57 
58  if (nonOrth == nNonOrthCorr)
59  {
60  // Calculate the conservative fluxes
61  phi = phiHbyA - p_rghEqn.flux();
62 
63  // Explicitly relax pressure for momentum corrector
64  p_rgh.relax();
65 
66  // Correct the momentum source with the pressure gradient flux
67  // calculated from the relaxed pressure
68  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rhorAUf);
69  U.correctBoundaryConditions();
70  fvOptions.correct(U);
71  }
72  }
73 
74  p = p_rgh + rho*gh;
75 
76  #include "continuityErrs.H"
77 
78  // For closed-volume cases adjust the pressure level
79  // to obey overall mass continuity
80  if (closedVolume && compressible)
81  {
84  p_rgh = p - rho*gh;
85  }
86 
87  rho = thermo.rho();
88  rho = max(rho, rhoMin[i]);
89  rho = min(rho, rhoMax[i]);
90  rho.relax();
91 
92  Info<< "Min/max rho:" << min(rho).value() << ' '
93  << max(rho).value() << endl;
94 }
phi
Definition: pEqn.H:20
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryField(),( phiHbyA.boundaryField() -MRF.relative(mesh.Sf().boundaryField()&U.boundaryField()) *rho.boundaryField() )/(mesh.magSf().boundaryField()*rhorAUf.boundaryField()))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedScalar compressibility
Definition: pEqn.H:39
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
volScalarField & p_rgh
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
dynamicFvMesh & mesh
surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
p
Definition: pEqn.H:59
const int nNonOrthCorr
fv::IOoptionList & fvOptions
const scalar pRefValue
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOMRFZoneList & MRF
fvVectorMatrix UEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
const volScalarField & psi
Definition: createFields.H:24
phiHbyA
Definition: pEqn.H:21
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
adjustPhi(phiHbyA, U, p_rgh)
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:154
volScalarField rAU(1.0/UEqn.A())
dimensionedScalar initialMass
Definition: createFields.H:83
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
const label pRefCell
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
rho
Definition: pEqn.H:1
psiReactionThermo & thermo
Definition: createFields.H:32
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
U
Definition: pEqn.H:82
const surfaceScalarField & ghf
HbyA
Definition: pEqn.H:7
const volScalarField & gh