phrghEqn.H
Go to the documentation of this file.
1 if (pimple.dict().lookupOrDefault<bool>("hydrostaticInitialization", false))
2 {
3  volScalarField& ph_rgh = regIOobject::store
4  (
6  (
7  IOobject
8  (
9  "ph_rgh",
10  "0",
11  mesh,
12  IOobject::MUST_READ,
13  IOobject::NO_WRITE
14  ),
15  mesh
16  )
17  );
18 
19  if (equal(runTime.value(), 0))
20  {
21  p = ph_rgh + rho*gh + pRef;
22  thermo.correct();
23  rho = thermo.rho();
24 
25  label nCorr
26  (
27  pimple.dict().lookupOrDefault<label>("nHydrostaticCorrectors", 5)
28  );
29 
30  for (label i=0; i<nCorr; i++)
31  {
33 
35  (
36  "phig",
37  -rhof*ghf*fvc::snGrad(rho)*mesh.magSf()
38  );
39 
40  // Update the pressure BCs to ensure flux consistency
41  constrainPressure(ph_rgh, rho, U, phig, rhof);
42 
43  fvScalarMatrix ph_rghEqn
44  (
45  fvm::laplacian(rhof, ph_rgh) == fvc::div(phig)
46  );
47 
48  ph_rghEqn.solve();
49 
50  p = ph_rgh + rho*gh + pRef;
51  thermo.correct();
52  rho = thermo.rho();
53 
54  Info<< "Hydrostatic pressure variation "
55  << (max(ph_rgh) - min(ph_rgh)).value() << endl;
56  }
57 
58  ph_rgh.write();
59 
60  p_rgh = ph_rgh;
61  }
62 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
scalar pRef
Definition: createFields.H:19
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
const volScalarField & gh
messageStream Info
virtual Ostream & write(const token &)=0
Write next token to stream.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
volScalarField & p_rgh
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45