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

localEulerDdtScheme.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::localEulerDdtScheme
00026 
00027 Description
00028     Local time-step first-order Euler implicit/explicit ddt.
00029     The reciprocal of the local time-step field is looked-up from the
00030     database with the name provided.
00031 
00032     This scheme should only be used for steady-state computations
00033     using transient codes where local time-stepping is preferably to
00034     under-relaxation for transport consistency reasons.
00035 
00036     See also CoEulerDdtScheme.
00037 
00038 SourceFiles
00039     localEulerDdtScheme.C
00040     localEulerDdtSchemes.C
00041 
00042 \*---------------------------------------------------------------------------*/
00043 
00044 #ifndef localEulerDdtScheme_H
00045 #define localEulerDdtScheme_H
00046 
00047 #include <finiteVolume/ddtScheme.H>
00048 
00049 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00050 
00051 namespace Foam
00052 {
00053 
00054 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00055 
00056 namespace fv
00057 {
00058 
00059 /*---------------------------------------------------------------------------*\
00060                        Class localEulerDdtScheme Declaration
00061 \*---------------------------------------------------------------------------*/
00062 
00063 template<class Type>
00064 class localEulerDdtScheme
00065 :
00066     public fv::ddtScheme<Type>
00067 {
00068     // Private Data
00069 
00070         //- Name of the reciprocal local time-step field
00071         word rDeltaTName_;
00072 
00073 
00074     // Private Member Functions
00075 
00076         //- Disallow default bitwise copy construct
00077         localEulerDdtScheme(const localEulerDdtScheme&);
00078 
00079         //- Disallow default bitwise assignment
00080         void operator=(const localEulerDdtScheme&);
00081 
00082         //- Return the reciprocal of the local time-step
00083         const volScalarField& localRDeltaT() const;
00084 
00085 
00086 public:
00087 
00088     //- Runtime type information
00089     TypeName("localEuler");
00090 
00091 
00092     // Constructors
00093 
00094         //- Construct from mesh and Istream
00095         localEulerDdtScheme(const fvMesh& mesh, Istream& is)
00096         :
00097             ddtScheme<Type>(mesh, is),
00098             rDeltaTName_(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         tmp<surfaceScalarField> meshPhi
00167         (
00168             const GeometricField<Type, fvPatchField, volMesh>&
00169         );
00170 };
00171 
00172 
00173 template<>
00174 tmp<surfaceScalarField> localEulerDdtScheme<scalar>::fvcDdtPhiCorr
00175 (
00176     const volScalarField& rA,
00177     const volScalarField& U,
00178     const surfaceScalarField& phi
00179 );
00180 
00181 
00182 template<>
00183 tmp<surfaceScalarField> localEulerDdtScheme<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 "localEulerDdtScheme.C"
00204 #endif
00205 
00206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00207 
00208 #endif
00209 
00210 // ************************************************************************* //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines