initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / ddtSchemes / CoEulerDdtScheme / CoEulerDdtScheme.H
blobf2b54169ed5f31e4c02a7dea2a1d41efa7c177a5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::fv::CoEulerDdtScheme
28 Description
29     Courant number limited first-order Euler implicit/explicit ddt.
31     The time-step is adjusted locally so that the local Courant number
32     does not exceed the specified value.
34     This scheme should only be used for steady-state computations
35     using transient codes where local time-stepping is preferably to
36     under-relaxation for transport consistency reasons.
38 SourceFiles
39     CoEulerDdtScheme.C
41 \*---------------------------------------------------------------------------*/
43 #ifndef CoEulerDdtScheme_H
44 #define CoEulerDdtScheme_H
46 #include "ddtScheme.H"
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 namespace Foam
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 namespace fv
58 /*---------------------------------------------------------------------------*\
59                        Class CoEulerDdtScheme Declaration
60 \*---------------------------------------------------------------------------*/
62 template<class Type>
63 class CoEulerDdtScheme
65     public fv::ddtScheme<Type>
67     // Private Data
69         //- Name of the flux field used to calculate the local time-step
70         word phiName_;
72         //- Name of the density field used to obtain the volumetric flux
73         //  from the mass flux if required
74         word rhoName_;
76         //- Maximum local Courant number
77         scalar maxCo_;
80     // Private Member Functions
82         //- Disallow default bitwise copy construct
83         CoEulerDdtScheme(const CoEulerDdtScheme&);
85         //- Disallow default bitwise assignment
86         void operator=(const CoEulerDdtScheme&);
88         //- Return the reciprocal of the Courant-number limited time-step
89         tmp<volScalarField> CorDeltaT() const;
91         //- Return the reciprocal of the face-Courant-number limited time-step
92         tmp<surfaceScalarField> CofrDeltaT() const;
95 public:
97     //- Runtime type information
98     TypeName("CoEuler");
101     // Constructors
103         //- Construct from mesh and Istream
104         CoEulerDdtScheme(const fvMesh& mesh, Istream& is)
105         :
106             ddtScheme<Type>(mesh, is),
107             phiName_(is),
108             rhoName_(is),
109             maxCo_(readScalar(is))
110         {}
113     // Member Functions
115         //- Return mesh reference
116         const fvMesh& mesh() const
117         {
118             return fv::ddtScheme<Type>::mesh();
119         }
121         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
122         (
123             const dimensioned<Type>&
124         );
126         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
127         (
128             const GeometricField<Type, fvPatchField, volMesh>&
129         );
131         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
132         (
133             const dimensionedScalar&,
134             const GeometricField<Type, fvPatchField, volMesh>&
135         );
137         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
138         (
139             const volScalarField&,
140             const GeometricField<Type, fvPatchField, volMesh>&
141         );
143         tmp<fvMatrix<Type> > fvmDdt
144         (
145             GeometricField<Type, fvPatchField, volMesh>&
146         );
148         tmp<fvMatrix<Type> > fvmDdt
149         (
150             const dimensionedScalar&,
151             GeometricField<Type, fvPatchField, volMesh>&
152         );
154         tmp<fvMatrix<Type> > fvmDdt
155         (
156             const volScalarField&,
157             GeometricField<Type, fvPatchField, volMesh>&
158         );
160         typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
162         tmp<fluxFieldType> fvcDdtPhiCorr
163         (
164             const volScalarField& rA,
165             const GeometricField<Type, fvPatchField, volMesh>& U,
166             const fluxFieldType& phi
167         );
169         tmp<fluxFieldType> fvcDdtPhiCorr
170         (
171             const volScalarField& rA,
172             const volScalarField& rho,
173             const GeometricField<Type, fvPatchField, volMesh>& U,
174             const fluxFieldType& phi
175         );
177         tmp<surfaceScalarField> meshPhi
178         (
179             const GeometricField<Type, fvPatchField, volMesh>&
180         );
184 template<>
185 tmp<surfaceScalarField> CoEulerDdtScheme<scalar>::fvcDdtPhiCorr
187     const volScalarField& rA,
188     const volScalarField& U,
189     const surfaceScalarField& phi
193 template<>
194 tmp<surfaceScalarField> CoEulerDdtScheme<scalar>::fvcDdtPhiCorr
196     const volScalarField& rA,
197     const volScalarField& rho,
198     const volScalarField& U,
199     const surfaceScalarField& phi
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 } // End namespace fv
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 } // End namespace Foam
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 #ifdef NoRepository
214 #   include "CoEulerDdtScheme.C"
215 #endif
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 #endif
221 // ************************************************************************* //