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