initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / finiteVolume / ddtSchemes / backwardDdtScheme / backwardDdtScheme.H
blob374d05031aa5cbca89f75c93e3600d0914996549
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::backwardDdtScheme
28 Description
29     Second-order backward-differencing ddt using the current and
30     two previous time-step values.
32 SourceFiles
33     backwardDdtScheme.C
35 \*---------------------------------------------------------------------------*/
37 #ifndef backwardDdtScheme_H
38 #define backwardDdtScheme_H
40 #include "ddtScheme.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace fv
52 /*---------------------------------------------------------------------------*\
53                        Class backwardDdtScheme Declaration
54 \*---------------------------------------------------------------------------*/
56 template<class Type>
57 class backwardDdtScheme
59     public fv::ddtScheme<Type>
61     // Private Member Functions
63         //- Return the current time-step
64         scalar deltaT_() const;
66         //- Return the previous time-step
67         scalar deltaT0_() const;
69         //- Return the previous time-step or GREAT if the old timestep field
70         //  wasn't available in which case Euler ddt is used
71         template<class GeoField>
72         scalar deltaT0_(const GeoField&) const;
74         //- Disallow default bitwise copy construct
75         backwardDdtScheme(const backwardDdtScheme&);
77         //- Disallow default bitwise assignment
78         void operator=(const backwardDdtScheme&);
81 public:
83     //- Runtime type information
84     TypeName("backward");
87     // Constructors
89         //- Construct from mesh
90         backwardDdtScheme(const fvMesh& mesh)
91         :
92             ddtScheme<Type>(mesh)
93         {}
95         //- Construct from mesh and Istream
96         backwardDdtScheme(const fvMesh& mesh, Istream& is)
97         :
98             ddtScheme<Type>(mesh, is)
99         {}
102     // Member Functions
104         //- Return mesh reference
105         const fvMesh& mesh() const
106         {
107             return fv::ddtScheme<Type>::mesh();
108         }
110         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
111         (
112             const dimensioned<Type>&
113         );
115         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
116         (
117             const GeometricField<Type, fvPatchField, volMesh>&
118         );
120         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
121         (
122             const dimensionedScalar&,
123             const GeometricField<Type, fvPatchField, volMesh>&
124         );
126         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
127         (
128             const volScalarField&,
129             const GeometricField<Type, fvPatchField, volMesh>&
130         );
132         tmp<fvMatrix<Type> > fvmDdt
133         (
134             GeometricField<Type, fvPatchField, volMesh>&
135         );
137         tmp<fvMatrix<Type> > fvmDdt
138         (
139             const dimensionedScalar&,
140             GeometricField<Type, fvPatchField, volMesh>&
141         );
143         tmp<fvMatrix<Type> > fvmDdt
144         (
145             const volScalarField&,
146             GeometricField<Type, fvPatchField, volMesh>&
147         );
149         typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
151         tmp<fluxFieldType> fvcDdtPhiCorr
152         (
153             const volScalarField& rA,
154             const GeometricField<Type, fvPatchField, volMesh>& U,
155             const fluxFieldType& phi
156         );
158         tmp<fluxFieldType> fvcDdtPhiCorr
159         (
160             const volScalarField& rA,
161             const volScalarField& rho,
162             const GeometricField<Type, fvPatchField, volMesh>& U,
163             const fluxFieldType& phi
164         );
167         tmp<surfaceScalarField> meshPhi
168         (
169             const GeometricField<Type, fvPatchField, volMesh>&
170         );
174 template<>
175 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
177     const volScalarField& rA,
178     const volScalarField& U,
179     const surfaceScalarField& phi
183 template<>
184 tmp<surfaceScalarField> backwardDdtScheme<scalar>::fvcDdtPhiCorr
186     const volScalarField& rA,
187     const volScalarField& rho,
188     const volScalarField& U,
189     const surfaceScalarField& phi
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 } // End namespace fv
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 } // End namespace Foam
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 #ifdef NoRepository
204 #   include "backwardDdtScheme.C"
205 #endif
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 #endif
211 // ************************************************************************* //