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 
10 (
11  "phid",
13  *(
14  (mesh.Sf() & fvc::interpolate(HbyA))
16  )
17 );
18 
19 MRF.makeRelative(fvc::interpolate(psi), phid);
20 
21 // Non-orthogonal pressure corrector loop
22 while (pimple.correctNonOrthogonal())
23 {
24  fvScalarMatrix pEqn
25  (
26  fvm::ddt(psi, p)
27  + fvm::div(phid, p)
29  ==
30  fvOptions(psi, p, rho.name())
31  );
32 
33  pEqn.solve();
34 
35  if (pimple.finalNonOrthogonalIter())
36  {
37  phi = pEqn.flux();
38  }
39 }
40 
41 #include "rhoEqn.H"
42 #include "compressibleContinuityErrs.H"
43 
45 U.correctBoundaryConditions();
46 fvOptions.correct(U);
47 K = 0.5*magSqr(U);
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
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