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-2021 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 
79 
80 const int Foam::Time::maxPrecision_(3 - log10(SMALL));
81 
83 
85 (
86  Foam::debug::infoSwitch("printExecutionFormat", 0)
87 );
88 
90 (
91  "printExecutionFormat",
92  int,
94 );
95 
96 
97 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
98 
100 {
101  bool adjustTime = false;
102  scalar timeToNextWrite = VGREAT;
103 
105  {
106  adjustTime = true;
107  timeToNextWrite = max
108  (
109  0.0,
111  );
112  }
113 
114  if (adjustTime)
115  {
116  scalar nSteps = timeToNextWrite/deltaT_;
117 
118  // For tiny deltaT the label can overflow!
119  if (nSteps < scalar(labelMax))
120  {
121  // nSteps can be < 1 so make sure at least 1
122  label nStepsToNextWrite = max(1, round(nSteps));
123 
124  scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
125 
126  // Control the increase of the time step to within a factor of 2
127  // and the decrease within a factor of 5.
128  if (newDeltaT >= deltaT_)
129  {
130  deltaT_ = min(newDeltaT, 2.0*deltaT_);
131  }
132  else
133  {
134  deltaT_ = max(newDeltaT, 0.2*deltaT_);
135  }
136  }
137  }
138 
139  functionObjects_.adjustTimeStep();
140 }
141 
142 
144 {
145  // default is to resume calculation from "latestTime"
146  const word startFrom = controlDict_.getOrDefault<word>
147  (
148  "startFrom",
149  "latestTime"
150  );
151 
152  if (startFrom == "startTime")
153  {
154  controlDict_.readEntry("startTime", startTime_);
155  }
156  else
157  {
158  // Search directory for valid time directories
159  instantList timeDirs = findTimes(path(), constant());
160 
161  const label nTimes = timeDirs.size();
162 
163  if (startFrom == "firstTime")
164  {
165  if (nTimes > 1 && timeDirs.first().name() == constant())
166  {
167  startTime_ = timeDirs[1].value();
168  }
169  else if (nTimes)
170  {
171  startTime_ = timeDirs.first().value();
172  }
173  }
174  else if (startFrom == "latestTime")
175  {
176  if (nTimes)
177  {
178  startTime_ = timeDirs.last().value();
179  }
180  }
181  else
182  {
183  FatalIOErrorInFunction(controlDict_)
184  << "expected startTime, firstTime or latestTime"
185  << " found '" << startFrom << "'"
186  << exit(FatalIOError);
187  }
188  }
189 
190  setTime(startTime_, 0);
191 
192  readDict();
193  deltaTSave_ = deltaT_;
194  deltaT0_ = deltaT_;
195 
196  // Check if time directory exists
197  // If not increase time precision to see if it is formatted differently.
198  // Note: do not use raw fileHandler().exists(...) since does not check
199  // alternative processorsDDD directories naming
200  if (fileHandler().filePath(timePath()).empty())
201  {
202  int oldPrecision = precision_;
203  int requiredPrecision = -1;
204  bool found = false;
205  word oldTime(timeName());
206  for
207  (
208  precision_ = maxPrecision_;
209  precision_ > oldPrecision;
210  --precision_
211  )
212  {
213  // Update the time formatting
214  setTime(startTime_, 0);
215 
216  word newTime(timeName());
217  if (newTime == oldTime)
218  {
219  break;
220  }
221  oldTime = newTime;
222 
223  // Check the existence of the time directory with the new format
224  //found = fileHandler().exists(timePath(), false);
225  const fileName dirName(fileHandler().filePath(timePath()));
226  found = !dirName.empty();
227 
228  if (found)
229  {
230  requiredPrecision = precision_;
231  }
232  }
233 
234  if (requiredPrecision > 0)
235  {
236  // Update the time precision
237  precision_ = requiredPrecision;
238 
239  // Update the time formatting
240  setTime(startTime_, 0);
241 
243  << "Increasing the timePrecision from " << oldPrecision
244  << " to " << precision_
245  << " to support the formatting of the current time directory "
246  << timeName() << nl << endl;
247  }
248  else
249  {
250  // Could not find time directory so assume it is not present
251  precision_ = oldPrecision;
252 
253  // Revert the time formatting
254  setTime(startTime_, 0);
255  }
256  }
257 
258  if (Pstream::parRun())
259  {
260  scalar sumStartTime = startTime_;
261  reduce(sumStartTime, sumOp<scalar>());
262  if
263  (
264  mag(Pstream::nProcs()*startTime_ - sumStartTime)
265  > Pstream::nProcs()*deltaT_/10.0
266  )
267  {
268  FatalIOErrorInFunction(controlDict_)
269  << "Start time is not the same for all processors" << nl
270  << "processor " << Pstream::myProcNo() << " has startTime "
271  << startTime_ << exit(FatalIOError);
272  }
273  }
274 
275  IOdictionary timeDict
276  (
277  IOobject
278  (
279  "time",
280  timeName(),
281  "uniform",
282  *this,
285  false
286  )
287  );
288 
289  // Read and set the deltaT only if time-step adjustment is active
290  // otherwise use the deltaT from the controlDict
291  if (controlDict_.getOrDefault("adjustTimeStep", false))
292  {
293  if (timeDict.readIfPresent("deltaT", deltaT_))
294  {
295  deltaTSave_ = deltaT_;
296  deltaT0_ = deltaT_;
297  }
298  }
299 
300  timeDict.readIfPresent("deltaT0", deltaT0_);
301 
302  if (timeDict.readIfPresent("index", startTimeIndex_))
303  {
304  timeIndex_ = startTimeIndex_;
305  }
306 
307 
308  // Check if values stored in time dictionary are consistent
309 
310  // 1. Based on time name
311  bool checkValue = true;
312 
313  string storedTimeName;
314  if (timeDict.readIfPresent("name", storedTimeName))
315  {
316  if (storedTimeName == timeName())
317  {
318  // Same time. No need to check stored value
319  checkValue = false;
320  }
321  }
322 
323  // 2. Based on time value
324  // (consistent up to the current time writing precision so it won't
325  // trigger if we just change the write precision)
326  if (checkValue)
327  {
328  scalar storedTimeValue;
329  if (timeDict.readIfPresent("value", storedTimeValue))
330  {
331  word storedTimeName(timeName(storedTimeValue));
332 
333  if (storedTimeName != timeName())
334  {
335  IOWarningInFunction(timeDict)
336  << "Time read from time dictionary " << storedTimeName
337  << " differs from actual time " << timeName() << '.' << nl
338  << " This may cause unexpected database behaviour."
339  << " If you are not interested" << nl
340  << " in preserving time state delete"
341  << " the time dictionary."
342  << endl;
343  }
344  }
345  }
346 }
347 
348 
349 void Foam::Time::setMonitoring(const bool forceProfiling)
350 {
351  const dictionary* profilingDict = controlDict_.findDict("profiling");
352  if (!profilingDict)
353  {
354  // ... or from etc/controlDict
355  profilingDict = debug::controlDict().findDict("profiling");
356  }
357 
358  // initialize profiling on request
359  // otherwise rely on profiling entry within controlDict
360  // and skip if 'active' keyword is explicitly set to false
361  if (forceProfiling)
362  {
364  (
365  IOobject
366  (
367  "profiling",
368  timeName(),
369  "uniform",
370  *this,
373  ),
374  *this
375  );
376  }
377  else if
378  (
379  profilingDict
380  && profilingDict->getOrDefault("active", true)
381  )
382  {
384  (
385  *profilingDict,
386  IOobject
387  (
388  "profiling",
389  timeName(),
390  "uniform",
391  *this,
394  ),
395  *this
396  );
397  }
398 
399  // Time objects not registered so do like objectRegistry::checkIn ourselves.
400  if (runTimeModifiable_)
401  {
402  // Monitor all files that controlDict depends on
403  fileHandler().addWatches(controlDict_, controlDict_.files());
404  }
405 
406  // Clear dependent files - not needed now
407  controlDict_.files().clear();
408 }
409 
410 
411 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
412 
414 (
415  const word& ctrlDictName,
416  const fileName& rootPath,
417  const fileName& caseName,
418  const word& systemName,
419  const word& constantName,
420  const bool enableFunctionObjects,
421  const bool enableLibs
422 )
423 :
424  TimePaths
425  (
426  rootPath,
427  caseName,
428  systemName,
429  constantName
430  ),
431 
432  objectRegistry(*this),
433 
434  loopProfiling_(nullptr),
435  libs_(),
436 
437  controlDict_
438  (
439  IOobject
440  (
441  ctrlDictName,
442  system(),
443  *this,
446  false
447  )
448  ),
449 
450  startTimeIndex_(0),
451  startTime_(0),
452  endTime_(0),
453 
454  stopAt_(saEndTime),
455  writeControl_(wcTimeStep),
456  writeInterval_(GREAT),
457  purgeWrite_(0),
458  subCycling_(0),
459  writeOnce_(false),
460  sigWriteNow_(*this, true),
461  sigStopAtWriteNow_(*this, true),
462  writeStreamOption_(IOstream::ASCII),
463  graphFormat_("raw"),
464  runTimeModifiable_(false),
465  functionObjects_(*this, false)
466 {
467  if (enableFunctionObjects)
468  {
469  functionObjects_.on();
470  }
471 
472  if (enableLibs)
473  {
474  libs_.open("libs", controlDict_);
475  }
476 
477  // Explicitly set read flags on objectRegistry so anything constructed
478  // from it reads as well (e.g. fvSolution).
480 
481  setControls();
482  setMonitoring();
483 }
484 
485 
487 (
488  const word& ctrlDictName,
489  const argList& args,
490  const word& systemName,
491  const word& constantName,
492  const bool enableFunctionObjects,
493  const bool enableLibs
494 )
495 :
496  TimePaths(args, systemName, constantName),
497 
498  objectRegistry(*this),
499 
500  loopProfiling_(nullptr),
501  libs_(),
502 
503  controlDict_
504  (
505  IOobject
506  (
507  ctrlDictName,
508  system(),
509  *this,
512  false
513  )
514  ),
515 
516  startTimeIndex_(0),
517  startTime_(0),
518  endTime_(0),
519 
520  stopAt_(saEndTime),
521  writeControl_(wcTimeStep),
522  writeInterval_(GREAT),
523  purgeWrite_(0),
524  subCycling_(0),
525  writeOnce_(false),
526  sigWriteNow_(*this, true),
527  sigStopAtWriteNow_(*this, true),
528  writeStreamOption_(IOstream::ASCII),
529  graphFormat_("raw"),
530  runTimeModifiable_(false),
531  functionObjects_(*this, false)
532 {
533  // Functions
534  //
535  // * '-withFunctionObjects' exists and used = enable
536  // * '-noFunctionObjects' exists and used = disable
537  // * default: no functions if there is no way to enable/disable them
538  if
539  (
540  argList::validOptions.found("withFunctionObjects")
541  ? args.found("withFunctionObjects")
542  : argList::validOptions.found("noFunctionObjects")
543  ? !args.found("noFunctionObjects")
544  : false
545  )
546  {
547  if (enableFunctionObjects)
548  {
549  functionObjects_.on();
550  }
551  }
552 
553  // Libraries
554  //
555  // * enable by default unless '-no-libs' option was used
556  if (enableLibs && !args.found("no-libs"))
557  {
558  libs_.open("libs", controlDict_);
559  }
560 
561  // Explicitly set read flags on objectRegistry so anything constructed
562  // from it reads as well (e.g. fvSolution).
564 
565  setControls();
566 
567  // '-profiling' = force profiling, ignore controlDict entry
568  setMonitoring(args.found("profiling"));
569 }
570 
571 
573 (
574  const dictionary& dict,
575  const fileName& rootPath,
576  const fileName& caseName,
577  const word& systemName,
578  const word& constantName,
579  const bool enableFunctionObjects,
580  const bool enableLibs
581 )
582 :
583  TimePaths
584  (
585  rootPath,
586  caseName,
587  systemName,
588  constantName
589  ),
590 
591  objectRegistry(*this),
592 
593  loopProfiling_(nullptr),
594  libs_(),
595 
596  controlDict_
597  (
598  IOobject
599  (
600  controlDictName,
601  system(),
602  *this,
605  false
606  ),
607  dict
608  ),
609 
610  startTimeIndex_(0),
611  startTime_(0),
612  endTime_(0),
613 
614  stopAt_(saEndTime),
615  writeControl_(wcTimeStep),
616  writeInterval_(GREAT),
617  purgeWrite_(0),
618  subCycling_(0),
619  writeOnce_(false),
620  sigWriteNow_(*this, true),
621  sigStopAtWriteNow_(*this, true),
622  writeStreamOption_(IOstream::ASCII),
623  graphFormat_("raw"),
624  runTimeModifiable_(false),
625  functionObjects_(*this, false)
626 {
627  if (enableFunctionObjects)
628  {
629  functionObjects_.on();
630  }
631 
632  if (enableLibs)
633  {
634  libs_.open("libs", controlDict_);
635  }
636 
637 
638  // Explicitly set read flags on objectRegistry so anything constructed
639  // from it reads as well (e.g. fvSolution).
641 
642  // Since could not construct regIOobject with setting:
643  controlDict_.readOpt(IOobject::MUST_READ_IF_MODIFIED);
644 
645  setControls();
646  setMonitoring();
647 }
648 
649 
651 (
652  const fileName& rootPath,
653  const fileName& caseName,
654  const word& systemName,
655  const word& constantName,
656  const bool enableFunctionObjects,
657  const bool enableLibs
658 )
659 :
660  TimePaths
661  (
662  rootPath,
663  caseName,
664  systemName,
665  constantName
666  ),
667 
668  objectRegistry(*this),
669 
670  loopProfiling_(nullptr),
671  libs_(),
672 
673  controlDict_
674  (
675  IOobject
676  (
677  controlDictName,
678  system(),
679  *this,
682  false
683  )
684  ),
685 
686  startTimeIndex_(0),
687  startTime_(0),
688  endTime_(0),
689 
690  stopAt_(saEndTime),
691  writeControl_(wcTimeStep),
692  writeInterval_(GREAT),
693  purgeWrite_(0),
694  subCycling_(0),
695  writeOnce_(false),
696  writeStreamOption_(IOstream::ASCII),
697  graphFormat_("raw"),
698  runTimeModifiable_(false),
699  functionObjects_(*this, false)
700 {
701  if (enableFunctionObjects)
702  {
703  functionObjects_.on();
704  }
705 
706  if (enableLibs)
707  {
708  libs_.open("libs", controlDict_);
709  }
710 
711  setMonitoring(); // for profiling etc
712 }
713 
714 
715 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
716 
718 {
719  return
721  (
722  fileName("."), // root-path
723  fileName("."), // case-name
724  false, // No enableFunctionObjects
725  false // No enableLibs
726  );
727 }
728 
729 
731 {
732  return
734  (
735  caseDir.path(), // root-path
736  caseDir.name(), // case-name
737  false, // No enableFunctionObjects
738  false // No enableLibs
739  );
740 }
741 
742 
744 {
745  return
747  (
749  args,
750  false, // No enableFunctionObjects
751  false // No enableLibs
752  );
753 }
754 
755 
756 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
757 
759 {
760  loopProfiling_.reset(nullptr);
761 
762  forAllReverse(controlDict_.watchIndices(), i)
763  {
764  fileHandler().removeWatch(controlDict_.watchIndices()[i]);
765  }
766 
767  // Destroy function objects first
768  functionObjects_.clear();
769 
770  // Clean up profiling
771  profiling::stop(*this);
772 
773  // Ensure all owned objects are also cleaned up now
775 }
776 
777 
778 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
779 
780 Foam::word Foam::Time::timeName(const scalar t, const int precision)
781 {
782  std::ostringstream buf;
783  buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
784  buf.precision(precision);
785  buf << t;
786  return buf.str();
787 }
788 
789 
791 {
792  return dimensionedScalar::name();
793 }
794 
795 
797 (
798  const fileName& dir,
799  const word& name,
800  const IOobject::readOption rOpt,
801  const word& stopInstance
802 ) const
803 {
804  IOobject startIO
805  (
806  name, // name might be empty!
807  timeName(),
808  dir,
809  *this,
810  rOpt
811  );
812 
813  IOobject io
814  (
815  fileHandler().findInstance
816  (
817  startIO,
818  timeOutputValue(),
819  stopInstance
820  )
821  );
822  return io.instance();
823 }
824 
825 
827 (
828  const fileName& directory,
829  const instant& t
830 ) const
831 {
832  // Simplified version: use findTimes (readDir + sort). The expensive
833  // bit is the readDir, not the sorting. Tbd: avoid calling findInstancePath
834  // from filePath.
835 
836  instantList timeDirs = findTimes(path(), constant());
837  // Note:
838  // - times will include constant (with value 0) as first element.
839  // For backwards compatibility make sure to find 0 in preference
840  // to constant.
841  // - list is sorted so could use binary search
842 
843  forAllReverse(timeDirs, i)
844  {
845  if (t.equal(timeDirs[i].value()))
846  {
847  return timeDirs[i].name();
848  }
849  }
850 
851  return word::null;
852 }
853 
854 
856 {
857  return findInstancePath(path(), t);
858 }
859 
860 
861 Foam::label Foam::Time::startTimeIndex() const
862 {
863  return startTimeIndex_;
864 }
865 
866 
868 {
869  return dimensionedScalar("startTime", dimTime, startTime_);
870 }
871 
872 
874 {
875  return dimensionedScalar("endTime", dimTime, endTime_);
876 }
877 
878 
880 {
881  return stopAt_;
882 }
883 
884 
885 bool Foam::Time::run() const
886 {
887  loopProfiling_.reset(nullptr);
888 
889  bool isRunning = value() < (endTime_ - 0.5*deltaT_);
890 
891  if (!subCycling_)
892  {
893  // Only execute when the condition is no longer true
894  // ie, when exiting the control loop
895  if (!isRunning && timeIndex_ != startTimeIndex_)
896  {
897  // Ensure functionObjects execute on last time step
898  // (and hence write uptodate functionObjectProperties)
899  {
900  addProfiling(fo, "functionObjects.execute()");
901  functionObjects_.execute();
902  }
903  {
904  addProfiling(fo, "functionObjects.end()");
905  functionObjects_.end();
906  }
907  }
908  }
909 
910  if (isRunning)
911  {
912  if (!subCycling_)
913  {
914  const_cast<Time&>(*this).readModifiedObjects();
915 
916  if (timeIndex_ == startTimeIndex_)
917  {
918  addProfiling(functionObjects, "functionObjects.start()");
919  functionObjects_.start();
920  }
921  else
922  {
923  addProfiling(functionObjects, "functionObjects.execute()");
924  functionObjects_.execute();
925  }
926 
927  // Check if the execution of functionObjects require re-reading
928  // any files. This moves effect of e.g. 'timeActivatedFileUpdate'
929  // one time step forward. Note that we cannot call
930  // readModifiedObjects from within timeActivatedFileUpdate since
931  // it might re-read the functionObjects themselves (and delete
932  // the timeActivatedFileUpdate one)
933  if (functionObjects_.filesModified())
934  {
935  const_cast<Time&>(*this).readModifiedObjects();
936  }
937  }
938 
939  // Update the "is-running" status following the
940  // possible side-effects from functionObjects
941  isRunning = value() < (endTime_ - 0.5*deltaT_);
942 
943  // (re)trigger profiling
944  if (profiling::active())
945  {
946  loopProfiling_.reset
947  (
948  new profilingTrigger("time.run() " + objectRegistry::name())
949  );
950  }
951  }
952 
953  return isRunning;
954 }
955 
956 
958 {
959  const bool isRunning = run();
960 
961  if (isRunning)
962  {
963  operator++();
964  }
965 
966  return isRunning;
967 }
968 
969 
970 bool Foam::Time::end() const
971 {
972  return value() > (endTime_ + 0.5*deltaT_);
973 }
974 
975 
976 bool Foam::Time::stopAt(const stopAtControls stopCtrl) const
977 {
978  if (stopCtrl == stopAtControls::saUnknown)
979  {
980  return false;
981  }
982 
983  const bool changed = (stopAt_ != stopCtrl);
984  stopAt_ = stopCtrl;
985  endTime_ = GREAT;
986 
987  // Adjust endTime
988  if (stopCtrl == stopAtControls::saEndTime)
989  {
990  controlDict_.readEntry("endTime", endTime_);
991  }
992 
993  return changed;
994 }
995 
996 
998 {
999  return controlDict_.getOrDefault("adjustTimeStep", false);
1000 }
1001 
1002 
1004 {
1005  value() = t.value();
1006  dimensionedScalar::name() = t.dimensionedScalar::name();
1007  timeIndex_ = t.timeIndex_;
1008  fileHandler().setTime(*this);
1009 }
1010 
1011 
1012 void Foam::Time::setTime(const instant& inst, const label newIndex)
1013 {
1014  value() = inst.value();
1015  dimensionedScalar::name() = inst.name();
1016  timeIndex_ = newIndex;
1017 
1018  IOdictionary timeDict
1019  (
1020  IOobject
1021  (
1022  "time",
1023  timeName(),
1024  "uniform",
1025  *this,
1028  false
1029  )
1030  );
1031 
1032  timeDict.readIfPresent("deltaT", deltaT_);
1033  timeDict.readIfPresent("deltaT0", deltaT0_);
1034  timeDict.readIfPresent("index", timeIndex_);
1035  fileHandler().setTime(*this);
1036 }
1037 
1038 
1039 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
1040 {
1041  setTime(newTime.value(), newIndex);
1042 }
1043 
1044 
1045 void Foam::Time::setTime(const scalar newTime, const label newIndex)
1046 {
1047  value() = newTime;
1048  dimensionedScalar::name() = timeName(timeToUserTime(newTime));
1049  timeIndex_ = newIndex;
1050  fileHandler().setTime(*this);
1051 }
1052 
1053 
1055 {
1056  setEndTime(endTime.value());
1057 }
1058 
1059 
1060 void Foam::Time::setEndTime(const scalar endTime)
1061 {
1062  endTime_ = endTime;
1063 }
1064 
1065 
1068  const dimensionedScalar& deltaT,
1069  const bool adjust
1070 )
1071 {
1072  setDeltaT(deltaT.value(), adjust);
1073 }
1074 
1075 
1076 void Foam::Time::setDeltaT(const scalar deltaT, const bool adjust)
1077 {
1078  deltaT_ = deltaT;
1079  deltaTchanged_ = true;
1080 
1081  if (adjust)
1082  {
1083  adjustDeltaT();
1084  }
1085 }
1086 
1087 
1088 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
1089 {
1090  #ifdef FULLDEBUG
1091  if (prevTimeState_)
1092  {
1094  << "previous time state already set" << nl
1095  << exit(FatalError);
1096  }
1097  #endif
1098 
1099  prevTimeState_.reset(new TimeState(*this));
1100 
1101  setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
1102  deltaT_ /= nSubCycles;
1103  deltaT0_ /= nSubCycles;
1104  deltaTSave_ = deltaT0_;
1105 
1106  subCycling_ = nSubCycles;
1107 
1108  return prevTimeState();
1109 }
1110 
1111 
1112 void Foam::Time::subCycleIndex(const label index)
1113 {
1114  // Only permit adjustment if sub-cycling was already active
1115  // and if the index is valid (positive, non-zero).
1116  // This avoids potential mixups for deleting.
1117 
1118  if (subCycling_ && index > 0)
1119  {
1120  subCycling_ = index;
1121  }
1122 }
1123 
1124 
1126 {
1127  if (subCycling_)
1128  {
1129  TimeState::operator=(prevTimeState());
1130  prevTimeState_.reset(nullptr);
1131  }
1132 
1133  subCycling_ = 0;
1134 }
1135 
1136 
1137 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1138 
1140 {
1141  return operator+=(deltaT.value());
1142 }
1143 
1144 
1145 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
1146 {
1147  setDeltaT(deltaT);
1148  return operator++();
1149 }
1150 
1151 
1153 {
1154  deltaT0_ = deltaTSave_;
1155  deltaTSave_ = deltaT_;
1156 
1157  // Save old time value and name
1158  const scalar oldTimeValue = timeToUserTime(value());
1159  const word oldTimeName = dimensionedScalar::name();
1160 
1161  // Increment time
1162  setTime(value() + deltaT_, timeIndex_ + 1);
1163 
1164  if (!subCycling_)
1165  {
1166  // If the time is very close to zero reset to zero
1167  if (mag(value()) < 10*SMALL*deltaT_)
1168  {
1169  setTime(0.0, timeIndex_);
1170  }
1171 
1172  if (sigStopAtWriteNow_.active() || sigWriteNow_.active())
1173  {
1174  // A signal might have been sent on one processor only
1175  // Reduce so all decide the same.
1176 
1177  label flag = 0;
1178  if (sigStopAtWriteNow_.active() && stopAt_ == saWriteNow)
1179  {
1180  flag += 1;
1181  }
1182  if (sigWriteNow_.active() && writeOnce_)
1183  {
1184  flag += 2;
1185  }
1186  reduce(flag, maxOp<label>());
1187 
1188  if (flag & 1)
1189  {
1190  stopAt_ = saWriteNow;
1191  }
1192  if (flag & 2)
1193  {
1194  writeOnce_ = true;
1195  }
1196  }
1197 
1198  writeTime_ = false;
1199 
1200  switch (writeControl_)
1201  {
1202  case wcNone:
1203  case wcUnknown:
1204  break;
1205 
1206  case wcTimeStep:
1207  writeTime_ = !(timeIndex_ % label(writeInterval_));
1208  break;
1209 
1210  case wcRunTime:
1211  case wcAdjustableRunTime:
1212  {
1213  const label writeIndex = label
1214  (
1215  ((value() - startTime_) + 0.5*deltaT_)
1216  / writeInterval_
1217  );
1218 
1219  if (writeIndex > writeTimeIndex_)
1220  {
1221  writeTime_ = true;
1222  writeTimeIndex_ = writeIndex;
1223  }
1224  }
1225  break;
1226 
1227  case wcCpuTime:
1228  {
1229  const label writeIndex = label
1230  (
1231  returnReduce(elapsedCpuTime(), maxOp<double>())
1232  / writeInterval_
1233  );
1234  if (writeIndex > writeTimeIndex_)
1235  {
1236  writeTime_ = true;
1237  writeTimeIndex_ = writeIndex;
1238  }
1239  }
1240  break;
1241 
1242  case wcClockTime:
1243  {
1244  const label writeIndex = label
1245  (
1246  returnReduce(elapsedClockTime(), maxOp<double>())
1247  / writeInterval_
1248  );
1249  if (writeIndex > writeTimeIndex_)
1250  {
1251  writeTime_ = true;
1252  writeTimeIndex_ = writeIndex;
1253  }
1254  }
1255  break;
1256  }
1257 
1258 
1259  // Check if endTime needs adjustment to stop at the next run()/end()
1260  if (!end())
1261  {
1262  if (stopAt_ == saNoWriteNow)
1263  {
1264  endTime_ = value();
1265  }
1266  else if (stopAt_ == saWriteNow)
1267  {
1268  endTime_ = value();
1269  writeTime_ = true;
1270  }
1271  else if (stopAt_ == saNextWrite && writeTime_ == true)
1272  {
1273  endTime_ = value();
1274  }
1275  }
1276 
1277  // Override writeTime if one-shot writing
1278  if (writeOnce_)
1279  {
1280  writeTime_ = true;
1281  writeOnce_ = false;
1282  }
1283 
1284  // Adjust the precision of the time directory name if necessary
1285  if (writeTime_)
1286  {
1287  // Tolerance used when testing time equivalence
1288  const scalar timeTol =
1289  max(min(pow(10.0, -precision_), 0.1*deltaT_), SMALL);
1290 
1291  // User-time equivalent of deltaT
1292  const scalar userDeltaT = timeToUserTime(deltaT_);
1293 
1294  // Time value obtained by reading timeName
1295  scalar timeNameValue = -VGREAT;
1296 
1297  // Check that new time representation differs from old one
1298  // reinterpretation of the word
1299  if
1300  (
1301  readScalar(dimensionedScalar::name(), timeNameValue)
1302  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1303  )
1304  {
1305  int oldPrecision = precision_;
1306  while
1307  (
1308  precision_ < maxPrecision_
1309  && readScalar(dimensionedScalar::name(), timeNameValue)
1310  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1311  )
1312  {
1313  precision_++;
1314  setTime(value(), timeIndex());
1315  }
1316 
1317  if (precision_ != oldPrecision)
1318  {
1320  << "Increased the timePrecision from " << oldPrecision
1321  << " to " << precision_
1322  << " to distinguish between timeNames at time "
1324  << endl;
1325 
1326  if (precision_ == maxPrecision_)
1327  {
1328  // Reached maxPrecision limit
1330  << "Current time name " << dimensionedScalar::name()
1331  << nl
1332  << " The maximum time precision has been reached"
1333  " which might result in overwriting previous"
1334  " results."
1335  << endl;
1336  }
1337 
1338  // Check if round-off error caused time-reversal
1339  scalar oldTimeNameValue = -VGREAT;
1340  if
1341  (
1342  readScalar(oldTimeName, oldTimeNameValue)
1343  && (
1344  sign(timeNameValue - oldTimeNameValue)
1345  != sign(deltaT_)
1346  )
1347  )
1348  {
1350  << "Current time name " << dimensionedScalar::name()
1351  << " is set to an instance prior to the "
1352  "previous one "
1353  << oldTimeName << nl
1354  << " This might result in temporal "
1355  "discontinuities."
1356  << endl;
1357  }
1358  }
1359  }
1360  }
1361  }
1362 
1363  return *this;
1364 }
1365 
1366 
1368 {
1369  return operator++();
1370 }
1371 
1372 
1373 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::Time::general
default float notation
Definition: Time.H:109
profiling.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::Time::Time
Time(const word &ctrlDictName, const argList &args, const bool enableFunctionObjects=true, const bool enableLibs=true)
Construct from dictionary name to read and argument list.
Definition: TimeI.H:31
Foam::maxOp
Definition: ops.H:223
Foam::Time::writeControls
writeControls
Write control options.
Definition: Time.H:84
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::functionObjectList::adjustTimeStep
bool adjustTimeStep()
Called at the end of Time::adjustDeltaT() if adjustTime is true.
Definition: functionObjectList.C:902
Foam::Time::end
virtual bool end() const
Return true if end of run,.
Definition: Time.C:970
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::IOobject::instance
const fileName & instance() const noexcept
Definition: IOobjectI.H:196
Foam::fileOperation::addWatches
virtual void addWatches(regIOobject &, const fileNameList &) const
Helper: add watches for list of regIOobjects.
Definition: fileOperation.C:881
Foam::fileName::path
static std::string path(const std::string &str)
Return directory path name (part before last /)
Definition: fileNameI.H:176
Foam::Time::writeControlNames
static const Enum< writeControls > writeControlNames
Names for writeControls.
Definition: Time.H:116
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:199
Foam::system
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition: MSwindows.C:1150
oldTime
Info<< "Creating field kinetic energy K\n"<< 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
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
Foam::Time::stopAt
virtual stopAtControls stopAt() const
Return the stop control information.
Definition: Time.C:879
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:123
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1485
Foam::FatalIOError
IOerror FatalIOError
Foam::objectRegistry::clear
void clear()
Clear all entries from the registry.
Definition: objectRegistry.C:318
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::Instant::equal
bool equal(scalar val) const
True if values are equal (includes SMALL for rounding)
Definition: Instant.C:59
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::Time::controlDictName
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:226
setTime
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
Foam::Time::fmtflags
fmtflags
Supported time directory name formats.
Definition: Time.H:107
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::sumOp
Definition: ops.H:213
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::debug::infoSwitch
int infoSwitch(const char *name, const int deflt=0)
Lookup info switch or add default value.
Definition: debug.C:231
Foam::Time::stopAtControlNames
static const Enum< stopAtControls > stopAtControlNames
Names for stopAtControls.
Definition: Time.H:119
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Time::precision_
static int precision_
Time directory name precision.
Definition: Time.H:186
Foam::Time::timeName
virtual word timeName() const
Return current time name.
Definition: Time.C:790
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Time::adjustDeltaT
void adjustDeltaT()
Adjust the time step so that writing occurs at the specified time.
Definition: Time.C:99
Foam::TimeState
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:51
Foam::Time::~Time
virtual ~Time()
Destructor.
Definition: Time.C:758
registerInfoSwitch
registerInfoSwitch("printExecutionFormat", int, Foam::Time::printExecutionFormat_)
Foam::Time::subCycle
virtual TimeState subCycle(const label nSubCycles)
Set time to sub-cycle for the given number of steps.
Definition: Time.C:1088
Foam::profiling::stop
static void stop(const Time &owner)
Stop profiling, cleanup pool if possible.
Definition: profiling.C:172
Foam::Time::setEndTime
virtual void setEndTime(const dimensionedScalar &endTime)
Reset end time.
Definition: Time.C:1054
Foam::profiling::active
static bool active()
True if profiling is allowed and is active.
Definition: profiling.C:111
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::fileOperation::setTime
virtual void setTime(const Time &) const
Callback for time change.
Definition: fileOperation.H:548
argList.H
Foam::TimeState::writeTimeIndex_
label writeTimeIndex_
Definition: TimeState.H:58
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::Time::operator+=
virtual Time & operator+=(const dimensionedScalar &deltaT)
Set deltaT to that specified and increment time via operator++()
Definition: Time.C:1139
timeName
word timeName
Definition: getTimeIndex.H:3
HashSet.H
dict
dictionary dict
Definition: searchingEngine.H:14
addProfiling
#define addProfiling(name, descr)
Define profiling trigger with specified name and description string.
Definition: profilingTrigger.H:113
Foam::argList::validOptions
static HashTable< string > validOptions
A list of valid options.
Definition: argList.H:212
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::Time::loop
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:957
Foam::Time::operator++
virtual Time & operator++()
Prefix increment,.
Definition: Time.C:1152
Foam::dimensioned< scalar >
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam::Time::wcAdjustableRunTime
"adjustable" / "adjustableRunTime"
Definition: Time.H:89
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::Time::New
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
PstreamReduceOps.H
Inter-processor communication reduction functions.
Foam::TimePaths
Address the time paths without using the Time class.
Definition: TimePaths.H:56
Foam::Time::findInstance
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Definition: Time.C:797
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fileOperation::removeWatch
virtual bool removeWatch(const label) const
Remove watch on a file (using handle)
Definition: fileOperation.C:857
Foam::profiling::initialize
static void initialize(const IOobject &ioObj, const Time &owner)
Singleton to initialize profiling pool, everything enabled.
Definition: profiling.C:146
Foam::Instant::name
const T & name() const
The name/key (const access)
Definition: Instant.H:111
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::Time::setDeltaT
virtual void setDeltaT(const dimensionedScalar &deltaT, const bool adjust=true)
Reset time step, normally also calling adjustDeltaT()
Definition: Time.C:1067
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
IOdictionary.H
Time.H
Foam::autoPtr< Foam::Time >
Foam::Time::findInstancePath
word findInstancePath(const fileName &directory, const instant &t) const
Definition: Time.C:827
Foam::Time::readModifiedObjects
void readModifiedObjects()
Read the objects that have been modified.
Definition: TimeIO.C:467
Foam::IOstreamOption::ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::profilingTrigger
Triggers for starting/stopping code profiling.
Definition: profilingTrigger.H:54
Foam::Time::printExecutionFormat_
static int printExecutionFormat_
Style for "ExecutionTime = " output.
Definition: Time.H:127
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< instant >
Foam::Time::maxPrecision_
static const int maxPrecision_
Maximum time directory name precision.
Definition: Time.H:189
Foam::Time::setTime
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:1003
Foam::Time::writeInterval_
scalar writeInterval_
Definition: Time.H:159
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::Instant::value
scalar value() const
The value (const access)
Definition: Instant.H:99
Foam::Time::format_
static fmtflags format_
Time directory name format.
Definition: Time.H:183
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::Time::startTime_
scalar startTime_
Definition: Time.H:151
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:186
Foam::TimeState::timeIndex_
label timeIndex_
Definition: TimeState.H:57
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:309
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::Time::isAdjustTimeStep
virtual bool isAdjustTimeStep() const
Return true if adjustTimeStep is true.
Definition: Time.C:997
registerSwitch.H
Foam::Time::stopAtControls
stopAtControls
Definition: Time.H:97
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
Foam::Time::endTime
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:873
Foam::instant
An instant of time. Contains the time value and name.
Definition: instant.H:52
Foam::Time::setMonitoring
void setMonitoring(const bool forceProfiling=false)
Set file monitoring, profiling, etc.
Definition: Time.C:349
Foam::debug::controlDict
dictionary & controlDict()
Definition: debug.C:143
Foam::IOobject::readOption
readOption
Enumeration defining the read options.
Definition: IOobject.H:183
Foam::Time::endSubCycle
virtual void endSubCycle()
Reset time after sub-cycling back to previous TimeState.
Definition: Time.C:1125
Foam::Time::setControls
void setControls()
Set the controls from the current controlDict.
Definition: Time.C:143
IOWarningInFunction
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
Definition: messageStream.H:340
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::IOobject::NO_READ
Definition: IOobject.H:188
args
Foam::argList args(argc, argv)
constant
constant condensation/saturation model.
Foam::Time::writeControl_
writeControls writeControl_
Definition: Time.H:157
Foam::Time::run
virtual bool run() const
Return true if run should continue,.
Definition: Time.C:885
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::TimeState::deltaT_
scalar deltaT_
Definition: TimeState.H:60
Foam::Time::startTime
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:867
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::Time::startTimeIndex
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:861
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::Time::subCycleIndex
virtual void subCycleIndex(const label index)
Adjust the reported sub-cycle index.
Definition: Time.C:1112
Foam::argList::found
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178