Loading...
Searching...
No Matches
pEqn.H
Go to the documentation of this file.
1{
2 rAU = 1.0/UEqn.A();
3 surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
4 volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p_rgh));
5 surfaceScalarField phiHbyA
6 (
7 "phiHbyA",
8 fvc::flux(HbyA)
9 );
10
11 if (p_rgh.needReference())
12 {
13 fvc::makeRelative(phiHbyA, U);
15 fvc::makeAbsolute(phiHbyA, U);
16 }
17
18 surfaceScalarField phig
19 (
20 (
21 interface.surfaceTensionForce()
22 - ghf*fvc::snGrad(rho)
23 )*faceMask*rAUf*mesh.magSf()
24 );
25
27
28 // Update the pressure BCs to ensure flux consistency
30
31 Pair<tmp<volScalarField>> vDotP = mixture->vDotP();
32 const volScalarField& vDotcP = vDotP[0]();
33 const volScalarField& vDotvP = vDotP[1]();
34
35 while (pimple.correctNonOrthogonal())
36 {
37 fvScalarMatrix p_rghEqn
38 (
39 fvc::div(phiHbyA) - fvm::laplacian(rAUf, p_rgh)
40 - (vDotvP - vDotcP)*(mixture->pSat() - rho*gh)
41 + fvm::Sp(vDotvP - vDotcP, p_rgh)
42 );
43
44
45 //p_rghEqn.setReference(pRefCell, pRefValue);
46 p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
47
48 p_rghEqn.solve(p_rgh.select(pimple.finalInnerIter()));
49
50 if (pimple.finalNonOrthogonalIter())
51 {
52 phi = phiHbyA + p_rghEqn.flux();
53
54 p_rgh.relax();
55
56 U =
57 cellMask
58 *(HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf));
59
60 U.correctBoundaryConditions();
61 fvOptions.correct(U);
62 }
63 }
64
65 #include "continuityErrs.H"
66
67 {
68 Uf = fvc::interpolate(U);
69 surfaceVectorField n(mesh.Sf()/mesh.magSf());
70 Uf += n*(phi/mesh.magSf() - (n & Uf));
71 }
72
73 // Make the fluxes relative to the mesh motion
74 fvc::makeRelative(phi, U);
75
76 // Zero faces H-I for transport Eq after pEq
77 phi *= faceMask;
78
79 p == p_rgh + rho*gh;
80
81 if (p_rgh.needReference())
82 {
83 p += dimensionedScalar
84 (
85 "p",
86 p.dimensions(),
87 pRefValue - getRefCellValue(p, pRefCell)
88 );
89 p_rgh = p - rho*gh;
90 }
91}
label n
volScalarField & p_rgh
fv::options & fvOptions
const scalar pRefValue
const surfaceScalarField & ghf
const label pRefCell
const volScalarField & gh
pimpleControl & pimple
U
Definition pEqn.H:72
volScalarField & p
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
fvVectorMatrix & UEqn
Definition UEqn.H:13
phiHbyA
Definition pcEqn.H:73
HbyA
Definition pcEqn.H:74
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
adjustPhi(phiHbyA, U, p_rgh)
tmp< volScalarField > rAU
const volScalarField & vDotcP
Definition pEqn.H:33
Pair< tmp< volScalarField > > vDotP
Definition pEqn.H:32
const volScalarField & vDotvP
Definition pEqn.H:34
faceMask
Definition setCellMask.H:43
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture