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

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