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-2024 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 "PstreamReduceOps.H"
29 #include "argList.H"
30 
31 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
36 
37  template<>
38  const char* Foam::NamedEnum
39  <
41  4
42  >::names[] =
43  {
44  "endTime",
45  "noWriteNow",
46  "writeNow",
47  "nextWrite"
48  };
49 
50  template<>
51  const char* Foam::NamedEnum
52  <
54  5
55  >::names[] =
56  {
57  "timeStep",
58  "runTime",
59  "adjustableRunTime",
60  "clockTime",
61  "cpuTime"
62  };
63 }
64 
67 
70 
72 
74 
76 
77 const int Foam::Time::maxPrecision_(3 - log10(small));
78 
80 
81 
82 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
83 
85 {
86  const scalar timeToNextAction = min
87  (
88  max
89  (
90  0,
91  (writeTimeIndex_ + 1)*writeInterval_ - (value() - beginTime_)
92  ),
93  functionObjects_.timeToNextAction()
94  );
95 
96  const scalar nSteps = timeToNextAction/deltaT_;
97 
98  // Ensure nStepsToNextWrite does not overflow
99  if (nSteps < labelMax)
100  {
101  // Allow the time-step to increase by up to 1%
102  // to accommodate the next write time before splitting
103  const label nStepsToNextWrite = label(max(nSteps, 1) + 0.99);
104 
105  deltaT_ = timeToNextAction/nStepsToNextWrite;
106  }
107 }
108 
109 
111 {
112  // default is to resume calculation from "latestTime"
113  const word startFrom = controlDict_.lookupOrDefault<word>
114  (
115  "startFrom",
116  "latestTime"
117  );
118 
119  if (startFrom == "startTime")
120  {
121  startTime_ = controlDict_.lookup<scalar>("startTime", userUnits());
122  }
123  else
124  {
125  // Search directory for valid time directories
126  const 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 
170  for
171  (
172  curPrecision_ = maxPrecision_;
173  curPrecision_ > oldPrecision;
174  curPrecision_--
175  )
176  {
177  // Update the time formatting
178  setTime(startTime_, 0);
179 
180  // Check the existence of the time directory with the new format
181  if (fileHandler().exists(timePath(), false, false))
182  {
183  requiredPrecision = curPrecision_;
184  }
185  }
186 
187  if (requiredPrecision > 0)
188  {
189  // Update the time precision
190  curPrecision_ = requiredPrecision;
191 
192  // Update the time formatting
193  setTime(startTime_, 0);
194 
196  << "Increasing the timePrecision from " << oldPrecision
197  << " to " << curPrecision_
198  << " to support the formatting of the current time directory "
199  << name() << nl << endl;
200  }
201  else
202  {
203  // Could not find time directory so assume it is not present
204  curPrecision_ = oldPrecision;
205 
206  // Revert the time formatting
207  setTime(startTime_, 0);
208  }
209  }
210 
211  if (Pstream::parRun())
212  {
213  scalar sumStartTime = startTime_;
214  reduce(sumStartTime, sumOp<scalar>());
215  if
216  (
217  mag(Pstream::nProcs()*startTime_ - sumStartTime)
218  > Pstream::nProcs()*deltaT_/10.0
219  )
220  {
221  FatalIOErrorInFunction(controlDict_)
222  << "Start time is not the same for all processors" << nl
223  << "processor " << Pstream::myProcNo() << " has startTime "
224  << startTime_ << exit(FatalIOError);
225  }
226  }
227 
228  timeIOdictionary timeDict
229  (
230  IOobject
231  (
232  "time",
233  name(),
234  "uniform",
235  *this,
238  false
239  )
240  );
241 
242  if (controlDict_.found("beginTime"))
243  {
244  beginTime_ = controlDict_.lookup<scalar>("beginTime", userUnits());
245  }
246  else if (timeDict.found("beginTime"))
247  {
248  beginTime_ = timeDict.lookup<scalar>("beginTime", userUnits());
249  }
250  else
251  {
252  beginTime_ = startTime_;
253  }
254 
255  // Read and set the deltaT only if time-step adjustment is active
256  // otherwise use the deltaT from the controlDict
257  if (controlDict_.lookupOrDefault<Switch>("adjustTimeStep", false))
258  {
259  if (timeDict.found("deltaT"))
260  {
261  deltaT_ = timeDict.lookup<scalar>("deltaT", userUnits());
262  deltaTSave_ = deltaT_;
263  deltaT0_ = deltaT_;
264  }
265  }
266 
267  if (timeDict.found("deltaT0"))
268  {
269  deltaT0_ = timeDict.lookup<scalar>("deltaT0", userUnits());
270  }
271 
272  if (timeDict.readIfPresent("index", startTimeIndex_))
273  {
274  timeIndex_ = startTimeIndex_;
275  }
276 
277  // Set writeTimeIndex_ to correspond to beginTime_
278  if
279  (
280  (
281  writeControl_ == writeControl::runTime
282  || writeControl_ == writeControl::adjustableRunTime
283  )
284  )
285  {
286  writeTimeIndex_ = label
287  (
288  ((value() - beginTime_) + 0.5*deltaT_)/writeInterval_
289  );
290  }
291 
292 
293  // Check if values stored in time dictionary are consistent
294 
295  // 1. Based on time name
296  bool checkValue = true;
297 
298  string storedTimeName;
299  if (timeDict.readIfPresent("name", storedTimeName))
300  {
301  if (storedTimeName == name())
302  {
303  // Same time. No need to check stored value
304  checkValue = false;
305  }
306  }
307 
308  // 2. Based on time value
309  // (consistent up to the current time writing precision so it won't
310  // trigger if we just change the write precision)
311  if (checkValue)
312  {
313  scalar storedTimeValue;
314  if (timeDict.readIfPresent("value", storedTimeValue))
315  {
316  word storedTimeName(timeName(storedTimeValue));
317 
318  if (storedTimeName != name())
319  {
320  IOWarningInFunction(timeDict)
321  << "Time read from time dictionary " << storedTimeName
322  << " differs from actual time " << name() << '.' << nl
323  << " This may cause unexpected database behaviour."
324  << " If you are not interested" << nl
325  << " in preserving time state delete"
326  << " the time dictionary."
327  << endl;
328  }
329  }
330  }
331 }
332 
333 
334 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
335 
337 (
338  const word& controlDictName,
339  const argList& args,
340  const bool enableFunctionObjects
341 )
342 :
343  TimePaths
344  (
345  args.parRunControl().parRun(),
346  args.rootPath(),
347  args.globalCaseName(),
348  args.caseName()
349  ),
350 
351  objectRegistry(*this),
352 
353  runTimeModifiable_(false),
354 
355  controlDict_
356  (
357  IOobject
358  (
359  controlDictName,
360  system(),
361  *this,
362  IOobject::MUST_READ_IF_MODIFIED,
363  IOobject::NO_WRITE
364  ),
365  *this
366  ),
367 
368  startTimeIndex_(0),
369  startTime_(0),
370  endTime_(0),
371  beginTime_(startTime_),
372 
373  userTime_(userTimes::userTime::New(controlDict_)),
374 
375  stopAt_(stopAtControl::endTime),
376  writeControl_(writeControl::timeStep),
377  writeInterval_(great),
378  purgeWrite_(0),
379  writeOnce_(false),
380 
381  subCycling_(false),
382 
383  sigWriteNow_(writeInfoHeader, *this),
384  sigStopAtWriteNow_(writeInfoHeader, *this),
385 
386  writeFormat_(IOstream::ASCII),
387  writeVersion_(IOstream::currentVersion),
388  writeCompression_(IOstream::UNCOMPRESSED),
389  cacheTemporaryObjects_(true),
390 
391  functionObjects_
392  (
393  *this,
394  enableFunctionObjects
395  ? argList::validOptions.found("functionObjects")
396  ? args.optionFound("functionObjects")
397  : !args.optionFound("noFunctionObjects")
398  : false
399  )
400 {
401  libs.open(controlDict_, "libs");
402 
403  // Explicitly set read flags on objectRegistry so anything constructed
404  // from it reads as well (e.g. fvSolution).
406 
407  if (args.options().found("case"))
408  {
409  const wordList switchSets
410  (
411  {
412  "InfoSwitches",
413  "OptimisationSwitches",
414  "DebugSwitches",
415  "DimensionedConstants",
416  "DimensionSets",
417  "UnitConversions"
418  }
419  );
420 
421  forAll(switchSets, i)
422  {
423  if (controlDict_.found(switchSets[i]))
424  {
425  IOWarningInFunction(controlDict_)
426  << switchSets[i]
427  << " in system/controlDict are only processed if "
428  << args.executable() << " is run in the "
429  << args.path() << " directory" << endl;
430  }
431  }
432  }
433 
434  setControls();
435 
436  // Add a watch on the controlDict and functions files
437  // after runTimeModifiable_ is set
438  controlDict_.addWatch();
439  functionObjects_.addWatch();
440 }
441 
442 
444 (
445  const word& controlDictName,
446  const fileName& rootPath,
447  const fileName& caseName,
448  const bool enableFunctionObjects
449 )
450 :
451  TimePaths(rootPath, caseName),
452 
453  objectRegistry(*this),
454 
455  runTimeModifiable_(false),
456 
457  controlDict_
458  (
459  IOobject
460  (
461  controlDictName,
462  system(),
463  *this,
464  IOobject::MUST_READ_IF_MODIFIED,
465  IOobject::NO_WRITE
466  ),
467  *this
468  ),
469 
470  startTimeIndex_(0),
471  startTime_(0),
472  endTime_(0),
473  beginTime_(startTime_),
474 
475  userTime_(userTimes::userTime::New(controlDict_)),
476 
477  stopAt_(stopAtControl::endTime),
478  writeControl_(writeControl::timeStep),
479  writeInterval_(great),
480  purgeWrite_(0),
481  writeOnce_(false),
482 
483  subCycling_(false),
484 
485  sigWriteNow_(writeInfoHeader, *this),
486  sigStopAtWriteNow_(writeInfoHeader, *this),
487 
488  writeFormat_(IOstream::ASCII),
489  writeVersion_(IOstream::currentVersion),
490  writeCompression_(IOstream::UNCOMPRESSED),
491  cacheTemporaryObjects_(true),
492 
493  functionObjects_(*this, enableFunctionObjects)
494 {
495  libs.open(controlDict_, "libs");
496 
497  // Explicitly set read flags on objectRegistry so anything constructed
498  // from it reads as well (e.g. fvSolution).
500 
501  setControls();
502 
503  // Add a watch on the controlDict and functions files
504  // after runTimeModifiable_ is set
505  controlDict_.addWatch();
506  functionObjects_.addWatch();
507 }
508 
509 
511 (
512  const dictionary& dict,
513  const fileName& rootPath,
514  const fileName& caseName,
515  const bool enableFunctionObjects
516 )
517 :
518  TimePaths(rootPath, caseName),
519 
520  objectRegistry(*this),
521 
522  runTimeModifiable_(false),
523 
524  controlDict_
525  (
526  IOobject
527  (
528  controlDictName,
529  system(),
530  *this,
531  IOobject::MUST_READ_IF_MODIFIED,
532  IOobject::NO_WRITE
533  ),
534  dict,
535  *this
536  ),
537 
538  startTimeIndex_(0),
539  startTime_(0),
540  endTime_(0),
541  beginTime_(startTime_),
542 
543  userTime_(userTimes::userTime::New(controlDict_)),
544 
545  stopAt_(stopAtControl::endTime),
546  writeControl_(writeControl::timeStep),
547  writeInterval_(great),
548  purgeWrite_(0),
549  writeOnce_(false),
550 
551  subCycling_(false),
552 
553  sigWriteNow_(writeInfoHeader, *this),
554  sigStopAtWriteNow_(writeInfoHeader, *this),
555 
556  writeFormat_(IOstream::ASCII),
557  writeVersion_(IOstream::currentVersion),
558  writeCompression_(IOstream::UNCOMPRESSED),
559  cacheTemporaryObjects_(true),
560 
561  functionObjects_(*this, enableFunctionObjects)
562 {
563  libs.open(controlDict_, "libs");
564 
565  // Explicitly set read flags on objectRegistry so anything constructed
566  // from it reads as well (e.g. fvSolution).
568 
569  setControls();
570 
571  // Add a watch on the controlDict and functions files
572  // after runTimeModifiable_ is set
573  controlDict_.addWatch();
574  functionObjects_.addWatch();
575 }
576 
577 
579 (
580  const fileName& rootPath,
581  const fileName& caseName,
582  const bool enableFunctionObjects
583 )
584 :
585  TimePaths(rootPath, caseName),
586 
587  objectRegistry(*this),
588 
589  runTimeModifiable_(false),
590 
591  controlDict_
592  (
593  IOobject
594  (
595  controlDictName,
596  system(),
597  *this,
598  IOobject::NO_READ,
599  IOobject::NO_WRITE
600  ),
601  *this
602  ),
603 
604  startTimeIndex_(0),
605  startTime_(0),
606  endTime_(0),
607  beginTime_(startTime_),
608 
609  userTime_(userTimes::userTime::New(controlDict_)),
610 
611  stopAt_(stopAtControl::endTime),
612  writeControl_(writeControl::timeStep),
613  writeInterval_(great),
614  purgeWrite_(0),
615  writeOnce_(false),
616 
617  subCycling_(false),
618 
619  writeFormat_(IOstream::ASCII),
620  writeVersion_(IOstream::currentVersion),
621  writeCompression_(IOstream::UNCOMPRESSED),
622  cacheTemporaryObjects_(true),
623 
624  functionObjects_(*this, enableFunctionObjects)
625 {
626  libs.open(controlDict_, "libs");
627 }
628 
629 
630 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
631 
633 {
634  // Destroy function objects first
635  functionObjects_.clear();
636 }
637 
638 
639 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
640 
641 Foam::word Foam::Time::timeName(const scalar t, const int precision)
642 {
643  std::ostringstream buf;
644  buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
645  buf.precision(precision);
646  buf << t;
647  return buf.str();
648 }
649 
650 
652 {
653  return findTimes(path(), constant());
654 }
655 
656 
658 (
659  const fileName& dir,
660  const word& name,
661  const IOobject::readOption rOpt,
662  const word& stopInstance
663 ) const
664 {
665  IOobject startIO
666  (
667  name, // name might be empty!
669  dir,
670  *this,
671  rOpt
672  );
673 
674  IOobject io
675  (
676  fileHandler().findInstance
677  (
678  startIO,
679  userTimeValue(),
680  stopInstance
681  )
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  const 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 
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  const 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
794  (
795  timeName(timeToUserTime(beginTime_)),
796  dimTime,
797  beginTime_
798  );
799 }
800 
801 
803 {
804  return dimensionedScalar
805  (
806  timeName(timeToUserTime(startTime_)),
807  dimTime,
808  startTime_
809  );
810 }
811 
812 
814 {
815  return dimensionedScalar
816  (
817  timeName(timeToUserTime(endTime_)),
818  dimTime,
819  endTime_
820  );
821 }
822 
823 
825 {
826  return *userTime_;
827 }
828 
829 
830 Foam::scalar Foam::Time::userTimeValue() const
831 {
832  return userTime_->timeToUserTime(value());
833 }
834 
835 
836 Foam::scalar Foam::Time::userDeltaTValue() const
837 {
838  return
839  userTime_->timeToUserTime(value())
840  - userTime_->timeToUserTime(value() - deltaT_);
841 }
842 
843 
844 Foam::scalar Foam::Time::userTimeToTime(const scalar tau) const
845 {
846  return userTime_->userTimeToTime(tau);
847 }
848 
849 
850 Foam::scalar Foam::Time::timeToUserTime(const scalar t) const
851 {
852  return userTime_->timeToUserTime(t);
853 }
854 
855 
857 {
858  if (name() == constant())
859  {
860  return constant();
861  }
862  else
863  {
864  return timeName(userTimeValue()) + userTime_->unitName();
865  }
866 }
867 
868 
870 {
871  return userTime_->units();
872 }
873 
874 
876 {
877  static const unitConversion unitSeconds(dimTime);
878 
879  switch (writeControl_)
880  {
881  case writeControl::timeStep:
882  return unitless;
883  case writeControl::runTime:
884  case writeControl::adjustableRunTime:
885  return userUnits();
886  case writeControl::cpuTime:
887  case writeControl::clockTime:
888  return unitSeconds;
889  }
890 
891  return unitNone;
892 }
893 
894 
896 {
897  return value() < (endTime_ - 0.5*deltaT_);
898 }
899 
900 
901 bool Foam::Time::run() const
902 {
903  bool running = this->running();
904 
905  if (!subCycling_)
906  {
907  if (!running && timeIndex_ != startTimeIndex_)
908  {
909  functionObjects_.execute();
910  functionObjects_.end();
911 
912  if (cacheTemporaryObjects_)
913  {
914  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
915  }
916  }
917  }
918 
919  if (running)
920  {
921  if (!subCycling_)
922  {
923  const_cast<Time&>(*this).readModifiedObjects();
924 
925  if (timeIndex_ == startTimeIndex_)
926  {
927  functionObjects_.start();
928  }
929  else
930  {
931  functionObjects_.execute();
932 
933  if (cacheTemporaryObjects_)
934  {
935  cacheTemporaryObjects_ = checkCacheTemporaryObjects();
936  }
937  }
938  }
939 
940  // Re-evaluate if running in case a function object has changed things
941  running = this->running();
942  }
943 
944  return running;
945 }
946 
947 
949 {
950  bool running = run();
951 
952  if (running)
953  {
954  operator++();
955  }
956 
957  return running;
958 }
959 
960 
961 bool Foam::Time::end() const
962 {
963  return value() > (endTime_ + 0.5*deltaT_);
964 }
965 
966 
967 bool Foam::Time::stopAt(const stopAtControl sa) const
968 {
969  const bool changed = (stopAt_ != sa);
970  stopAt_ = sa;
971 
972  // adjust endTime
973  if (sa == stopAtControl::endTime)
974  {
975  controlDict_.lookup("endTime") >> endTime_;
976  }
977  else
978  {
979  endTime_ = great;
980  }
981  return changed;
982 }
983 
984 
985 void Foam::Time::setTime(const Time& t)
986 {
987  value() = t.value();
988  dimensionedScalar::name() = t.dimensionedScalar::name();
989  timeIndex_ = t.timeIndex_;
990  fileHandler().setTime(*this);
991 }
992 
993 
994 void Foam::Time::setTime(const instant& inst, const label newIndex)
995 {
996  value() = userTimeToTime(inst.value());
997  dimensionedScalar::name() = inst.name();
998  timeIndex_ = newIndex;
999 
1000  timeIOdictionary timeDict
1001  (
1002  IOobject
1003  (
1004  "time",
1005  name(),
1006  "uniform",
1007  *this,
1010  false
1011  )
1012  );
1013 
1014  if (timeDict.found("deltaT"))
1015  {
1016  deltaT_ = timeDict.lookup<scalar>("deltaT", userUnits());
1017  }
1018 
1019  if (timeDict.found("deltaT0"))
1020  {
1021  deltaT0_ = timeDict.lookup<scalar>("deltaT0", userUnits());
1022  }
1023 
1024  timeDict.readIfPresent("index", timeIndex_);
1025 
1026  fileHandler().setTime(*this);
1027 }
1028 
1029 
1030 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
1031 {
1032  setTime(newTime.value(), newIndex);
1033 }
1034 
1035 
1036 void Foam::Time::setTime(const scalar newTime, const label newIndex)
1037 {
1038  value() = newTime;
1039  dimensionedScalar::name() = timeName(timeToUserTime(newTime));
1040  timeIndex_ = newIndex;
1041  fileHandler().setTime(*this);
1042 }
1043 
1044 
1046 {
1047  setEndTime(endTime.value());
1048 }
1049 
1050 
1051 void Foam::Time::setEndTime(const scalar endTime)
1052 {
1053  endTime_ = endTime;
1054 }
1055 
1056 
1058 {
1059  setDeltaT(deltaT.value());
1060 }
1061 
1062 
1063 void Foam::Time::setDeltaT(const scalar deltaT)
1064 {
1065  setDeltaTNoAdjust(deltaT);
1066 
1067  if (writeControl_ == writeControl::adjustableRunTime)
1068  {
1069  adjustDeltaT();
1070  }
1071 }
1072 
1073 
1074 void Foam::Time::setDeltaTNoAdjust(const scalar deltaT)
1075 {
1076  deltaT_ = deltaT;
1077  deltaTchanged_ = true;
1078 }
1079 
1080 
1081 void Foam::Time::setWriteInterval(const scalar writeInterval)
1082 {
1083  if (writeInterval_ == great || !equal(writeInterval, writeInterval_))
1084  {
1085  writeInterval_ = writeInterval;
1086 
1087  if
1088  (
1089  writeControl_ == writeControl::runTime
1090  || writeControl_ == writeControl::adjustableRunTime
1091  )
1092  {
1093  // Recalculate writeTimeIndex_ for consistency with the new
1094  // writeInterval
1095  writeTimeIndex_ = label
1096  (
1097  ((value() - beginTime_) + 0.5*deltaT_)/writeInterval_
1098  );
1099  }
1100  else if (writeControl_ == writeControl::timeStep)
1101  {
1102  // Set to the nearest integer
1103  writeInterval_ = label(writeInterval + 0.5);
1104  }
1105  }
1106 }
1107 
1108 
1110 {
1111  subCycling_ = true;
1112  prevTimeState_.set(new TimeState(*this));
1113 
1114  setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
1115  deltaT_ /= nSubCycles;
1116  deltaT0_ /= nSubCycles;
1117  deltaTSave_ = deltaT0_;
1118 
1119  return prevTimeState();
1120 }
1121 
1122 
1124 {
1125  if (subCycling_)
1126  {
1127  subCycling_ = false;
1128  TimeState::operator=(prevTimeState());
1129  prevTimeState_.clear();
1130  }
1131 }
1132 
1133 
1134 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1135 
1137 {
1138  return operator+=(deltaT.value());
1139 }
1140 
1141 
1142 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
1143 {
1144  setDeltaT(deltaT);
1145  return operator++();
1146 }
1147 
1148 
1150 {
1151  deltaT0_ = deltaTSave_;
1152  deltaTSave_ = deltaT_;
1153 
1154  // Save old time value and name
1155  const scalar oldTimeValue = timeToUserTime(value());
1156  const word oldTimeName = dimensionedScalar::name();
1157 
1158  // Increment time
1159  setTime(value() + deltaT_, timeIndex_ + 1);
1160 
1161  if (!subCycling_)
1162  {
1163  // If the time is very close to zero reset to zero
1164  if (mag(value()) < 10*small*deltaT_)
1165  {
1166  setTime(0, timeIndex_);
1167  }
1168 
1169  if (sigStopAtWriteNow_.active() || sigWriteNow_.active())
1170  {
1171  // A signal might have been sent on one processor only
1172  // Reduce so all decide the same.
1173 
1174  label flag = 0;
1175  if
1176  (
1177  sigStopAtWriteNow_.active()
1178  && stopAt_ == stopAtControl::writeNow
1179  )
1180  {
1181  flag += 1;
1182  }
1183  if (sigWriteNow_.active() && writeOnce_)
1184  {
1185  flag += 2;
1186  }
1187  reduce(flag, maxOp<label>());
1188 
1189  if (flag & 1)
1190  {
1191  stopAt_ = stopAtControl::writeNow;
1192  }
1193  if (flag & 2)
1194  {
1195  writeOnce_ = true;
1196  }
1197  }
1198 
1199  writeTime_ = false;
1200 
1201  switch (writeControl_)
1202  {
1203  case writeControl::timeStep:
1204  writeTime_ = !(timeIndex_ % label(writeInterval_));
1205  break;
1206 
1207  case writeControl::runTime:
1208  case writeControl::adjustableRunTime:
1209  {
1210  label writeIndex = label
1211  (
1212  ((value() - beginTime_) + 0.5*deltaT_)
1213  / writeInterval_
1214  );
1215 
1216  if (writeIndex > writeTimeIndex_)
1217  {
1218  writeTime_ = true;
1219  writeTimeIndex_ = writeIndex;
1220  }
1221  }
1222  break;
1223 
1224  case writeControl::cpuTime:
1225  {
1226  label writeIndex = label
1227  (
1228  returnReduce(elapsedCpuTime(), maxOp<double>())
1229  / writeInterval_
1230  );
1231  if (writeIndex > writeTimeIndex_)
1232  {
1233  writeTime_ = true;
1234  writeTimeIndex_ = writeIndex;
1235  }
1236  }
1237  break;
1238 
1239  case writeControl::clockTime:
1240  {
1241  label writeIndex = label
1242  (
1243  returnReduce(label(elapsedClockTime()), maxOp<label>())
1244  / writeInterval_
1245  );
1246  if (writeIndex > writeTimeIndex_)
1247  {
1248  writeTime_ = true;
1249  writeTimeIndex_ = writeIndex;
1250  }
1251  }
1252  break;
1253  }
1254 
1255 
1256  // Check if endTime needs adjustment to stop at the next run()/end()
1257  if (!end())
1258  {
1259  if (stopAt_ == stopAtControl::noWriteNow)
1260  {
1261  endTime_ = value();
1262  }
1263  else if (stopAt_ == stopAtControl::writeNow)
1264  {
1265  endTime_ = value();
1266  writeTime_ = true;
1267  }
1268  else if (stopAt_ == stopAtControl::nextWrite && writeTime_ == true)
1269  {
1270  endTime_ = value();
1271  }
1272  }
1273 
1274  // Override writeTime if one-shot writing
1275  if (writeOnce_)
1276  {
1277  writeTime_ = true;
1278  writeOnce_ = false;
1279  }
1280 
1281  // Adjust the precision of the time name if necessary
1282  {
1283  // User-time equivalent of deltaT
1284  const scalar userDeltaT = userDeltaTValue();
1285 
1286  // Tolerance used when testing time equivalence
1287  const scalar timeTol =
1288  max(min(pow(10.0, -precision_), 0.1*userDeltaT), small);
1289 
1290  // Time value obtained by reading timeName
1291  scalar timeNameValue = -vGreat;
1292 
1293  // Check that new time representation differs from old one
1294  // reinterpretation of the word
1295  if
1296  (
1297  readScalar(dimensionedScalar::name().c_str(), timeNameValue)
1298  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1299  )
1300  {
1301  int oldPrecision = curPrecision_;
1302  while
1303  (
1304  curPrecision_ < maxPrecision_
1305  && readScalar(dimensionedScalar::name().c_str(), timeNameValue)
1306  && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol)
1307  )
1308  {
1309  curPrecision_++;
1310  setTime(value(), timeIndex());
1311  }
1312 
1313  if (curPrecision_ != oldPrecision)
1314  {
1316  << "Increased the timePrecision from " << oldPrecision
1317  << " to " << curPrecision_
1318  << " to distinguish between timeNames at time "
1320  << endl;
1321 
1322  if (curPrecision_ == maxPrecision_)
1323  {
1324  // Reached maxPrecision limit
1326  << "Current time name " << dimensionedScalar::name()
1327  << nl
1328  << " The maximum time precision has been reached"
1329  " which might result in overwriting previous"
1330  " results."
1331  << exit(FatalError);
1332  }
1333 
1334  // Check if round-off error caused time-reversal
1335  scalar oldTimeNameValue = -vGreat;
1336  if
1337  (
1338  readScalar(oldTimeName.c_str(), oldTimeNameValue)
1339  && (
1340  sign(timeNameValue - oldTimeNameValue)
1341  != sign(deltaT_)
1342  )
1343  )
1344  {
1346  << "Current time name " << dimensionedScalar::name()
1347  << " is set to an instance prior to the "
1348  "previous one "
1349  << oldTimeName << nl
1350  << " This might result in temporal "
1351  "discontinuities."
1352  << endl;
1353  }
1354  }
1355  }
1356  }
1357  }
1358 
1359  return *this;
1360 }
1361 
1362 
1364 {
1365  return operator++();
1366 }
1367 
1368 
1369 // ************************************************************************* //
Inter-processor communication reduction functions.
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
readOption readOpt() const
Definition: IOobject.H:360
An IOstream is an abstract base class for all input/output systems; be they streams,...
Definition: IOstream.H:72
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A class for addressing time paths without using the Time class.
Definition: TimePaths.H:49
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:55
friend class Time
Declare friendship with the Time class.
Definition: TimeState.H:70
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
virtual bool run() const
Return true if run should continue,.
Definition: Time.C:901
static int precision_
Time directory name precision.
Definition: Time.H:165
scalar timeToUserTime(const scalar t) const
Convert the real-time (s) into user-time (e.g. CA deg)
Definition: Time.C:850
instantList times() const
Search the case for valid time directories.
Definition: Time.C:651
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:985
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:869
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:802
static const NamedEnum< stopAtControl, 4 > stopAtControlNames
Definition: Time.H:117
const unitConversion & writeIntervalUnits() const
Return the write interval units.
Definition: Time.C:875
static int curPrecision_
Current time directory name precision adjusted as necessary.
Definition: Time.H:170
void adjustDeltaT()
Adjust the time step so that writing occurs at the specified time.
Definition: Time.C:84
format
Supported time directory name formats.
Definition: Time.H:111
scalar userTimeToTime(const scalar tau) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: Time.C:844
scalar writeInterval_
Definition: Time.H:138
virtual dimensionedScalar beginTime() const
Return begin time (initial start time)
Definition: Time.C:791
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:830
scalar userDeltaTValue() const
Return user time step value.
Definition: Time.C:836
static const int maxPrecision_
Maximum time directory name precision.
Definition: Time.H:173
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:658
virtual void setWriteInterval(const scalar writeInterval)
Reset the write interval.
Definition: Time.C:1081
virtual void setEndTime(const dimensionedScalar &)
Reset end time.
Definition: Time.C:1045
word userTimeName() const
Return current user time name with units.
Definition: Time.C:856
virtual bool running() const
Return true if run should continue without any side effects.
Definition: Time.C:895
void setControls()
Set the controls from the current controlDict.
Definition: Time.C:110
word findInstancePath(const fileName &path, const instant &) const
Search the case for the time directory path.
Definition: Time.C:689
virtual void setDeltaT(const dimensionedScalar &)
Reset time step.
Definition: Time.C:1057
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:208
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:641
stopAtControl
Stop-run control options.
Definition: Time.H:102
static format format_
Time directory name format.
Definition: Time.H:162
virtual bool stopAt(const stopAtControl) const
Adjust the current stopAtControl. Note that this value.
Definition: Time.C:967
const userTimes::userTime & userTime() const
Return the userTime.
Definition: Time.C:824
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:785
virtual void endSubCycle()
Reset time after sub-cycling back to previous TimeState.
Definition: Time.C:1123
virtual void setDeltaTNoAdjust(const scalar)
Reset time step without additional adjustment or modification.
Definition: Time.C:1074
virtual Time & operator++()
Prefix increment,.
Definition: Time.C:1149
static const NamedEnum< writeControl, 5 > writeControlNames
Definition: Time.H:119
void readModifiedObjects()
Read the objects that have been modified.
Definition: TimeIO.C:200
scalar beginTime_
Definition: Time.H:129
virtual TimeState subCycle(const label nSubCycles)
Set time to sub-cycle for the given number of steps.
Definition: Time.C:1109
writeControl
Write control options.
Definition: Time.H:92
virtual bool loop()
Return true if run should continue and if so increment time.
Definition: Time.C:948
virtual dimensionedScalar endTime() const
Return end time.
Definition: Time.C:813
virtual ~Time()
Destructor.
Definition: Time.C:632
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
virtual Time & operator+=(const dimensionedScalar &)
Set deltaT to that specified and increment time via operator++()
Definition: Time.C:1136
virtual bool end() const
Return true if end of run,.
Definition: Time.C:961
instant findClosestTime(const scalar) const
Search the case for the time closest to the given time.
Definition: Time.C:723
T & last()
Return the last element of the list.
Definition: UListI.H:128
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:103
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
const Foam::HashTable< string > & options() const
Return options.
Definition: argListI.H:96
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
const scalar & value() const
Return const reference to value.
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.
A class for handling file names.
Definition: fileName.H:82
virtual void setTime(const Time &) const
Callback for time change.
scalar timeToNextAction()
Return the time to the next write.
An instant of time. Contains the time value and name.
Definition: instant.H:67
scalar value() const
Value (const access)
Definition: instant.H:114
const word & name() const
Name (const access)
Definition: instant.H:126
bool equal(const scalar) const
Comparison used for instants to be equal.
Definition: instant.C:59
Registry of regIOobjects.
void addWatch()
Add file watch on object (if registered and READ_IF_MODIFIED)
Definition: regIOobject.C:259
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
An abstract class for the user time description.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
static instantList timeDirs
Definition: globalFoam.H:44
runTimeSource setTime(sourceTimes[sourceTimeIndex], sourceTimeIndex)
word timeName
Definition: getTimeIndex.H:3
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dlLibraryTable libs
Table of loaded dynamic libraries.
void adjustDeltaT(Time &runTime, const PtrList< solver > &solvers)
Adjust the time-step according to the solver maxDeltaT.
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1230
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
void setDeltaT(Time &runTime, const PtrList< solver > &solvers)
Set the initial time-step according to the solver maxDeltaT.
dimensionedScalar sign(const dimensionedScalar &ds)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
bool writeInfoHeader
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
dimensionedScalar log10(const dimensionedScalar &ds)
const unitConversion unitNone
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:409
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
static const label labelMax
Definition: label.H:62
const unitConversion unitless
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
Foam::argList args(argc, argv)