Time.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "Time.H"
30 #include "PstreamReduceOps.H"
31 #include "argList.H"
32 #include "HashSet.H"
33 #include "profiling.H"
34 #include "IOdictionary.H"
35 #include "registerSwitch.H"
36 #include <sstream>
37 
38 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(Time, 0);
43 }
44 
45 const Foam::Enum
46 <
48 >
50 ({
51  { stopAtControls::saEndTime, "endTime" },
52  { stopAtControls::saNoWriteNow, "noWriteNow" },
53  { stopAtControls::saWriteNow, "writeNow" },
54  { stopAtControls::saNextWrite, "nextWrite" },
55  // Leave saUnknown untabulated - fallback to flag unknown settings
56 });
57 
58 
59 const Foam::Enum
60 <
62 >
64 ({
65  { writeControls::wcNone, "none" },
66  { writeControls::wcTimeStep, "timeStep" },
67  { writeControls::wcRunTime, "runTime" },
68  { writeControls::wcAdjustableRunTime, "adjustable" },
69  { writeControls::wcAdjustableRunTime, "adjustableRunTime" },
70  { writeControls::wcClockTime, "clockTime" },
71  { writeControls::wcCpuTime, "cpuTime" },
72  // Leave wcUnknown untabulated - fallback to flag unknown settings
73 });
74 
75 
77 (
78  IOstreamOption::floatFormat::general
79 );
80 
82 
83 const int Foam::Time::maxPrecision_(3 - log10(SMALL));
84 
86 
88 (
89  Foam::debug::infoSwitch("printExecutionFormat", 0)
90 );
91 
93 (
94  "printExecutionFormat",
95  int,
97 );
98 
99 
100 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
101 
103 {
104  bool adjustTime = false;
105  scalar timeToNextWrite = VGREAT;
106 
108  {
109  adjustTime = true;
110  timeToNextWrite = max
111  (
112  0.0,
114  );
115  }
116 
117  if (adjustTime)
118  {
119  scalar nSteps = timeToNextWrite/deltaT_;
120 
121  // For tiny deltaT the label can overflow!
122  if (nSteps < scalar(labelMax))
123  {
124  // nSteps can be < 1 so make sure at least 1
125  label nStepsToNextWrite = max(1, round(nSteps));
126 
127  scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
128 
129  // Control the increase of the time step to within a factor of 2
130  // and the decrease within a factor of 5.
131  if (newDeltaT >= deltaT_)
132  {
133  deltaT_ = min(newDeltaT, 2.0*deltaT_);
134  }
135  else
136  {
137  deltaT_ = max(newDeltaT, 0.2*deltaT_);
138  }
139  }
140  }
141 
142  functionObjects_.adjustTimeStep();
143 }
144 
145 
147 {
148  // default is to resume calculation from "latestTime"
149  const word startFrom = controlDict_.getOrDefault<word>
150  (
151  "startFrom",
152  "latestTime"
153  );
154 
155  if (startFrom == "startTime")
156  {
157  controlDict_.readEntry("startTime", startTime_);
158  }
159  else
160  {
161  // Search directory for valid time directories
162  instantList timeDirs = findTimes(path(), constant());
163 
164  const label nTimes = timeDirs.size();
165 
166  if (startFrom == "firstTime")
167  {
168  if (nTimes > 1 && timeDirs.front().name() == constant())
169  {
170  startTime_ = timeDirs[1].value();
171  }
172  else if (nTimes)
173  {
174  startTime_ = timeDirs.front().value();
175  }
176  }
177  else if (startFrom == "latestTime")
178  {
179  if (nTimes)
180  {
181  startTime_ = timeDirs.back().value();
182  }
183  }
184  else
185  {
186  FatalIOErrorInFunction(controlDict_)
187  << "expected startTime, firstTime or latestTime"
188  << " found '" << startFrom << "'"
189  << exit(FatalIOError);
190  }
191  }
192 
193  setTime(startTime_, 0);
194 
195  readDict();
196  deltaTSave_ = deltaT_;
197  deltaT0_ = deltaT_;
198 
199  // Check if time directory exists
200  // If not increase time precision to see if it is formatted differently.
201  // Note: - do not use raw fileHandler().exists(...) since does not check
202  // alternative processorsDDD directories naming
203  // - do not look for compressed (is directory)
204 
205  const fileName timeDir(fileHandler().filePath(timePath(), false));
206  if (returnReduceAnd(timeDir.empty()))
207  {
208  const int oldPrecision = precision_;
209  int requiredPrecision = -1;
210  word oldTime(timeName());
211  for
212  (
213  precision_ = maxPrecision_;
214  precision_ > oldPrecision;
215  --precision_
216  )
217  {
218  // Update the time formatting
219  setTime(startTime_, 0);
220 
221  word newTime(timeName());
222  if (newTime == oldTime)
223  {
224  break;
225  }
226  oldTime = std::move(newTime);
227 
228  // Does a time directory exist with the new format?
229  const fileName timeDir(fileHandler().filePath(timePath(), false));
230  if (returnReduceOr(!timeDir.empty()))
231  {
232  requiredPrecision = precision_;
233  }
234  }
235 
236  if (requiredPrecision > 0)
237  {
238  // Update the time precision
239  precision_ = requiredPrecision;
240 
241  // Update the time formatting
242  setTime(startTime_, 0);
243 
245  << "Increasing the timePrecision from " << oldPrecision
246  << " to " << precision_
247  << " to support the formatting of the current time directory "
248  << timeName() << nl << endl;
249  }
250  else
251  {
252  // Could not find time directory so assume it is not present
253  precision_ = oldPrecision;
254 
255  // Revert the time formatting
256  setTime(startTime_, 0);
257  }
258  }
259 
260  if (Pstream::parRun())
261  {
262  scalar sumStartTime = startTime_;
263  reduce(sumStartTime, sumOp<scalar>());
264  if
265  (
266  mag(Pstream::nProcs()*startTime_ - sumStartTime)
267  > Pstream::nProcs()*deltaT_/10.0
268  )
269  {
270  FatalIOErrorInFunction(controlDict_)
271  << "Start time is not the same for all processors" << nl
272  << "processor " << Pstream::myProcNo() << " has startTime "
273  << startTime_ << exit(FatalIOError);
274  }
275  }
276 
277  IOdictionary timeDict
278  (
279  IOobject
280  (
281  "time",
282  timeName(),
283  "uniform",
284  *this,
288  )
289  );
290 
291  // Read and set the deltaT only if time-step adjustment is active
292  // otherwise use the deltaT from the controlDict
293  if (controlDict_.getOrDefault("adjustTimeStep", false))
294  {
295  if (timeDict.readIfPresent("deltaT", deltaT_))
296  {
297  deltaTSave_ = deltaT_;
298  deltaT0_ = deltaT_;
299  }
300  }
301 
302  timeDict.readIfPresent("deltaT0", deltaT0_);
303 
304  if (timeDict.readIfPresent("index", startTimeIndex_))
305  {
306  timeIndex_ = startTimeIndex_;
307  }
308 
309 
310  // Check if values stored in time dictionary are consistent
311 
312  // 1. Based on time name
313  bool checkValue = true;
314 
315  string storedTimeName;
316  if (timeDict.readIfPresent("name", storedTimeName))
317  {
318  if (storedTimeName == timeName())
319  {
320  // Same time. No need to check stored value
321  checkValue = false;
322  }
323  }
324 
325  // 2. Based on time value
326  // (consistent up to the current time writing precision so it won't
327  // trigger if we just change the write precision)
328  if (checkValue)
329  {
330  scalar storedTimeValue;
331  if (timeDict.readIfPresent("value", storedTimeValue))
332  {
333  word storedTimeName(timeName(storedTimeValue));
334 
335  if (storedTimeName != timeName())
336  {
337  IOWarningInFunction(timeDict)
338  << "Time read from time dictionary " << storedTimeName
339  << " differs from actual time " << timeName() << '.' << nl
340  << " This may cause unexpected database behaviour."
341  << " If you are not interested" << nl
342  << " in preserving time state delete"
343  << " the time dictionary."
344  << endl;
345  }
346  }
347  }
348 }
349 
350 
351 void Foam::Time::setMonitoring(const bool forceProfiling)
352 {
353  const dictionary* profilingDict = controlDict_.findDict("profiling");
354  if (!profilingDict)
355  {
356  // ... or from etc/controlDict
357  profilingDict = debug::controlDict().findDict("profiling");
358  }
359 
360  // initialize profiling on request
361  // otherwise rely on profiling entry within controlDict
362  // and skip if 'active' keyword is explicitly set to false
363  if (forceProfiling)
364  {
366  (
367  IOobject
368  (
369  "profiling",
370  timeName(),
371  "uniform",
372  *this,
376  ),
377  *this
378  );
379  }
380  else if
381  (
382  profilingDict
383  && profilingDict->getOrDefault("active", true)
384  )
385  {
387  (
388  *profilingDict,
389  IOobject
390  (
391  "profiling",
392  timeName(),
393  "uniform",
394  *this,
398  ),
399  *this
400  );
401  }
402 
403  // Time objects not registered so do like objectRegistry::checkIn ourselves.
404  if (runTimeModifiable_)
405  {
406  // Monitor all files that controlDict depends on
407  // Files might have been set during token reading so only on master
408  // processor. Broadcast names to all processors
409  // (even although they are only checked on master)
410  // so that the watched states are synchronised
411 
412  Pstream::broadcast(controlDict_.files(), UPstream::worldComm);
413  fileHandler().addWatches(controlDict_, controlDict_.files());
414  }
415 
416  // Clear dependent files - not needed now
417  controlDict_.files().clear();
418 }
419 
420 
421 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
422 
424 (
425  const word& ctrlDictName,
426  const fileName& rootPath,
427  const fileName& caseName,
428  const word& systemDirName,
429  const word& constantDirName,
430  const bool enableFunctionObjects,
431  const bool enableLibs,
432  IOobjectOption::readOption rOpt // (default: READ_MODIFIED)
433 )
434 :
435  TimePaths(rootPath, caseName, systemDirName, constantDirName),
436  objectRegistry(*this),
437  loopProfiling_(nullptr),
438  libs_(),
439 
440  controlDict_
441  (
442  IOobject
443  (
444  ctrlDictName,
445  system(),
446  *this,
447  rOpt,
448  IOobjectOption::NO_WRITE,
449  IOobjectOption::NO_REGISTER
450  )
451  ),
452 
453  startTimeIndex_(0),
454  startTime_(0),
455  endTime_(0),
456 
457  stopAt_(saEndTime),
458  writeControl_(wcTimeStep),
459  writeInterval_(GREAT),
460  purgeWrite_(0),
461  subCycling_(0),
462  writeOnce_(false),
463  sigWriteNow_(*this, true),
464  sigStopAtWriteNow_(*this, true),
465  writeStreamOption_(IOstreamOption::ASCII),
466  graphFormat_("raw"),
467  runTimeModifiable_(false),
468  cacheTemporaryObjects_(true),
469  functionObjects_(*this, false)
470 {
471  if (enableFunctionObjects)
472  {
473  functionObjects_.on();
474  }
475 
476  if (enableLibs)
477  {
478  libs_.open("libs", controlDict_);
479  }
480 
481  // Explicitly set read flags on objectRegistry so anything constructed
482  // from it reads as well (e.g. fvSolution).
485  setControls();
486  setMonitoring();
487 }
488 
489 
491 (
492  const word& ctrlDictName,
493  const argList& args,
494  const word& systemDirName,
495  const word& constantDirName,
496  const bool enableFunctionObjects,
497  const bool enableLibs,
498  IOobjectOption::readOption rOpt // (default: READ_MODIFIED)
499 )
500 :
501  TimePaths(args, systemDirName, constantDirName),
502  objectRegistry(*this),
503  loopProfiling_(nullptr),
504  libs_(),
505 
506  controlDict_
507  (
508  IOobject
509  (
510  ctrlDictName,
511  system(),
512  *this,
513  rOpt,
514  IOobjectOption::NO_WRITE,
515  IOobjectOption::NO_REGISTER
516  )
517  ),
518 
519  startTimeIndex_(0),
520  startTime_(0),
521  endTime_(0),
522 
523  stopAt_(saEndTime),
524  writeControl_(wcTimeStep),
525  writeInterval_(GREAT),
526  purgeWrite_(0),
527  subCycling_(0),
528  writeOnce_(false),
529  sigWriteNow_(*this, true),
530  sigStopAtWriteNow_(*this, true),
531  writeStreamOption_(IOstreamOption::ASCII),
532  graphFormat_("raw"),
533  runTimeModifiable_(false),
534  cacheTemporaryObjects_(true),
535  functionObjects_(*this, false)
536 {
537  // Functions
538  //
539  // * '-withFunctionObjects' exists and used = enable
540  // * '-noFunctionObjects' exists and used = disable
541  // * default: no functions if there is no way to enable/disable them
542  if (enableFunctionObjects && args.allowFunctionObjects())
543  {
544  functionObjects_.on();
545  }
546 
547  // Libraries
548  //
549  // * enable by default unless '-no-libs' option was used
550  if (enableLibs && args.allowLibs())
551  {
552  libs_.open("libs", controlDict_);
553  }
554 
555  // Explicitly set read flags on objectRegistry so anything constructed
556  // from it reads as well (e.g. fvSolution).
558 
559  setControls();
561  // '-profiling' = force profiling, ignore controlDict entry
562  setMonitoring(args.found("profiling"));
563 }
564 
565 
567 (
568  const dictionary& dict,
569  const fileName& rootPath,
570  const fileName& caseName,
571  const word& systemDirName,
572  const word& constantDirName,
573  const bool enableFunctionObjects,
574  const bool enableLibs,
575  IOobjectOption::readOption rOpt // (default: READ_MODIFIED)
576 )
577 :
578  TimePaths(rootPath, caseName, systemDirName, constantDirName),
579  objectRegistry(*this),
580  loopProfiling_(nullptr),
581  libs_(),
582 
583  controlDict_
584  (
585  // Construct no-read, copying initial contents
586  IOobject
587  (
588  Time::controlDictName,
589  system(),
590  *this,
591  IOobjectOption::NO_READ,
592  IOobjectOption::NO_WRITE,
593  IOobjectOption::NO_REGISTER
594  ),
595  dict
596  ),
597 
598  startTimeIndex_(0),
599  startTime_(0),
600  endTime_(0),
601 
602  stopAt_(saEndTime),
603  writeControl_(wcTimeStep),
604  writeInterval_(GREAT),
605  purgeWrite_(0),
606  subCycling_(0),
607  writeOnce_(false),
608  sigWriteNow_(*this, true),
609  sigStopAtWriteNow_(*this, true),
610  writeStreamOption_(IOstreamOption::ASCII),
611  graphFormat_("raw"),
612  runTimeModifiable_(false),
613  cacheTemporaryObjects_(true),
614  functionObjects_(*this, false)
615 {
616  if (enableFunctionObjects)
617  {
618  functionObjects_.on();
619  }
620 
621  if (enableLibs)
622  {
623  libs_.open("libs", controlDict_);
624  }
625 
626 
627  // Explicitly set read flags on objectRegistry so anything constructed
628  // from it reads as well (e.g. fvSolution).
630 
631  // Since could not construct controlDict with setting:
632  controlDict_.readOpt(rOpt);
634  setControls();
635  setMonitoring();
636 }
637 
638 
640 (
641  const fileName& rootPath,
642  const fileName& caseName,
643  const word& systemDirName,
644  const word& constantDirName,
645  const bool enableFunctionObjects,
646  const bool enableLibs
647 )
648 :
649  TimePaths(rootPath, caseName, systemDirName, constantDirName),
650  objectRegistry(*this),
651  loopProfiling_(nullptr),
652  libs_(),
653 
654  controlDict_
655  (
656  IOobject
657  (
658  Time::controlDictName,
659  system(),
660  *this,
661  IOobjectOption::NO_READ,
662  IOobjectOption::NO_WRITE,
663  IOobjectOption::NO_REGISTER
664  )
665  ),
666 
667  startTimeIndex_(0),
668  startTime_(0),
669  endTime_(0),
670 
671  stopAt_(saEndTime),
672  writeControl_(wcTimeStep),
673  writeInterval_(GREAT),
674  purgeWrite_(0),
675  subCycling_(0),
676  writeOnce_(false),
677  writeStreamOption_(IOstreamOption::ASCII),
678  graphFormat_("raw"),
679  runTimeModifiable_(false),
680  cacheTemporaryObjects_(true),
681  functionObjects_(*this, false)
682 {
683  if (enableFunctionObjects)
684  {
685  functionObjects_.on();
686  }
687 
688  if (enableLibs)
689  {
690  libs_.open("libs", controlDict_);
691  }
693  setMonitoring(); // for profiling etc
694 }
695 
696 
697 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
698 
700 {
701  loopProfiling_.reset(nullptr);
702 
703  forAllReverse(controlDict_.watchIndices(), i)
704  {
705  fileHandler().removeWatch(controlDict_.watchIndices()[i]);
706  }
707 
708  // Destroy function objects first
709  functionObjects_.clear();
710 
711  // Clean up profiling
712  profiling::stop(*this);
713 
714  // Ensure all owned objects are also cleaned up now
716 }
717 
718 
719 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
720 
721 Foam::word Foam::Time::timeName(const scalar t, const int precision)
722 {
723  std::ostringstream buf;
724  buf.setf(std::ios_base::fmtflags(format_), std::ios_base::floatfield);
725  buf.precision(precision);
726  buf << t;
727  return buf.str();
728 }
729 
730 
732 (
733  const fileName& directory,
734  const word& name,
736  const word& stopInstance,
737  const bool constant_fallback
738 ) const
739 {
740  // Note: name can empty (ie, search for directory only)
741  IOobject startIO(name, timeName(), directory, *this, rOpt);
742 
743  IOobject io
744  (
745  fileHandler().findInstance
746  (
747  startIO,
748  timeOutputValue(),
749  stopInstance,
750  constant_fallback
751  )
752  );
753  return io.instance();
754 }
755 
757 Foam::label Foam::Time::startTimeIndex() const
758 {
759  return startTimeIndex_;
760 }
761 
764 {
765  return dimensionedScalar("startTime", dimTime, startTime_);
766 }
767 
770 {
771  return dimensionedScalar("endTime", dimTime, endTime_);
772 }
773 
776 {
777  return stopAt_;
778 }
779 
780 
781 bool Foam::Time::run() const
782 {
783  loopProfiling_.reset(nullptr);
784 
785  bool isRunning = value() < (endTime_ - 0.5*deltaT_);
786 
787  if (!subCycling_)
788  {
789  // Only execute when the condition is no longer true
790  // ie, when exiting the control loop
791  if (!isRunning && timeIndex_ != startTimeIndex_)
792  {
793  // Ensure functionObjects execute on last time step
794  // (and hence write uptodate functionObjectProperties)
795  {
796  addProfiling(fo, "functionObjects.execute()");
797  functionObjects_.execute();
798  }
799  {
800  addProfiling(fo, "functionObjects.end()");
801  functionObjects_.end();
802  }
803 
804  if (cacheTemporaryObjects_)
805  {
806  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
807  }
808  }
809  }
810 
811  if (isRunning)
812  {
813  if (!subCycling_)
814  {
815  const_cast<Time&>(*this).readModifiedObjects();
816 
817  if (timeIndex_ == startTimeIndex_)
818  {
819  addProfiling(functionObjects, "functionObjects.start()");
820  functionObjects_.start();
821  }
822  else
823  {
824  addProfiling(functionObjects, "functionObjects.execute()");
825  functionObjects_.execute();
826  }
827 
828  // Check if the execution of functionObjects require re-reading
829  // any files. This moves effect of e.g. 'timeActivatedFileUpdate'
830  // one time step forward. Note that we cannot call
831  // readModifiedObjects from within timeActivatedFileUpdate since
832  // it might re-read the functionObjects themselves (and delete
833  // the timeActivatedFileUpdate one)
834  if (functionObjects_.filesModified())
835  {
836  const_cast<Time&>(*this).readModifiedObjects();
837  }
838 
839  if (cacheTemporaryObjects_)
840  {
841  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
842  }
843  }
844 
845  // Update the "is-running" status following the
846  // possible side-effects from functionObjects
847  isRunning = value() < (endTime_ - 0.5*deltaT_);
848 
849  // (re)trigger profiling
850  if (profiling::active())
851  {
852  loopProfiling_.reset
853  (
854  new profilingTrigger("time.run() " + objectRegistry::name())
855  );
856  }
857  }
858 
859  return isRunning;
860 }
861 
862 
863 bool Foam::Time::loop()
864 {
865  const bool isRunning = run();
866 
867  if (isRunning)
868  {
869  operator++();
870  }
871 
872  return isRunning;
873 }
874 
876 bool Foam::Time::end() const
877 {
878  return value() > (endTime_ + 0.5*deltaT_);
879 }
880 
881 
882 bool Foam::Time::stopAt(const stopAtControls stopCtrl) const
883 {
884  if (stopCtrl == stopAtControls::saUnknown)
885  {
886  return false;
887  }
888 
889  const bool changed = (stopAt_ != stopCtrl);
890  stopAt_ = stopCtrl;
891  endTime_ = GREAT;
892 
893  // Adjust endTime
894  if (stopCtrl == stopAtControls::saEndTime)
895  {
896  controlDict_.readEntry("endTime", endTime_);
897  }
898 
899  return changed;
900 }
901 
903 bool Foam::Time::isAdjustTimeStep() const
904 {
905  return controlDict_.getOrDefault("adjustTimeStep", false);
906 }
907 
908 
909 void Foam::Time::setTime(const Time& t)
910 {
911  resetTimeState(t.timeName(), t.value(), t.timeIndex());
912  fileHandler().setTime(*this);
913 }
914 
915 
916 void Foam::Time::setTime(const instant& inst, const label newIndex)
917 {
918  resetTimeState(inst.name(), inst.value(), newIndex);
919 
920  IOdictionary timeDict
921  (
922  IOobject
923  (
924  "time",
925  timeName(),
926  "uniform",
927  *this,
931  )
932  );
933 
934  timeDict.readIfPresent("deltaT", deltaT_);
935  timeDict.readIfPresent("deltaT0", deltaT0_);
936  timeDict.readIfPresent("index", timeIndex_);
937  fileHandler().setTime(*this);
938 }
939 
941 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
942 {
943  setTime(newTime.value(), newIndex);
944 }
945 
946 
947 void Foam::Time::setTime(const scalar newTime, const label newIndex)
948 {
949  resetTimeState
950  (
951  timeName(timeToUserTime(newTime)),
952  newTime,
953  newIndex
954  );
955  fileHandler().setTime(*this);
956 }
957 
959 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
960 {
961  setEndTime(endTime.value());
962 }
963 
964 
965 void Foam::Time::setEndTime(const scalar endTime)
966 {
967  endTime_ = endTime;
968 }
969 
970 
972 (
973  const dimensionedScalar& deltaT,
974  const bool adjust
975 )
976 {
977  setDeltaT(deltaT.value(), adjust);
978 }
979 
980 
981 void Foam::Time::setDeltaT(const scalar deltaT, const bool adjust)
982 {
983  deltaT_ = deltaT;
984  deltaTchanged_ = true;
985 
986  if (adjust)
987  {
988  adjustDeltaT();
989  }
990 }
991 
992 
993 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
994 {
995  #ifdef FULLDEBUG
996  if (prevTimeState_)
997  {
999  << "previous time state already set" << nl
1000  << exit(FatalError);
1001  }
1002  #endif
1003 
1004  prevTimeState_.reset(new TimeState(*this));
1005 
1006  setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
1007  deltaT_ /= nSubCycles;
1008  deltaT0_ /= nSubCycles;
1009  deltaTSave_ = deltaT0_;
1011  subCycling_ = nSubCycles;
1012 
1013  return prevTimeState();
1014 }
1015 
1016 
1017 void Foam::Time::subCycleIndex(const label index)
1018 {
1019  // Only permit adjustment if sub-cycling was already active
1020  // and if the index is valid (positive, non-zero).
1021  // This avoids potential mixups for deleting.
1022 
1023  if (subCycling_ && index > 0)
1024  {
1025  subCycling_ = index;
1026  }
1027 }
1028 
1029 
1031 {
1032  if (subCycling_)
1033  {
1034  TimeState::operator=(prevTimeState());
1035  prevTimeState_.reset(nullptr);
1036  }
1038  subCycling_ = 0;
1039 }
1040 
1041 
1042 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1045 {
1046  return operator+=(deltaT.value());
1047 }
1048 
1049 
1050 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
1051 {
1052  setDeltaT(deltaT);
1053  return operator++();
1054 }
1055 
1056 
1058 {
1059  deltaT0_ = deltaTSave_;
1060  deltaTSave_ = deltaT_;
1061 
1062  // Save old time value and name
1063  const scalar oldTimeValue = timeToUserTime(value());
1064  const word oldTimeName = dimensionedScalar::name();
1065 
1066  // Increment time
1067  setTime(value() + deltaT_, timeIndex_ + 1);
1068 
1069  if (!subCycling_)
1070  {
1071  // If the time is very close to zero reset to zero
1072  if (mag(value()) < 10*SMALL*deltaT_)
1073  {
1074  setTime(0.0, timeIndex_);
1075  }
1076 
1077  if (sigStopAtWriteNow_.active() || sigWriteNow_.active())
1078  {
1079  // A signal might have been sent on one processor only
1080  // Reduce so all decide the same.
1081 
1082  label flag = 0;
1083  if (sigStopAtWriteNow_.active() && stopAt_ == saWriteNow)
1084  {
1085  flag += 1;
1086  }
1087  if (sigWriteNow_.active() && writeOnce_)
1088  {
1089  flag += 2;
1090  }
1091  reduce(flag, maxOp<label>());
1092 
1093  if (flag & 1)
1094  {
1095  stopAt_ = saWriteNow;
1096  }
1097  if (flag & 2)
1098  {
1099  writeOnce_ = true;
1100  }
1101  }
1102 
1103  writeTime_ = false;
1104 
1105  switch (writeControl_)
1106  {
1107  case wcNone:
1108  case wcUnknown:
1109  break;
1110 
1111  case wcTimeStep:
1112  writeTime_ = !(timeIndex_ % label(writeInterval_));
1113  break;
1114 
1115  case wcRunTime:
1116  case wcAdjustableRunTime:
1117  {
1118  const label writeIndex = label
1119  (
1120  ((value() - startTime_) + 0.5*deltaT_)
1121  / writeInterval_
1122  );
1123 
1124  if (writeIndex > writeTimeIndex_)
1125  {
1126  writeTime_ = true;
1127  writeTimeIndex_ = writeIndex;
1128  }
1129  }
1130  break;
1131 
1132  case wcCpuTime:
1133  {
1134  const label writeIndex = label
1135  (
1136  returnReduce(elapsedCpuTime(), maxOp<double>())
1137  / writeInterval_
1138  );
1139  if (writeIndex > writeTimeIndex_)
1140  {
1141  writeTime_ = true;
1142  writeTimeIndex_ = writeIndex;
1143  }
1144  }
1145  break;
1146 
1147  case wcClockTime:
1148  {
1149  const label writeIndex = label
1150  (
1151  returnReduce(elapsedClockTime(), maxOp<double>())
1152  / writeInterval_
1153  );
1154  if (writeIndex > writeTimeIndex_)
1155  {
1156  writeTime_ = true;
1157  writeTimeIndex_ = writeIndex;
1158  }
1159  }
1160  break;
1161  }
1162 
1163 
1164  // Check if endTime needs adjustment to stop at the next run()/end()
1165  if (!end())
1166  {
1167  if (stopAt_ == saNoWriteNow)
1168  {
1169  endTime_ = value();
1170  }
1171  else if (stopAt_ == saWriteNow)
1172  {
1173  endTime_ = value();
1174  writeTime_ = true;
1175  }
1176  else if (stopAt_ == saNextWrite && writeTime_ == true)
1177  {
1178  endTime_ = value();
1179  }
1180  }
1181 
1182  // Override writeTime if one-shot writing
1183  if (writeOnce_)
1184  {
1185  writeTime_ = true;
1186  writeOnce_ = false;
1187  }
1188 
1189  // Adjust the precision of the time directory name if necessary
1190  if (writeTime_)
1191  {
1192  // User-time equivalent of deltaT
1193  const scalar userDeltaT =
1194  timeToUserTime(value()) - timeToUserTime(value() - deltaT_);
1195 
1196  // Tolerance used when testing time equivalence
1197  const scalar timeTol =
1198  max(min(pow(scalar(10), -precision_), 0.1*userDeltaT), SMALL);
1199 
1200  // Time value obtained by reading timeName
1201  scalar timeNameValue = -VGREAT;
1202 
1203  // Check that new time representation differs from old one
1204  // reinterpretation of the word
1205  if
1206  (
1207  readScalar(dimensionedScalar::name(), timeNameValue)
1208  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1209  )
1210  {
1211  int oldPrecision = precision_;
1212  while
1213  (
1214  precision_ < maxPrecision_
1215  && readScalar(dimensionedScalar::name(), timeNameValue)
1216  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1217  )
1218  {
1219  precision_++;
1220  setTime(value(), timeIndex());
1221  }
1222 
1223  if (precision_ != oldPrecision)
1224  {
1226  << "Increased the timePrecision from " << oldPrecision
1227  << " to " << precision_
1228  << " to distinguish between timeNames at time "
1230  << endl;
1231 
1232  if (precision_ == maxPrecision_)
1233  {
1234  // Reached maxPrecision limit
1236  << "Current time name " << dimensionedScalar::name()
1237  << nl
1238  << " The maximum time precision has been reached"
1239  " which might result in overwriting previous"
1240  " results."
1241  << endl;
1242  }
1243 
1244  // Check if round-off error caused time-reversal
1245  scalar oldTimeNameValue = -VGREAT;
1246  if
1247  (
1248  readScalar(oldTimeName, oldTimeNameValue)
1249  && (
1250  sign(timeNameValue - oldTimeNameValue)
1251  != sign(deltaT_)
1252  )
1253  )
1254  {
1256  << "Current time name " << dimensionedScalar::name()
1257  << " is set to an instance prior to the "
1258  "previous one "
1259  << oldTimeName << nl
1260  << " This might result in temporal "
1261  "discontinuities."
1262  << endl;
1263  }
1264  }
1265  }
1266  }
1267  }
1268 
1269  return *this;
1270 }
1271 
1272 
1274 {
1275  return operator++();
1276 }
1277 
1278 
1279 // ************************************************************************* //
virtual Time & operator++()
Prefix increment,.
Definition: Time.C:1050
dimensionedScalar sign(const dimensionedScalar &ds)
const scalar & value() const noexcept
Return const reference to value.
dictionary dict
registerInfoSwitch("printExecutionFormat", int, Foam::Time::printExecutionFormat_)
Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Definition: createK.H:12
Address the time paths without using the Time class.
Definition: TimePaths.H:51
scalar deltaT_
Definition: TimeState.H:67
A class for handling file names.
Definition: fileName.H:72
readOption readOpt() const noexcept
Get the read option.
Inter-processor communication reduction functions.
static int printExecutionFormat_
Style for "ExecutionTime = " output.
Definition: Time.H:125
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void adjustDeltaT()
Adjust the time step so that writing occurs at the specified time.
Definition: Time.C:95
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
virtual ~Time()
Destructor.
Definition: Time.C:692
virtual bool isAdjustTimeStep() const
Return true if adjustTimeStep is true.
Definition: Time.C:896
scalar writeInterval_
Definition: Time.H:164
virtual stopAtControls stopAt() const
Return the stop control information.
Definition: Time.C:768
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual void setEndTime(const dimensionedScalar &endTime)
Reset end time.
Definition: Time.C:952
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
floatFormat
Float formats (eg, time directory name formats)
writeControls
Write control options.
Definition: Time.H:82
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:856
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:756
static bool active() noexcept
True if profiling is allowed and is active.
Definition: profiling.C:107
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
void on()
Switch the function objects on.
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
static const Enum< stopAtControls > stopAtControlNames
Names for stopAtControls.
Definition: Time.H:114
A simple container for options an IOstream can normally have.
Time(const word &ctrlDictName, const argList &args, const bool enableFunctionObjects=true, const bool enableLibs=true, IOobjectOption::readOption rOpt=IOobjectOption::READ_MODIFIED)
Construct from argument list, reading from specified control dictionary name.
Definition: TimeI.H:24
bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:762
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
scalar startTime_
Definition: Time.H:156
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:421
writeControls writeControl_
Definition: Time.H:162
void setMonitoring(const bool forceProfiling=false)
Set file monitoring, profiling, etc.
Definition: Time.C:344
const word & timeName() const noexcept
Return the current time name.
Definition: TimeStateI.H:30
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
int infoSwitch(const char *name, const int deflt=0)
Lookup info switch or add default value.
Definition: debug.C:228
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
label writeTimeIndex_
Definition: TimeState.H:65
bool allowFunctionObjects() const
The controlDict &#39;functions&#39; entry is allowed to be used.
Definition: argList.C:2159
virtual bool run() const
Return true if run should continue,.
Definition: Time.C:774
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
static const int maxPrecision_
Maximum time directory name precision.
Definition: Time.H:207
scalar value() const noexcept
The value (const access)
Definition: Instant.H:134
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
word timeName
Definition: getTimeIndex.H:3
virtual TimeState subCycle(const label nSubCycles)
Set time to sub-cycle for the given number of steps.
Definition: Time.C:986
"adjustable" / "adjustableRunTime"
Definition: Time.H:87
bool allowLibs() const
The controlDict &#39;libs&#39; entry is allowed to be used. (eg, has not been disabled by the -no-libs option...
Definition: argList.C:2177
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1077
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const T & name() const noexcept
The name/key (const access)
Definition: Instant.H:144
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
word findInstance(const fileName &directory, const word &name=word::null, IOobjectOption::readOption rOpt=IOobjectOption::MUST_READ, const word &stopInstance=word::null, const bool constant_fallback=true) const
Return time instance (location) of directory containing the file name (eg, used in reading mesh data)...
Definition: Time.C:725
static int precision_
Time directory name precision.
Definition: Time.H:202
A class for handling words, derived from Foam::string.
Definition: word.H:63
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:118
static void initialize(const IOobject &ioObj, const Time &owner)
Singleton to initialize profiling pool, everything enabled.
Definition: profiling.C:142
static const Enum< writeControls > writeControlNames
Names for writeControls.
Definition: Time.H:109
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static void stop(const Time &owner)
Stop profiling, cleanup pool if possible.
Definition: profiling.C:168
virtual void subCycleIndex(const label index)
Adjust the reported sub-cycle index.
Definition: Time.C:1010
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
static IOstreamOption::floatFormat format_
Format for time directory names (general | fixed | scientific)
Definition: Time.H:197
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
const word & name() const noexcept
Return const reference to name.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
An instant of time. Contains the time value and name. Uses Foam::Time when formatting the name...
Definition: instant.H:53
virtual bool end() const
Return true if end of run,.
Definition: Time.C:869
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual Time & operator+=(const dimensionedScalar &deltaT)
Set deltaT to that specified and increment time via operator++()
Definition: Time.C:1037
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
virtual void setDeltaT(const dimensionedScalar &deltaT, const bool adjust=true)
Reset time step, normally also calling adjustDeltaT()
Definition: Time.C:965
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void clear()
Clear all entries from the registry.
constexpr label labelMax
Definition: label.H:55
Reading is optional [identical to READ_IF_PRESENT].
virtual void endSubCycle()
Reset time after sub-cycling back to previous TimeState.
Definition: Time.C:1023
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:437
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition: POSIX.C:1702
stopAtControls
Stop-run control options, which are primarily used when altering the stopAt condition.
Definition: Time.H:97
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
List< instant > instantList
List of instants.
Definition: instantList.H:41
Registry of regIOobjects.
bool open(bool verbose=true)
Open named, but unopened libraries. These names will normally have been added with push_back()...
Foam::argList args(argc, argv)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
dictionary & controlDict()
The central control dictionary, the contents of which are either taken directly from the FOAM_CONTROL...
Definition: debug.C:142
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:750
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Request registration (bool: true)
dimensionedScalar log10(const dimensionedScalar &ds)
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
void setControls()
Set the controls from the current controlDict.
Definition: Time.C:139
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
readOption
Enumeration defining read preferences.