pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 volScalarField rAU(1.0/UEqn.A());
6 
8 
10 (
11  "phiHbyA",
12  (
13  fvc::flux(rho*HbyA)
15  )
16  + phig
17 );
18 
19 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
20 
21 // Update the pressure BCs to ensure flux consistency
23 
24 while (pimple.correctNonOrthogonal())
25 {
26  fvScalarMatrix p_rghEqn
27  (
29  + fvc::ddt(psi, rho)*gh
30  + fvc::ddt(psi)*pRef
31  + fvc::div(phiHbyA)
33  ==
34  parcels.Srho()
35  + surfaceFilm.Srho()
36  + fvOptions(psi, p_rgh, rho.name())
37  );
38 
39  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
40 
41  if (pimple.finalNonOrthogonalIter())
42  {
43  phi = phiHbyA + p_rghEqn.flux();
44  U = HbyA + rAU*fvc::reconstruct((p_rghEqn.flux() + phig)/rhorAUf);
45  U.correctBoundaryConditions();
46  fvOptions.correct(U);
47  }
48 }
49 
50 p = p_rgh + rho*gh + pRef;
51 
52 #include "rhoEqn.H"
53 #include "compressibleContinuityErrs.H"
54 
55 K = 0.5*magSqr(U);
56 
57 if (thermo.dpdt())
58 {
59  dpdt = fvc::ddt(p);
60 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
p
Definition: pEqn.H:50
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
const dictionary & pimple
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
CGAL::Exact_predicates_exact_constructions_kernel K
fv::options & fvOptions
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
const surfaceScalarField & ghf
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
IOMRFZoneList & MRF
volScalarField & dpdt
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
Solve the continuity for density.
p_rgh
Definition: pEqn.H:134
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
scalar pRef
Definition: createFields.H:19
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
phi
Definition: pEqn.H:18
const volScalarField & psi
filmModelType & surfaceFilm
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
rho
Definition: pEqn.H:1
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
volScalarField rAU(1.0/UEqn.A())