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