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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2015-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "localEulerDdtScheme.H"
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
40{
41 return
42 word(mesh.ddtScheme("default"))
44}
45
46
48(
49 const fvMesh& mesh
50)
51{
52 return mesh.objectRegistry::lookupObject<volScalarField>
53 (
54 mesh.time().subCycling() ? rSubDeltaTName : rDeltaTName
55 );
56}
57
58
60(
61 const fvMesh& mesh
62)
63{
64 return mesh.objectRegistry::lookupObject<surfaceScalarField>
65 (
66 rDeltaTfName
67 );
68}
69
70
72(
73 const fvMesh& mesh,
74 const label nAlphaSubCycles
75)
76{
78 (
80 (
81 rSubDeltaTName,
83 *mesh.objectRegistry::lookupObject<volScalarField>
84 (
85 rDeltaTName
86 )
87 )
88 );
89}
90
91
92// ************************************************************************* //
virtual bool enabled() const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
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:60
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:72
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:48
ITstream & ddtScheme(const word &name) const
Get ddt scheme for given name, or default.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
label nAlphaSubCycles(alphaControls.get< label >("nAlphaSubCycles"))