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::flux(HbyA)
9 + fvc::interpolate(rho*rAU)*fvc::ddtCorr(U, phi)
10 );
11
12 surfaceScalarField phig
13 (
14 (
15 mixture.surfaceTensionForce()
16 - ghf*fvc::snGrad(rho)
17 )*rAUf*mesh.magSf()
18 );
19
21
22 // Update the pressure BCs to ensure flux consistency
24
25 PtrList<fvScalarMatrix> p_rghEqnComps(mixture.phases().size());
26
27 label phasei = 0;
28 forAllConstIters(mixture.phases(), phase)
29 {
30 const rhoThermo& thermo = phase().thermo();
31 const volScalarField& rho = thermo.rho()();
32
34 (
35 phasei,
36 (
37 fvc::ddt(rho) + thermo.psi()*correction(fvm::ddt(p_rgh))
38 + fvc::div(phi, rho) - fvc::Sp(fvc::div(phi), rho)
39 ).ptr()
40 );
41
42 ++phasei;
43 }
44
45 // Cache p_rgh prior to solve for density update
46 volScalarField p_rgh_0(p_rgh);
47
48 while (pimple.correctNonOrthogonal())
49 {
50 fvScalarMatrix p_rghEqnIncomp
51 (
52 fvc::div(phiHbyA)
53 - fvm::laplacian(rAUf, p_rgh)
54 );
55
56 tmp<fvScalarMatrix> p_rghEqnComp;
57
58 phasei = 0;
59 forAllConstIters(mixture.phases(), phase)
60 {
61 tmp<fvScalarMatrix> hmm
62 (
63 (max(phase(), scalar(0))/phase().thermo().rho())
65 );
66
67 if (phasei == 0)
68 {
69 p_rghEqnComp = hmm;
70 }
71 else
72 {
73 p_rghEqnComp.ref() += hmm;
74 }
75
76 ++phasei;
77 }
78
79 solve
80 (
81 p_rghEqnComp
82 + p_rghEqnIncomp,
83 mesh.solver(p_rgh.select(pimple.finalInnerIter()))
84 );
85
86 if (pimple.finalNonOrthogonalIter())
87 {
88 phasei = 0;
89 for (phaseModel& phase : mixture.phases())
90 {
91 phase.dgdt() =
92 pos0(phase)
93 *(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
94
95 ++phasei;
96 }
97
98 phi = phiHbyA + p_rghEqnIncomp.flux();
99
100 U = HbyA
101 + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
102 U.correctBoundaryConditions();
103 }
104 }
105
106 p = max(p_rgh + mixture.rho()*gh, pMin);
107
108 // Update densities from change in p_rgh
109 mixture.correctRho(p_rgh - p_rgh_0);
110 rho = mixture.rho();
111
112 // Correct p_rgh for consistency with p and the updated densities
114 p_rgh.correctBoundaryConditions();
115
116 K = 0.5*magSqr(U);
117
118 Info<< "max(U) " << max(mag(U)).value() << endl;
119 Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
120}
CGAL::Exact_predicates_exact_constructions_kernel K
Y[inertIndex] max(0.0)
volScalarField & p_rgh
surfaceScalarField & phi
const surfaceScalarField & ghf
const volScalarField & gh
pimpleControl & pimple
const dimensionedScalar & pMin
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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))
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
volScalarField p_rgh_0(p_rgh)
PtrList< fvScalarMatrix > p_rghEqnComps(mixture.phases().size())
label phasei
Definition: pEqn.H:27
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
CEqn solve()
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:39
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278