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