pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 volScalarField rAU(1.0/UEqn.A());
5 
6 volVectorField HbyA("HbyA", U);
7 HbyA = rAU*UEqn.H();
8 
9 if (pimple.transonic())
10 {
12  (
13  "phid",
15  *(
16  (fvc::interpolate(HbyA) & mesh.Sf())
18  )
19  );
20 
21  MRF.makeRelative(fvc::interpolate(psi), phid);
22 
23  while (pimple.correctNonOrthogonal())
24  {
25  fvScalarMatrix pEqn
26  (
27  fvm::ddt(psi, p)
28  + fvm::div(phid, p)
30  ==
31  coalParcels.Srho()
32  + fvOptions(psi, p, rho.name())
33  );
34 
35  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
36 
37  if (pimple.finalNonOrthogonalIter())
38  {
39  phi == pEqn.flux();
40  }
41  }
42 }
43 else
44 {
46  (
47  "phiHbyA",
48  (
49  (fvc::interpolate(rho*HbyA) & mesh.Sf())
51  )
52  );
53 
54  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
55 
56  while (pimple.correctNonOrthogonal())
57  {
58  fvScalarMatrix pEqn
59  (
60  fvm::ddt(psi, p)
61  + fvc::div(phiHbyA)
63  ==
64  coalParcels.Srho()
65  + fvOptions(psi, p, rho.name())
66  );
67 
68  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
69 
70  if (pimple.finalNonOrthogonalIter())
71  {
72  phi = phiHbyA + pEqn.flux();
73  }
74  }
75 }
76 
77 #include "rhoEqn.H"
78 #include "compressibleContinuityErrs.H"
79 
80 U = HbyA - rAU*fvc::grad(p);
81 U.correctBoundaryConditions();
82 fvOptions.correct(U);
83 
84 K = 0.5*magSqr(U);
85 
86 if (thermo.dpdt())
87 {
88  dpdt = fvc::ddt(p);
89 }
volScalarField & dpdt
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())
Solve the continuity for density.
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