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 OpenCFD Ltd.
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
30
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace functionObjects
37{
40 (
44 );
45}
46}
47
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51Foam::functionObjects::
52setTimeStepFaRegionsFunctionObject::
53setTimeStepFaRegionsFunctionObject
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().lookupOrDefault<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
114Foam::scalar Foam::functionObjects::setTimeStepFaRegionsFunctionObject::
115regionDeltaT() const
116{
117 const wordList names(time_.sortedNames<regionFaModel>());
118
119 scalar Co = 0.0;
120
121 forAll (names, i)
122 {
123 const auto* regionFa = time_.cfindObject<regionFaModel>(names[i]);
124
125 if (regionFa)
126 {
127 const scalar regionCo = regionFa->CourantNumber();
128 if (regionCo > Co)
129 {
130 Co = regionCo;
131 }
132 }
133 }
134
135 if (names.size() > 0)
136 {
137 const scalar regionFaMaxCo =
138 time_.controlDict().get<scalar>("regionFaMaxCo");
139
140 const scalar maxDeltaTFact = regionFaMaxCo/(Co + SMALL);
141 const scalar deltaTFact =
142 min(min(maxDeltaTFact, 1.0 + 0.1*maxDeltaTFact), 1.2);
143
144 return deltaTFact*time_.deltaTValue();
145 }
146
147 return time_.deltaTValue();
148}
149
150
152{
153 return true;
154}
155
156
158{
159 return true;
160}
161
162
163// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
This function object controls the time step for classes of the type regionFaModel....
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
virtual bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
Virtual base class for function objects with a reference to Time.
Base class for area region models.
virtual scalar CourantNumber() const
Courant number of the region.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Namespace for OpenFOAM.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333