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  UEqn.clear();
8 
10 
12  (
13  "phiHbyA",
14  (fvc::interpolate(HbyA) & mesh.Sf())
15  );
16 
17  MRF.makeRelative(phiHbyA);
18 
20 
22 
23  // Update the fixedFluxPressure BCs to ensure flux consistency
25  (
26  p_rgh.boundaryField(),
27  (
28  phiHbyA.boundaryField()
29  - MRF.relative(mesh.Sf().boundaryField() & U.boundaryField())
30  )/(mesh.magSf().boundaryField()*rAUf.boundaryField())
31  );
32 
33  while (simple.correctNonOrthogonal())
34  {
35  fvScalarMatrix p_rghEqn
36  (
38  );
39 
40  p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
41 
42  p_rghEqn.solve();
43 
44  if (simple.finalNonOrthogonalIter())
45  {
46  // Calculate the conservative fluxes
47  phi = phiHbyA - p_rghEqn.flux();
48 
49  // Explicitly relax pressure for momentum corrector
50  p_rgh.relax();
51 
52  // Correct the momentum source with the pressure gradient flux
53  // calculated from the relaxed pressure
54  U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rAUf);
55  U.correctBoundaryConditions();
56  fvOptions.correct(U);
57  }
58  }
59 
60  #include "continuityErrs.H"
61 
62  p = p_rgh + rhok*gh;
63 
64  if (p_rgh.needReference())
65  {
67  (
68  "p",
69  p.dimensions(),
71  );
72  p_rgh = p - rhok*gh;
73  }
74 }
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()))
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
volScalarField & p_rgh
tmp< surfaceScalarField > interpolate(const RhoType &rho)
dynamicFvMesh & mesh
surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
p
Definition: pEqn.H:59
rhok
Definition: TEqn.H:27
fv::IOoptionList & fvOptions
const scalar pRefValue
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
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
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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())
const label pRefCell
const dictionary & simple
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
U
Definition: pEqn.H:82
const surfaceScalarField & ghf
HbyA
Definition: pEqn.H:7
const volScalarField & gh