pEqn.H
Go to the documentation of this file.
1 {
2  rho = thermo.rho();
3  rho.relax();
4 
5  volScalarField rAU("rAU", 1.0/UEqn().A());
7 
8  volVectorField HbyA("HbyA", U);
9  HbyA = rAU*UEqn().H();
10  UEqn.clear();
11 
13 
15  (
16  "phiHbyA",
17  (fvc::interpolate(rho*HbyA) & mesh.Sf())
18  );
19 
20  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
21 
23 
25 
26  // Update the fixedFluxPressure BCs to ensure flux consistency
28  (
29  p_rgh.boundaryField(),
30  (
31  phiHbyA.boundaryField()
32  - MRF.relative(mesh.Sf().boundaryField() & U.boundaryField())
33  *rho.boundaryField()
34  )/(mesh.magSf().boundaryField()*rhorAUf.boundaryField())
35  );
36 
37 
38  while (simple.correctNonOrthogonal())
39  {
40  fvScalarMatrix p_rghEqn
41  (
43  );
44 
45  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
46  p_rghEqn.solve();
47 
48  if (simple.finalNonOrthogonalIter())
49  {
50  // Calculate the conservative fluxes
51  phi = phiHbyA - p_rghEqn.flux();
52 
53  // Explicitly relax pressure for momentum corrector
54  p_rgh.relax();
55 
56  // Correct the momentum source with the pressure gradient flux
57  // calculated from the relaxed pressure
58  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rhorAUf);
59  U.correctBoundaryConditions();
60  fvOptions.correct(U);
61  }
62  }
63 
64  #include "continuityErrs.H"
65 
66  p = p_rgh + rho*gh;
67 
69  bool compressible = (compressibility.value() > SMALL);
70 
71  // For closed-volume cases adjust the pressure level
72  // to obey overall mass continuity
74  {
75  if(!compressible)
76  {
78  (
79  "p",
80  p.dimensions(),
82  );
83  }
84  else
85  {
88  }
89  p_rgh = p - rho*gh;
90  }
91 
92  rho = thermo.rho();
93  rho.relax();
94  Info<< "rho max/min : " << max(rho).value() << " " << min(rho).value()
95  << endl;
96 }
Foam::fvc::reconstruct
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcReconstruct.C:54
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:27
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
p
p
Definition: pEqn.H:62
ghf
const surfaceScalarField & ghf
Definition: setRegionFluidFields.H:35
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:86
phig
surfaceScalarField phig(-rAUf *ghf *fvc::snGrad(rhok) *mesh.magSf())
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
setSnGrad< fixedFluxPressureFvPatchScalarField >
setSnGrad< fixedFluxPressureFvPatchScalarField >(p_rgh.boundaryField(),(phiHbyA.boundaryField() - MRF.relative(mesh.Sf().boundaryField() &U.boundaryField()))/(mesh.magSf().boundaryField() *rAUf.boundaryField()))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
phi
phi
Definition: pEqn.H:20
A
simpleMatrix< scalar > A(Nc)
simple
Simple relative velocity model.
U
U
Definition: pEqn.H:46
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Info
messageStream Info
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:15
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:34
rhorAUf
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
adjustPhi
adjustPhi(phiHbyA, U, p)
Foam::getRefCellValue
scalar getRefCellValue(const volScalarField &field, const label refCelli)
Return the current value of field in the reference cell.
Definition: findRefCell.C:130
HbyA
HbyA
Definition: pEqn.H:4
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
compressible
bool compressible
Definition: pEqn.H:40
closedVolume
bool closedVolume
Definition: pEqn.H:24
initialMass
dimensionedScalar initialMass
Definition: createFields.H:83
phiHbyA
phiHbyA
Definition: pEqn.H:21
compressibility
dimensionedScalar compressibility
Definition: pEqn.H:39
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:31
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16