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 + fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, Uf)
10 );
11
12 if (p_rgh.needReference())
13 {
14 fvc::makeRelative(phiHbyA, U);
16 fvc::makeAbsolute(phiHbyA, U);
17 }
18
19 surfaceScalarField phig
20 (
21 (
22 interface.surfaceTensionForce()
23 - ghf*fvc::snGrad(rho)
24 )*rAUf*mesh.magSf()
25 );
26
28
29 // Update the pressure BCs to ensure flux consistency
31
32 Pair<tmp<volScalarField>> vDotP = mixture->vDotP();
33 const volScalarField& vDotcP = vDotP[0]();
34 const volScalarField& vDotvP = vDotP[1]();
35
36 while (pimple.correctNonOrthogonal())
37 {
38 fvScalarMatrix p_rghEqn
39 (
40 fvc::div(phiHbyA) - fvm::laplacian(rAUf, p_rgh)
41 - (vDotvP - vDotcP)*(mixture->pSat() - rho*gh)
42 + fvm::Sp(vDotvP - vDotcP, p_rgh)
43 );
44
45 p_rghEqn.setReference(pRefCell, pRefValue);
46
47 p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
48
49 if (pimple.finalNonOrthogonalIter())
50 {
51 phi = phiHbyA + p_rghEqn.flux();
52
53 U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf);
54 U.correctBoundaryConditions();
55 fvOptions.correct(U);
56 }
57 }
58
59 {
60 Uf = fvc::interpolate(U);
61 surfaceVectorField n(mesh.Sf()/mesh.magSf());
62 Uf += n*(phi/mesh.magSf() - (n & Uf));
63 }
64
65 // Make the fluxes relative to the mesh motion
66 fvc::makeRelative(phi, U);
67
68 p == p_rgh + rho*gh;
69
70 if (p_rgh.needReference())
71 {
72 p += dimensionedScalar
73 (
74 "p",
75 p.dimensions(),
76 pRefValue - getRefCellValue(p, pRefCell)
77 );
78 p_rgh = p - rho*gh;
79 }
80}
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
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