Time.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "Time.H"
27 #include "timeIOdictionary.H"
28 #include "argList.H"
29 
30 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(Time, 0);
35 
36  template<>
37  const char* Foam::NamedEnum
38  <
40  4
41  >::names[] =
42  {
43  "endTime",
44  "noWriteNow",
45  "writeNow",
46  "nextWrite"
47  };
48 
49  template<>
50  const char* Foam::NamedEnum
51  <
53  5
54  >::names[] =
55  {
56  "timeStep",
57  "runTime",
58  "adjustableRunTime",
59  "clockTime",
60  "cpuTime"
61  };
62 }
63 
66 
69 
71 
73 
74 int Foam::Time::curPrecision_(Foam::Time::precision_);
75 
76 const int Foam::Time::maxPrecision_(3 - log10(small));
77 
79 
80 
81 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
82 
84 {
85  const scalar timeToNextWrite = min
86  (
87  max
88  (
89  0,
90  (writeTimeIndex_ + 1)*writeInterval_ - (value() - beginTime_)
91  ),
92  functionObjects_.timeToNextWrite()
93  );
94 
95  const scalar nSteps = timeToNextWrite/deltaT_;
96 
97  // Ensure nStepsToNextWrite does not overflow
98  if (nSteps < labelMax)
99  {
100  // Allow the time-step to increase by up to 1%
101  // to accommodate the next write time before splitting
102  const label nStepsToNextWrite = label(max(nSteps, 1) + 0.99);
103 
104  deltaT_ = timeToNextWrite/nStepsToNextWrite;
105  }
106 }
107 
108 
110 {
111  // default is to resume calculation from "latestTime"
112  const word startFrom = controlDict_.lookupOrDefault<word>
113  (
114  "startFrom",
115  "latestTime"
116  );
117 
118  if (startFrom == "startTime")
119  {
120  controlDict_.lookup("startTime") >> startTime_;
121  startTime_ = userTimeToTime(startTime_);
122  }
123  else
124  {
125  // Search directory for valid time directories
126  instantList timeDirs = findTimes(path(), constant());
127 
128  if (startFrom == "firstTime")
129  {
130  if (timeDirs.size())
131  {
132  if (timeDirs[0].name() == constant() && timeDirs.size() >= 2)
133  {
134  startTime_ = userTimeToTime(timeDirs[1].value());
135  }
136  else
137  {
138  startTime_ = userTimeToTime(timeDirs[0].value());
139  }
140  }
141  }
142  else if (startFrom == "latestTime")
143  {
144  if (timeDirs.size())
145  {
146  startTime_ = userTimeToTime(timeDirs.last().value());
147  }
148  }
149  else
150  {
151  FatalIOErrorInFunction(controlDict_)
152  << "expected startTime, firstTime or latestTime"
153  << " found '" << startFrom << "'"
154  << exit(FatalIOError);
155  }
156  }
157 
158  setTime(startTime_, 0);
159  readDict();
160  deltaTSave_ = deltaT_;
161  deltaT0_ = deltaT_;
162 
163  // Check if time directory exists
164  // If not increase time precision to see if it is formatted differently.
165  if (!fileHandler().exists(timePath(), false, false))
166  {
167  int oldPrecision = curPrecision_;
168  int requiredPrecision = -1;
169  bool found = false;
170  word oldTime(timeName());
171  for
172  (
173  curPrecision_ = maxPrecision_;
174  curPrecision_ > oldPrecision;
175  curPrecision_--
176  )
177  {
178  // Update the time formatting
179  setTime(startTime_, 0);
180 
181  word newTime(timeName());
182 
183  if (newTime == oldTime)
184  {
185  break;
186  }
187  oldTime = newTime;
188 
189  // Check the existence of the time directory with the new format
190  found = fileHandler().exists(timePath(), false, false);
191 
192  if (found)
193  {
194  requiredPrecision = curPrecision_;
195  }
196  }
197 
198  if (requiredPrecision > 0)
199  {
200  // Update the time precision
201  curPrecision_ = requiredPrecision;
202 
203  // Update the time formatting
204  setTime(startTime_, 0);
205 
207  << "Increasing the timePrecision from " << oldPrecision
208  << " to " << curPrecision_
209  << " to support the formatting of the current time directory "
210  << timeName() << nl << endl;
211  }
212  else
213  {
214  // Could not find time directory so assume it is not present
215  curPrecision_ = oldPrecision;
216 
217  // Revert the time formatting
218  setTime(startTime_, 0);
219  }
220  }
221 
222  if (Pstream::parRun())
223  {
224  scalar sumStartTime = startTime_;
225  reduce(sumStartTime, sumOp<scalar>());
226  if
227  (
228  mag(Pstream::nProcs()*startTime_ - sumStartTime)
229  > Pstream::nProcs()*deltaT_/10.0
230  )
231  {
232  FatalIOErrorInFunction(controlDict_)
233  << "Start time is not the same for all processors" << nl
234  << "processor " << Pstream::myProcNo() << " has startTime "
235  << startTime_ << exit(FatalIOError);
236  }
237  }
238 
239  timeIOdictionary timeDict
240  (
241  IOobject
242  (
243  "time",
244  timeName(),
245  "uniform",
246  *this,
249  false
250  )
251  );
252 
253  if (timeDict.found("beginTime"))
254  {
255  beginTime_ = userTimeToTime(timeDict.lookup<scalar>("beginTime"));
256  }
257  else
258  {
259  beginTime_ = startTime_;
260  }
261 
262  // Read and set the deltaT only if time-step adjustment is active
263  // otherwise use the deltaT from the controlDict
264  if (controlDict_.lookupOrDefault<Switch>("adjustTimeStep", false))
265  {
266  if (timeDict.found("deltaT"))
267  {
268  deltaT_ = userTimeToTime(timeDict.lookup<scalar>("deltaT"));
269  deltaTSave_ = deltaT_;
270  deltaT0_ = deltaT_;
271  }
272  }
273 
274  if (timeDict.found("deltaT0"))
275  {
276  deltaT0_ = userTimeToTime(timeDict.lookup<scalar>("deltaT0"));
277  }
278 
279  if (timeDict.readIfPresent("index", startTimeIndex_))
280  {
281  timeIndex_ = startTimeIndex_;
282  }
283 
284  // Set writeTimeIndex_ to correspond to beginTime_
285  if
286  (
287  (
288  writeControl_ == writeControl::runTime
289  || writeControl_ == writeControl::adjustableRunTime
290  )
291  )
292  {
293  writeTimeIndex_ = label
294  (
295  ((value() - beginTime_) + 0.5*deltaT_)/writeInterval_
296  );
297  }
298 
299 
300  // Check if values stored in time dictionary are consistent
301 
302  // 1. Based on time name
303  bool checkValue = true;
304 
305  string storedTimeName;
306  if (timeDict.readIfPresent("name", storedTimeName))
307  {
308  if (storedTimeName == timeName())
309  {
310  // Same time. No need to check stored value
311  checkValue = false;
312  }
313  }
314 
315  // 2. Based on time value
316  // (consistent up to the current time writing precision so it won't
317  // trigger if we just change the write precision)
318  if (checkValue)
319  {
320  scalar storedTimeValue;
321  if (timeDict.readIfPresent("value", storedTimeValue))
322  {
323  word storedTimeName(timeName(storedTimeValue));
324 
325  if (storedTimeName != timeName())
326  {
327  IOWarningInFunction(timeDict)
328  << "Time read from time dictionary " << storedTimeName
329  << " differs from actual time " << timeName() << '.' << nl
330  << " This may cause unexpected database behaviour."
331  << " If you are not interested" << nl
332  << " in preserving time state delete"
333  << " the time dictionary."
334  << endl;
335  }
336  }
337  }
338 }
339 
340 
341 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
342 
344 (
345  const word& controlDictName,
346  const fileName& rootPath,
347  const fileName& caseName,
348  const word& systemName,
349  const word& constantName,
350  const bool enableFunctionObjects
351 )
352 :
353  TimePaths
354  (
355  rootPath,
356  caseName,
357  systemName,
358  constantName
359  ),
360 
361  objectRegistry(*this),
362 
363  runTimeModifiable_(false),
364 
365  controlDict_
366  (
367  IOobject
368  (
369  controlDictName,
370  system(),
371  *this,
374  )
375  ),
376 
377  startTimeIndex_(0),
378  startTime_(0),
379  endTime_(0),
380  beginTime_(startTime_),
381 
382  userTime_(userTimes::userTime::New(controlDict_)),
383 
384  stopAt_(stopAtControl::endTime),
385  writeControl_(writeControl::timeStep),
386  writeInterval_(great),
387  purgeWrite_(0),
388  writeOnce_(false),
389 
390  subCycling_(false),
391 
392  sigWriteNow_(writeInfoHeader, *this),
393  sigStopAtWriteNow_(writeInfoHeader, *this),
394 
395  writeFormat_(IOstream::ASCII),
396  writeVersion_(IOstream::currentVersion),
397  writeCompression_(IOstream::UNCOMPRESSED),
398  graphFormat_("raw"),
399  cacheTemporaryObjects_(true),
400 
401  functionObjects_(*this, enableFunctionObjects)
402 {
403  libs.open(controlDict_, "libs");
404 
405  // Explicitly set read flags on objectRegistry so anything constructed
406  // from it reads as well (e.g. fvSolution).
408 
409  setControls();
410 
411  // Add a watch on the controlDict file after runTimeModifiable_ is set
412  controlDict_.addWatch();
413 }
414 
415 
417 (
418  const word& controlDictName,
419  const argList& args,
420  const word& systemName,
421  const word& constantName
422 )
423 :
424  TimePaths
425  (
426  args.parRunControl().parRun(),
427  args.rootPath(),
428  args.globalCaseName(),
429  args.caseName(),
430  systemName,
431  constantName
432  ),
433 
434  objectRegistry(*this),
435 
436  runTimeModifiable_(false),
437 
438  controlDict_
439  (
440  IOobject
441  (
442  controlDictName,
443  system(),
444  *this,
447  )
448  ),
449 
450  startTimeIndex_(0),
451  startTime_(0),
452  endTime_(0),
453  beginTime_(startTime_),
454 
455  userTime_(userTimes::userTime::New(controlDict_)),
456 
457  stopAt_(stopAtControl::endTime),
458  writeControl_(writeControl::timeStep),
459  writeInterval_(great),
460  purgeWrite_(0),
461  writeOnce_(false),
462 
463  subCycling_(false),
464 
465  sigWriteNow_(writeInfoHeader, *this),
466  sigStopAtWriteNow_(writeInfoHeader, *this),
467 
468  writeFormat_(IOstream::ASCII),
469  writeVersion_(IOstream::currentVersion),
470  writeCompression_(IOstream::UNCOMPRESSED),
471  graphFormat_("raw"),
472  cacheTemporaryObjects_(true),
473 
474  functionObjects_
475  (
476  *this,
477  argList::validOptions.found("withFunctionObjects")
478  ? args.optionFound("withFunctionObjects")
479  : !args.optionFound("noFunctionObjects")
480  )
481 {
482  libs.open(controlDict_, "libs");
483 
484  // Explicitly set read flags on objectRegistry so anything constructed
485  // from it reads as well (e.g. fvSolution).
487 
488  if (args.options().found("case"))
489  {
490  const wordList switchSets
491  (
492  {
493  "InfoSwitches",
494  "OptimisationSwitches",
495  "DebugSwitches",
496  "DimensionedConstants",
497  "DimensionSets"
498  }
499  );
500 
501  forAll(switchSets, i)
502  {
503  if (controlDict_.found(switchSets[i]))
504  {
505  IOWarningInFunction(controlDict_)
506  << switchSets[i]
507  << " in system/controlDict are only processed if "
508  << args.executable() << " is run in the "
509  << args.path() << " directory" << endl;
510  }
511  }
512  }
513 
514  setControls();
515 
516  // Add a watch on the controlDict file after runTimeModifiable_ is set
517  controlDict_.addWatch();
518 }
519 
520 
522 (
523  const dictionary& dict,
524  const fileName& rootPath,
525  const fileName& caseName,
526  const word& systemName,
527  const word& constantName,
528  const bool enableFunctionObjects
529 )
530 :
531  TimePaths
532  (
533  rootPath,
534  caseName,
535  systemName,
536  constantName
537  ),
538 
539  objectRegistry(*this),
540 
541  runTimeModifiable_(false),
542 
543  controlDict_
544  (
545  IOobject
546  (
547  controlDictName,
548  system(),
549  *this,
552  ),
553  dict
554  ),
555 
556  startTimeIndex_(0),
557  startTime_(0),
558  endTime_(0),
559  beginTime_(startTime_),
560 
561  userTime_(userTimes::userTime::New(controlDict_)),
562 
563  stopAt_(stopAtControl::endTime),
564  writeControl_(writeControl::timeStep),
565  writeInterval_(great),
566  purgeWrite_(0),
567  writeOnce_(false),
568 
569  subCycling_(false),
570 
571  sigWriteNow_(writeInfoHeader, *this),
572  sigStopAtWriteNow_(writeInfoHeader, *this),
573 
574  writeFormat_(IOstream::ASCII),
575  writeVersion_(IOstream::currentVersion),
576  writeCompression_(IOstream::UNCOMPRESSED),
577  graphFormat_("raw"),
578  cacheTemporaryObjects_(true),
579 
580  functionObjects_(*this, enableFunctionObjects)
581 {
582  libs.open(controlDict_, "libs");
583 
584  // Explicitly set read flags on objectRegistry so anything constructed
585  // from it reads as well (e.g. fvSolution).
587 
588  setControls();
589 
590  // Add a watch on the controlDict file after runTimeModifiable_ is set
591  controlDict_.addWatch();
592 }
593 
594 
596 (
597  const fileName& rootPath,
598  const fileName& caseName,
599  const word& systemName,
600  const word& constantName,
601  const bool enableFunctionObjects
602 )
603 :
604  TimePaths
605  (
606  rootPath,
607  caseName,
608  systemName,
609  constantName
610  ),
611 
612  objectRegistry(*this),
613 
614  runTimeModifiable_(false),
615 
616  controlDict_
617  (
618  IOobject
619  (
620  controlDictName,
621  system(),
622  *this,
625  )
626  ),
627 
628  startTimeIndex_(0),
629  startTime_(0),
630  endTime_(0),
631  beginTime_(startTime_),
632 
633  userTime_(userTimes::userTime::New(controlDict_)),
634 
635  stopAt_(stopAtControl::endTime),
636  writeControl_(writeControl::timeStep),
637  writeInterval_(great),
638  purgeWrite_(0),
639  writeOnce_(false),
640 
641  subCycling_(false),
642 
643  writeFormat_(IOstream::ASCII),
644  writeVersion_(IOstream::currentVersion),
645  writeCompression_(IOstream::UNCOMPRESSED),
646  graphFormat_("raw"),
647  cacheTemporaryObjects_(true),
648 
649  functionObjects_(*this, enableFunctionObjects)
650 {
651  libs.open(controlDict_, "libs");
652 }
653 
654 
655 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
656 
658 {
659  // Destroy function objects first
660  functionObjects_.clear();
661 }
662 
663 
664 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
665 
666 Foam::word Foam::Time::timeName(const scalar t, const int precision)
667 {
668  std::ostringstream buf;
669  buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
670  buf.precision(precision);
671  buf << t;
672  return buf.str();
673 }
674 
675 
677 {
678  return dimensionedScalar::name();
679 }
680 
681 
683 {
684  return findTimes(path(), constant());
685 }
686 
687 
689 (
690  const fileName& dir,
691  const word& name,
692  const IOobject::readOption rOpt,
693  const word& stopInstance
694 ) const
695 {
696  IOobject startIO
697  (
698  name, // name might be empty!
699  timeName(),
700  dir,
701  *this,
702  rOpt
703  );
704 
705  IOobject io
706  (
707  fileHandler().findInstance
708  (
709  startIO,
710  userTimeValue(),
711  stopInstance
712  )
713  );
714  return io.instance();
715 }
716 
717 
719 (
720  const fileName& directory,
721  const instant& t
722 ) const
723 {
724  // Simplified version: use findTimes (readDir + sort). The expensive
725  // bit is the readDir, not the sorting. Tbd: avoid calling findInstancePath
726  // from filePath.
727 
728  instantList timeDirs = findTimes(path(), constant());
729  // Note:
730  // - times will include constant (with value 0) as first element.
731  // For backwards compatibility make sure to find 0 in preference
732  // to constant.
733  // - list is sorted so could use binary search
734 
735  forAllReverse(timeDirs, i)
736  {
737  if (t.equal(timeDirs[i].value()))
738  {
739  return timeDirs[i].name();
740  }
741  }
742 
743  return word::null;
744 }
745 
746 
748 {
749  return findInstancePath(path(), t);
750 }
751 
752 
754 {
755  instantList timeDirs = findTimes(path(), constant());
756 
757  // There is only one time (likely "constant") so return it
758  if (timeDirs.size() == 1)
759  {
760  return timeDirs[0];
761  }
762 
763  if (t < timeDirs[1].value())
764  {
765  return timeDirs[1];
766  }
767  else if (t > timeDirs.last().value())
768  {
769  return timeDirs.last();
770  }
771 
772  label nearestIndex = -1;
773  scalar deltaT = great;
774 
775  for (label timei=1; timei < timeDirs.size(); ++timei)
776  {
777  scalar diff = mag(timeDirs[timei].value() - t);
778  if (diff < deltaT)
779  {
780  deltaT = diff;
781  nearestIndex = timei;
782  }
783  }
784 
785  return timeDirs[nearestIndex];
786 }
787 
788 
790 (
791  const instantList& timeDirs,
792  const scalar t,
793  const word& constantName
794 )
795 {
796  label nearestIndex = -1;
797  scalar deltaT = great;
798 
799  forAll(timeDirs, timei)
800  {
801  if (timeDirs[timei].name() == constantName) continue;
802 
803  scalar diff = mag(timeDirs[timei].value() - t);
804  if (diff < deltaT)
805  {
806  deltaT = diff;
807  nearestIndex = timei;
808  }
809  }
810 
811  return nearestIndex;
812 }
813 
814 
816 {
817  return startTimeIndex_;
818 }
819 
820 
822 {
823  return dimensionedScalar("beginTime", dimTime, beginTime_);
824 }
825 
826 
828 {
829  return dimensionedScalar("startTime", dimTime, startTime_);
830 }
831 
832 
834 {
835  return dimensionedScalar("endTime", dimTime, endTime_);
836 }
837 
838 
840 {
841  return *userTime_;
842 }
843 
844 
845 Foam::scalar Foam::Time::userTimeValue() const
846 {
847  return userTime_->timeToUserTime(value());
848 }
849 
850 
851 Foam::scalar Foam::Time::userTimeToTime(const scalar tau) const
852 {
853  return userTime_->userTimeToTime(tau);
854 }
855 
856 
857 Foam::scalar Foam::Time::timeToUserTime(const scalar t) const
858 {
859  return userTime_->timeToUserTime(t);
860 }
861 
862 
864 {
865  return timeName(userTimeValue()) + userTime_->unit();
866 }
867 
868 
870 {
871  return value() < (endTime_ - 0.5*deltaT_);
872 }
873 
874 
875 bool Foam::Time::run() const
876 {
877  bool running = this->running();
878 
879  if (!subCycling_)
880  {
881  if (!running && timeIndex_ != startTimeIndex_)
882  {
883  functionObjects_.execute();
884  functionObjects_.end();
885 
886  if (cacheTemporaryObjects_)
887  {
888  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
889  }
890  }
891  }
892 
893  if (running)
894  {
895  if (!subCycling_)
896  {
897  const_cast<Time&>(*this).readModifiedObjects();
898 
899  if (timeIndex_ == startTimeIndex_)
900  {
901  functionObjects_.start();
902  }
903  else
904  {
905  functionObjects_.execute();
906 
907  if (cacheTemporaryObjects_)
908  {
909  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
910  }
911  }
912  }
913 
914  // Re-evaluate if running in case a function object has changed things
915  running = this->running();
916  }
917 
918  return running;
919 }
920 
921 
923 {
924  bool running = run();
925 
926  if (running)
927  {
928  operator++();
929  }
930 
931  return running;
932 }
933 
934 
935 bool Foam::Time::end() const
936 {
937  return value() > (endTime_ + 0.5*deltaT_);
938 }
939 
940 
941 bool Foam::Time::stopAt(const stopAtControl sa) const
942 {
943  const bool changed = (stopAt_ != sa);
944  stopAt_ = sa;
945 
946  // adjust endTime
947  if (sa == stopAtControl::endTime)
948  {
949  controlDict_.lookup("endTime") >> endTime_;
950  }
951  else
952  {
953  endTime_ = great;
954  }
955  return changed;
956 }
957 
958 
959 void Foam::Time::setTime(const Time& t)
960 {
961  value() = t.value();
962  dimensionedScalar::name() = t.dimensionedScalar::name();
963  timeIndex_ = t.timeIndex_;
964  fileHandler().setTime(*this);
965 }
966 
967 
968 void Foam::Time::setTime(const instant& inst, const label newIndex)
969 {
970  value() = userTimeToTime(inst.value());
971  dimensionedScalar::name() = inst.name();
972  timeIndex_ = newIndex;
973 
974  timeIOdictionary timeDict
975  (
976  IOobject
977  (
978  "time",
979  timeName(),
980  "uniform",
981  *this,
984  false
985  )
986  );
987 
988  if (timeDict.found("deltaT"))
989  {
990  deltaT_ = userTimeToTime(timeDict.lookup<scalar>("deltaT"));
991  }
992 
993  if (timeDict.found("deltaT0"))
994  {
995  deltaT0_ = userTimeToTime(timeDict.lookup<scalar>("deltaT0"));
996  }
997 
998  timeDict.readIfPresent("index", timeIndex_);
999 
1000  fileHandler().setTime(*this);
1001 }
1002 
1003 
1004 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
1005 {
1006  setTime(newTime.value(), newIndex);
1007 }
1008 
1009 
1010 void Foam::Time::setTime(const scalar newTime, const label newIndex)
1011 {
1012  value() = newTime;
1013  dimensionedScalar::name() = timeName(timeToUserTime(newTime));
1014  timeIndex_ = newIndex;
1015  fileHandler().setTime(*this);
1016 }
1017 
1018 
1020 {
1021  setEndTime(endTime.value());
1022 }
1023 
1024 
1025 void Foam::Time::setEndTime(const scalar endTime)
1026 {
1027  endTime_ = endTime;
1028 }
1029 
1030 
1032 {
1033  setDeltaT(deltaT.value());
1034 }
1035 
1036 
1037 void Foam::Time::setDeltaT(const scalar deltaT)
1038 {
1039  setDeltaTNoAdjust(deltaT);
1040 
1041  if (writeControl_ == writeControl::adjustableRunTime)
1042  {
1043  adjustDeltaT();
1044  }
1045 }
1046 
1047 
1048 void Foam::Time::setDeltaTNoAdjust(const scalar deltaT)
1049 {
1050  deltaT_ = deltaT;
1051  deltaTchanged_ = true;
1052 }
1053 
1054 
1055 void Foam::Time::setWriteInterval(const scalar writeInterval)
1056 {
1057  if (writeInterval_ == great || !equal(writeInterval, writeInterval_))
1058  {
1059  writeInterval_ = writeInterval;
1060 
1061  if
1062  (
1063  writeControl_ == writeControl::runTime
1064  || writeControl_ == writeControl::adjustableRunTime
1065  )
1066  {
1067  // Recalculate writeTimeIndex_ for consistency with the new
1068  // writeInterval
1069  writeTimeIndex_ = label
1070  (
1071  ((value() - beginTime_) + 0.5*deltaT_)/writeInterval_
1072  );
1073  }
1074  else if (writeControl_ == writeControl::timeStep)
1075  {
1076  // Set to the nearest integer
1077  writeInterval_ = label(writeInterval + 0.5);
1078  }
1079  }
1080 }
1081 
1082 
1084 {
1085  subCycling_ = true;
1086  prevTimeState_.set(new TimeState(*this));
1087 
1088  setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
1089  deltaT_ /= nSubCycles;
1090  deltaT0_ /= nSubCycles;
1091  deltaTSave_ = deltaT0_;
1092 
1093  return prevTimeState();
1094 }
1095 
1096 
1098 {
1099  if (subCycling_)
1100  {
1101  subCycling_ = false;
1102  TimeState::operator=(prevTimeState());
1103  prevTimeState_.clear();
1104  }
1105 }
1106 
1107 
1108 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1109 
1111 {
1112  return operator+=(deltaT.value());
1113 }
1114 
1115 
1116 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
1117 {
1118  setDeltaT(deltaT);
1119  return operator++();
1120 }
1121 
1122 
1124 {
1125  deltaT0_ = deltaTSave_;
1126  deltaTSave_ = deltaT_;
1127 
1128  // Save old time value and name
1129  const scalar oldTimeValue = timeToUserTime(value());
1130  const word oldTimeName = dimensionedScalar::name();
1131 
1132  // Increment time
1133  setTime(value() + deltaT_, timeIndex_ + 1);
1134 
1135  if (!subCycling_)
1136  {
1137  // If the time is very close to zero reset to zero
1138  if (mag(value()) < 10*small*deltaT_)
1139  {
1140  setTime(0, timeIndex_);
1141  }
1142 
1143  if (sigStopAtWriteNow_.active() || sigWriteNow_.active())
1144  {
1145  // A signal might have been sent on one processor only
1146  // Reduce so all decide the same.
1147 
1148  label flag = 0;
1149  if
1150  (
1151  sigStopAtWriteNow_.active()
1152  && stopAt_ == stopAtControl::writeNow
1153  )
1154  {
1155  flag += 1;
1156  }
1157  if (sigWriteNow_.active() && writeOnce_)
1158  {
1159  flag += 2;
1160  }
1161  reduce(flag, maxOp<label>());
1162 
1163  if (flag & 1)
1164  {
1165  stopAt_ = stopAtControl::writeNow;
1166  }
1167  if (flag & 2)
1168  {
1169  writeOnce_ = true;
1170  }
1171  }
1172 
1173  writeTime_ = false;
1174 
1175  switch (writeControl_)
1176  {
1177  case writeControl::timeStep:
1178  writeTime_ = !(timeIndex_ % label(writeInterval_));
1179  break;
1180 
1181  case writeControl::runTime:
1182  case writeControl::adjustableRunTime:
1183  {
1184  label writeIndex = label
1185  (
1186  ((value() - beginTime_) + 0.5*deltaT_)
1187  / writeInterval_
1188  );
1189 
1190  if (writeIndex > writeTimeIndex_)
1191  {
1192  writeTime_ = true;
1193  writeTimeIndex_ = writeIndex;
1194  }
1195  }
1196  break;
1197 
1198  case writeControl::cpuTime:
1199  {
1200  label writeIndex = label
1201  (
1202  returnReduce(elapsedCpuTime(), maxOp<double>())
1203  / writeInterval_
1204  );
1205  if (writeIndex > writeTimeIndex_)
1206  {
1207  writeTime_ = true;
1208  writeTimeIndex_ = writeIndex;
1209  }
1210  }
1211  break;
1212 
1213  case writeControl::clockTime:
1214  {
1215  label writeIndex = label
1216  (
1217  returnReduce(label(elapsedClockTime()), maxOp<label>())
1218  / writeInterval_
1219  );
1220  if (writeIndex > writeTimeIndex_)
1221  {
1222  writeTime_ = true;
1223  writeTimeIndex_ = writeIndex;
1224  }
1225  }
1226  break;
1227  }
1228 
1229 
1230  // Check if endTime needs adjustment to stop at the next run()/end()
1231  if (!end())
1232  {
1233  if (stopAt_ == stopAtControl::noWriteNow)
1234  {
1235  endTime_ = value();
1236  }
1237  else if (stopAt_ == stopAtControl::writeNow)
1238  {
1239  endTime_ = value();
1240  writeTime_ = true;
1241  }
1242  else if (stopAt_ == stopAtControl::nextWrite && writeTime_ == true)
1243  {
1244  endTime_ = value();
1245  }
1246  }
1247 
1248  // Override writeTime if one-shot writing
1249  if (writeOnce_)
1250  {
1251  writeTime_ = true;
1252  writeOnce_ = false;
1253  }
1254 
1255  // Adjust the precision of the time directory name if necessary
1256  if (writeTime_)
1257  {
1258  // User-time equivalent of deltaT
1259  const scalar userDeltaT =
1260  timeToUserTime(value()) - timeToUserTime(value() - deltaT_);
1261 
1262  // Tolerance used when testing time equivalence
1263  const scalar timeTol =
1264  max(min(pow(10.0, -precision_), 0.1*userDeltaT), small);
1265 
1266  // Time value obtained by reading timeName
1267  scalar timeNameValue = -vGreat;
1268 
1269  // Check that new time representation differs from old one
1270  // reinterpretation of the word
1271  if
1272  (
1273  readScalar(dimensionedScalar::name().c_str(), timeNameValue)
1274  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1275  )
1276  {
1277  int oldPrecision = curPrecision_;
1278  while
1279  (
1280  curPrecision_ < maxPrecision_
1281  && readScalar(dimensionedScalar::name().c_str(), timeNameValue)
1282  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1283  )
1284  {
1285  curPrecision_++;
1286  setTime(value(), timeIndex());
1287  }
1288 
1289  if (curPrecision_ != oldPrecision)
1290  {
1292  << "Increased the timePrecision from " << oldPrecision
1293  << " to " << curPrecision_
1294  << " to distinguish between timeNames at time "
1296  << endl;
1297 
1298  if (curPrecision_ == maxPrecision_)
1299  {
1300  // Reached maxPrecision limit
1302  << "Current time name " << dimensionedScalar::name()
1303  << nl
1304  << " The maximum time precision has been reached"
1305  " which might result in overwriting previous"
1306  " results."
1307  << endl;
1308  }
1309 
1310  // Check if round-off error caused time-reversal
1311  scalar oldTimeNameValue = -vGreat;
1312  if
1313  (
1314  readScalar(oldTimeName.c_str(), oldTimeNameValue)
1315  && (
1316  sign(timeNameValue - oldTimeNameValue)
1317  != sign(deltaT_)
1318  )
1319  )
1320  {
1322  << "Current time name " << dimensionedScalar::name()
1323  << " is set to an instance prior to the "
1324  "previous one "
1325  << oldTimeName << nl
1326  << " This might result in temporal "
1327  "discontinuities."
1328  << endl;
1329  }
1330  }
1331  }
1332  }
1333  }
1334 
1335  return *this;
1336 }
1337 
1338 
1340 {
1341  return operator++();
1342 }
1343 
1344 
1345 // ************************************************************************* //
virtual Time & operator++()
Prefix increment,.
Definition: Time.C:1123
dimensionedScalar sign(const dimensionedScalar &ds)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:403
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
A class for addressing time paths without using the Time class.
Definition: TimePaths.H:48
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A class for handling file names.
Definition: fileName.H:79
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:52
const fileName & globalCaseName() const
Return case name.
Definition: argListI.H:54
virtual dimensionedScalar beginTime() const
Return begin time (initial start time)
Definition: Time.C:821
const word & name() const
Name (const access)
Definition: instant.H:126
void adjustDeltaT()
Adjust the time step so that writing occurs at the specified time.
Definition: Time.C:83
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual ~Time()
Destructor.
Definition: Time.C:657
bool writeInfoHeader
const fileName & rootPath() const
Return root path.
Definition: argListI.H:42
static autoPtr< userTime > New(const dictionary &controlDict)
Definition: userTimeNew.C:32
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:922
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:827
bool exists(IOobject &io) const
Does ioobject exist. Is either a directory (empty name()) or.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: Time.C:689
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:833
stopAtControl
Stop-run control options.
Definition: Time.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:116
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
virtual word timeName() const
Return current time name.
Definition: Time.C:676
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
dlLibraryTable libs
Table of loaded dynamic libraries.
virtual bool run() const
Return true if run should continue,.
Definition: Time.C:875
const ParRunControl & parRunControl() const
Return parRunControl.
Definition: argListI.H:60
static const int maxPrecision_
Maximum time directory name precision.
Definition: Time.H:169
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:845
const userTimes::userTime & userTime() const
Return the userTime.
Definition: Time.C:839
virtual TimeState subCycle(const label nSubCycles)
Set time to sub-cycle for the given number of steps.
Definition: Time.C:1083
const dimensionSet dimTime
virtual bool stopAt(const stopAtControl) const
Adjust the current stopAtControl. Note that this value.
Definition: Time.C:941
const Foam::HashTable< string > & options() const
Return options.
Definition: argListI.H:96
static int precision_
Time directory name precision.
Definition: Time.H:161
A class for handling words, derived from string.
Definition: word.H:59
virtual bool running() const
Return true if run should continue without any side effects.
Definition: Time.C:869
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
static int curPrecision_
Current time directory name precision adjusted as necessary.
Definition: Time.H:166
const Type & value() const
Return const reference to value.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
static const word null
An empty word.
Definition: word.H:77
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static const label labelMax
Definition: label.H:62
word timeName
Definition: getTimeIndex.H:3
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:959
format
Supported time directory name formats.
Definition: Time.H:108
const fileOperation & fileHandler()
Get current file handler.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
Time(const word &name, const argList &args, const word &systemName="system", const word &constantName="constant")
Construct given name of dictionary to read and argument list.
Definition: Time.C:417
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:207
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:48
static format format_
Time directory name format.
Definition: Time.H:158
static const char nl
Definition: Ostream.H:260
virtual void setDeltaTNoAdjust(const scalar)
Reset time step without additional adjustment or modification.
Definition: Time.C:1048
defineTypeNameAndDebug(combustionModel, 0)
scalar value() const
Value (const access)
Definition: instant.H:114
static const NamedEnum< stopAtControl, 4 > stopAtControlNames_
Definition: Time.H:128
const word & name() const
Return const reference to name.
bool open(const fileName &libName, const bool verbose=true)
Open the named library, optionally with warnings if problems occur.
void readModifiedObjects()
Read the objects that have been modified.
Definition: TimeIO.C:227
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
rho oldTime()
virtual void setEndTime(const dimensionedScalar &)
Reset end time.
Definition: Time.C:1019
static const NamedEnum< writeControl, 5 > writeControlNames_
Definition: Time.H:131
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
An instant of time. Contains the time value and name.
Definition: instant.H:66
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
virtual bool end() const
Return true if end of run,.
Definition: Time.C:935
virtual void setTime(const Time &) const
Callback for time change.
#define WarningInFunction
Report a warning using Foam::Warning.
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
word userTimeName() const
Return current user time name with units.
Definition: Time.C:863
label timeIndex
Definition: getTimeIndex.H:4
instantList times() const
Search the case for valid time directories.
Definition: Time.C:682
instant findClosestTime(const scalar) const
Search the case for the time closest to the given time.
Definition: Time.C:753
virtual void setDeltaT(const dimensionedScalar &)
Reset time step.
Definition: Time.C:1031
writeControl
Write control options.
Definition: Time.H:89
scalar userTimeToTime(const scalar tau) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: Time.C:851
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void endSubCycle()
Reset time after sub-cycling back to previous TimeState.
Definition: Time.C:1097
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
virtual void setWriteInterval(const scalar writeInterval)
Reset the write interval.
Definition: Time.C:1055
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (e.g. CA deg)
Definition: Time.C:857
bool parRun() const
Definition: parRun.H:77
Registry of regIOobjects.
T & last()
Return the last element of the list.
Definition: UListI.H:128
word findInstancePath(const fileName &path, const instant &) const
Search the case for the time directory path.
Definition: Time.C:719
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:815
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
bool found
dimensionedScalar log10(const dimensionedScalar &ds)
static HashTable< string > validOptions
A list of valid options.
Definition: argList.H:156
virtual Time & operator+=(const dimensionedScalar &)
Set deltaT to that specified and increment time via operator++()
Definition: Time.C:1110
bool equal(const scalar) const
Comparison used for instants to be equal.
Definition: instant.C:59
void setControls()
Set the controls from the current controlDict.
Definition: Time.C:109
Namespace for OpenFOAM.
static label findClosestTimeIndex(const instantList &, const scalar, const word &constantName="constant")
Search instantList for the time index closest to the given time.
Definition: Time.C:790
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1230
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
IOerror FatalIOError