Skip to content

Commit 2e2c72c

Browse files
committed
Merge branch 'dev_doury' of https://gitlab.com/phorgue/porousMultiphaseFoam into dev_doury
2 parents 0c056e2 + ea33999 commit 2e2c72c

File tree

7 files changed

+50
-11
lines changed

7 files changed

+50
-11
lines changed

solvers/impesEvapFoam/CoatsNo.H

+8-3
Original file line numberDiff line numberDiff line change
@@ -20,9 +20,14 @@
2020

2121
CFLCoats /= mesh.V();
2222
CFLUse = gMax(CFLCoats);
23-
if (runTime.value() < 500.0)
24-
maxDeltaTFact = 0.005/(CFLUse + SMALL);
25-
else
23+
if (runTime.value() < 100000000.0)
24+
maxDeltaTFact = (maxCo/3.0)/(CFLUse + SMALL);
25+
else if (runTime.value() < 1000000000.0)
26+
{
27+
scalar max = (maxCo/3.0) + (runTime.value()-100000000.0) * (maxCo - (maxCo/3.0))/(900000000.0);
28+
maxDeltaTFact = max/(CFLUse + SMALL);
29+
}
30+
else
2631
maxDeltaTFact = maxCo/(CFLUse + SMALL);
2732

2833
Info<< "Coats CFL Number mean: " << gAverage(CFLCoats) << " max: " << CFLUse << endl;

solvers/impesEvapFoam/SEqn.H

+1-1
Original file line numberDiff line numberDiff line change
@@ -40,7 +40,7 @@
4040
}
4141

4242
Info << " Lost Volume of Water" << "=" << lost_water_vol << "m3" << endl;
43-
Info << " Lost Volume of Water" << "=" << gain_water_vol << "m3" << endl;
43+
Info << " Gain Volume of Water" << "=" << gain_water_vol << "m3" << endl;
4444
Info << "Saturation b: " << " Min(Sb) = " << gMin(Sb) << " Max(Sb) = " << gMax(Sb) << endl;
4545

4646
}

solvers/impesEvapFoam/TEqn.H

+3
Original file line numberDiff line numberDiff line change
@@ -13,6 +13,8 @@
1313
// );
1414
//TEqn.solve();
1515
//Info << T << endl;
16+
17+
// Info << (hwa-hwb) * fvc::interpolate(m) << endl;
1618
fvScalarMatrix TEqn
1719
(
1820
rhoa * eps * (fvc::ddt(u_a) - fvc::ddt(u_a,Sb))
@@ -22,5 +24,6 @@
2224
== fvm::laplacian(lambda_pm,T) + (hwa-hwb) * m - wallSourceTerm
2325
);
2426
TEqn.solve();
27+
//Info << "gradT" << rhoa * fvc::div(phia,(Xw*hwa+(1-Xw)*hair)) + rhob * fvc::div(phib*(hwb)) - fvm::laplacian(lambda_pm,T) << endl;
2528

2629
}

solvers/impesEvapFoam/createFields.H

+27-2
Original file line numberDiff line numberDiff line change
@@ -133,14 +133,38 @@ dimensionedScalar hair(transportProperties.lookupOrDefault("hair",dimensionedSca
133133
// non-wetting phase properties
134134
autoPtr<incompressiblePhase> phasea = incompressiblePhase::New(mesh,transportProperties,"a");
135135
volVectorField& Ua = phasea->U();
136-
surfaceScalarField& phia = phasea->phi();
136+
//surfaceScalarField& phia = phasea->phi();
137+
surfaceScalarField phia
138+
(
139+
IOobject
140+
(
141+
"phia",
142+
runTime.timeName(),
143+
mesh,
144+
IOobject::NO_READ,
145+
IOobject::AUTO_WRITE
146+
),
147+
phasea->phi()
148+
);
137149
const dimensionedScalar& rhoa = phasea->rho();
138150
const dimensionedScalar& mua = phasea->mu();
139151

140152
// wetting phase properties
141153
autoPtr<incompressiblePhase> phaseb = incompressiblePhase::New(mesh,transportProperties,"b");
142154
volVectorField& Ub = phaseb->U();
143-
surfaceScalarField& phib = phaseb->phi();
155+
//surfaceScalarField& phib = phaseb->phi();
156+
surfaceScalarField phib
157+
(
158+
IOobject
159+
(
160+
"phib",
161+
runTime.timeName(),
162+
mesh,
163+
IOobject::NO_READ,
164+
IOobject::AUTO_WRITE
165+
),
166+
phaseb->phi()
167+
);
144168
const dimensionedScalar& rhob = phaseb->rho();
145169
const dimensionedScalar& mub = phaseb->mu();
146170

@@ -491,6 +515,7 @@ forAll(mesh.boundaryMesh(),patchi)
491515
waterMassBalanceCSV << " Xw";
492516
waterMassBalanceCSV << " p";
493517
}
518+
waterMassBalanceCSV << "Heat Flux(" << T.boundaryField()[patchi].patch().name() << ")";
494519
}
495520
waterMassBalanceCSV << " Sb_var" ;
496521
waterMassBalanceCSV << endl;

solvers/impesEvapFoam/pEqn.H

+5-4
Original file line numberDiff line numberDiff line change
@@ -1,21 +1,22 @@
11
{
2+
23
fvScalarMatrix pEqn
34
(
4-
fvm::laplacian(-Mf, p) + fvc::div(phiG)
5+
fvm::laplacian(-Mf, p) + fvc::div(phiG)
56
// capillary term
67
+ fvc::div(phiPc)*activateCapillarity
78
==
89
// event source terms
910
- sourceTerm
1011
);
11-
//pEqn.setReference(860,0);
12+
// pEqn.setReference(860,0);
1213
pEqn.solve();
1314

1415
phiP = pEqn.flux();
1516

16-
phi = phiP+phiG+phiPc*activateCapillarity;//
17+
phi = phiP+phiG+phiPc*activateCapillarity;
1718

18-
phib = Fbf*phiP + (Lbf/Lf)*phiG + phiPc*activateCapillarity ; // ;
19+
phib = Fbf*phiP + (Lbf/Lf)*phiG + phiPc*activateCapillarity;
1920
phia = phi - phib;
2021

2122
U = fvc::reconstruct(phi);

solvers/impesEvapFoam/updateSbProperties.H

+1
Original file line numberDiff line numberDiff line change
@@ -1,3 +1,4 @@
1+
12
//- relative permeability computation
23
krModel->correct(Sb);
34
kraf = fvc::interpolate(kra,"kra");

solvers/impesEvapFoam/waterMassBalance.H

+5-1
Original file line numberDiff line numberDiff line change
@@ -8,9 +8,12 @@
88
Info << phib.boundaryField()[patchi].patch().name() << " = " << gSum(phib.boundaryField()[patchi]) << " ; ";
99
}
1010
}
11-
11+
1212
Info << endl;
1313

14+
surfaceScalarField heatFlux = fvc::interpolate(-lambda_pm)*(fvc::snGrad(T)) ;//+ (hwa-hwb) * fvc::interpolate(m);
15+
16+
1417
//- write fields using usual openfoam rules
1518
runTime.write();
1619

@@ -26,6 +29,7 @@
2629
waterMassBalanceCSV << " " << gSum(Xw.boundaryField()[patchi]);
2730
waterMassBalanceCSV << " " << gSum(p.boundaryField()[patchi]);
2831
}
32+
waterMassBalanceCSV << " " << gSum(heatFlux.boundaryField()[patchi]);
2933
}
3034
waterMassBalanceCSV << " " << fvc::domainIntegrate(fvc::ddt(Sb)*eps).value();
3135
waterMassBalanceCSV << endl;

0 commit comments

Comments
 (0)