setTimeStepFaRegionsFunctionObject.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) 2020-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 
30 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace functionObjects
37 {
38  defineTypeNameAndDebug(setTimeStepFaRegionsFunctionObject, 0);
40  (
41  functionObject,
42  setTimeStepFaRegionsFunctionObject,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::functionObjects::
52 setTimeStepFaRegionsFunctionObject::
53 setTimeStepFaRegionsFunctionObject
54 (
55  const word& name,
56  const Time& runTime,
57  const dictionary& dict
58 )
59 :
61 {
62  read(dict);
63 }
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
70  // Wanted timestep
71  scalar newDeltaT = regionDeltaT();
72 
73  static label index = -1;
74 
75  if ((time_.timeIndex() != index) && (newDeltaT < time_.deltaTValue()))
76  {
77  // Store current time so we don't get infinite recursion (since
78  // setDeltaT calls adjustTimeStep() again)
79  index = time_.timeIndex();
80 
81  // Set time, allow deltaT to be adjusted for writeInterval purposes
82  const_cast<Time&>(time_).setDeltaT(newDeltaT, false);
83 
84  return true;
85  }
86 
87  return false;
88 }
89 
90 
92 (
93  const dictionary& dict
94 )
95 {
97  {
98  // Ensure that adjustTimeStep is active
99  if (!time_.controlDict().getOrDefault<bool>("adjustTimeStep", false))
100  {
102  << "Need to set 'adjustTimeStep' true to allow timestep control"
103  << nl
104  << exit(FatalIOError);
105  }
106 
107  return true;
108  }
109 
110  return false;
111 }
112 
113 
114 Foam::scalar Foam::functionObjects::setTimeStepFaRegionsFunctionObject::
115 regionDeltaT() const
116 {
117  bool adjust = false;
118  scalar Co = 0;
119 
120  for (const regionFaModel& reg : time_.cobjects<regionFaModel>())
121  {
122  const scalar regionCo = reg.CourantNumber();
123  if (Co < regionCo)
124  {
125  Co = regionCo;
126  adjust = true;
127  }
128  }
129 
130  if (adjust)
131  {
132  const scalar regionFaMaxCo =
133  time_.controlDict().get<scalar>("regionFaMaxCo");
134 
135  const scalar maxDeltaTFact = regionFaMaxCo/(Co + SMALL);
136  const scalar deltaTFact =
137  min(min(maxDeltaTFact, 1.0 + 0.1*maxDeltaTFact), 1.2);
138 
139  return deltaTFact*time_.deltaTValue();
140  }
141 
142  return time_.deltaTValue();
143 }
144 
147 {
148  return true;
149 }
150 
151 
153 {
154  return true;
155 }
156 
157 
158 // ************************************************************************* //
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
engineTime & runTime
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Base class for area region models.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Virtual base class for function objects with a reference to Time.
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
const Time & time_
Reference to the time database.