Loading...
Searching...
No Matches
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 tmp<volScalarField> trho(thermo.rho());
32 const volScalarField& rho = trho();
33
35 (
36 phasei,
37 (
38 fvc::ddt(rho) + thermo.psi()*correction(fvm::ddt(p_rgh))
39 + fvc::div(phi, rho) - fvc::Sp(fvc::div(phi), rho)
40 ).ptr()
41 );
42
43 ++phasei;
44 }
45
46 // Cache p_rgh prior to solve for density update
47 volScalarField p_rgh_0(p_rgh);
48
49 while (pimple.correctNonOrthogonal())
50 {
51 fvScalarMatrix p_rghEqnIncomp
52 (
53 fvc::div(phiHbyA)
54 - fvm::laplacian(rAUf, p_rgh)
55 );
56
57 tmp<fvScalarMatrix> p_rghEqnComp;
58
59 phasei = 0;
60 forAllConstIters(mixture.phases(), phase)
61 {
62 tmp<fvScalarMatrix> hmm
63 (
64 (max(phase(), scalar(0))/phase().thermo().rho())
66 );
67
68 if (phasei == 0)
69 {
70 p_rghEqnComp = hmm;
71 }
72 else
73 {
74 p_rghEqnComp.ref() += hmm;
75 }
76
77 ++phasei;
78 }
79
80 solve
81 (
82 p_rghEqnComp
83 + p_rghEqnIncomp,
84 p_rgh.select(pimple.finalInnerIter())
85 );
86
87 if (pimple.finalNonOrthogonalIter())
88 {
89 phasei = 0;
90 for (phaseModel& phase : mixture.phases())
91 {
92 phase.dgdt() =
93 pos0(phase)
94 *(p_rghEqnComps[phasei] & p_rgh)/phase.thermo().rho();
95
96 ++phasei;
97 }
98
99 phi = phiHbyA + p_rghEqnIncomp.flux();
100
101 U = HbyA
102 + rAU*fvc::reconstruct((phig + p_rghEqnIncomp.flux())/rAUf);
103 U.correctBoundaryConditions();
104 }
105 }
106
107 p = max(p_rgh + mixture.rho()*gh, pMin);
108
109 // Update densities from change in p_rgh
110 mixture.correctRho(p_rgh - p_rgh_0);
111 rho = mixture.rho();
112
113 // Correct p_rgh for consistency with p and the updated densities
115 p_rgh.correctBoundaryConditions();
116
117 K = 0.5*magSqr(U);
118
119 Info<< "max(U) " << max(mag(U)).value() << endl;
120 Info<< "min(p_rgh) " << min(p_rgh).value() << endl;
121}
CGAL::Exact_predicates_exact_constructions_kernel K
volScalarField & p_rgh
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
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
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
tmp< volScalarField > rAU
volScalarField p_rgh_0(p_rgh)
PtrList< fvScalarMatrix > p_rghEqnComps(mixture.phases().size())
label phasei
Definition pEqn.H:27
rho
Definition pEqn.H:111
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tmp< volScalarField > trho
CEqn solve()
psiReactionThermo & thermo
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235