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 "SLTSDdtScheme.H"
00027 #include <finiteVolume/surfaceInterpolate.H>
00028 #include <finiteVolume/fvMatrices.H>
00029
00030
00031
00032 namespace Foam
00033 {
00034
00035
00036
00037 namespace fv
00038 {
00039
00040
00041
00042 template<class Type>
00043 void SLTSDdtScheme<Type>::relaxedDiag
00044 (
00045 scalarField& rD,
00046 const surfaceScalarField& phi
00047 ) const
00048 {
00049 const unallocLabelList& owner = mesh().owner();
00050 const unallocLabelList& neighbour = mesh().neighbour();
00051 scalarField diag(rD.size(), 0.0);
00052
00053 forAll(owner, faceI)
00054 {
00055 if (phi[faceI] > 0.0)
00056 {
00057 diag[owner[faceI]] += phi[faceI];
00058 rD[neighbour[faceI]] += phi[faceI];
00059 }
00060 else
00061 {
00062 diag[neighbour[faceI]] -= phi[faceI];
00063 rD[owner[faceI]] -= phi[faceI];
00064 }
00065 }
00066
00067 forAll(phi.boundaryField(), patchi)
00068 {
00069 const fvsPatchScalarField& pphi = phi.boundaryField()[patchi];
00070 const unallocLabelList& faceCells = pphi.patch().patch().faceCells();
00071
00072 forAll(pphi, patchFacei)
00073 {
00074 if (pphi[patchFacei] > 0.0)
00075 {
00076 diag[faceCells[patchFacei]] += pphi[patchFacei];
00077 }
00078 else
00079 {
00080 rD[faceCells[patchFacei]] -= pphi[patchFacei];
00081 }
00082 }
00083 }
00084
00085 rD += (1.0/alpha_ - 2.0)*diag;
00086 }
00087
00088
00089 template<class Type>
00090 tmp<volScalarField> SLTSDdtScheme<Type>::SLrDeltaT() const
00091 {
00092 const surfaceScalarField& phi =
00093 mesh().objectRegistry::lookupObject<surfaceScalarField>(phiName_);
00094
00095 const dimensionedScalar& deltaT = mesh().time().deltaT();
00096
00097 tmp<volScalarField> trDeltaT
00098 (
00099 new volScalarField
00100 (
00101 IOobject
00102 (
00103 "rDeltaT",
00104 phi.instance(),
00105 mesh()
00106 ),
00107 mesh(),
00108 dimensionedScalar("rDeltaT", dimless/dimTime, 0.0),
00109 zeroGradientFvPatchScalarField::typeName
00110 )
00111 );
00112
00113 volScalarField& rDeltaT = trDeltaT();
00114
00115 relaxedDiag(rDeltaT, phi);
00116
00117 if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
00118 {
00119 rDeltaT.internalField() = max
00120 (
00121 rDeltaT.internalField()/mesh().V(),
00122 scalar(1)/deltaT.value()
00123 );
00124 }
00125 else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
00126 {
00127 const volScalarField& rho =
00128 mesh().objectRegistry::lookupObject<volScalarField>(rhoName_)
00129 .oldTime();
00130
00131 rDeltaT.internalField() = max
00132 (
00133 rDeltaT.internalField()/(rho.internalField()*mesh().V()),
00134 scalar(1)/deltaT.value()
00135 );
00136 }
00137 else
00138 {
00139 FatalErrorIn("SLTSDdtScheme<Type>::CofrDeltaT() const")
00140 << "Incorrect dimensions of phi: " << phi.dimensions()
00141 << abort(FatalError);
00142 }
00143
00144 rDeltaT.correctBoundaryConditions();
00145
00146 return trDeltaT;
00147 }
00148
00149
00150 template<class Type>
00151 tmp<GeometricField<Type, fvPatchField, volMesh> >
00152 SLTSDdtScheme<Type>::fvcDdt
00153 (
00154 const dimensioned<Type>& dt
00155 )
00156 {
00157 volScalarField rDeltaT = SLrDeltaT();
00158
00159 IOobject ddtIOobject
00160 (
00161 "ddt("+dt.name()+')',
00162 mesh().time().timeName(),
00163 mesh()
00164 );
00165
00166 if (mesh().moving())
00167 {
00168 tmp<GeometricField<Type, fvPatchField, volMesh> > tdtdt
00169 (
00170 new GeometricField<Type, fvPatchField, volMesh>
00171 (
00172 ddtIOobject,
00173 mesh(),
00174 dimensioned<Type>
00175 (
00176 "0",
00177 dt.dimensions()/dimTime,
00178 pTraits<Type>::zero
00179 )
00180 )
00181 );
00182
00183 tdtdt().internalField() =
00184 rDeltaT.internalField()*dt.value()*(1.0 - mesh().V0()/mesh().V());
00185
00186 return tdtdt;
00187 }
00188 else
00189 {
00190 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00191 (
00192 new GeometricField<Type, fvPatchField, volMesh>
00193 (
00194 ddtIOobject,
00195 mesh(),
00196 dimensioned<Type>
00197 (
00198 "0",
00199 dt.dimensions()/dimTime,
00200 pTraits<Type>::zero
00201 ),
00202 calculatedFvPatchField<Type>::typeName
00203 )
00204 );
00205 }
00206 }
00207
00208
00209 template<class Type>
00210 tmp<GeometricField<Type, fvPatchField, volMesh> >
00211 SLTSDdtScheme<Type>::fvcDdt
00212 (
00213 const GeometricField<Type, fvPatchField, volMesh>& vf
00214 )
00215 {
00216 volScalarField rDeltaT = SLrDeltaT();
00217
00218 IOobject ddtIOobject
00219 (
00220 "ddt("+vf.name()+')',
00221 mesh().time().timeName(),
00222 mesh()
00223 );
00224
00225 if (mesh().moving())
00226 {
00227 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00228 (
00229 new GeometricField<Type, fvPatchField, volMesh>
00230 (
00231 ddtIOobject,
00232 mesh(),
00233 rDeltaT.dimensions()*vf.dimensions(),
00234 rDeltaT.internalField()*
00235 (
00236 vf.internalField()
00237 - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00238 ),
00239 rDeltaT.boundaryField()*
00240 (
00241 vf.boundaryField() - vf.oldTime().boundaryField()
00242 )
00243 )
00244 );
00245 }
00246 else
00247 {
00248 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00249 (
00250 new GeometricField<Type, fvPatchField, volMesh>
00251 (
00252 ddtIOobject,
00253 rDeltaT*(vf - vf.oldTime())
00254 )
00255 );
00256 }
00257 }
00258
00259
00260 template<class Type>
00261 tmp<GeometricField<Type, fvPatchField, volMesh> >
00262 SLTSDdtScheme<Type>::fvcDdt
00263 (
00264 const dimensionedScalar& rho,
00265 const GeometricField<Type, fvPatchField, volMesh>& vf
00266 )
00267 {
00268 volScalarField rDeltaT = SLrDeltaT();
00269
00270 IOobject ddtIOobject
00271 (
00272 "ddt("+rho.name()+','+vf.name()+')',
00273 mesh().time().timeName(),
00274 mesh()
00275 );
00276
00277 if (mesh().moving())
00278 {
00279 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00280 (
00281 new GeometricField<Type, fvPatchField, volMesh>
00282 (
00283 ddtIOobject,
00284 mesh(),
00285 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00286 rDeltaT.internalField()*rho.value()*
00287 (
00288 vf.internalField()
00289 - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00290 ),
00291 rDeltaT.boundaryField()*rho.value()*
00292 (
00293 vf.boundaryField() - vf.oldTime().boundaryField()
00294 )
00295 )
00296 );
00297 }
00298 else
00299 {
00300 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00301 (
00302 new GeometricField<Type, fvPatchField, volMesh>
00303 (
00304 ddtIOobject,
00305 rDeltaT*rho*(vf - vf.oldTime())
00306 )
00307 );
00308 }
00309 }
00310
00311
00312 template<class Type>
00313 tmp<GeometricField<Type, fvPatchField, volMesh> >
00314 SLTSDdtScheme<Type>::fvcDdt
00315 (
00316 const volScalarField& rho,
00317 const GeometricField<Type, fvPatchField, volMesh>& vf
00318 )
00319 {
00320 volScalarField rDeltaT = SLrDeltaT();
00321
00322 IOobject ddtIOobject
00323 (
00324 "ddt("+rho.name()+','+vf.name()+')',
00325 mesh().time().timeName(),
00326 mesh()
00327 );
00328
00329 if (mesh().moving())
00330 {
00331 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00332 (
00333 new GeometricField<Type, fvPatchField, volMesh>
00334 (
00335 ddtIOobject,
00336 mesh(),
00337 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00338 rDeltaT.internalField()*
00339 (
00340 rho.internalField()*vf.internalField()
00341 - rho.oldTime().internalField()
00342 *vf.oldTime().internalField()*mesh().V0()/mesh().V()
00343 ),
00344 rDeltaT.boundaryField()*
00345 (
00346 rho.boundaryField()*vf.boundaryField()
00347 - rho.oldTime().boundaryField()
00348 *vf.oldTime().boundaryField()
00349 )
00350 )
00351 );
00352 }
00353 else
00354 {
00355 return tmp<GeometricField<Type, fvPatchField, volMesh> >
00356 (
00357 new GeometricField<Type, fvPatchField, volMesh>
00358 (
00359 ddtIOobject,
00360 rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
00361 )
00362 );
00363 }
00364 }
00365
00366
00367 template<class Type>
00368 tmp<fvMatrix<Type> >
00369 SLTSDdtScheme<Type>::fvmDdt
00370 (
00371 GeometricField<Type, fvPatchField, volMesh>& vf
00372 )
00373 {
00374 tmp<fvMatrix<Type> > tfvm
00375 (
00376 new fvMatrix<Type>
00377 (
00378 vf,
00379 vf.dimensions()*dimVol/dimTime
00380 )
00381 );
00382
00383 fvMatrix<Type>& fvm = tfvm();
00384
00385 scalarField rDeltaT = SLrDeltaT()().internalField();
00386
00387 Info<< "max/min rDeltaT " << max(rDeltaT) << " " << min(rDeltaT) << endl;
00388
00389 fvm.diag() = rDeltaT*mesh().V();
00390
00391 if (mesh().moving())
00392 {
00393 fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V0();
00394 }
00395 else
00396 {
00397 fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V();
00398 }
00399
00400 return tfvm;
00401 }
00402
00403
00404 template<class Type>
00405 tmp<fvMatrix<Type> >
00406 SLTSDdtScheme<Type>::fvmDdt
00407 (
00408 const dimensionedScalar& rho,
00409 GeometricField<Type, fvPatchField, volMesh>& vf
00410 )
00411 {
00412 tmp<fvMatrix<Type> > tfvm
00413 (
00414 new fvMatrix<Type>
00415 (
00416 vf,
00417 rho.dimensions()*vf.dimensions()*dimVol/dimTime
00418 )
00419 );
00420 fvMatrix<Type>& fvm = tfvm();
00421
00422 scalarField rDeltaT = SLrDeltaT()().internalField();
00423
00424 fvm.diag() = rDeltaT*rho.value()*mesh().V();
00425
00426 if (mesh().moving())
00427 {
00428 fvm.source() = rDeltaT
00429 *rho.value()*vf.oldTime().internalField()*mesh().V0();
00430 }
00431 else
00432 {
00433 fvm.source() = rDeltaT
00434 *rho.value()*vf.oldTime().internalField()*mesh().V();
00435 }
00436
00437 return tfvm;
00438 }
00439
00440
00441 template<class Type>
00442 tmp<fvMatrix<Type> >
00443 SLTSDdtScheme<Type>::fvmDdt
00444 (
00445 const volScalarField& rho,
00446 GeometricField<Type, fvPatchField, volMesh>& vf
00447 )
00448 {
00449 tmp<fvMatrix<Type> > tfvm
00450 (
00451 new fvMatrix<Type>
00452 (
00453 vf,
00454 rho.dimensions()*vf.dimensions()*dimVol/dimTime
00455 )
00456 );
00457 fvMatrix<Type>& fvm = tfvm();
00458
00459 scalarField rDeltaT = SLrDeltaT()().internalField();
00460
00461 fvm.diag() = rDeltaT*rho.internalField()*mesh().V();
00462
00463 if (mesh().moving())
00464 {
00465 fvm.source() = rDeltaT
00466 *rho.oldTime().internalField()
00467 *vf.oldTime().internalField()*mesh().V0();
00468 }
00469 else
00470 {
00471 fvm.source() = rDeltaT
00472 *rho.oldTime().internalField()
00473 *vf.oldTime().internalField()*mesh().V();
00474 }
00475
00476 return tfvm;
00477 }
00478
00479
00480 template<class Type>
00481 tmp<typename SLTSDdtScheme<Type>::fluxFieldType>
00482 SLTSDdtScheme<Type>::fvcDdtPhiCorr
00483 (
00484 const volScalarField& rA,
00485 const GeometricField<Type, fvPatchField, volMesh>& U,
00486 const fluxFieldType& phi
00487 )
00488 {
00489 IOobject ddtIOobject
00490 (
00491 "ddtPhiCorr(" + rA.name() + ',' + U.name() + ',' + phi.name() + ')',
00492 mesh().time().timeName(),
00493 mesh()
00494 );
00495
00496 if (mesh().moving())
00497 {
00498 return tmp<fluxFieldType>
00499 (
00500 new fluxFieldType
00501 (
00502 ddtIOobject,
00503 mesh(),
00504 dimensioned<typename flux<Type>::type>
00505 (
00506 "0",
00507 rA.dimensions()*phi.dimensions()/dimTime,
00508 pTraits<typename flux<Type>::type>::zero
00509 )
00510 )
00511 );
00512 }
00513 else
00514 {
00515 volScalarField rDeltaT = SLrDeltaT();
00516
00517 return tmp<fluxFieldType>
00518 (
00519 new fluxFieldType
00520 (
00521 ddtIOobject,
00522 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())*
00523 (
00524 fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00525 - (fvc::interpolate(rDeltaT*rA*U.oldTime()) & mesh().Sf())
00526 )
00527 )
00528 );
00529 }
00530 }
00531
00532
00533 template<class Type>
00534 tmp<typename SLTSDdtScheme<Type>::fluxFieldType>
00535 SLTSDdtScheme<Type>::fvcDdtPhiCorr
00536 (
00537 const volScalarField& rA,
00538 const volScalarField& rho,
00539 const GeometricField<Type, fvPatchField, volMesh>& U,
00540 const fluxFieldType& phi
00541 )
00542 {
00543 IOobject ddtIOobject
00544 (
00545 "ddtPhiCorr("
00546 + rA.name() + ',' + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
00547 mesh().time().timeName(),
00548 mesh()
00549 );
00550
00551 if (mesh().moving())
00552 {
00553 return tmp<fluxFieldType>
00554 (
00555 new fluxFieldType
00556 (
00557 ddtIOobject,
00558 mesh(),
00559 dimensioned<typename flux<Type>::type>
00560 (
00561 "0",
00562 rA.dimensions()*rho.dimensions()*phi.dimensions()/dimTime,
00563 pTraits<typename flux<Type>::type>::zero
00564 )
00565 )
00566 );
00567 }
00568 else
00569 {
00570 volScalarField rDeltaT = SLrDeltaT();
00571
00572 if
00573 (
00574 U.dimensions() == dimVelocity
00575 && phi.dimensions() == dimVelocity*dimArea
00576 )
00577 {
00578 return tmp<fluxFieldType>
00579 (
00580 new fluxFieldType
00581 (
00582 ddtIOobject,
00583 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
00584 *(
00585 fvc::interpolate(rDeltaT*rA*rho.oldTime())*phi.oldTime()
00586 - (fvc::interpolate(rDeltaT*rA*rho.oldTime()*U.oldTime())
00587 & mesh().Sf())
00588 )
00589 )
00590 );
00591 }
00592 else if
00593 (
00594 U.dimensions() == dimVelocity
00595 && phi.dimensions() == dimDensity*dimVelocity*dimArea
00596 )
00597 {
00598 return tmp<fluxFieldType>
00599 (
00600 new fluxFieldType
00601 (
00602 ddtIOobject,
00603 this->fvcDdtPhiCoeff
00604 (
00605 U.oldTime(),
00606 phi.oldTime()/fvc::interpolate(rho.oldTime())
00607 )
00608 *(
00609 fvc::interpolate(rDeltaT*rA*rho.oldTime())
00610 *phi.oldTime()/fvc::interpolate(rho.oldTime())
00611 - (
00612 fvc::interpolate
00613 (
00614 rDeltaT*rA*rho.oldTime()*U.oldTime()
00615 ) & mesh().Sf()
00616 )
00617 )
00618 )
00619 );
00620 }
00621 else if
00622 (
00623 U.dimensions() == dimDensity*dimVelocity
00624 && phi.dimensions() == dimDensity*dimVelocity*dimArea
00625 )
00626 {
00627 return tmp<fluxFieldType>
00628 (
00629 new fluxFieldType
00630 (
00631 ddtIOobject,
00632 this->fvcDdtPhiCoeff
00633 (
00634 rho.oldTime(),
00635 U.oldTime(),
00636 phi.oldTime()
00637 )
00638 *(
00639 fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00640 - (
00641 fvc::interpolate(rDeltaT*rA*U.oldTime())&mesh().Sf()
00642 )
00643 )
00644 )
00645 );
00646 }
00647 else
00648 {
00649 FatalErrorIn
00650 (
00651 "SLTSDdtScheme<Type>::fvcDdtPhiCorr"
00652 ) << "dimensions of phi are not correct"
00653 << abort(FatalError);
00654
00655 return fluxFieldType::null();
00656 }
00657 }
00658 }
00659
00660
00661 template<class Type>
00662 tmp<surfaceScalarField> SLTSDdtScheme<Type>::meshPhi
00663 (
00664 const GeometricField<Type, fvPatchField, volMesh>&
00665 )
00666 {
00667 return tmp<surfaceScalarField>
00668 (
00669 new surfaceScalarField
00670 (
00671 IOobject
00672 (
00673 "meshPhi",
00674 mesh().time().timeName(),
00675 mesh()
00676 ),
00677 mesh(),
00678 dimensionedScalar("0", dimVolume/dimTime, 0.0)
00679 )
00680 );
00681 }
00682
00683
00684
00685
00686 }
00687
00688
00689
00690 }
00691
00692