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

backwardDdtScheme.H

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 Class
00025     Foam::fv::backwardDdtScheme
00026 
00027 Description
00028     Second-order backward-differencing ddt using the current and
00029     two previous time-step values.
00030 
00031 SourceFiles
00032     backwardDdtScheme.C
00033 
00034 \*---------------------------------------------------------------------------*/
00035 
00036 #ifndef backwardDdtScheme_H
00037 #define backwardDdtScheme_H
00038 
00039 #include <finiteVolume/ddtScheme.H>
00040 
00041 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00042 
00043 namespace Foam
00044 {
00045 
00046 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00047 
00048 namespace fv
00049 {
00050 
00051 /*---------------------------------------------------------------------------*\
00052                        Class backwardDdtScheme Declaration
00053 \*---------------------------------------------------------------------------*/
00054 
00055 template<class Type>
00056 class backwardDdtScheme
00057 :
00058     public fv::ddtScheme<Type>
00059 {
00060     // Private Member Functions
00061 
00062         //- Return the current time-step
00063         scalar deltaT_() const;
00064 
00065         //- Return the previous time-step
00066         scalar deltaT0_() const;
00067 
00068         //- Return the previous time-step or GREAT if the old timestep field
00069         //  wasn't available in which case Euler ddt is used
00070         template<class GeoField>
00071         scalar deltaT0_(const GeoField&) const;
00072 
00073         //- Disallow default bitwise copy construct
00074         backwardDdtScheme(const backwardDdtScheme&);
00075 
00076         //- Disallow default bitwise assignment
00077         void operator=(const backwardDdtScheme&);
00078 
00079 
00080 public:
00081 
00082     //- Runtime type information
00083     TypeName("backward");
00084 
00085 
00086     // Constructors
00087 
00088         //- Construct from mesh
00089         backwardDdtScheme(const fvMesh& mesh)
00090         :
00091             ddtScheme<Type>(mesh)
00092         {}
00093 
00094         //- Construct from mesh and Istream
00095         backwardDdtScheme(const fvMesh& mesh, Istream& is)
00096         :
00097             ddtScheme<Type>(mesh, is)
00098         {}
00099 
00100 
00101     // Member Functions
00102 
00103         //- Return mesh reference
00104         const fvMesh& mesh() const
00105         {
00106             return fv::ddtScheme<Type>::mesh();
00107         }
00108 
00109         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00110         (
00111             const dimensioned<Type>&
00112         );
00113 
00114         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00115         (
00116             const GeometricField<Type, fvPatchField, volMesh>&
00117         );
00118 
00119         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00120         (
00121             const dimensionedScalar&,
00122             const GeometricField<Type, fvPatchField, volMesh>&
00123         );
00124 
00125         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00126         (
00127             const volScalarField&,
00128             const GeometricField<Type, fvPatchField, volMesh>&
00129         );
00130 
00131         tmp<fvMatrix<Type> > fvmDdt
00132         (
00133             GeometricField<Type, fvPatchField, volMesh>&
00134         );
00135 
00136         tmp<fvMatrix<Type> > fvmDdt
00137         (
00138             const dimensionedScalar&,
00139             GeometricField<Type, fvPatchField, volMesh>&
00140         );
00141 
00142         tmp<fvMatrix<Type> > fvmDdt
00143         (
00144             const volScalarField&,
00145             GeometricField<Type, fvPatchField, volMesh>&
00146         );
00147 
00148         typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
00149 
00150         tmp<fluxFieldType> fvcDdtPhiCorr
00151         (
00152             const volScalarField& rA,
00153             const GeometricField<Type, fvPatchField, volMesh>& U,
00154             const fluxFieldType& phi
00155         );
00156 
00157         tmp<fluxFieldType> fvcDdtPhiCorr
00158         (
00159             const volScalarField& rA,
00160             const volScalarField& rho,
00161             const GeometricField<Type, fvPatchField, volMesh>& U,
00162             const fluxFieldType& phi
00163         );
00164 
00165 
00166         tmp<surfaceScalarField> meshPhi
00167         (
00168             const GeometricField<Type, fvPatchField, volMesh>&
00169         );
00170 };
00171 
00172 
00173 template<>
00174 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
00175 (
00176     const volScalarField& rA,
00177     const volScalarField& U,
00178     const surfaceScalarField& phi
00179 );
00180 
00181 
00182 template<>
00183 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
00184 (
00185     const volScalarField& rA,
00186     const volScalarField& rho,
00187     const volScalarField& U,
00188     const surfaceScalarField& phi
00189 );
00190 
00191 
00192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00193 
00194 } // End namespace fv
00195 
00196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00197 
00198 } // End namespace Foam
00199 
00200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00201 
00202 #ifdef NoRepository
00203 #   include <finiteVolume/backwardDdtScheme.C>
00204 #endif
00205 
00206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00207 
00208 #endif
00209 
00210 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines