localEulerDdt.C
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) 2015-2022 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 
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
38 {
39  return
40  word(mesh.schemes().ddt("default"))
42 }
43 
44 
46 (
47  const fvMesh& mesh
48 )
49 {
50  return mesh.objectRegistry::lookupObject<volScalarField>
51  (
52  mesh.time().subCycling() ? rSubDeltaTName : rDeltaTName
53  );
54 }
55 
56 
58 (
59  const fvMesh& mesh
60 )
61 {
62  return mesh.objectRegistry::lookupObject<surfaceScalarField>
63  (
64  rDeltaTfName
65  );
66 }
67 
68 
70 (
71  const fvMesh& mesh,
73 )
74 {
75  return tmp<volScalarField>
76  (
77  new volScalarField
78  (
79  rSubDeltaTName,
81  *mesh.objectRegistry::lookupObject<volScalarField>
82  (
83  rDeltaTName
84  )
85  )
86  );
87 }
88 
89 
90 // ************************************************************************* //
const label nAlphaSubCycles(alphaControls.lookup< label >("nAlphaSubCycles"))
Generic GeometricField class.
bool subCycling() const
Return true if time currently being sub-cycled, otherwise false.
Definition: Time.H:438
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1751
ITstream & ddt(const word &name) const
Definition: fvSchemes.C:456
Local time-step first-order Euler implicit/explicit ddt.
static word rSubDeltaTName
Name of the reciprocal local sub-cycling time-step field.
static const surfaceScalarField & localRDeltaTf(const fvMesh &mesh)
Return the reciprocal of the local face time-step.
Definition: localEulerDdt.C:58
static word rDeltaTfName
Name of the reciprocal local face time-step field.
static tmp< volScalarField > localRSubDeltaT(const fvMesh &mesh, const label nAlphaSubCycles)
Calculate and return the reciprocal of the local sub-cycling.
Definition: localEulerDdt.C:70
static word rDeltaTName
Name of the reciprocal local time-step field.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
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