fluidMaxDeltaT.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) 2023-2024 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 "fluidMaxDeltaT.H"
27 #include "fluidSolver.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
37 
39  (
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& name,
53  const Time& runTime,
54  const dictionary& dict
55 )
56 :
57  fvMeshFunctionObject(name, runTime, dict)
58 {
59  read(dict);
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 {
74 
75  maxCoPtr_ =
77  (
78  "maxCo",
79  time_.userUnits(),
80  dimless,
81  dict
82  );
83  maxDeltaTPtr_ =
85  (
86  "maxDeltaT",
87  time_.userUnits(),
88  time_.userUnits(),
89  dict
90  );
91 
92  return true;
93 }
94 
95 
97 {
98  return true;
99 }
100 
101 
103 {
104  return true;
105 }
106 
107 
109 {
110  scalar deltaT = maxDeltaTPtr_().value(time_.value());
111 
112  const scalar CoNum =
113  mesh_.lookupObject<solvers::fluidSolver>(solver::typeName).CoNum;
114 
115  if (CoNum > small)
116  {
117  const scalar maxCo = maxCoPtr_().value(time_.value());
118 
119  deltaT = min(deltaT, maxCo/CoNum*time_.deltaTValue());
120  }
121 
122  return deltaT;
123 }
124 
125 
126 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Returns the maximum time-step evaluated from time-dependent maximum Courant number and maximum time-s...
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
fluidMaxDeltaT(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
virtual bool execute()
Execute, currently does nothing.
virtual bool write()
Write the time step value.
virtual bool read(const dictionary &)
Read the controls.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool read(const dictionary &)
Read optional controls.
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
A class for handling words, derived from string.
Definition: word.H:62
scalar CoNum
scalar maxCo
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict