38namespace functionObjects
47 { controlMode::TIME,
"time" },
48 { controlMode::TRIGGER,
"trigger" },
49 { controlMode::TIME_OR_TRIGGER,
"timeOrTrigger" },
50 { controlMode::TIME_AND_TRIGGER,
"timeAndTrigger" }
56void Foam::functionObjects::timeControl::readControls()
79 nStepsToStartTimeChange_ = 3;
82 "nStepsToStartTimeChange",
83 nStepsToStartTimeChange_
89bool Foam::functionObjects::timeControl::active()
const
91 label triggeri = time_.functionObjects().triggerIndex();
94 time_.value() >= (timeStart_ - 0.5*time_.deltaTValue())
95 && time_.value() <= (timeEnd_ + 0.5*time_.deltaTValue());
97 bool inTrigger = triggeri >= triggerStart_ && triggeri <= triggerEnd_;
100 <<
name() <<
" mode:" << controlModeNames_[controlMode_] <<
nl
101 <<
" - time:" << time_.value()
102 <<
" timeStart:" << timeStart_
103 <<
" timeEnd:" << timeEnd_
104 <<
" inTime:" << inTime <<
nl
105 <<
" - triggeri:" << triggeri
106 <<
" triggerStart:" << triggerStart_
107 <<
" triggerEnd:" << triggerEnd_
108 <<
" inTrigger:" << inTrigger
111 switch (controlMode_)
113 case controlMode::TIME:
117 case controlMode::TRIGGER:
121 case controlMode::TIME_OR_TRIGGER:
123 return inTime || inTrigger;
125 case controlMode::TIME_AND_TRIGGER:
127 return inTime && inTrigger;
132 <<
"Unhandled enumeration: " << controlModeNames_[controlMode_]
141Foam::scalar Foam::functionObjects::timeControl::calcExpansion
143 const scalar startRatio,
148 scalar ratio = startRatio;
167 for (label iter = 0; iter < 100; iter++)
170 scalar
f = (
y-1)*
pow(ratio,
n)+1-
y*
pow(ratio,
n-1);
171 scalar dfdratio = (
y-1)*
n*
pow(ratio,
n-1)-
y*(
n-1)*
pow(ratio,
n-2);
172 scalar dratio =
f/(dfdratio + SMALL);
184 <<
"Did not converge to find new timestep growth factor given "
185 <<
"overall factor " <<
y <<
" and " <<
n <<
" steps to do it in."
186 <<
endl <<
" Returning current best guess " << ratio <<
endl;
192void Foam::functionObjects::timeControl::calcDeltaTCoeff
194 scalar& requiredDeltaTCoeff,
195 const scalar wantedDT,
197 const scalar presentTime,
198 const scalar timeToNextWrite,
199 const bool rampDirectionUp
202 const scalar writeInterval = writeControl_.interval();
207 scalar newDeltaT = deltaT0_;
209 if (seriesDTCoeff_ != GREAT)
211 newDeltaT *= seriesDTCoeff_;
219 scalar requiredTimeInterval = newDeltaT;
224 if (requiredDeltaTCoeff != 1.0)
226 requiredTimeInterval *=
227 (
pow(requiredDeltaTCoeff, nSteps) - 1)
228 /(requiredDeltaTCoeff - 1);
233 scalar timeToNextMultiple = -presentTime;
237 timeToNextMultiple +=
240 (presentTime + requiredTimeInterval)
247 timeToNextMultiple +=
250 (presentTime - requiredTimeInterval)
257 if (timeToNextWrite <= timeToNextMultiple)
262 newDeltaT = deltaT0_;
263 if (timeToNextWrite < newDeltaT)
265 newDeltaT = timeToNextWrite;
269 if (requiredTimeInterval > writeInterval)
272 <<
"With given ratio needed time span "
273 << requiredTimeInterval
274 <<
" exceeds available writeInterval "
275 << writeInterval <<
nl
276 <<
"Disabling all future time step ramping"
278 deltaTCoeff_ = GREAT;
279 newDeltaT = wantedDT;
284 seriesDTCoeff_ = GREAT;
287 Info<<
"Disabling ramping until time "
288 << presentTime+timeToNextWrite <<
endl;
291 requiredDeltaTCoeff = newDeltaT/deltaT0_;
300 scalar
y = timeToNextMultiple/wantedDT;
301 label requiredSteps = nSteps;
303 scalar ratioEstimate = deltaTCoeff_;
304 scalar ratioMax = deltaTCoeff_;
306 if (seriesDTCoeff_ != GREAT)
308 ratioEstimate = seriesDTCoeff_;
311 if (!rampDirectionUp)
313 ratioEstimate = 1/ratioEstimate;
314 ratioMax = 1/ratioMax;
318 bool searchConverged =
false;
319 for (label iter = 0; iter < 100; iter++)
322 const scalar newRatio = calcExpansion
328 const scalar deltaTLoop =
330 /
pow(newRatio,
mag(requiredSteps)-1);
331 scalar firstDeltaRatio = deltaTLoop/deltaT0_;
335 *(
pow(newRatio,
mag(requiredSteps))-1)
340 Info<<
" nSteps " << requiredSteps
341 <<
" ratioEstimate " << ratioEstimate
342 <<
" a0 " << deltaTLoop
343 <<
" firstDeltaRatio " << firstDeltaRatio
344 <<
" Sn " <<
Sn <<
" Sn+Time " <<
Sn+presentTime
345 <<
" seriesRatio "<< newRatio <<
endl;
352 (firstDeltaRatio < 1.0 && rampDirectionUp)
353 || (firstDeltaRatio > 1.0 && !rampDirectionUp)
358 requiredSteps =
mag(nSteps);
361 Info<<
"firstDeltaRatio " << firstDeltaRatio <<
" rampDir"
362 << rampDirectionUp <<
" newRatio " << newRatio
363 <<
" y " <<
y <<
" steps " << requiredSteps <<
endl;
368 if (firstDeltaRatio > ratioMax && rampDirectionUp)
373 Info<<
"First ratio " << firstDeltaRatio
374 <<
" exceeds threshold " << ratioMax <<
nl
375 <<
"Increasing required steps " << requiredSteps
379 else if (firstDeltaRatio < ratioMax && !rampDirectionUp)
384 Info<<
"First ratio " << firstDeltaRatio
385 <<
" exceeds threshold " << ratioMax <<
nl
386 <<
"Decreasing required steps " << requiredSteps
392 (newRatio > ratioMax && rampDirectionUp)
393 || (newRatio < ratioMax && !rampDirectionUp)
397 requiredSteps = nSteps;
400 Info<<
"Series ratio " << newRatio <<
"exceeds threshold "
401 << ratioMax <<
nl <<
"Consider next deltaT multiple "
402 <<
y*wantedDT + presentTime <<
endl;
407 seriesDTCoeff_ = newRatio;
409 if (requiredSteps == 1)
411 Sn = deltaT0_*firstDeltaRatio;
412 seriesDTCoeff_ = GREAT;
417 remainder(
Sn + presentTime, wantedDT) / wantedDT;
419 if (
mag(jumpError) > ROOTSMALL)
421 requiredSteps = label(timeToNextWrite/wantedDT);
422 firstDeltaRatio = timeToNextWrite/requiredSteps/deltaT0_;
426 Info<<
"All conditions satisfied deltaT0_:" << deltaT0_
427 <<
" calculated deltaTCoeff:" << firstDeltaRatio
428 <<
" series ratio set to:" << seriesDTCoeff_ <<
endl;
430 searchConverged =
true;
431 requiredDeltaTCoeff = firstDeltaRatio;
436 if (!searchConverged)
442 <<
"Did not converge to find new timestep growth factor"
443 <<
" given overall factor " <<
y
444 <<
" and " << requiredSteps <<
" steps to do it in."
445 <<
nl <<
" Falling back to non-adjusted deltatT "
448 requiredDeltaTCoeff = 1.0;
474 executeTimeIndex_(-1),
476 seriesDTCoeff_(GREAT)
491 ||
dict.found(
"timeStart")
492 ||
dict.found(
"timeEnd")
493 ||
dict.found(
"triggerStart")
494 ||
dict.found(
"triggerEnd")
507 deltaT0_ = time_.deltaTValue();
509 if (active() && (postProcess || executeControl_.execute()))
511 executeTimeIndex_ = time_.timeIndex();
524 foPtr_->execute(subIndex);
533 if (active() && (postProcess || writeControl_.execute()))
536 if (executeTimeIndex_ != time_.timeIndex())
538 executeTimeIndex_ = time_.timeIndex();
551 if (active() && (executeControl_.execute() || writeControl_.execute()))
566 writeControl_.control()
570 const scalar presentTime = time_.value();
575 foPtr_->adjustTimeStep();
576 const scalar wantedDT = time_.deltaTValue();
581 const label writeTimeIndex = writeControl_.executionIndex();
582 const scalar writeInterval = writeControl_.interval();
585 scalar timeToNextWrite =
586 (writeTimeIndex+1)*writeInterval
587 - (presentTime - time_.startTime().value());
588 if (timeToNextWrite <= 0.0)
590 timeToNextWrite = writeInterval;
598 scalar intervalError =
599 remainder(writeInterval, wantedDT)/writeInterval;
602 mag(intervalError) > ROOTSMALL
603 || deltaTCoeff_ == GREAT
606 scalar deltaT = time_.deltaTValue();
607 scalar nSteps = timeToNextWrite/deltaT;
614 deltaTCoeff_ != GREAT
615 || nSteps < nStepsToStartTimeChange_
621 label nStepsToNextWrite =
max(1, round(nSteps));
622 scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
625 scalar clipThreshold = 2;
626 if (deltaTCoeff_ != GREAT)
628 clipThreshold = deltaTCoeff_;
631 if (newDeltaT >= deltaT)
633 deltaT =
min(newDeltaT, clipThreshold*deltaT);
637 clipThreshold = 1/clipThreshold;
638 deltaT =
max(newDeltaT, clipThreshold*deltaT);
641 const_cast<Time&
>(time_).setDeltaT(deltaT,
false);
647 scalar requiredDeltaTCoeff = deltaTCoeff_;
653 if (seriesDTCoeff_ != GREAT)
655 requiredDeltaTCoeff = seriesDTCoeff_;
660 requiredDeltaTCoeff = deltaTCoeff_;
671 if (wantedDT < deltaT0_)
699 if (wantedDT < deltaT0_)
701 requiredDeltaTCoeff = 1/requiredDeltaTCoeff;
714 if (timeToNextWrite > wantedDT)
716 requiredDeltaTCoeff = wantedDT/deltaT0_;
720 requiredDeltaTCoeff = timeToNextWrite/deltaT0_;
724 if (requiredDeltaTCoeff > deltaTCoeff_ && debug)
727 <<
"Required deltaTCoeff "<< requiredDeltaTCoeff
728 <<
" is larger than allowed value "<< deltaTCoeff_
746 scalar newDeltaT = deltaT0_*requiredDeltaTCoeff;
750 const_cast<Time&
>(time_).setDeltaT(newDeltaT,
false);
752 if (seriesDTCoeff_ < 1.0)
754 requiredDeltaTCoeff = 1/requiredDeltaTCoeff;
755 seriesDTCoeff_ = 1/seriesDTCoeff_;
761 foPtr_->adjustTimeStep();
764 if (deltaTCoeff_ != GREAT)
768 scalar requiredDeltaTCoeff =
773 min(deltaTCoeff_, wantedDT/deltaT0_)
777 scalar newDeltaT = deltaT0_*requiredDeltaTCoeff;
781 const_cast<Time&
>( time_).setDeltaT(newDeltaT,
false);
786 const_cast<Time&
>( time_).setDeltaT(wantedDT,
false);
801 writeControl_.read(
dict);
802 executeControl_.read(
dict);
806 return foPtr_->read(
dict);
815 return active() ? foPtr_->filesModified() :
false;
823 foPtr_->updateMesh(mpm);
832 foPtr_->movePoints(
mesh);
Macros for easy insertion into run-time selection tables.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
EnumType getOrDefault(const word &key, const dictionary &dict, const EnumType deflt, const bool failsafe=false) const
virtual bool read()
Re-read model coefficients if they have changed.
virtual scalar userTimeToTime(const scalar theta) const
Convert the user-time (e.g. CA deg) to real-time (s).
scalar deltaTValue() const noexcept
Return time step value.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base-class for Time/database function objects.
Computes the magnitude of an input field.
Wrapper around functionObjects to add time control.
static const Enum< controlMode > controlModeNames_
virtual bool filesModified() const
Did any file get changed during execution?
static bool entriesPresent(const dictionary &dict)
Helper function to identify if a timeControl object is present.
virtual bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
virtual bool execute()
Called at each ++ or += of the time-loop.
virtual bool write()
Called at each ++ or += of the time-loop.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
Virtual base class for function objects with a reference to Time.
const Time & time_
Reference to the time database.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void movePoints()
Update for new mesh geometry.
void updateMesh()
Update for new mesh topology.
Mesh consisting of general polyhedral cells.
General time dependent execution controller. The execution parameters are given by the "Control" and ...
static bool entriesPresent(const dictionary &dict, const word &prefix)
Identify if a timeControl object is present in the dictionary.
@ ocAdjustableRunTime
Currently identical to "runTime".
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
errorManip< error > abort(error &err)
static scalar Sn(const scalar a, const scalar x)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
constexpr char nl
The newline '\n' character (0x0a)