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;
7 
8  volScalarField rAU(1.0/UEqn.A());
10 
11  volVectorField HbyA("HbyA", U);
12  HbyA = rAU*UEqn.H();
13 
14  if (pimple.transonic())
15  {
17  (
18  "phiHbyA",
19  (
20  (fvc::interpolate(HbyA) & mesh.Sf())
22  )
23  );
24 
25  MRF.makeRelative(phiHbyA);
26 
27  surfaceScalarField phid("phid", fvc::interpolate(thermo.psi())*phiHbyA);
28 
29  phiHbyA *= fvc::interpolate(rho);
30 
32  (
33  fvc::ddt(rho) + fvc::div(phiHbyA)
34  + correction(psi*fvm::ddt(p) + fvm::div(phid, p))
35  );
36 
37  while (pimple.correctNonOrthogonal())
38  {
39  fvScalarMatrix pEqn
40  (
41  pDDtEqn
43  ==
44  fvOptions(psi, p, rho.name())
45  );
46 
47  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
48 
49  if (pimple.finalNonOrthogonalIter())
50  {
51  phi = phiHbyA + pEqn.flux();
52  }
53  }
54  }
55  else
56  {
58  (
59  "phiHbyA",
60  (
61  (fvc::interpolate(rho*HbyA) & mesh.Sf())
63  )
64  );
65 
66  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
67 
69  (
71  + fvc::div(phiHbyA)
72  ==
73  fvOptions(psi, p, rho.name())
74  );
75 
76  while (pimple.correctNonOrthogonal())
77  {
78  fvScalarMatrix pEqn
79  (
80  pDDtEqn
82  );
83 
84  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
85 
86  if (pimple.finalNonOrthogonalIter())
87  {
88  phi = phiHbyA + pEqn.flux();
89  }
90  }
91  }
92 
93  // Second part of thermodynamic density update
94  thermo.rho() += psi*p;
95 
96  #include "rhoEqn.H"
97  #include "compressibleContinuityErrs.H"
98 
99  U = HbyA - rAU*fvc::grad(p);
100  U.correctBoundaryConditions();
101  fvOptions.correct(U);
102  K = 0.5*magSqr(U);
103 
104  if (thermo.dpdt())
105  {
106  dpdt = fvc::ddt(p);
107  }
108 }
volScalarField & dpdt
fvScalarMatrix pDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p))+fvc::div(phiHbyA)==fvOptions(psi, p, rho.name()))
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
phi
Definition: pEqn.H:20
CGAL::Exact_predicates_exact_constructions_kernel K
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
const dictionary & pimple
p
Definition: pEqn.H:59
fv::IOoptionList & fvOptions
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))
surfaceScalarField phid("phid", fvc::interpolate(psi)*( (mesh.Sf()&fvc::interpolate(HbyA)) +rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho) ))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
const volScalarField & psi
Definition: createFields.H:24
phiHbyA
Definition: pEqn.H:21
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
volScalarField rAU(1.0/UEqn.A())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
rho
Definition: pEqn.H:1
psiReactionThermo & thermo
Definition: createFields.H:32
U
Definition: pEqn.H:82
HbyA
Definition: pEqn.H:7