Forward compatibility: flex
[foam-extend-3.2.git] / src / dbns / timeStepping / EulerLocalDdtScheme / EulerLocalDdtScheme.H
blob0ab283ca98b7e2765a71c4868aa791f81be9c4ee
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::fv::EulerLocalDdtScheme
27 Description
28     Local time-setp limited first-order Euler implicit/explicit ddt.
30     The local time-step field has to be registered and provided outside of this
31     class. If the physical time-step deltaT is already divided by a multi-stage
32     coefficient this scheme can be also used for multi-stage time integration.
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 Author
39     Oliver Borm
40     Aleksandar Jemcov
42 SourceFiles
43     EulerLocalDdtScheme.C
45 \*---------------------------------------------------------------------------*/
47 #ifndef EulerLocalDdtScheme_H
48 #define EulerLocalDdtScheme_H
50 #include "ddtScheme.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 namespace Foam
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 namespace fv
62 /*---------------------------------------------------------------------------*\
63                        Class EulerLocalDdtScheme Declaration
64 \*---------------------------------------------------------------------------*/
66 template<class Type>
67 class EulerLocalDdtScheme
69     public fv::ddtScheme<Type>
71     // Private Data
73         //- Name of the physical time-step (field)
74         word deltaTName_;
76         //- Name of the local time-step field
77         word deltaTauName_;
79     // Private Member Functions
81         //- Disallow default bitwise copy construct
82         EulerLocalDdtScheme(const EulerLocalDdtScheme&);
84         //- Disallow default bitwise assignment
85         void operator=(const EulerLocalDdtScheme&);
87 public:
89     //- Runtime type information
90     TypeName("EulerLocal");
93     // Constructors
95         //- Construct from mesh and Istream
96         EulerLocalDdtScheme(const fvMesh& mesh, Istream& is)
97         :
98             ddtScheme<Type>(mesh, is),
99             deltaTName_(is),
100             deltaTauName_(is)
101         {}
104     // Member Functions
106         //- Return mesh reference
107         const fvMesh& mesh() const
108         {
109             return fv::ddtScheme<Type>::mesh();
110         }
112         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
113         (
114             const dimensioned<Type>&
115         );
117         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
118         (
119             const GeometricField<Type, fvPatchField, volMesh>&
120         );
122         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
123         (
124             const dimensionedScalar&,
125             const GeometricField<Type, fvPatchField, volMesh>&
126         );
128         tmp<GeometricField<Type, fvPatchField, volMesh> > fvcDdt
129         (
130             const volScalarField&,
131             const GeometricField<Type, fvPatchField, volMesh>&
132         );
134         tmp<fvMatrix<Type> > fvmDdt
135         (
136             GeometricField<Type, fvPatchField, volMesh>&
137         );
139         tmp<fvMatrix<Type> > fvmDdt
140         (
141             const dimensionedScalar&,
142             GeometricField<Type, fvPatchField, volMesh>&
143         );
145         tmp<fvMatrix<Type> > fvmDdt
146         (
147             const volScalarField&,
148             GeometricField<Type, fvPatchField, volMesh>&
149         );
151         typedef typename ddtScheme<Type>::fluxFieldType fluxFieldType;
153         tmp<fluxFieldType> fvcDdtPhiCorr
154         (
155             const volScalarField& rA,
156             const GeometricField<Type, fvPatchField, volMesh>& U,
157             const fluxFieldType& phi
158         );
160         tmp<fluxFieldType> fvcDdtPhiCorr
161         (
162             const volScalarField& rA,
163             const volScalarField& rho,
164             const GeometricField<Type, fvPatchField, volMesh>& U,
165             const fluxFieldType& phi
166         );
168         tmp<surfaceScalarField> meshPhi
169         (
170             const GeometricField<Type, fvPatchField, volMesh>&
171         );
175 template<>
176 tmp<surfaceScalarField> EulerLocalDdtScheme<scalar>::fvcDdtPhiCorr
178     const volScalarField& rA,
179     const volScalarField& U,
180     const surfaceScalarField& phi
184 template<>
185 tmp<surfaceScalarField> EulerLocalDdtScheme<scalar>::fvcDdtPhiCorr
187     const volScalarField& rA,
188     const volScalarField& rho,
189     const volScalarField& U,
190     const surfaceScalarField& phi
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 } // End namespace fv
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 } // End namespace Foam
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 #ifdef NoRepository
205 #   include "EulerLocalDdtScheme.C"
206 #endif
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 #endif
212 // ************************************************************************* //