FreeFOAM The Cross-Platform CFD Toolkit
localEulerDdtScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
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.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fv::localEulerDdtScheme
26 
27 Description
28  Local time-step first-order Euler implicit/explicit ddt.
29  The reciprocal of the local time-step field is looked-up from the
30  database with the name provided.
31 
32  This scheme should only be used for steady-state computations
33  using transient codes where local time-stepping is preferably to
34  under-relaxation for transport consistency reasons.
35 
36  See also CoEulerDdtScheme.
37 
38 SourceFiles
39  localEulerDdtScheme.C
40  localEulerDdtSchemes.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef localEulerDdtScheme_H
45 #define localEulerDdtScheme_H
46 
47 #include <finiteVolume/ddtScheme.H>
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace fv
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  Class localEulerDdtScheme Declaration
61 \*---------------------------------------------------------------------------*/
62 
63 template<class Type>
65 :
66  public fv::ddtScheme<Type>
67 {
68  // Private Data
69 
70  //- Name of the reciprocal local time-step field
71  word rDeltaTName_;
72 
73 
74  // Private Member Functions
75 
76  //- Disallow default bitwise copy construct
78 
79  //- Disallow default bitwise assignment
80  void operator=(const localEulerDdtScheme&);
81 
82  //- Return the reciprocal of the local time-step
83  const volScalarField& localRDeltaT() const;
84 
85 
86 public:
87 
88  //- Runtime type information
89  TypeName("localEuler");
90 
91 
92  // Constructors
93 
94  //- Construct from mesh and Istream
96  :
97  ddtScheme<Type>(mesh, is),
98  rDeltaTName_(is)
99  {}
100 
101 
102  // Member Functions
103 
104  //- Return mesh reference
105  const fvMesh& mesh() const
106  {
107  return fv::ddtScheme<Type>::mesh();
108  }
109 
111  (
112  const dimensioned<Type>&
113  );
114 
116  (
118  );
119 
121  (
122  const dimensionedScalar&,
124  );
125 
127  (
128  const volScalarField&,
130  );
131 
133  (
135  );
136 
138  (
139  const dimensionedScalar&,
141  );
142 
144  (
145  const volScalarField&,
147  );
148 
150 
152  (
153  const volScalarField& rA,
155  const fluxFieldType& phi
156  );
157 
159  (
160  const volScalarField& rA,
161  const volScalarField& rho,
163  const fluxFieldType& phi
164  );
165 
167  (
169  );
170 };
171 
172 
173 template<>
175 (
176  const volScalarField& rA,
177  const volScalarField& U,
178  const surfaceScalarField& phi
179 );
180 
181 
182 template<>
184 (
185  const volScalarField& rA,
186  const volScalarField& rho,
187  const volScalarField& U,
188  const surfaceScalarField& phi
189 );
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace fv
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 } // End namespace Foam
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 #ifdef NoRepository
203 # include "localEulerDdtScheme.C"
204 #endif
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 #endif
209 
210 // ************************************************************************* //