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) 2011-2016 OpenFOAM Foundation
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 
30  The reciprocal of the local time-step field is looked-up from the database.
31 
32  This scheme should only be used for steady-state computations using
33  transient codes where local time-stepping is preferably to under-relaxation
34  for transport consistency reasons.
35 
36 See also
37  Foam::fv::CoEulerDdtScheme
38 
39 SourceFiles
40  localEulerDdt.C
41  localEulerDdtScheme.C
42  localEulerDdtSchemes.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef localEulerDdtScheme_H
47 #define localEulerDdtScheme_H
48 
49 #include "ddtScheme.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace fv
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class localEulerDdt Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 class localEulerDdt
66 {
67 public:
68 
69  //- Name of the reciprocal local time-step field
70  static word rDeltaTName;
71 
72  //- Name of the reciprocal local face time-step field
73  static word rDeltaTfName;
74 
75  //- Name of the reciprocal local sub-cycling time-step field
76  static word rSubDeltaTName;
77 
78 
79  // Constructors
80 
82  {}
83 
84 
85  // Member Functions
86 
87  //- Return true if LTS is enabled
88  static bool enabled(const fvMesh& mesh);
89 
90  //- Return the reciprocal of the local time-step
91  // looked-up from the objectRegistry
92  static const volScalarField& localRDeltaT(const fvMesh& mesh);
93 
94  //- Return the reciprocal of the local face time-step
95  // looked-up from the objectRegistry
96  static const surfaceScalarField& localRDeltaTf(const fvMesh& mesh);
97 
98  //- Calculate and return the reciprocal of the local sub-cycling
99  // time-step
101  (
102  const fvMesh& mesh,
103  const label nAlphaSubCycles
104  );
105 };
106 
107 
108 /*---------------------------------------------------------------------------*\
109  Class localEulerDdtScheme Declaration
110 \*---------------------------------------------------------------------------*/
111 
112 template<class Type>
114 :
115  public localEulerDdt,
116  public fv::ddtScheme<Type>
117 {
118  // Private Member Functions
119 
120  //- Disallow default bitwise copy construct
122 
123  //- Disallow default bitwise assignment
124  void operator=(const localEulerDdtScheme&);
125 
126  //- Return the reciprocal of the local time-step
127  const volScalarField& localRDeltaT() const;
128 
129  //- Return the reciprocal of the local face time-step
130  const surfaceScalarField& localRDeltaTf() const;
131 
132 
133 public:
134 
135  //- Runtime type information
136  TypeName("localEuler");
137 
138 
139  // Constructors
140 
141  //- Construct from mesh
143  :
144  ddtScheme<Type>(mesh)
145  {}
146 
147  //- Construct from mesh and Istream
149  :
150  ddtScheme<Type>(mesh, is)
151  {}
152 
153 
154  // Member Functions
155 
156  //- Return mesh reference
157  const fvMesh& mesh() const
158  {
159  return fv::ddtScheme<Type>::mesh();
160  }
161 
163  (
164  const dimensioned<Type>&
165  );
166 
168  (
170  );
171 
173  (
174  const dimensionedScalar&,
176  );
177 
179  (
180  const volScalarField&,
182  );
183 
185  (
186  const volScalarField& alpha,
187  const volScalarField& rho,
189  );
190 
192  (
194  );
195 
196  tmp<fvMatrix<Type>> fvmDdt
197  (
199  );
200 
201  tmp<fvMatrix<Type>> fvmDdt
202  (
203  const dimensionedScalar&,
205  );
206 
207  tmp<fvMatrix<Type>> fvmDdt
208  (
209  const volScalarField&,
211  );
212 
213  tmp<fvMatrix<Type>> fvmDdt
214  (
215  const volScalarField& alpha,
216  const volScalarField& rho,
218  );
221 
222  tmp<fluxFieldType> fvcDdtUfCorr
223  (
226  );
227 
228  tmp<fluxFieldType> fvcDdtPhiCorr
229  (
231  const fluxFieldType& phi
232  );
233 
234  tmp<fluxFieldType> fvcDdtUfCorr
235  (
236  const volScalarField& rho,
239  );
240 
241  tmp<fluxFieldType> fvcDdtPhiCorr
242  (
243  const volScalarField& rho,
245  const fluxFieldType& phi
246  );
247 
249  (
251  );
252 };
253 
254 
255 template<>
257 (
260 );
261 
262 template<>
264 (
265  const volScalarField& U,
266  const surfaceScalarField& phi
267 );
268 
269 template<>
271 (
272  const volScalarField& rho,
273  const volScalarField& U,
274  const surfaceScalarField& Uf
275 );
276 
277 template<>
279 (
280  const volScalarField& rho,
281  const volScalarField& U,
282  const surfaceScalarField& phi
283 );
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace fv
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 } // End namespace Foam
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 #ifdef NoRepository
297  #include "localEulerDdtScheme.C"
298 #endif
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 #endif
303 
304 // ************************************************************************* //
Local time-step first-order Euler implicit/explicit ddt.
surfaceScalarField & phi
U
Definition: pEqn.H:83
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Uf
Definition: pEqn.H:67
Generic dimensioned Type class.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
const fvMesh & mesh() const
Return mesh reference.
Definition: ddtScheme.H:135
dynamicFvMesh & mesh
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
static word rDeltaTfName
Name of the reciprocal local face time-step field.
const volScalarField & psi
A class for managing temporary objects.
Definition: PtrList.H:53
static word rDeltaTName
Name of the reciprocal local time-step field.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static word rSubDeltaTName
Name of the reciprocal local sub-cycling time-step field.
Namespace for OpenFOAM.
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))