00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "multiphaseMixture.H"
00027 #include <multiphaseMixture/alphaContactAngleFvPatchScalarField.H>
00028 #include <OpenFOAM/Time.H>
00029 #include <OpenFOAM/subCycle.H>
00030 #include <finiteVolume/fvCFD.H>
00031
00032
00033
00034 const scalar Foam::multiphaseMixture::convertToRad =
00035 Foam::mathematicalConstant::pi/180.0;
00036
00037
00038
00039
00040 void Foam::multiphaseMixture::calcAlphas()
00041 {
00042 scalar level = 0.0;
00043 alphas_ == 0.0;
00044
00045 forAllIter(PtrDictionary<phase>, phases_, iter)
00046 {
00047 alphas_ += level*iter();
00048 level += 1.0;
00049 }
00050
00051 alphas_.correctBoundaryConditions();
00052 }
00053
00054
00055
00056
00057 Foam::multiphaseMixture::multiphaseMixture
00058 (
00059 const volVectorField& U,
00060 const surfaceScalarField& phi
00061 )
00062 :
00063 transportModel(U, phi),
00064 phases_(lookup("phases"), phase::iNew(U, phi)),
00065 refPhase_(*phases_.lookup(word(lookup("refPhase")))),
00066
00067 mesh_(U.mesh()),
00068 U_(U),
00069 phi_(phi),
00070
00071 rhoPhi_
00072 (
00073 IOobject
00074 (
00075 "rho*phi",
00076 mesh_.time().timeName(),
00077 mesh_,
00078 IOobject::NO_READ,
00079 IOobject::NO_WRITE
00080 ),
00081 mesh_,
00082 dimensionedScalar("rho*phi", dimMass/dimTime, 0.0)
00083 ),
00084
00085 alphas_
00086 (
00087 IOobject
00088 (
00089 "alphas",
00090 mesh_.time().timeName(),
00091 mesh_,
00092 IOobject::NO_READ,
00093 IOobject::AUTO_WRITE
00094 ),
00095 mesh_,
00096 dimensionedScalar("alphas", dimless, 0.0),
00097 zeroGradientFvPatchScalarField::typeName
00098 ),
00099
00100 sigmas_(lookup("sigmas")),
00101 dimSigma_(1, 0, -2, 0, 0),
00102 deltaN_
00103 (
00104 "deltaN",
00105 1e-8/pow(average(mesh_.V()), 1.0/3.0)
00106 )
00107 {
00108 calcAlphas();
00109 alphas_.write();
00110
00111 forAllIter(PtrDictionary<phase>, phases_, iter)
00112 {
00113 alphaTable_.add(iter());
00114 }
00115 }
00116
00117
00118
00119
00120 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::rho() const
00121 {
00122 PtrDictionary<phase>::const_iterator iter = phases_.begin();
00123
00124 tmp<volScalarField> trho = iter()*iter().rho();
00125
00126 for(++iter; iter != phases_.end(); ++iter)
00127 {
00128 trho() += iter()*iter().rho();
00129 }
00130
00131 return trho;
00132 }
00133
00134
00135 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::mu() const
00136 {
00137 PtrDictionary<phase>::const_iterator iter = phases_.begin();
00138
00139 tmp<volScalarField> tmu = iter()*iter().rho()*iter().nu();
00140
00141 for(++iter; iter != phases_.end(); ++iter)
00142 {
00143 tmu() += iter()*iter().rho()*iter().nu();
00144 }
00145
00146 return tmu;
00147 }
00148
00149
00150 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::muf() const
00151 {
00152 PtrDictionary<phase>::const_iterator iter = phases_.begin();
00153
00154 tmp<surfaceScalarField> tmuf =
00155 fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
00156
00157 for(++iter; iter != phases_.end(); ++iter)
00158 {
00159 tmuf() +=
00160 fvc::interpolate(iter())*iter().rho()*fvc::interpolate(iter().nu());
00161 }
00162
00163 return tmuf;
00164 }
00165
00166
00167 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::nu() const
00168 {
00169 return mu()/rho();
00170 }
00171
00172
00173 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::nuf() const
00174 {
00175 return muf()/fvc::interpolate(rho());
00176 }
00177
00178
00179 Foam::tmp<Foam::surfaceScalarField>
00180 Foam::multiphaseMixture::surfaceTensionForce() const
00181 {
00182 tmp<surfaceScalarField> tstf
00183 (
00184 new surfaceScalarField
00185 (
00186 IOobject
00187 (
00188 "surfaceTensionForce",
00189 mesh_.time().timeName(),
00190 mesh_
00191 ),
00192 mesh_,
00193 dimensionedScalar
00194 (
00195 "surfaceTensionForce",
00196 dimensionSet(1, -2, -2, 0, 0),
00197 0.0
00198 )
00199 )
00200 );
00201
00202 surfaceScalarField& stf = tstf();
00203
00204 forAllConstIter(PtrDictionary<phase>, phases_, iter1)
00205 {
00206 const phase& alpha1 = iter1();
00207
00208 PtrDictionary<phase>::const_iterator iter2 = iter1;
00209 ++iter2;
00210
00211 for(; iter2 != phases_.end(); ++iter2)
00212 {
00213 const phase& alpha2 = iter2();
00214
00215 sigmaTable::const_iterator sigma =
00216 sigmas_.find(interfacePair(alpha1, alpha2));
00217
00218 if (sigma == sigmas_.end())
00219 {
00220 FatalErrorIn("multiphaseMixture::surfaceTensionForce() const")
00221 << "Cannot find interface " << interfacePair(alpha1, alpha2)
00222 << " in list of sigma values"
00223 << exit(FatalError);
00224 }
00225
00226 stf += dimensionedScalar("sigma", dimSigma_, sigma())
00227 *fvc::interpolate(K(alpha1, alpha2))*
00228 (
00229 fvc::interpolate(alpha2)*fvc::snGrad(alpha1)
00230 - fvc::interpolate(alpha1)*fvc::snGrad(alpha2)
00231 );
00232 }
00233 }
00234
00235 return tstf;
00236 }
00237
00238
00239 void Foam::multiphaseMixture::solve()
00240 {
00241 forAllIter(PtrDictionary<phase>, phases_, iter)
00242 {
00243 iter().correct();
00244 }
00245
00246 const Time& runTime = mesh_.time();
00247
00248 label nAlphaSubCycles
00249 (
00250 readLabel
00251 (
00252 mesh_.solutionDict().subDict("PISO").lookup("nAlphaSubCycles")
00253 )
00254 );
00255
00256 label nAlphaCorr
00257 (
00258 readLabel(mesh_.solutionDict().subDict("PISO").lookup("nAlphaCorr"))
00259 );
00260
00261 bool cycleAlpha
00262 (
00263 Switch(mesh_.solutionDict().subDict("PISO").lookup("cycleAlpha"))
00264 );
00265
00266 scalar cAlpha
00267 (
00268 readScalar(mesh_.solutionDict().subDict("PISO").lookup("cAlpha"))
00269 );
00270
00271
00272 volScalarField& alpha = phases_.first();
00273
00274 if (nAlphaSubCycles > 1)
00275 {
00276 surfaceScalarField rhoPhiSum = 0.0*rhoPhi_;
00277 dimensionedScalar totalDeltaT = runTime.deltaT();
00278
00279 for
00280 (
00281 subCycle<volScalarField> alphaSubCycle(alpha, nAlphaSubCycles);
00282 !(++alphaSubCycle).end();
00283 )
00284 {
00285 solveAlphas(nAlphaCorr, cycleAlpha, cAlpha);
00286 rhoPhiSum += (runTime.deltaT()/totalDeltaT)*rhoPhi_;
00287 }
00288
00289 rhoPhi_ = rhoPhiSum;
00290 }
00291 else
00292 {
00293 solveAlphas(nAlphaCorr, cycleAlpha, cAlpha);
00294 }
00295 }
00296
00297
00298 void Foam::multiphaseMixture::correct()
00299 {}
00300
00301
00302 Foam::tmp<Foam::surfaceVectorField> Foam::multiphaseMixture::nHatfv
00303 (
00304 const volScalarField& alpha1,
00305 const volScalarField& alpha2
00306 ) const
00307 {
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317 surfaceVectorField gradAlphaf =
00318 fvc::interpolate(alpha2)*fvc::interpolate(fvc::grad(alpha1))
00319 - fvc::interpolate(alpha1)*fvc::interpolate(fvc::grad(alpha2));
00320
00321
00322 return gradAlphaf/(mag(gradAlphaf) + deltaN_);
00323 }
00324
00325
00326 Foam::tmp<Foam::surfaceScalarField> Foam::multiphaseMixture::nHatf
00327 (
00328 const volScalarField& alpha1,
00329 const volScalarField& alpha2
00330 ) const
00331 {
00332
00333 return nHatfv(alpha1, alpha2) & mesh_.Sf();
00334 }
00335
00336
00337
00338
00339
00340
00341
00342
00343 void Foam::multiphaseMixture::correctContactAngle
00344 (
00345 const phase& alpha1,
00346 const phase& alpha2,
00347 surfaceVectorField::GeometricBoundaryField& nHatb
00348 ) const
00349 {
00350 const volScalarField::GeometricBoundaryField& gbf
00351 = refPhase_.boundaryField();
00352
00353 const fvBoundaryMesh& boundary = mesh_.boundary();
00354
00355 forAll(boundary, patchi)
00356 {
00357 if (isA<alphaContactAngleFvPatchScalarField>(gbf[patchi]))
00358 {
00359 const alphaContactAngleFvPatchScalarField& acap =
00360 refCast<const alphaContactAngleFvPatchScalarField>(gbf[patchi]);
00361
00362 vectorField& nHatPatch = nHatb[patchi];
00363
00364 vectorField AfHatPatch =
00365 mesh_.Sf().boundaryField()[patchi]
00366 /mesh_.magSf().boundaryField()[patchi];
00367
00368 alphaContactAngleFvPatchScalarField::thetaPropsTable::
00369 const_iterator tp =
00370 acap.thetaProps().find(interfacePair(alpha1, alpha2));
00371
00372 if (tp == acap.thetaProps().end())
00373 {
00374 FatalErrorIn
00375 (
00376 "multiphaseMixture::correctContactAngle"
00377 "(const phase& alpha1, const phase& alpha2, "
00378 "fvPatchVectorFieldField& nHatb) const"
00379 ) << "Cannot find interface " << interfacePair(alpha1, alpha2)
00380 << "\n in table of theta properties for patch "
00381 << acap.patch().name()
00382 << exit(FatalError);
00383 }
00384
00385 bool matched = (tp.key().first() == alpha1.name());
00386
00387 scalar theta0 = convertToRad*tp().theta0(matched);
00388 scalarField theta(boundary[patchi].size(), theta0);
00389
00390 scalar uTheta = tp().uTheta();
00391
00392
00393 if (uTheta > SMALL)
00394 {
00395 scalar thetaA = convertToRad*tp().thetaA(matched);
00396 scalar thetaR = convertToRad*tp().thetaR(matched);
00397
00398
00399 vectorField Uwall =
00400 U_.boundaryField()[patchi].patchInternalField()
00401 - U_.boundaryField()[patchi];
00402 Uwall -= (AfHatPatch & Uwall)*AfHatPatch;
00403
00404
00405 vectorField nWall =
00406 nHatPatch - (AfHatPatch & nHatPatch)*AfHatPatch;
00407
00408
00409 nWall /= (mag(nWall) + SMALL);
00410
00411
00412
00413 scalarField uwall = nWall & Uwall;
00414
00415 theta += (thetaA - thetaR)*tanh(uwall/uTheta);
00416 }
00417
00418
00419
00420
00421 scalarField a12 = nHatPatch & AfHatPatch;
00422
00423 scalarField b1 = cos(theta);
00424
00425 scalarField b2(nHatPatch.size());
00426
00427 forAll(b2, facei)
00428 {
00429 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
00430 }
00431
00432 scalarField det = 1.0 - a12*a12;
00433
00434 scalarField a = (b1 - a12*b2)/det;
00435 scalarField b = (b2 - a12*b1)/det;
00436
00437 nHatPatch = a*AfHatPatch + b*nHatPatch;
00438
00439 nHatPatch /= (mag(nHatPatch) + deltaN_.value());
00440 }
00441 }
00442 }
00443
00444
00445 Foam::tmp<Foam::volScalarField> Foam::multiphaseMixture::K
00446 (
00447 const phase& alpha1,
00448 const phase& alpha2
00449 ) const
00450 {
00451 tmp<surfaceVectorField> tnHatfv = nHatfv(alpha1, alpha2);
00452
00453 correctContactAngle(alpha1, alpha2, tnHatfv().boundaryField());
00454
00455
00456 return -fvc::div(tnHatfv & mesh_.Sf());
00457 }
00458
00459
00460 Foam::tmp<Foam::surfaceScalarField>
00461 Foam::multiphaseMixture::nearInterface() const
00462 {
00463 tmp<surfaceScalarField> tnearInt
00464 (
00465 new surfaceScalarField
00466 (
00467 IOobject
00468 (
00469 "nearInterface",
00470 mesh_.time().timeName(),
00471 mesh_
00472 ),
00473 mesh_,
00474 dimensionedScalar("nearInterface", dimless, 0.0)
00475 )
00476 );
00477
00478 forAllConstIter(PtrDictionary<phase>, phases_, iter)
00479 {
00480 surfaceScalarField alphaf = fvc::interpolate(iter());
00481 tnearInt() = max(tnearInt(), pos(alphaf - 0.01)*pos(0.99 - alphaf));
00482 }
00483
00484 return tnearInt;
00485 }
00486
00487
00488 void Foam::multiphaseMixture::solveAlphas
00489 (
00490 const label nAlphaCorr,
00491 const bool cycleAlpha,
00492 const scalar cAlpha
00493 )
00494 {
00495 static label nSolves=-1;
00496 nSolves++;
00497
00498 word alphaScheme("div(phi,alpha)");
00499 word alphacScheme("div(phirb,alpha)");
00500
00501 tmp<fv::convectionScheme<scalar> > mvConvection
00502 (
00503 fv::convectionScheme<scalar>::New
00504 (
00505 mesh_,
00506 alphaTable_,
00507 phi_,
00508 mesh_.divScheme(alphaScheme)
00509 )
00510 );
00511
00512 surfaceScalarField phic = mag(phi_/mesh_.magSf());
00513 phic = min(cAlpha*phic, max(phic));
00514
00515 for (int gCorr=0; gCorr<nAlphaCorr; gCorr++)
00516 {
00517 phase* refPhasePtr = &refPhase_;
00518
00519 if (cycleAlpha)
00520 {
00521 PtrDictionary<phase>::iterator refPhaseIter = phases_.begin();
00522 for(label i=0; i<nSolves%phases_.size(); i++)
00523 {
00524 ++refPhaseIter;
00525 }
00526 refPhasePtr = &refPhaseIter();
00527 }
00528
00529 phase& refPhase = *refPhasePtr;
00530
00531 volScalarField refPhaseNew = refPhase;
00532 refPhaseNew == 1.0;
00533
00534 rhoPhi_ = phi_*refPhase.rho();
00535
00536 forAllIter(PtrDictionary<phase>, phases_, iter)
00537 {
00538 phase& alpha = iter();
00539
00540 if (&alpha == &refPhase) continue;
00541
00542 fvScalarMatrix alphaEqn
00543 (
00544 fvm::ddt(alpha)
00545 + mvConvection->fvmDiv(phi_, alpha)
00546 );
00547
00548 forAllIter(PtrDictionary<phase>, phases_, iter2)
00549 {
00550 phase& alpha2 = iter2();
00551
00552 if (&alpha2 == &alpha) continue;
00553
00554 surfaceScalarField phir = phic*nHatf(alpha, alpha2);
00555 surfaceScalarField phirb12 =
00556 -fvc::flux(-phir, alpha2, alphacScheme);
00557
00558 alphaEqn += fvm::div(phirb12, alpha, alphacScheme);
00559 }
00560
00561 alphaEqn.solve(mesh_.solver("alpha"));
00562
00563 rhoPhi_ += alphaEqn.flux()*(alpha.rho() - refPhase.rho());
00564
00565 Info<< alpha.name() << " volume fraction, min, max = "
00566 << alpha.weightedAverage(mesh_.V()).value()
00567 << ' ' << min(alpha).value()
00568 << ' ' << max(alpha).value()
00569 << endl;
00570
00571 refPhaseNew == refPhaseNew - alpha;
00572 }
00573
00574 refPhase == refPhaseNew;
00575
00576 Info<< refPhase.name() << " volume fraction, min, max = "
00577 << refPhase.weightedAverage(mesh_.V()).value()
00578 << ' ' << min(refPhase).value()
00579 << ' ' << max(refPhase).value()
00580 << endl;
00581 }
00582
00583 calcAlphas();
00584 }
00585
00586
00587 bool Foam::multiphaseMixture::read()
00588 {
00589 if (transportModel::read())
00590 {
00591 bool readOK = true;
00592
00593 PtrList<entry> phaseData(lookup("phases"));
00594 label phasei = 0;
00595
00596 forAllIter(PtrDictionary<phase>, phases_, iter)
00597 {
00598 readOK &= iter().read(phaseData[phasei++].dict());
00599 }
00600
00601 lookup("sigmas") >> sigmas_;
00602
00603 return readOK;
00604 }
00605 else
00606 {
00607 return false;
00608 }
00609 }
00610
00611
00612