FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

boundedBackwardDdtScheme.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023     
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "boundedBackwardDdtScheme.H"
00027 #include <finiteVolume/surfaceInterpolate.H>
00028 #include <finiteVolume/fvcDiv.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034 
00035 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00036 
00037 namespace fv
00038 {
00039 
00040 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00041 
00042 scalar boundedBackwardDdtScheme::deltaT_() const
00043 {
00044     return mesh().time().deltaT().value();
00045 }
00046 
00047 
00048 scalar boundedBackwardDdtScheme::deltaT0_() const
00049 {
00050     return mesh().time().deltaT0().value();
00051 }
00052 
00053 
00054 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00055 
00056 tmp<volScalarField>
00057 boundedBackwardDdtScheme::fvcDdt
00058 (
00059     const dimensionedScalar& dt
00060 )
00061 {
00062     // No change compared to backward
00063 
00064     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00065 
00066     IOobject ddtIOobject
00067     (
00068         "ddt("+dt.name()+')',
00069         mesh().time().timeName(),
00070         mesh()
00071     );
00072 
00073     scalar deltaT = deltaT_();
00074     scalar deltaT0 = deltaT0_();
00075 
00076     scalar coefft   = 1 + deltaT/(deltaT + deltaT0);
00077     scalar coefft00 = deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
00078     scalar coefft0  = coefft + coefft00;
00079 
00080     if (mesh().moving())
00081     {
00082         tmp<volScalarField> tdtdt
00083         (
00084             new volScalarField
00085             (
00086                 ddtIOobject,
00087                 mesh(),
00088                 dimensionedScalar
00089                 (
00090                     "0",
00091                     dt.dimensions()/dimTime,
00092                     0.0
00093                 )
00094             )
00095         );
00096 
00097         tdtdt().internalField() = rDeltaT.value()*dt.value()*
00098         (
00099             coefft - (coefft0*mesh().V0() - coefft00*mesh().V00())/mesh().V()
00100         );
00101 
00102         return tdtdt;
00103     }
00104     else
00105     {
00106         return tmp<volScalarField>
00107         (
00108             new volScalarField
00109             (
00110                 ddtIOobject,
00111                 mesh(),
00112                 dimensionedScalar
00113                 (
00114                     "0",
00115                     dt.dimensions()/dimTime,
00116                     0.0
00117                 ),
00118                 calculatedFvPatchScalarField::typeName
00119             )
00120         );
00121     }
00122 }
00123 
00124 
00125 tmp<volScalarField>
00126 boundedBackwardDdtScheme::fvcDdt
00127 (
00128     const volScalarField& vf
00129 )
00130 {
00131     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00132 
00133     IOobject ddtIOobject
00134     (
00135         "ddt("+vf.name()+')',
00136         mesh().time().timeName(),
00137         mesh()
00138     );
00139 
00140     scalar deltaT = deltaT_();
00141     scalar deltaT0 = deltaT0_(vf);
00142 
00143     // Calculate unboundedness indicator
00144     // Note: all times moved by one because access to internal field
00145     // copies current field into the old-time level.
00146     volScalarField phict =
00147         mag
00148         (
00149             vf.oldTime().oldTime()
00150           - vf.oldTime().oldTime().oldTime()
00151         )/
00152         (
00153             mag
00154             (
00155                 vf.oldTime()
00156               - vf.oldTime().oldTime()
00157             )
00158           + dimensionedScalar("small", vf.dimensions(), VSMALL)
00159         );
00160 
00161     volScalarField limiter(pos(phict) - pos(phict - scalar(1)));
00162 
00163     volScalarField coefft   = scalar(1) + limiter*deltaT/(deltaT + deltaT0);
00164     volScalarField coefft00 = limiter*sqr(deltaT)/(deltaT0*(deltaT + deltaT0));
00165     volScalarField coefft0  = coefft + coefft00;
00166 
00167     if (mesh().moving())
00168     {
00169         return tmp<volScalarField>
00170         (
00171             new volScalarField
00172             (
00173                 ddtIOobject,
00174                 mesh(),
00175                 rDeltaT.dimensions()*vf.dimensions(),
00176                 rDeltaT.value()*
00177                 (
00178                     coefft*vf.internalField() -
00179                     (
00180                         coefft0.internalField()
00181                         *vf.oldTime().internalField()*mesh().V0()
00182                       - coefft00.internalField()
00183                         *vf.oldTime().oldTime().internalField()
00184                        *mesh().V00()
00185                     )/mesh().V()
00186                 ),
00187                 rDeltaT.value()*
00188                 (
00189                     coefft.boundaryField()*vf.boundaryField() -
00190                     (
00191                         coefft0.boundaryField()*
00192                         vf.oldTime().boundaryField()
00193                       - coefft00.boundaryField()*
00194                         vf.oldTime().oldTime().boundaryField()
00195                     )
00196                 )
00197             )
00198         );
00199     }
00200     else
00201     {
00202         return tmp<volScalarField>
00203         (
00204             new volScalarField
00205             (
00206                 ddtIOobject,
00207                 rDeltaT*
00208                 (
00209                     coefft*vf
00210                   - coefft0*vf.oldTime()
00211                   + coefft00*vf.oldTime().oldTime()
00212                 )
00213             )
00214         );
00215     }
00216 }
00217 
00218 
00219 tmp<volScalarField>
00220 boundedBackwardDdtScheme::fvcDdt
00221 (
00222     const dimensionedScalar& rho,
00223     const volScalarField& vf
00224 )
00225 {
00226     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00227 
00228     IOobject ddtIOobject
00229     (
00230         "ddt("+rho.name()+','+vf.name()+')',
00231         mesh().time().timeName(),
00232         mesh()
00233     );
00234 
00235     scalar deltaT = deltaT_();
00236     scalar deltaT0 = deltaT0_(vf);
00237 
00238     // Calculate unboundedness indicator
00239     // Note: all times moved by one because access to internal field
00240     // copies current field into the old-time level.
00241     volScalarField phict =
00242         mag
00243         (
00244             vf.oldTime().oldTime()
00245           - vf.oldTime().oldTime().oldTime()
00246         )/
00247         (
00248             mag
00249             (
00250                 vf.oldTime()
00251               - vf.oldTime().oldTime()
00252             )
00253           + dimensionedScalar("small", vf.dimensions(), VSMALL)
00254         );
00255 
00256     volScalarField limiter(pos(phict) - pos(phict - scalar(1)));
00257 
00258     volScalarField coefft   = scalar(1) + limiter*deltaT/(deltaT + deltaT0);
00259     volScalarField coefft00 = limiter*sqr(deltaT)/(deltaT0*(deltaT + deltaT0));
00260     volScalarField coefft0  = coefft + coefft00;
00261 
00262     if (mesh().moving())
00263     {
00264         return tmp<volScalarField>
00265         (
00266             new volScalarField
00267             (
00268                 ddtIOobject,
00269                 mesh(),
00270                 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00271                 rDeltaT.value()*rho.value()*
00272                 (
00273                     coefft*vf.internalField() -
00274                     (
00275                         coefft0.internalField()*
00276                         vf.oldTime().internalField()*mesh().V0()
00277                       - coefft00.internalField()*
00278                         vf.oldTime().oldTime().internalField()
00279                        *mesh().V00()
00280                     )/mesh().V()
00281                 ),
00282                 rDeltaT.value()*rho.value()*
00283                 (
00284                     coefft.boundaryField()*vf.boundaryField() -
00285                     (
00286                         coefft0.boundaryField()*
00287                         vf.oldTime().boundaryField()
00288                       - coefft00.boundaryField()*
00289                         vf.oldTime().oldTime().boundaryField()
00290                     )
00291                 )
00292             )
00293         );
00294     }
00295     else
00296     {
00297         return tmp<volScalarField>
00298         (
00299             new volScalarField
00300             (
00301                 ddtIOobject,
00302                 rDeltaT*rho*
00303                 (
00304                     coefft*vf
00305                   - coefft0*vf.oldTime()
00306                  + coefft00*vf.oldTime().oldTime()
00307                 )
00308             )
00309         );
00310     }
00311 }
00312 
00313 
00314 tmp<volScalarField>
00315 boundedBackwardDdtScheme::fvcDdt
00316 (
00317     const volScalarField& rho,
00318     const volScalarField& vf
00319 )
00320 {
00321     dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();
00322 
00323     IOobject ddtIOobject
00324     (
00325         "ddt("+rho.name()+','+vf.name()+')',
00326         mesh().time().timeName(),
00327         mesh()
00328     );
00329 
00330     scalar deltaT = deltaT_();
00331     scalar deltaT0 = deltaT0_(vf);
00332 
00333     // Calculate unboundedness indicator
00334     // Note: all times moved by one because access to internal field
00335     // copies current field into the old-time level.
00336     volScalarField phict =
00337         mag
00338         (
00339             rho.oldTime().oldTime()*vf.oldTime().oldTime()
00340           - rho.oldTime().oldTime().oldTime()*vf.oldTime().oldTime().oldTime()
00341         )/
00342         (
00343             mag
00344             (
00345                 rho.oldTime()*vf.oldTime()
00346               - rho.oldTime().oldTime()*vf.oldTime().oldTime()
00347             )
00348           + dimensionedScalar("small", rho.dimensions()*vf.dimensions(), VSMALL)
00349         );
00350 
00351     volScalarField limiter(pos(phict) - pos(phict - scalar(1)));
00352 
00353     volScalarField coefft   = scalar(1) + limiter*deltaT/(deltaT + deltaT0);
00354     volScalarField coefft00 = limiter*sqr(deltaT)/(deltaT0*(deltaT + deltaT0));
00355     volScalarField coefft0  = coefft + coefft00;
00356 
00357     if (mesh().moving())
00358     {
00359         return tmp<volScalarField>
00360         (
00361             new volScalarField
00362             (
00363                 ddtIOobject,
00364                 mesh(),
00365                 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00366                 rDeltaT.value()*
00367                 (
00368                     coefft*rho.internalField()*vf.internalField() -
00369                     (
00370                         coefft0.internalField()*
00371                         rho.oldTime().internalField()*
00372                         vf.oldTime().internalField()*mesh().V0()
00373                       - coefft00.internalField()*
00374                         rho.oldTime().oldTime().internalField()
00375                        *vf.oldTime().oldTime().internalField()*mesh().V00()
00376                     )/mesh().V()
00377                 ),
00378                 rDeltaT.value()*
00379                 (
00380                     coefft.boundaryField()*vf.boundaryField() -
00381                     (
00382                         coefft0.boundaryField()*
00383                         rho.oldTime().boundaryField()*
00384                         vf.oldTime().boundaryField()
00385                       - coefft00.boundaryField()*
00386                         rho.oldTime().oldTime().boundaryField()*
00387                         vf.oldTime().oldTime().boundaryField()
00388                     )
00389                 )
00390             )
00391         );
00392     }
00393     else
00394     {
00395         return tmp<volScalarField>
00396         (
00397             new volScalarField
00398             (
00399                 ddtIOobject,
00400                 rDeltaT*
00401                 (
00402                     coefft*rho*vf
00403                   - coefft0*rho.oldTime()*vf.oldTime()
00404                   + coefft00*rho.oldTime().oldTime()*vf.oldTime().oldTime()
00405                 )
00406             )
00407         );
00408     }
00409 }
00410 
00411 
00412 tmp<fvScalarMatrix>
00413 boundedBackwardDdtScheme::fvmDdt
00414 (
00415     volScalarField& vf
00416 )
00417 {
00418     tmp<fvScalarMatrix> tfvm
00419     (
00420         new fvScalarMatrix
00421         (
00422             vf,
00423             vf.dimensions()*dimVol/dimTime
00424         )
00425     );
00426 
00427     fvScalarMatrix& fvm = tfvm();
00428 
00429     scalar rDeltaT = 1.0/deltaT_();
00430 
00431     scalar deltaT = deltaT_();
00432     scalar deltaT0 = deltaT0_(vf);
00433 
00434     // Calculate unboundedness indicator
00435     // Note: all times moved by one because access to internal field
00436     // copies current field into the old-time level.
00437     scalarField phict =
00438         mag
00439         (
00440             vf.oldTime().oldTime().internalField()
00441           - vf.oldTime().oldTime().oldTime().internalField()
00442         )/
00443         (
00444             mag
00445             (
00446                 vf.oldTime().internalField()
00447               - vf.oldTime().oldTime().internalField()
00448             )
00449             + VSMALL
00450         );
00451 
00452     scalarField limiter(pos(phict) - pos(phict - 1.0));
00453 
00454     scalarField coefft   = 1.0 + limiter*deltaT/(deltaT + deltaT0);
00455     scalarField coefft00 = limiter*deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
00456     scalarField coefft0  = coefft + coefft00;
00457 
00458     fvm.diag() = (coefft*rDeltaT)*mesh().V();
00459 
00460     if (mesh().moving())
00461     {
00462         fvm.source() = rDeltaT*
00463         (
00464             coefft0*vf.oldTime().internalField()*mesh().V0()
00465           - coefft00*vf.oldTime().oldTime().internalField()
00466            *mesh().V00()
00467         );
00468     }
00469     else
00470     {
00471         fvm.source() = rDeltaT*mesh().V()*
00472         (
00473             coefft0*vf.oldTime().internalField()
00474           - coefft00*vf.oldTime().oldTime().internalField()
00475         );
00476     }
00477 
00478     return tfvm;
00479 }
00480 
00481 
00482 tmp<fvScalarMatrix>
00483 boundedBackwardDdtScheme::fvmDdt
00484 (
00485     const dimensionedScalar& rho,
00486     volScalarField& vf
00487 )
00488 {
00489     tmp<fvScalarMatrix> tfvm
00490     (
00491         new fvScalarMatrix
00492         (
00493             vf,
00494             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00495         )
00496     );
00497     fvScalarMatrix& fvm = tfvm();
00498 
00499     scalar rDeltaT = 1.0/deltaT_();
00500 
00501     scalar deltaT = deltaT_();
00502     scalar deltaT0 = deltaT0_(vf);
00503 
00504     // Calculate unboundedness indicator
00505     // Note: all times moved by one because access to internal field
00506     // copies current field into the old-time level.
00507     scalarField phict =
00508         mag
00509         (
00510             vf.oldTime().oldTime().internalField()
00511           - vf.oldTime().oldTime().oldTime().internalField()
00512         )/
00513         (
00514             mag
00515             (
00516                 vf.oldTime().internalField()
00517               - vf.oldTime().oldTime().internalField()
00518             )
00519             + VSMALL
00520         );
00521 
00522     scalarField limiter(pos(phict) - pos(phict - 1.0));
00523 
00524     scalarField coefft   = 1.0 + limiter*deltaT/(deltaT + deltaT0);
00525     scalarField coefft00 = limiter*deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
00526     scalarField coefft0  = coefft + coefft00;
00527 
00528     fvm.diag() = (coefft*rDeltaT*rho.value())*mesh().V();
00529 
00530     if (mesh().moving())
00531     {
00532         fvm.source() = rDeltaT*rho.value()*
00533         (
00534             coefft0*vf.oldTime().internalField()*mesh().V0()
00535           - coefft00*vf.oldTime().oldTime().internalField()
00536            *mesh().V00()
00537         );
00538     }
00539     else
00540     {
00541         fvm.source() = rDeltaT*mesh().V()*rho.value()*
00542         (
00543             coefft0*vf.oldTime().internalField()
00544           - coefft00*vf.oldTime().oldTime().internalField()
00545         );
00546     }
00547 
00548     return tfvm;
00549 }
00550 
00551 
00552 tmp<fvScalarMatrix>
00553 boundedBackwardDdtScheme::fvmDdt
00554 (
00555     const volScalarField& rho,
00556     volScalarField& vf
00557 )
00558 {
00559     tmp<fvScalarMatrix> tfvm
00560     (
00561         new fvScalarMatrix
00562         (
00563             vf,
00564             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00565         )
00566     );
00567     fvScalarMatrix& fvm = tfvm();
00568 
00569     scalar rDeltaT = 1.0/deltaT_();
00570 
00571     scalar deltaT = deltaT_();
00572     scalar deltaT0 = deltaT0_(vf);
00573 
00574     // Calculate unboundedness indicator
00575     // Note: all times moved by one because access to internal field
00576     // copies current field into the old-time level.
00577     scalarField phict =
00578         mag
00579         (
00580             rho.oldTime().oldTime().internalField()*
00581             vf.oldTime().oldTime().internalField()
00582           - rho.oldTime().oldTime().oldTime().internalField()*
00583             vf.oldTime().oldTime().oldTime().internalField()
00584         )/
00585         (
00586             mag
00587             (
00588                 rho.oldTime().internalField()*
00589                 vf.oldTime().internalField()
00590               - rho.oldTime().oldTime().internalField()*
00591                 vf.oldTime().oldTime().internalField()
00592             )
00593             + VSMALL
00594         );
00595 
00596     scalarField limiter(pos(phict) - pos(phict - 1.0));
00597 
00598     scalarField coefft   = 1.0 + limiter*deltaT/(deltaT + deltaT0);
00599     scalarField coefft00 = limiter*deltaT*deltaT/(deltaT0*(deltaT + deltaT0));
00600     scalarField coefft0  = coefft + coefft00;
00601 
00602     fvm.diag() = (coefft*rDeltaT)*rho.internalField()*mesh().V();
00603 
00604     if (mesh().moving())
00605     {
00606         fvm.source() = rDeltaT*
00607         (
00608             coefft0*rho.oldTime().internalField()
00609            *vf.oldTime().internalField()*mesh().V0()
00610           - coefft00*rho.oldTime().oldTime().internalField()
00611            *vf.oldTime().oldTime().internalField()*mesh().V00()
00612         );
00613     }
00614     else
00615     {
00616         fvm.source() = rDeltaT*mesh().V()*
00617         (
00618             coefft0*rho.oldTime().internalField()
00619            *vf.oldTime().internalField()
00620           - coefft00*rho.oldTime().oldTime().internalField()
00621            *vf.oldTime().oldTime().internalField()
00622         );
00623     }
00624 
00625     return tfvm;
00626 }
00627 
00628 
00629 tmp<surfaceScalarField> boundedBackwardDdtScheme::fvcDdtPhiCorr
00630 (
00631     const volScalarField& rA,
00632     const volScalarField& U,
00633     const surfaceScalarField& phi
00634 )
00635 {
00636     notImplemented
00637     (
00638         "boundedBackwardDdtScheme::fvcDdtPhiCorr"
00639     );
00640 
00641     return surfaceScalarField::null();
00642 }
00643 
00644 
00645 tmp<surfaceScalarField> boundedBackwardDdtScheme::fvcDdtPhiCorr
00646 (
00647     const volScalarField& rA,
00648     const volScalarField& rho,
00649     const volScalarField& U,
00650     const surfaceScalarField& phi
00651 )
00652 {
00653     notImplemented
00654     (
00655         "boundedBackwardDdtScheme::fvcDdtPhiCorr"
00656     );
00657 
00658     return surfaceScalarField::null();
00659 }
00660 
00661 
00662 tmp<surfaceScalarField> boundedBackwardDdtScheme::meshPhi
00663 (
00664     const volScalarField& vf
00665 )
00666 {
00667     notImplemented
00668     (
00669         "boundedBackwardDdtScheme::meshPhi(const volScalarField& vf)"
00670     );
00671 
00672     return surfaceScalarField::null();
00673 }
00674 
00675 
00676 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00677 
00678 } // End namespace fv
00679 
00680 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00681 
00682 } // End namespace Foam
00683 
00684 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines