pEqn.H
Go to the documentation of this file.
1 rho = thermo.rho();
2 
3 // Thermodynamic density needs to be updated by psi*d(p) after the
4 // pressure solution
6 
7 volScalarField rAU(1.0/UEqn.A());
11 (
12  "phiHbyA",
13  (
14  fvc::flux(rho*HbyA)
16  )
17 );
18 
19 MRF.makeRelative(fvc::interpolate(rho), phiHbyA);
20 
21 // Update the pressure BCs to ensure flux consistency
23 
25 (
27  + fvc::div(phiHbyA)
28  ==
29  parcels.Srho()
30  + fvOptions(psi, p, rho.name())
31 );
32 
33 while (pimple.correctNonOrthogonal())
34 {
35  fvScalarMatrix pEqn
36  (
37  pDDtEqn
39  );
40 
41  pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
42 
43  if (pimple.finalNonOrthogonalIter())
44  {
45  phi = phiHbyA + pEqn.flux();
46  }
47 }
48 
49 p.relax();
50 
51 // Thermodynamic density update
52 thermo.correctRho(psi*p - psip0);
53 
54 #include "rhoEqn.H" // NOTE: flux and time scales now inconsistent
55 #include "compressibleContinuityErrs.H"
56 
57 U = HbyA - rAU*fvc::grad(p);
58 U.correctBoundaryConditions();
59 fvOptions.correct(U);
60 K = 0.5*magSqr(U);
61 
62 if (thermo.dpdt())
63 {
64  dpdt = fvc::ddt(p);
65 }
66 
67 rho = thermo.rho();
68 rho = max(rho, rhoMin);
69 rho = min(rho, rhoMax);
70 
71 Info<< "p min/max = " << min(p).value() << ", " << max(p).value() << endl;
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
PtrList< dimensionedScalar > rhoMax(fluidRegions.size())
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
Solve the continuity for density.
p
Definition: pEqn.H:50
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
phiHbyA
Definition: pEqn.H:20
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
PtrList< dimensionedScalar > rhoMin(fluidRegions.size())
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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))
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
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
fvScalarMatrix pDDtEqn(fvc::ddt(rho)+psi *correction(fvm::ddt(p))+fvc::div(phiHbyA)==parcels.Srho()+fvOptions(psi, p, rho.name()))
dimensioned< scalar > magSqr(const dimensioned< Type > &)
constrainPressure(p, rho, U, phiHbyA, rhorAUf, MRF)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
volVectorField & HbyA
Definition: pEqn.H:13
fvVectorMatrix & UEqn
Definition: UEqn.H:13
messageStream Info
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
const volScalarField psip0(psi *p)
volScalarField rAU(1.0/UEqn.A())