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
13 {
14 surfaceScalarField faceMaskOld
15 (
16 localMin<scalar>(mesh).interpolate(cellMask.oldTime())
17 );
18 phiHbyA += faceMaskOld*fvc::ddtCorr(U, Uf);
19 }
20
21 if (p_rgh.needReference())
22 {
23 fvc::makeRelative(phiHbyA, U);
25 fvc::makeAbsolute(phiHbyA, U);
26 }
27
28 surfaceScalarField phig
29 (
30 (
31 interface.surfaceTensionForce()
32 - ghf*fvc::snGrad(rho)
33 )*faceMask*rAUf*mesh.magSf()
34 );
35
37
38 // Update the pressure BCs to ensure flux consistency
40
41 Pair<tmp<volScalarField>> vDotP = mixture->vDotP();
42 const volScalarField& vDotcP = vDotP[0]();
43 const volScalarField& vDotvP = vDotP[1]();
44
45 while (pimple.correctNonOrthogonal())
46 {
47 fvScalarMatrix p_rghEqn
48 (
49 fvc::div(phiHbyA) - fvm::laplacian(rAUf, p_rgh)
50 - (vDotvP - vDotcP)*(mixture->pSat() - rho*gh)
51 + fvm::Sp(vDotvP - vDotcP, p_rgh)
52 );
53
54
55 //p_rghEqn.setReference(pRefCell, pRefValue);
56 p_rghEqn.setReference(pRefCell, getRefCellValue(p_rgh, pRefCell));
57
58 p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
59
60 if (pimple.finalNonOrthogonalIter())
61 {
62 phi = phiHbyA + p_rghEqn.flux();
63
64 p_rgh.relax();
65
66 U =
67 cellMask
68 *(HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf));
69
70 U.correctBoundaryConditions();
71 fvOptions.correct(U);
72 }
73 }
74
75 #include "continuityErrs.H"
76
77 {
78 Uf = fvc::interpolate(U);
79 surfaceVectorField n(mesh.Sf()/mesh.magSf());
80 Uf += n*(phi/mesh.magSf() - (n & Uf));
81 }
82
83 // Make the fluxes relative to the mesh motion
84 fvc::makeRelative(phi, U);
85
86 // Zero faces H-I for transport Eq after pEq
87 phi *= faceMask;
88
89 p == p_rgh + rho*gh;
90
91 if (p_rgh.needReference())
92 {
93 p += dimensionedScalar
94 (
95 "p",
96 p.dimensions(),
97 pRefValue - getRefCellValue(p, pRefCell)
98 );
99 p_rgh = p - rho*gh;
100 }
101}
label n
volScalarField & p_rgh
fv::options & fvOptions
surfaceScalarField & phi
const scalar pRefValue
const surfaceScalarField & ghf
const label pRefCell
const volScalarField & gh
pimpleControl & pimple
U
Definition: pEqn.H:72
volScalarField & p
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
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
mesh interpolate(rAU)
ddtCorr
Definition: readControls.H:9
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
adjustPhi(phiHbyA, U, p_rgh)
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
const volScalarField & vDotcP
Definition: pEqn.H:33
Pair< tmp< volScalarField > > vDotP
Definition: pEqn.H:32
const volScalarField & vDotvP
Definition: pEqn.H:34
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39