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