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

CoEulerDdtScheme.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::CoEulerDdtScheme
00026 
00027 Description
00028     Courant number limited first-order Euler implicit/explicit ddt.
00029 
00030     The time-step is adjusted locally so that the local Courant number
00031     does not exceed the specified value.
00032 
00033     This scheme should only be used for steady-state computations
00034     using transient codes where local time-stepping is preferably to
00035     under-relaxation for transport consistency reasons.
00036 
00037 SourceFiles
00038     CoEulerDdtScheme.C
00039 
00040 \*---------------------------------------------------------------------------*/
00041 
00042 #ifndef CoEulerDdtScheme_H
00043 #define CoEulerDdtScheme_H
00044 
00045 #include <finiteVolume/ddtScheme.H>
00046 
00047 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00048 
00049 namespace Foam
00050 {
00051 
00052 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00053 
00054 namespace fv
00055 {
00056 
00057 /*---------------------------------------------------------------------------*\
00058                        Class CoEulerDdtScheme Declaration
00059 \*---------------------------------------------------------------------------*/
00060 
00061 template<class Type>
00062 class CoEulerDdtScheme
00063 :
00064     public fv::ddtScheme<Type>
00065 {
00066     // Private Data
00067 
00068         //- Name of the flux field used to calculate the local time-step
00069         word phiName_;
00070 
00071         //- Name of the density field used to obtain the volumetric flux
00072         //  from the mass flux if required
00073         word rhoName_;
00074 
00075         //- Maximum local Courant number
00076         scalar maxCo_;
00077 
00078 
00079     // Private Member Functions
00080 
00081         //- Disallow default bitwise copy construct
00082         CoEulerDdtScheme(const CoEulerDdtScheme&);
00083 
00084         //- Disallow default bitwise assignment
00085         void operator=(const CoEulerDdtScheme&);
00086 
00087         //- Return the reciprocal of the Courant-number limited time-step
00088         tmp<volScalarField> CorDeltaT() const;
00089 
00090         //- Return the reciprocal of the face-Courant-number limited time-step
00091         tmp<surfaceScalarField> CofrDeltaT() const;
00092 
00093 
00094 public:
00095 
00096     //- Runtime type information
00097     TypeName("CoEuler");
00098 
00099 
00100     // Constructors
00101 
00102         //- Construct from mesh and Istream
00103         CoEulerDdtScheme(const fvMesh& mesh, Istream& is)
00104         :
00105             ddtScheme<Type>(mesh, is),
00106             phiName_(is),
00107             rhoName_(is),
00108             maxCo_(readScalar(is))
00109         {}
00110 
00111 
00112     // Member Functions
00113 
00114         //- Return mesh reference
00115         const fvMesh& mesh() const
00116         {
00117             return fv::ddtScheme<Type>::mesh();
00118         }
00119 
00120         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00121         (
00122             const dimensioned<Type>&
00123         );
00124 
00125         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00126         (
00127             const GeometricField<Type, fvPatchField, volMesh>&
00128         );
00129 
00130         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00131         (
00132             const dimensionedScalar&,
00133             const GeometricField<Type, fvPatchField, volMesh>&
00134         );
00135 
00136         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
00137         (
00138             const volScalarField&,
00139             const GeometricField<Type, fvPatchField, volMesh>&
00140         );
00141 
00142         tmp<fvMatrix<Type> > fvmDdt
00143         (
00144             GeometricField<Type, fvPatchField, volMesh>&
00145         );
00146 
00147         tmp<fvMatrix<Type> > fvmDdt
00148         (
00149             const dimensionedScalar&,
00150             GeometricField<Type, fvPatchField, volMesh>&
00151         );
00152 
00153         tmp<fvMatrix<Type> > fvmDdt
00154         (
00155             const volScalarField&,
00156             GeometricField<Type, fvPatchField, volMesh>&
00157         );
00158 
00159         typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
00160 
00161         tmp<fluxFieldType> fvcDdtPhiCorr
00162         (
00163             const volScalarField& rA,
00164             const GeometricField<Type, fvPatchField, volMesh>& U,
00165             const fluxFieldType& phi
00166         );
00167 
00168         tmp<fluxFieldType> fvcDdtPhiCorr
00169         (
00170             const volScalarField& rA,
00171             const volScalarField& rho,
00172             const GeometricField<Type, fvPatchField, volMesh>& U,
00173             const fluxFieldType& phi
00174         );
00175 
00176         tmp<surfaceScalarField> meshPhi
00177         (
00178             const GeometricField<Type, fvPatchField, volMesh>&
00179         );
00180 };
00181 
00182 
00183 template<>
00184 tmp<surfaceScalarField> CoEulerDdtScheme<scalar>::fvcDdtPhiCorr
00185 (
00186     const volScalarField& rA,
00187     const volScalarField& U,
00188     const surfaceScalarField& phi
00189 );
00190 
00191 
00192 template<>
00193 tmp<surfaceScalarField> CoEulerDdtScheme<scalar>::fvcDdtPhiCorr
00194 (
00195     const volScalarField& rA,
00196     const volScalarField& rho,
00197     const volScalarField& U,
00198     const surfaceScalarField& phi
00199 );
00200 
00201 
00202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00203 
00204 } // End namespace fv
00205 
00206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00207 
00208 } // End namespace Foam
00209 
00210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00211 
00212 #ifdef NoRepository
00213 #   include <finiteVolume/CoEulerDdtScheme.C>
00214 #endif
00215 
00216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00217 
00218 #endif
00219 
00220 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines