OpenFOAM logo
Open Source CFD Toolkit

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