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 
15  (
16  "phiHbyA",
17  (
18  (fvc::interpolate(rho*HbyA) & mesh.Sf())
20  )
21  );
22 
23  MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
24 
26  (
28  + fvc::div(phiHbyA)
29  ==
30  parcels.Srho()
31  + fvOptions(psi, p, rho.name())
32  );
33 
34  while (pimple.correctNonOrthogonal())
35  {
36  fvScalarMatrix pEqn
37  (
38  pDDtEqn
40  );
41 
42  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
43 
44  if (pimple.finalNonOrthogonalIter())
45  {
46  phi = phiHbyA + pEqn.flux();
47  }
48  }
49 
50  p.relax();
51 
52  // Second part of thermodynamic density update
53  thermo.rho() += psi*p;
54 
55  #include "rhoEqn.H" // NOTE: flux and time scales now inconsistent
56  #include "compressibleContinuityErrs.H"
57 
58  U = HbyA - rAU*fvc::grad(p);
59  U.correctBoundaryConditions();
60  fvOptions.correct(U);
61  K = 0.5*magSqr(U);
62 
63  if (thermo.dpdt())
64  {
65  dpdt = fvc::ddt(p);
66  }
67 
68  rho = thermo.rho();
69  rho = max(rho, rhoMin);
70  rho = min(rho, rhoMax);
71 
72  Info<< "p min/max = " << min(p).value() << ", " << max(p).value() << endl;
73 }
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
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Solve the continuity for density.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOMRFZoneList & MRF
fvVectorMatrix UEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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())
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
U
Definition: pEqn.H:82
HbyA
Definition: pEqn.H:7