pEqn.H
Go to the documentation of this file.
1 {
2  rho = thermo.rho();
3 
4  // Thermodynamic density needs to be updated by psi*d(p) after the
5  // pressure solution - done in 2 parts. Part 1:
6  thermo.rho() -= psi*p_rgh;
7 
8  volScalarField rAU(1.0/UEqn.A());
11 
13 
15  (
16  "phiHbyA",
17  (
18  fvc::flux(rho*HbyA)
20  )
21  + phig
22  );
23 
24  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
25 
26  // Update the pressure BCs to ensure flux consistency
28 
30  (
32  + fvc::div(phiHbyA)
33  ==
34  fvOptions(psi, p_rgh, rho.name())
35  );
36 
37  while (pimple.correctNonOrthogonal())
38  {
39  fvScalarMatrix p_rghEqn
40  (
42  - fvm::laplacian(rhorAUf, p_rgh)
43  );
44 
45  p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
46 
47  if (pimple.finalNonOrthogonalIter())
48  {
49  // Calculate the conservative fluxes
50  phi = phiHbyA + p_rghEqn.flux();
51 
52  // Explicitly relax pressure for momentum corrector
53  p_rgh.relax();
54 
55  // Correct the momentum source with the pressure gradient flux
56  // calculated from the relaxed pressure
57  U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rhorAUf);
58  U.correctBoundaryConditions();
59  fvOptions.correct(U);
60  K = 0.5*magSqr(U);
61  }
62  }
63 
64  p = p_rgh + rho*gh;
65 
66  // Second part of thermodynamic density update
67  thermo.rho() += psi*p_rgh;
68 
69  if (thermo.dpdt())
70  {
71  dpdt = fvc::ddt(p);
72  }
73 
74  #include "rhoEqn.H"
75  #include "compressibleContinuityErrs.H"
76 }
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
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
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
const dictionary & pimple
tmp< surfaceScalarField > interpolate(const RhoType &rho)
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
surfaceScalarField phig("phig",-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
p_rgh
Definition: pEqn.H:120
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
fvScalarMatrix p_rghDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p_rgh))+fvc::div(phiHbyA)==fvOptions(psi, p_rgh, rho.name()))
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
const volScalarField & gh
phi
Definition: pEqn.H:18
const volScalarField & psi
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())