pEqn.H
Go to the documentation of this file.
1 {
2  volScalarField rAU("rAU", 1.0/UEqn.A());
4 
5  volVectorField HbyA("HbyA", U);
6  HbyA = rAU*UEqn.H();
7 
9  (
10  "phiHbyA",
11  (fvc::interpolate(HbyA) & mesh.Sf())
13  );
14  MRF.makeRelative(phiHbyA);
16 
18  (
19  (
21  )*rAUf*mesh.magSf()
22  );
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  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
34  );
35 
36  while (pimple.correctNonOrthogonal())
37  {
38  fvScalarMatrix p_rghEqn
39  (
41  );
42 
43  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
44 
45  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
46 
47  if (pimple.finalNonOrthogonalIter())
48  {
49  phi = phiHbyA - p_rghEqn.flux();
50 
51  p_rgh.relax();
52 
53  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rAUf);
54  U.correctBoundaryConditions();
55  fvOptions.correct(U);
56  }
57  }
58 
59  #include "continuityErrs.H"
60 
61  p == p_rgh + rho*gh;
62 
63  if (p_rgh.needReference())
64  {
66  (
67  "p",
68  p.dimensions(),
70  );
71  p_rgh = p - rho*gh;
72  }
73 }
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
phig
surfaceScalarField phig(-rAUf *ghf *fvc::snGrad(rhok) *mesh.magSf())
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
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()))
phi
phi
Definition: pEqn.H:20
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::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
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
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())
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
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
phiHbyA
phiHbyA
Definition: pEqn.H:21
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:31
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))