forked from EastEriq/phaseFieldFoam-6
-
Notifications
You must be signed in to change notification settings - Fork 1
/
pEqn.H
74 lines (58 loc) · 2.01 KB
/
pEqn.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
{
//-Calculate descretization coefficients based on the definition of Ueqn
volScalarField rAU("rAU",1.0/UEqn.A());
surfaceScalarField rAUf = fvc::interpolate(rAU);
surfaceScalarField rhof = twoPhaseProperties.rhoMixF(scalar(0.5)*(alpha1.oldTime() + alpha1));
//-Solve for the velocity based on the velocity terms described in UEqn.
volVectorField HbyA("HbyA",U);
HbyA = rAU*UEqn.H();
/*
// ddtPhiCorr //
This term accounts for the divergence of the face velocity field by taking out the difference between the interpolated velocity and the flux
*/
surfaceScalarField phiHbyA
(
"phiHbyA",
(fvc::interpolate(HbyA) & mesh.Sf())
+ fvc::interpolate(rho*rAU)*fvc::ddtCorr(U,phi)
);
adjustPhi(phiHbyA,U,p_rgh);
phi = phiHbyA;
surfaceScalarField phig
(
(KSaG)*rAUf
);
phiHbyA += phig;
//-Update the fixedFluxPressure BCs to ensure flux consistency
// from https://www.cfd-online.com/Forums/openfoam-programming-development/187041-problem-compiling-solver-made-different-version-v2-0-v4-1-a.html
constrainPressure(p_rgh, U, phiHbyA, rAUf);
while (pimple.correctNonOrthogonal())
{
fvScalarMatrix p_rghEqn
(
fvm::laplacian(rAUf,p_rgh) == fvc::div(phi)
);
p_rghEqn.setReference(pRefCell,getRefCellValue(p_rgh,pRefCell));
p_rghEqn.solve(mesh.solver(p_rgh.select(pimple.finalInnerIter())));
if (pimple.finalNonOrthogonalIter())
{
phi = phiHbyA - p_rghEqn.flux();
U = HbyA + rAU*fvc::reconstruct((phig - p_rghEqn.flux())/rAUf);
U.correctBoundaryConditions();
fvOptions.correct(U);
}
}
fvc::makeRelative(phi,U);
#include "continuityErrs.H"
p == p_rgh + rho*gh;
if (p_rgh.needReference())
{
p += dimensionedScalar
(
"p",
p.dimensions(),
pRefValue - getRefCellValue(p,pRefCell)
);
p_rgh = p - rho*gh;
}
}