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