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

localEulerDdtScheme.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 "localEulerDdtScheme.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 const volScalarField& localEulerDdtScheme<Type>::localRDeltaT() const
00044 {
00045     return mesh().objectRegistry::lookupObject<volScalarField>(rDeltaTName_);
00046 }
00047 
00048 
00049 template<class Type>
00050 tmp<GeometricField<Type, fvPatchField, volMesh> >
00051 localEulerDdtScheme<Type>::fvcDdt
00052 (
00053     const dimensioned<Type>& dt
00054 )
00055 {
00056     const volScalarField& rDeltaT = localRDeltaT();
00057 
00058     IOobject ddtIOobject
00059     (
00060         "ddt(" + dt.name() + ')',
00061         mesh().time().timeName(),
00062         mesh()
00063     );
00064 
00065     if (mesh().moving())
00066     {
00067         tmp<GeometricField<Type, fvPatchField, volMesh> > tdtdt
00068         (
00069             new GeometricField<Type, fvPatchField, volMesh>
00070             (
00071                 ddtIOobject,
00072                 mesh(),
00073                 dimensioned<Type>
00074                 (
00075                     "0",
00076                     dt.dimensions()/dimTime,
00077                     pTraits<Type>::zero
00078                 )
00079             )
00080         );
00081 
00082         tdtdt().internalField() =
00083             rDeltaT.internalField()*dt.value()*(1.0 - mesh().V0()/mesh().V());
00084 
00085         return tdtdt;
00086     }
00087     else
00088     {
00089         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00090         (
00091             new GeometricField<Type, fvPatchField, volMesh>
00092             (
00093                 ddtIOobject,
00094                 mesh(),
00095                 dimensioned<Type>
00096                 (
00097                     "0",
00098                     dt.dimensions()/dimTime,
00099                     pTraits<Type>::zero
00100                 ),
00101                 calculatedFvPatchField<Type>::typeName
00102             )
00103         );
00104     }
00105 }
00106 
00107 
00108 template<class Type>
00109 tmp<GeometricField<Type, fvPatchField, volMesh> >
00110 localEulerDdtScheme<Type>::fvcDdt
00111 (
00112     const GeometricField<Type, fvPatchField, volMesh>& vf
00113 )
00114 {
00115     const volScalarField& rDeltaT = localRDeltaT();
00116 
00117     IOobject ddtIOobject
00118     (
00119         "ddt(" + vf.name() + ')',
00120         mesh().time().timeName(),
00121         mesh()
00122     );
00123 
00124     if (mesh().moving())
00125     {
00126         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00127         (
00128             new GeometricField<Type, fvPatchField, volMesh>
00129             (
00130                 ddtIOobject,
00131                 mesh(),
00132                 rDeltaT.dimensions()*vf.dimensions(),
00133                 rDeltaT.internalField()*
00134                 (
00135                     vf.internalField()
00136                   - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00137                 ),
00138                 rDeltaT.boundaryField()*
00139                 (
00140                     vf.boundaryField() - vf.oldTime().boundaryField()
00141                 )
00142             )
00143         );
00144     }
00145     else
00146     {
00147         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00148         (
00149             new GeometricField<Type, fvPatchField, volMesh>
00150             (
00151                 ddtIOobject,
00152                 rDeltaT*(vf - vf.oldTime())
00153             )
00154         );
00155     }
00156 }
00157 
00158 
00159 template<class Type>
00160 tmp<GeometricField<Type, fvPatchField, volMesh> >
00161 localEulerDdtScheme<Type>::fvcDdt
00162 (
00163     const dimensionedScalar& rho,
00164     const GeometricField<Type, fvPatchField, volMesh>& vf
00165 )
00166 {
00167     const volScalarField& rDeltaT = localRDeltaT();
00168 
00169     IOobject ddtIOobject
00170     (
00171         "ddt(" + rho.name() + ',' + vf.name() + ')',
00172         mesh().time().timeName(),
00173         mesh()
00174     );
00175 
00176     if (mesh().moving())
00177     {
00178         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00179         (
00180             new GeometricField<Type, fvPatchField, volMesh>
00181             (
00182                 ddtIOobject,
00183                 mesh(),
00184                 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00185                 rDeltaT.internalField()*rho.value()*
00186                 (
00187                     vf.internalField()
00188                   - vf.oldTime().internalField()*mesh().V0()/mesh().V()
00189                 ),
00190                 rDeltaT.boundaryField()*rho.value()*
00191                 (
00192                     vf.boundaryField() - vf.oldTime().boundaryField()
00193                 )
00194             )
00195         );
00196     }
00197     else
00198     {
00199         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00200         (
00201             new GeometricField<Type, fvPatchField, volMesh>
00202             (
00203                 ddtIOobject,
00204                 rDeltaT*rho*(vf - vf.oldTime())
00205             )
00206         );
00207     }
00208 }
00209 
00210 
00211 template<class Type>
00212 tmp<GeometricField<Type, fvPatchField, volMesh> >
00213 localEulerDdtScheme<Type>::fvcDdt
00214 (
00215     const volScalarField& rho,
00216     const GeometricField<Type, fvPatchField, volMesh>& vf
00217 )
00218 {
00219     const volScalarField& rDeltaT = localRDeltaT();
00220 
00221     IOobject ddtIOobject
00222     (
00223         "ddt(" + rho.name() + ',' + vf.name() + ')',
00224         mesh().time().timeName(),
00225         mesh()
00226     );
00227 
00228     if (mesh().moving())
00229     {
00230         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00231         (
00232             new GeometricField<Type, fvPatchField, volMesh>
00233             (
00234                 ddtIOobject,
00235                 mesh(),
00236                 rDeltaT.dimensions()*rho.dimensions()*vf.dimensions(),
00237                 rDeltaT.internalField()*
00238                 (
00239                     rho.internalField()*vf.internalField()
00240                   - rho.oldTime().internalField()
00241                    *vf.oldTime().internalField()*mesh().V0()/mesh().V()
00242                 ),
00243                 rDeltaT.boundaryField()*
00244                 (
00245                     rho.boundaryField()*vf.boundaryField()
00246                   - rho.oldTime().boundaryField()
00247                    *vf.oldTime().boundaryField()
00248                 )
00249             )
00250         );
00251     }
00252     else
00253     {
00254         return tmp<GeometricField<Type, fvPatchField, volMesh> >
00255         (
00256             new GeometricField<Type, fvPatchField, volMesh>
00257             (
00258                 ddtIOobject,
00259                 rDeltaT*(rho*vf - rho.oldTime()*vf.oldTime())
00260             )
00261         );
00262     }
00263 }
00264 
00265 
00266 template<class Type>
00267 tmp<fvMatrix<Type> >
00268 localEulerDdtScheme<Type>::fvmDdt
00269 (
00270     GeometricField<Type, fvPatchField, volMesh>& vf
00271 )
00272 {
00273     tmp<fvMatrix<Type> > tfvm
00274     (
00275         new fvMatrix<Type>
00276         (
00277             vf,
00278             vf.dimensions()*dimVol/dimTime
00279         )
00280     );
00281 
00282     fvMatrix<Type>& fvm = tfvm();
00283 
00284     const scalarField& rDeltaT = localRDeltaT().internalField();
00285 
00286     fvm.diag() = rDeltaT*mesh().V();
00287 
00288     if (mesh().moving())
00289     {
00290         fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V0();
00291     }
00292     else
00293     {
00294         fvm.source() = rDeltaT*vf.oldTime().internalField()*mesh().V();
00295     }
00296 
00297     return tfvm;
00298 }
00299 
00300 
00301 template<class Type>
00302 tmp<fvMatrix<Type> >
00303 localEulerDdtScheme<Type>::fvmDdt
00304 (
00305     const dimensionedScalar& rho,
00306     GeometricField<Type, fvPatchField, volMesh>& vf
00307 )
00308 {
00309     tmp<fvMatrix<Type> > tfvm
00310     (
00311         new fvMatrix<Type>
00312         (
00313             vf,
00314             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00315         )
00316     );
00317     fvMatrix<Type>& fvm = tfvm();
00318 
00319     const scalarField& rDeltaT = localRDeltaT().internalField();
00320 
00321     fvm.diag() = rDeltaT*rho.value()*mesh().V();
00322 
00323     if (mesh().moving())
00324     {
00325         fvm.source() = rDeltaT
00326             *rho.value()*vf.oldTime().internalField()*mesh().V0();
00327     }
00328     else
00329     {
00330         fvm.source() = rDeltaT
00331             *rho.value()*vf.oldTime().internalField()*mesh().V();
00332     }
00333 
00334     return tfvm;
00335 }
00336 
00337 
00338 template<class Type>
00339 tmp<fvMatrix<Type> >
00340 localEulerDdtScheme<Type>::fvmDdt
00341 (
00342     const volScalarField& rho,
00343     GeometricField<Type, fvPatchField, volMesh>& vf
00344 )
00345 {
00346     tmp<fvMatrix<Type> > tfvm
00347     (
00348         new fvMatrix<Type>
00349         (
00350             vf,
00351             rho.dimensions()*vf.dimensions()*dimVol/dimTime
00352         )
00353     );
00354     fvMatrix<Type>& fvm = tfvm();
00355 
00356     const scalarField& rDeltaT = localRDeltaT().internalField();
00357 
00358     fvm.diag() = rDeltaT*rho.internalField()*mesh().V();
00359 
00360     if (mesh().moving())
00361     {
00362         fvm.source() = rDeltaT
00363             *rho.oldTime().internalField()
00364             *vf.oldTime().internalField()*mesh().V0();
00365     }
00366     else
00367     {
00368         fvm.source() = rDeltaT
00369             *rho.oldTime().internalField()
00370             *vf.oldTime().internalField()*mesh().V();
00371     }
00372 
00373     return tfvm;
00374 }
00375 
00376 
00377 template<class Type>
00378 tmp<typename localEulerDdtScheme<Type>::fluxFieldType>
00379 localEulerDdtScheme<Type>::fvcDdtPhiCorr
00380 (
00381     const volScalarField& rA,
00382     const GeometricField<Type, fvPatchField, volMesh>& U,
00383     const fluxFieldType& phi
00384 )
00385 {
00386     IOobject ddtIOobject
00387     (
00388         "ddtPhiCorr(" + rA.name() + ',' + U.name() + ',' + phi.name() + ')',
00389         mesh().time().timeName(),
00390         mesh()
00391     );
00392 
00393     if (mesh().moving())
00394     {
00395         return tmp<fluxFieldType>
00396         (
00397             new fluxFieldType
00398             (
00399                 ddtIOobject,
00400                 mesh(),
00401                 dimensioned<typename flux<Type>::type>
00402                 (
00403                     "0",
00404                     rA.dimensions()*phi.dimensions()/dimTime,
00405                     pTraits<typename flux<Type>::type>::zero
00406                 )
00407             )
00408         );
00409     }
00410     else
00411     {
00412         const volScalarField& rDeltaT = localRDeltaT();
00413 
00414         return tmp<fluxFieldType>
00415         (
00416             new fluxFieldType
00417             (
00418                 ddtIOobject,
00419                 this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())*
00420                 (
00421                     fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00422                   - (fvc::interpolate(rDeltaT*rA*U.oldTime()) & mesh().Sf())
00423                 )
00424             )
00425         );
00426     }
00427 }
00428 
00429 
00430 template<class Type>
00431 tmp<typename localEulerDdtScheme<Type>::fluxFieldType>
00432 localEulerDdtScheme<Type>::fvcDdtPhiCorr
00433 (
00434     const volScalarField& rA,
00435     const volScalarField& rho,
00436     const GeometricField<Type, fvPatchField, volMesh>& U,
00437     const fluxFieldType& phi
00438 )
00439 {
00440     IOobject ddtIOobject
00441     (
00442         "ddtPhiCorr("
00443       + rA.name() + ',' + rho.name() + ',' + U.name() + ',' + phi.name() + ')',
00444         mesh().time().timeName(),
00445         mesh()
00446     );
00447 
00448     if (mesh().moving())
00449     {
00450         return tmp<fluxFieldType>
00451         (
00452             new fluxFieldType
00453             (
00454                 ddtIOobject,
00455                 mesh(),
00456                 dimensioned<typename flux<Type>::type>
00457                 (
00458                     "0",
00459                     rA.dimensions()*rho.dimensions()*phi.dimensions()/dimTime,
00460                     pTraits<typename flux<Type>::type>::zero
00461                 )
00462             )
00463         );
00464     }
00465     else
00466     {
00467         const volScalarField& rDeltaT = localRDeltaT();
00468 
00469         if
00470         (
00471             U.dimensions() == dimVelocity
00472          && phi.dimensions() == dimVelocity*dimArea
00473         )
00474         {
00475             return tmp<fluxFieldType>
00476             (
00477                 new fluxFieldType
00478                 (
00479                     ddtIOobject,
00480 //                     (
00481 //                         scalar(1)
00482 //                       - (min(rDeltaT)/fvc::interpolate(rDeltaT))
00483 //                        *(scalar(1) - fvcDdtPhiCoeff(U.oldTime(), phi.oldTime()))
00484 //                     )
00485                     this->fvcDdtPhiCoeff(U.oldTime(), phi.oldTime())
00486                     //0.95
00487                    *(
00488                         fvc::interpolate(rDeltaT*rA*rho.oldTime())*phi.oldTime()
00489                       - (fvc::interpolate(rDeltaT*rA*rho.oldTime()*U.oldTime())
00490                       & mesh().Sf())
00491                     )
00492                 )
00493             );
00494         }
00495         else if
00496         (
00497             U.dimensions() == dimVelocity
00498          && phi.dimensions() == dimDensity*dimVelocity*dimArea
00499         )
00500         {
00501             return tmp<fluxFieldType>
00502             (
00503                 new fluxFieldType
00504                 (
00505                     ddtIOobject,
00506                     this->fvcDdtPhiCoeff
00507                     (
00508                         U.oldTime(),
00509                         phi.oldTime()/fvc::interpolate(rho.oldTime())
00510                     )
00511                    *(
00512                         fvc::interpolate(rDeltaT*rA*rho.oldTime())
00513                        *phi.oldTime()/fvc::interpolate(rho.oldTime())
00514                       - (
00515                             fvc::interpolate
00516                             (
00517                                 rDeltaT*rA*rho.oldTime()*U.oldTime()
00518                             ) & mesh().Sf()
00519                         )
00520                     )
00521                 )
00522             );
00523         }
00524         else if
00525         (
00526             U.dimensions() == dimDensity*dimVelocity
00527          && phi.dimensions() == dimDensity*dimVelocity*dimArea
00528         )
00529         {
00530             return tmp<fluxFieldType>
00531             (
00532                 new fluxFieldType
00533                 (
00534                     ddtIOobject,
00535                     this->fvcDdtPhiCoeff
00536                     (
00537                         rho.oldTime(),
00538                         U.oldTime(),
00539                         phi.oldTime()
00540                     )
00541                    *(
00542                         fvc::interpolate(rDeltaT*rA)*phi.oldTime()
00543                       - (
00544                             fvc::interpolate(rDeltaT*rA*U.oldTime())&mesh().Sf()
00545                         )
00546                     )
00547                 )
00548             );
00549         }
00550         else
00551         {
00552             FatalErrorIn
00553             (
00554                 "localEulerDdtScheme<Type>::fvcDdtPhiCorr"
00555             )   << "dimensions of phi are not correct"
00556                 << abort(FatalError);
00557 
00558             return fluxFieldType::null();
00559         }
00560     }
00561 }
00562 
00563 
00564 template<class Type>
00565 tmp<surfaceScalarField> localEulerDdtScheme<Type>::meshPhi
00566 (
00567     const GeometricField<Type, fvPatchField, volMesh>&
00568 )
00569 {
00570     return tmp<surfaceScalarField>
00571     (
00572         new surfaceScalarField
00573         (
00574             IOobject
00575             (
00576                 "meshPhi",
00577                 mesh().time().timeName(),
00578                 mesh()
00579             ),
00580             mesh(),
00581             dimensionedScalar("0", dimVolume/dimTime, 0.0)
00582         )
00583     );
00584 }
00585 
00586 
00587 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00588 
00589 } // End namespace fv
00590 
00591 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00592 
00593 } // End namespace Foam
00594 
00595 // ************************************************************************* //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines