localEulerDdt.C
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) 2015-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 \*---------------------------------------------------------------------------*/
25 
26 #include "localEulerDdtScheme.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
37 {
38  return
39  word(mesh.ddtScheme("default"))
41 }
42 
43 
45 (
46  const fvMesh& mesh
47 )
48 {
49  return mesh.objectRegistry::lookupObject<volScalarField>
50  (
52  );
53 }
54 
55 
57 (
58  const fvMesh& mesh,
59  const label nAlphaSubCycles
60 )
61 {
62  return tmp<volScalarField>
63  (
64  new volScalarField
65  (
67  nAlphaSubCycles
68  *mesh.objectRegistry::lookupObject<volScalarField>
69  (
71  )
72  )
73  );
74 }
75 
76 
77 // ************************************************************************* //
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 tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:57
bool subCycling() const
Return true if time currently being sub-cycled, otherwise false.
Definition: Time.H:444
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:36
A class for handling words, derived from string.
Definition: word.H:59
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:45
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:54
static word rDeltaTName
Name of the reciprocal local time-step field.
ITstream & ddtScheme(const word &name) const
Definition: fvSchemes.C:360
static word rSubDeltaTName
Name of the reciprocal local sub-cycling time-step field.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243