pEqn.H
Go to the documentation of this file.
1{
2 if (pimple.nCorrPIMPLE() == 1)
3 {
4 p =
5 (
6 rho
7 - alphal*rhol0
8 - ((alphav*psiv + alphal*psil) - psi)*pSat
9 )/psi;
10 }
11
12 surfaceScalarField rhof("rhof", fvc::interpolate(rho));
13
14 volScalarField rAU(1.0/UEqn.A());
15 surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho*rAU));
16 volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
17
18 phi = fvc::flux(HbyA)
19 + rhorAUf*fvc::ddtCorr(U, Uf);
20 fvc::makeRelative(phi, U);
21
22 surfaceScalarField phiGradp(rhorAUf*mesh.magSf()*fvc::snGrad(p));
23
24 phi -= phiGradp/rhof;
25
26 volScalarField rho0(rho - psi*p);
27
28 while (pimple.correctNonOrthogonal())
29 {
30 fvScalarMatrix pEqn
31 (
32 fvc::ddt(rho)
33 + psi*correction(fvm::ddt(p))
34 + fvc::div(phi, rho)
35 + fvc::div(phiGradp)
36 - fvm::laplacian(rhorAUf, p)
37 );
38
39 pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
40
41 if (pimple.finalNonOrthogonalIter())
42 {
43 phi += (phiGradp + pEqn.flux())/rhof;
44 }
45 }
46
47 Info<< "Predicted p max-min : " << max(p).value()
48 << " " << min(p).value() << endl;
49
50 rho == max(rho0 + psi*p, rhoMin);
51
52 #include "alphavPsi.H"
53
54 p =
55 (
56 rho
57 - alphal*rhol0
58 - ((alphav*psiv + alphal*psil) - psi)*pSat
59 )/psi;
60
61 p.correctBoundaryConditions();
62
63 Info<< "Phase-change corrected p max-min : " << max(p).value()
64 << " " << min(p).value() << endl;
65
66 // Correct velocity
67
68 U = HbyA - rAU*fvc::grad(p);
69
70 // Remove the swirl component of velocity for "wedge" cases
71 if (pimple.dict().found("removeSwirl"))
72 {
73 label swirlCmpt(pimple.dict().get<label>("removeSwirl"));
74
75 Info<< "Removing swirl component-" << swirlCmpt << " of U" << endl;
76 U.field().replace(swirlCmpt, Zero);
77 }
78
79 U.correctBoundaryConditions();
80
81 Info<< "max(U) " << max(mag(U)).value() << endl;
82
83 {
84 Uf = fvc::interpolate(U);
85 surfaceVectorField n(mesh.Sf()/mesh.magSf());
86 Uf += n*(phi/mesh.magSf() - (n & Uf));
87 }
88}
label n
Y[inertIndex] max(0.0)
alphal
Definition: alphavPsi.H:12
surfaceScalarField & phi
pimpleControl & pimple
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & psi
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
fvVectorMatrix & UEqn
Definition: UEqn.H:13
HbyA
Definition: pcEqn.H:74
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
const dimensionedScalar rhoMin
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
surfaceScalarField phiGradp(rhorAUf *mesh.magSf() *fvc::snGrad(p))
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
scalar rho0