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

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