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 "MULES.H"
00027 #include <finiteVolume/upwind.H>
00028 #include <finiteVolume/uncorrectedSnGrad.H>
00029 #include <finiteVolume/gaussConvectionScheme.H>
00030 #include <finiteVolume/gaussLaplacianScheme.H>
00031 #include <finiteVolume/uncorrectedSnGrad.H>
00032 #include <finiteVolume/surfaceInterpolate.H>
00033 #include <finiteVolume/fvcSurfaceIntegrate.H>
00034 #include <finiteVolume/slicedSurfaceFields.H>
00035 #include <OpenFOAM/syncTools.H>
00036
00037 #include <finiteVolume/fvCFD.H>
00038
00039
00040
00041 template<class RhoType, class SpType, class SuType>
00042 void Foam::MULES::explicitSolve
00043 (
00044 const RhoType& rho,
00045 volScalarField& psi,
00046 const surfaceScalarField& phi,
00047 surfaceScalarField& phiPsi,
00048 const SpType& Sp,
00049 const SuType& Su,
00050 const scalar psiMax,
00051 const scalar psiMin
00052 )
00053 {
00054 Info<< "MULES: Solving for " << psi.name() << endl;
00055
00056 const fvMesh& mesh = psi.mesh();
00057 psi.correctBoundaryConditions();
00058
00059 surfaceScalarField phiBD = upwind<scalar>(psi.mesh(), phi).flux(psi);
00060
00061 surfaceScalarField& phiCorr = phiPsi;
00062 phiCorr -= phiBD;
00063
00064 scalarField allLambda(mesh.nFaces(), 1.0);
00065
00066 slicedSurfaceScalarField lambda
00067 (
00068 IOobject
00069 (
00070 "lambda",
00071 mesh.time().timeName(),
00072 mesh,
00073 IOobject::NO_READ,
00074 IOobject::NO_WRITE,
00075 false
00076 ),
00077 mesh,
00078 dimless,
00079 allLambda,
00080 false
00081 );
00082
00083 limiter
00084 (
00085 allLambda,
00086 rho,
00087 psi,
00088 phiBD,
00089 phiCorr,
00090 Sp.field(),
00091 Su.field(),
00092 psiMax,
00093 psiMin,
00094 3
00095 );
00096
00097 phiPsi = phiBD + lambda*phiCorr;
00098
00099 scalarField& psiIf = psi;
00100 const scalarField& psi0 = psi.oldTime();
00101 const scalar deltaT = mesh.time().deltaT().value();
00102
00103 psiIf = 0.0;
00104 fvc::surfaceIntegrate(psiIf, phiPsi);
00105
00106 if (mesh.moving())
00107 {
00108 psiIf =
00109 (
00110 mesh.Vsc0()*rho.oldTime()*psi0/(deltaT*mesh.Vsc())
00111 + Su.field()
00112 - psiIf
00113 )/(rho/deltaT - Sp.field());
00114 }
00115 else
00116 {
00117 psiIf =
00118 (
00119 rho.oldTime()*psi0/deltaT
00120 + Su.field()
00121 - psiIf
00122 )/(rho/deltaT - Sp.field());
00123 }
00124
00125 psi.correctBoundaryConditions();
00126 }
00127
00128
00129 template<class RhoType, class SpType, class SuType>
00130 void Foam::MULES::implicitSolve
00131 (
00132 const RhoType& rho,
00133 volScalarField& psi,
00134 const surfaceScalarField& phi,
00135 surfaceScalarField& phiPsi,
00136 const SpType& Sp,
00137 const SuType& Su,
00138 const scalar psiMax,
00139 const scalar psiMin
00140 )
00141 {
00142 const fvMesh& mesh = psi.mesh();
00143
00144 const dictionary& MULEScontrols = mesh.solverDict(psi.name());
00145
00146 label maxIter
00147 (
00148 readLabel(MULEScontrols.lookup("maxIter"))
00149 );
00150
00151 label nLimiterIter
00152 (
00153 readLabel(MULEScontrols.lookup("nLimiterIter"))
00154 );
00155
00156 scalar maxUnboundedness
00157 (
00158 readScalar(MULEScontrols.lookup("maxUnboundedness"))
00159 );
00160
00161 scalar CoCoeff
00162 (
00163 readScalar(MULEScontrols.lookup("CoCoeff"))
00164 );
00165
00166 scalarField allCoLambda(mesh.nFaces());
00167
00168 {
00169 surfaceScalarField Cof =
00170 mesh.time().deltaT()*mesh.surfaceInterpolation::deltaCoeffs()
00171 *mag(phi)/mesh.magSf();
00172
00173 slicedSurfaceScalarField CoLambda
00174 (
00175 IOobject
00176 (
00177 "CoLambda",
00178 mesh.time().timeName(),
00179 mesh,
00180 IOobject::NO_READ,
00181 IOobject::NO_WRITE,
00182 false
00183 ),
00184 mesh,
00185 dimless,
00186 allCoLambda,
00187 false
00188 );
00189
00190 CoLambda == 1.0/max(CoCoeff*Cof, scalar(1));
00191 }
00192
00193 scalarField allLambda(allCoLambda);
00194
00195
00196 slicedSurfaceScalarField lambda
00197 (
00198 IOobject
00199 (
00200 "lambda",
00201 mesh.time().timeName(),
00202 mesh,
00203 IOobject::NO_READ,
00204 IOobject::NO_WRITE,
00205 false
00206 ),
00207 mesh,
00208 dimless,
00209 allLambda,
00210 false
00211 );
00212
00213 linear<scalar> CDs(mesh);
00214 upwind<scalar> UDs(mesh, phi);
00215
00216
00217 fvScalarMatrix psiConvectionDiffusion
00218 (
00219 fvm::ddt(rho, psi)
00220 + fv::gaussConvectionScheme<scalar>(mesh, phi, UDs).fvmDiv(phi, psi)
00221
00222
00223 - fvm::Sp(Sp, psi)
00224 - Su
00225 );
00226
00227 surfaceScalarField phiBD = psiConvectionDiffusion.flux();
00228
00229 surfaceScalarField& phiCorr = phiPsi;
00230 phiCorr -= phiBD;
00231
00232 for (label i=0; i<maxIter; i++)
00233 {
00234 if (i != 0 && i < 4)
00235 {
00236 allLambda = allCoLambda;
00237 }
00238
00239 limiter
00240 (
00241 allLambda,
00242 rho,
00243 psi,
00244 phiBD,
00245 phiCorr,
00246 Sp.field(),
00247 Su.field(),
00248 psiMax,
00249 psiMin,
00250 nLimiterIter
00251 );
00252
00253 solve
00254 (
00255 psiConvectionDiffusion + fvc::div(lambda*phiCorr),
00256 MULEScontrols
00257 );
00258
00259 scalar maxPsiM1 = gMax(psi.internalField()) - 1.0;
00260 scalar minPsi = gMin(psi.internalField());
00261
00262 scalar unboundedness = max(max(maxPsiM1, 0.0), -min(minPsi, 0.0));
00263
00264 if (unboundedness < maxUnboundedness)
00265 {
00266 break;
00267 }
00268 else
00269 {
00270 Info<< "MULES: max(" << psi.name() << " - 1) = " << maxPsiM1
00271 << " min(" << psi.name() << ") = " << minPsi << endl;
00272
00273 phiBD = psiConvectionDiffusion.flux();
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297 }
00298 }
00299
00300 phiPsi = psiConvectionDiffusion.flux() + lambda*phiCorr;
00301 }
00302
00303
00304 template<class RhoType, class SpType, class SuType>
00305 void Foam::MULES::limiter
00306 (
00307 scalarField& allLambda,
00308 const RhoType& rho,
00309 const volScalarField& psi,
00310 const surfaceScalarField& phiBD,
00311 const surfaceScalarField& phiCorr,
00312 const SpType& Sp,
00313 const SuType& Su,
00314 const scalar psiMax,
00315 const scalar psiMin,
00316 const label nLimiterIter
00317 )
00318 {
00319 const scalarField& psiIf = psi;
00320 const volScalarField::GeometricBoundaryField& psiBf = psi.boundaryField();
00321
00322 const scalarField& psi0 = psi.oldTime();
00323
00324 const fvMesh& mesh = psi.mesh();
00325
00326 const unallocLabelList& owner = mesh.owner();
00327 const unallocLabelList& neighb = mesh.neighbour();
00328 tmp<volScalarField::DimensionedInternalField> tVsc = mesh.Vsc();
00329 const scalarField& V = tVsc();
00330 const scalar deltaT = mesh.time().deltaT().value();
00331
00332 const scalarField& phiBDIf = phiBD;
00333 const surfaceScalarField::GeometricBoundaryField& phiBDBf =
00334 phiBD.boundaryField();
00335
00336 const scalarField& phiCorrIf = phiCorr;
00337 const surfaceScalarField::GeometricBoundaryField& phiCorrBf =
00338 phiCorr.boundaryField();
00339
00340 slicedSurfaceScalarField lambda
00341 (
00342 IOobject
00343 (
00344 "lambda",
00345 mesh.time().timeName(),
00346 mesh,
00347 IOobject::NO_READ,
00348 IOobject::NO_WRITE,
00349 false
00350 ),
00351 mesh,
00352 dimless,
00353 allLambda,
00354 false
00355 );
00356
00357 scalarField& lambdaIf = lambda;
00358 surfaceScalarField::GeometricBoundaryField& lambdaBf =
00359 lambda.boundaryField();
00360
00361 scalarField psiMaxn(psiIf.size(), psiMin);
00362 scalarField psiMinn(psiIf.size(), psiMax);
00363
00364 scalarField sumPhiBD(psiIf.size(), 0.0);
00365
00366 scalarField sumPhip(psiIf.size(), VSMALL);
00367 scalarField mSumPhim(psiIf.size(), VSMALL);
00368
00369 forAll(phiCorrIf, facei)
00370 {
00371 label own = owner[facei];
00372 label nei = neighb[facei];
00373
00374 psiMaxn[own] = max(psiMaxn[own], psiIf[nei]);
00375 psiMinn[own] = min(psiMinn[own], psiIf[nei]);
00376
00377 psiMaxn[nei] = max(psiMaxn[nei], psiIf[own]);
00378 psiMinn[nei] = min(psiMinn[nei], psiIf[own]);
00379
00380 sumPhiBD[own] += phiBDIf[facei];
00381 sumPhiBD[nei] -= phiBDIf[facei];
00382
00383 scalar phiCorrf = phiCorrIf[facei];
00384
00385 if (phiCorrf > 0.0)
00386 {
00387 sumPhip[own] += phiCorrf;
00388 mSumPhim[nei] += phiCorrf;
00389 }
00390 else
00391 {
00392 mSumPhim[own] -= phiCorrf;
00393 sumPhip[nei] -= phiCorrf;
00394 }
00395 }
00396
00397 forAll(phiCorrBf, patchi)
00398 {
00399 const fvPatchScalarField& psiPf = psiBf[patchi];
00400 const scalarField& phiBDPf = phiBDBf[patchi];
00401 const scalarField& phiCorrPf = phiCorrBf[patchi];
00402
00403 const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
00404
00405 if (psiPf.coupled())
00406 {
00407 scalarField psiPNf = psiPf.patchNeighbourField();
00408
00409 forAll(phiCorrPf, pFacei)
00410 {
00411 label pfCelli = pFaceCells[pFacei];
00412
00413 psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPNf[pFacei]);
00414 psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPNf[pFacei]);
00415 }
00416 }
00417 else
00418 {
00419 forAll(phiCorrPf, pFacei)
00420 {
00421 label pfCelli = pFaceCells[pFacei];
00422
00423 psiMaxn[pfCelli] = max(psiMaxn[pfCelli], psiPf[pFacei]);
00424 psiMinn[pfCelli] = min(psiMinn[pfCelli], psiPf[pFacei]);
00425 }
00426 }
00427
00428 forAll(phiCorrPf, pFacei)
00429 {
00430 label pfCelli = pFaceCells[pFacei];
00431
00432 sumPhiBD[pfCelli] += phiBDPf[pFacei];
00433
00434 scalar phiCorrf = phiCorrPf[pFacei];
00435
00436 if (phiCorrf > 0.0)
00437 {
00438 sumPhip[pfCelli] += phiCorrf;
00439 }
00440 else
00441 {
00442 mSumPhim[pfCelli] -= phiCorrf;
00443 }
00444 }
00445 }
00446
00447 psiMaxn = min(psiMaxn, psiMax);
00448 psiMinn = max(psiMinn, psiMin);
00449
00450
00451
00452
00453
00454 if (mesh.moving())
00455 {
00456 tmp<volScalarField::DimensionedInternalField> V0 = mesh.Vsc0();
00457
00458 psiMaxn =
00459 V*((rho/deltaT - Sp)*psiMaxn - Su)
00460 - (V0()/deltaT)*rho.oldTime()*psi0
00461 + sumPhiBD;
00462
00463 psiMinn =
00464 V*(Su - (rho/deltaT - Sp)*psiMinn)
00465 + (V0/deltaT)*rho.oldTime()*psi0
00466 - sumPhiBD;
00467 }
00468 else
00469 {
00470 psiMaxn =
00471 V*((rho/deltaT - Sp)*psiMaxn - Su - (rho.oldTime()/deltaT)*psi0)
00472 + sumPhiBD;
00473
00474 psiMinn =
00475 V*(Su - (rho/deltaT - Sp)*psiMinn + (rho.oldTime()/deltaT)*psi0)
00476 - sumPhiBD;
00477 }
00478
00479 scalarField sumlPhip(psiIf.size());
00480 scalarField mSumlPhim(psiIf.size());
00481
00482 for(int j=0; j<nLimiterIter; j++)
00483 {
00484 sumlPhip = 0.0;
00485 mSumlPhim = 0.0;
00486
00487 forAll(lambdaIf, facei)
00488 {
00489 label own = owner[facei];
00490 label nei = neighb[facei];
00491
00492 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
00493
00494 if (lambdaPhiCorrf > 0.0)
00495 {
00496 sumlPhip[own] += lambdaPhiCorrf;
00497 mSumlPhim[nei] += lambdaPhiCorrf;
00498 }
00499 else
00500 {
00501 mSumlPhim[own] -= lambdaPhiCorrf;
00502 sumlPhip[nei] -= lambdaPhiCorrf;
00503 }
00504 }
00505
00506 forAll(lambdaBf, patchi)
00507 {
00508 scalarField& lambdaPf = lambdaBf[patchi];
00509 const scalarField& phiCorrfPf = phiCorrBf[patchi];
00510
00511 const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
00512
00513 forAll(lambdaPf, pFacei)
00514 {
00515 label pfCelli = pFaceCells[pFacei];
00516
00517 scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
00518
00519 if (lambdaPhiCorrf > 0.0)
00520 {
00521 sumlPhip[pfCelli] += lambdaPhiCorrf;
00522 }
00523 else
00524 {
00525 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
00526 }
00527 }
00528 }
00529
00530 forAll (sumlPhip, celli)
00531 {
00532 sumlPhip[celli] =
00533 max(min
00534 (
00535 (sumlPhip[celli] + psiMaxn[celli])/mSumPhim[celli],
00536 1.0), 0.0
00537 );
00538
00539 mSumlPhim[celli] =
00540 max(min
00541 (
00542 (mSumlPhim[celli] + psiMinn[celli])/sumPhip[celli],
00543 1.0), 0.0
00544 );
00545 }
00546
00547 const scalarField& lambdam = sumlPhip;
00548 const scalarField& lambdap = mSumlPhim;
00549
00550 forAll(lambdaIf, facei)
00551 {
00552 if (phiCorrIf[facei] > 0.0)
00553 {
00554 lambdaIf[facei] = min
00555 (
00556 lambdaIf[facei],
00557 min(lambdap[owner[facei]], lambdam[neighb[facei]])
00558 );
00559 }
00560 else
00561 {
00562 lambdaIf[facei] = min
00563 (
00564 lambdaIf[facei],
00565 min(lambdam[owner[facei]], lambdap[neighb[facei]])
00566 );
00567 }
00568 }
00569
00570
00571 forAll(lambdaBf, patchi)
00572 {
00573 fvsPatchScalarField& lambdaPf = lambdaBf[patchi];
00574 const scalarField& phiCorrfPf = phiCorrBf[patchi];
00575
00576 const labelList& pFaceCells = mesh.boundary()[patchi].faceCells();
00577
00578 forAll(lambdaPf, pFacei)
00579 {
00580 label pfCelli = pFaceCells[pFacei];
00581
00582 if (phiCorrfPf[pFacei] > 0.0)
00583 {
00584 lambdaPf[pFacei] = min(lambdaPf[pFacei], lambdap[pfCelli]);
00585 }
00586 else
00587 {
00588 lambdaPf[pFacei] = min(lambdaPf[pFacei], lambdam[pfCelli]);
00589 }
00590 }
00591 }
00592
00593 syncTools::syncFaceList(mesh, allLambda, minEqOp<scalar>(), false);
00594 }
00595 }
00596
00597
00598