pEqn.H
Go to the documentation of this file.
1{
2 volScalarField rAU("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::interpolate(HbyA) & mesh.Sf())
9 + fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, phi)
10 );
12
13 surfaceScalarField phig
14 (
15 (
16 interface.surfaceTensionForce()
17 - ghf*fvc::snGrad(rho)
18 )*rAUf*mesh.magSf()
19 );
20
22
23 // Update the pressure BCs to ensure flux consistency
25
26 Pair<tmp<volScalarField>> vDot = mixture->vDot();
27 const volScalarField& vDotc = vDot[0]();
28 const volScalarField& vDotv = vDot[1]();
29
30 while (pimple.correctNonOrthogonal())
31 {
32 fvScalarMatrix p_rghEqn
33 (
34 fvc::div(phiHbyA)
35 - fvm::laplacian(rAUf, p_rgh)
36 ==
37 vDotv + vDotc
38 );
39
40 p_rghEqn.setReference(pRefCell, pRefValue);
41
42 p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
43
44 if (pimple.finalNonOrthogonalIter())
45 {
46 phi = phiHbyA + p_rghEqn.flux();
47
48 U = HbyA + rAU*fvc::reconstruct((phig + p_rghEqn.flux())/rAUf);
49 U.correctBoundaryConditions();
50 fvOptions.correct(U);
51 }
52 }
53
55
56 if (p_rgh.needReference())
57 {
58 p += dimensionedScalar
59 (
60 "p",
61 p.dimensions(),
62 pRefValue - getRefCellValue(p, pRefCell)
63 );
64 p_rgh = p - rho*gh;
65 }
66}
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
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
adjustPhi(phiHbyA, U, p_rgh)
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
const volScalarField & vDotv
Definition: pEqn.H:28
const volScalarField & vDotc
Definition: pEqn.H:27
Pair< tmp< volScalarField > > vDot
Definition: pEqn.H:26
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
interfaceProperties interface(alpha1, U, thermo->transportPropertiesDict())
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39