localEulerDdtScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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  //- Return the reciprocal of the local time-step
121  const volScalarField& localRDeltaT() const;
122 
123  //- Return the reciprocal of the local face time-step
124  const surfaceScalarField& localRDeltaTf() const;
125 
126 
127 public:
128 
129  //- Runtime type information
130  TypeName("localEuler");
131 
132 
133  // Constructors
134 
135  //- Construct from mesh
137  :
138  ddtScheme<Type>(mesh)
139  {}
140 
141  //- Construct from mesh and Istream
143  :
144  ddtScheme<Type>(mesh, is)
145  {}
146 
147  //- Disallow default bitwise copy construction
148  localEulerDdtScheme(const localEulerDdtScheme&) = delete;
149 
150 
151  // Member Functions
152 
153  using ddtScheme<Type>::mesh;
154 
156  (
157  const dimensioned<Type>&
158  );
159 
161  (
163  );
164 
166  (
167  const dimensionedScalar&,
169  );
170 
172  (
173  const volScalarField&,
175  );
176 
178  (
179  const volScalarField& alpha,
180  const volScalarField& rho,
182  );
183 
185  (
187  );
188 
189  virtual tmp<fvMatrix<Type>> fvmDdt
190  (
192  );
193 
194  virtual tmp<fvMatrix<Type>> fvmDdt
195  (
196  const dimensionedScalar&,
198  );
199 
200  virtual tmp<fvMatrix<Type>> fvmDdt
201  (
202  const volScalarField&,
204  );
205 
206  virtual tmp<fvMatrix<Type>> fvmDdt
207  (
208  const volScalarField& alpha,
209  const volScalarField& rho,
211  );
214 
215  // using ddtScheme<Type>::fvcDdtPhiCoeff;
216 
217  // virtual tmp<surfaceScalarField> fvcDdtPhiCoeff
218  // (
219  // const GeometricField<Type, fvPatchField, volMesh>& U,
220  // const fluxFieldType& phi,
221  // const fluxFieldType& phiCorr
222  // );
223 
224  virtual tmp<fluxFieldType> fvcDdtUfCorr
225  (
228  );
229 
230  virtual tmp<fluxFieldType> fvcDdtPhiCorr
231  (
233  const fluxFieldType& phi
234  );
235 
236  virtual tmp<fluxFieldType> fvcDdtUfCorr
237  (
238  const volScalarField& rho,
241  );
242 
243  virtual tmp<fluxFieldType> fvcDdtPhiCorr
244  (
245  const volScalarField& rho,
247  const fluxFieldType& phi
248  );
249 
251  (
253  );
254 
255 
256  // Member Operators
257 
258  //- Disallow default bitwise assignment
259  void operator=(const localEulerDdtScheme&) = delete;
260 };
261 
262 
263 template<>
265 (
268 );
269 
270 template<>
272 (
273  const volScalarField& U,
274  const surfaceScalarField& phi
275 );
276 
277 template<>
279 (
280  const volScalarField& rho,
281  const volScalarField& U,
282  const surfaceScalarField& Uf
283 );
284 
285 template<>
287 (
288  const volScalarField& rho,
289  const volScalarField& U,
290  const surfaceScalarField& phi
291 );
292 
293 
294 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
295 
296 } // End namespace fv
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 } // End namespace Foam
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 #ifdef NoRepository
305  #include "localEulerDdtScheme.C"
306 #endif
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 #endif
311 
312 // ************************************************************************* //
Local time-step first-order Euler implicit/explicit ddt.
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual tmp< fluxFieldType > fvcDdtUfCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Generic dimensioned Type class.
Abstract base class for ddt schemes.
Definition: ddtScheme.H:64
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)
virtual tmp< fluxFieldType > fvcDdtPhiCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const fluxFieldType &phi)
phi
Definition: correctPhi.H:3
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
autoPtr< surfaceVectorField > Uf
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:70
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
U
Definition: pEqn.H:72
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.
static word rSubDeltaTName
Name of the reciprocal local sub-cycling time-step field.
Namespace for OpenFOAM.