ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / OpenFOAM / db / Time / Time.C
blob4068357705a37717d2109b6a7d3a4969e6d38e29
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "Time.H"
27 #include "PstreamReduceOps.H"
28 #include "argList.H"
30 #include <sstream>
32 // * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(Foam::Time, 0);
36 namespace Foam
38     template<>
39     const char* Foam::NamedEnum
40     <
41         Foam::Time::stopAtControls,
42         4
43     >::names[] =
44     {
45         "endTime",
46         "noWriteNow",
47         "writeNow",
48         "nextWrite"
49     };
51     template<>
52     const char* Foam::NamedEnum
53     <
54         Foam::Time::writeControls,
55         5
56     >::names[] =
57     {
58         "timeStep",
59         "runTime",
60         "adjustableRunTime",
61         "clockTime",
62         "cpuTime"
63     };
66 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
67     Foam::Time::stopAtControlNames_;
69 const Foam::NamedEnum<Foam::Time::writeControls, 5>
70     Foam::Time::writeControlNames_;
72 Foam::Time::fmtflags Foam::Time::format_(Foam::Time::general);
73 int Foam::Time::precision_(6);
75 Foam::word Foam::Time::controlDictName("controlDict");
78 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
80 void Foam::Time::adjustDeltaT()
82     if (writeControl_ == wcAdjustableRunTime)
83     {
84         scalar timeToNextWrite = max
85         (
86             0.0,
87             (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
88         );
90         scalar nSteps = timeToNextWrite/deltaT_ - SMALL;
92         // For tiny deltaT the label can overflow!
93         if (nSteps < labelMax)
94         {
95             label nStepsToNextWrite = label(nSteps) + 1;
97             scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
99             // Control the increase of the time step to within a factor of 2
100             // and the decrease within a factor of 5.
101             if (newDeltaT >= deltaT_)
102             {
103                 deltaT_ = min(newDeltaT, 2.0*deltaT_);
104             }
105             else
106             {
107                 deltaT_ = max(newDeltaT, 0.2*deltaT_);
108             }
109         }
110     }
114 void Foam::Time::setControls()
116     // default is to resume calculation from "latestTime"
117     const word startFrom = controlDict_.lookupOrDefault<word>
118     (
119         "startFrom",
120         "latestTime"
121     );
123     if (startFrom == "startTime")
124     {
125         controlDict_.lookup("startTime") >> startTime_;
126     }
127     else
128     {
129         // Search directory for valid time directories
130         instantList timeDirs = findTimes(path());
132         if (startFrom == "firstTime")
133         {
134             if (timeDirs.size())
135             {
136                 if (timeDirs[0].name() == constant() && timeDirs.size() >= 2)
137                 {
138                     startTime_ = timeDirs[1].value();
139                 }
140                 else
141                 {
142                     startTime_ = timeDirs[0].value();
143                 }
144             }
145         }
146         else if (startFrom == "latestTime")
147         {
148             if (timeDirs.size())
149             {
150                 startTime_ = timeDirs.last().value();
151             }
152         }
153         else
154         {
155             FatalIOErrorIn("Time::setControls()", controlDict_)
156                 << "expected startTime, firstTime or latestTime"
157                 << " found '" << startFrom << "'"
158                 << exit(FatalIOError);
159         }
160     }
162     setTime(startTime_, 0);
164     readDict();
165     deltaTSave_ = deltaT_;
166     deltaT0_ = deltaT_;
168     if (Pstream::parRun())
169     {
170         scalar sumStartTime = startTime_;
171         reduce(sumStartTime, sumOp<scalar>());
172         if
173         (
174             mag(Pstream::nProcs()*startTime_ - sumStartTime)
175           > Pstream::nProcs()*deltaT_/10.0
176         )
177         {
178             FatalIOErrorIn("Time::setControls()", controlDict_)
179                 << "Start time is not the same for all processors" << nl
180                 << "processor " << Pstream::myProcNo() << " has startTime "
181                 << startTime_ << exit(FatalIOError);
182         }
183     }
185     IOdictionary timeDict
186     (
187         IOobject
188         (
189             "time",
190             timeName(),
191             "uniform",
192             *this,
193             IOobject::READ_IF_PRESENT,
194             IOobject::NO_WRITE,
195             false
196         )
197     );
199     if (timeDict.readIfPresent("deltaT", deltaT_))
200     {
201         deltaTSave_ = deltaT_;
202         deltaT0_ = deltaT_;
203     }
205     timeDict.readIfPresent("deltaT0", deltaT0_);
207     if (timeDict.readIfPresent("index", startTimeIndex_))
208     {
209         timeIndex_ = startTimeIndex_;
210     }
214 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
216 Foam::Time::Time
218     const word& controlDictName,
219     const fileName& rootPath,
220     const fileName& caseName,
221     const word& systemName,
222     const word& constantName,
223     const bool enableFunctionObjects
226     TimePaths
227     (
228         rootPath,
229         caseName,
230         systemName,
231         constantName
232     ),
234     objectRegistry(*this),
236     libs_(),
238     controlDict_
239     (
240         IOobject
241         (
242             controlDictName,
243             system(),
244             *this,
245             IOobject::MUST_READ_IF_MODIFIED,
246             IOobject::NO_WRITE,
247             false
248         )
249     ),
251     startTimeIndex_(0),
252     startTime_(0),
253     endTime_(0),
255     stopAt_(saEndTime),
256     writeControl_(wcTimeStep),
257     writeInterval_(GREAT),
258     purgeWrite_(0),
259     subCycling_(false),
261     writeFormat_(IOstream::ASCII),
262     writeVersion_(IOstream::currentVersion),
263     writeCompression_(IOstream::UNCOMPRESSED),
264     graphFormat_("raw"),
265     runTimeModifiable_(true),
267     functionObjects_(*this, enableFunctionObjects)
269     libs_.open(controlDict_, "libs");
271     // Explicitly set read flags on objectRegistry so anything constructed
272     // from it reads as well (e.g. fvSolution).
273     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
275     setControls();
277     // Time objects not registered so do like objectRegistry::checkIn ourselves.
278     if (runTimeModifiable_)
279     {
280         monitorPtr_.reset
281         (
282             new fileMonitor
283             (
284                 regIOobject::fileModificationChecking == inotify
285              || regIOobject::fileModificationChecking == inotifyMaster
286             )
287         );
289         // File might not exist yet.
290         fileName f(controlDict_.filePath());
292         if (!f.size())
293         {
294             // We don't have this file but would like to re-read it.
295             // Possibly if in master-only reading mode. Use a non-existing
296             // file to keep fileMonitor synced.
297             f = controlDict_.objectPath();
298         }
300         controlDict_.watchIndex() = addWatch(f);
301     }
305 Foam::Time::Time
307     const word& controlDictName,
308     const argList& args,
309     const word& systemName,
310     const word& constantName
313     TimePaths
314     (
315         args.rootPath(),
316         args.caseName(),
317         systemName,
318         constantName
319     ),
321     objectRegistry(*this),
323     libs_(),
325     controlDict_
326     (
327         IOobject
328         (
329             controlDictName,
330             system(),
331             *this,
332             IOobject::MUST_READ_IF_MODIFIED,
333             IOobject::NO_WRITE,
334             false
335         )
336     ),
338     startTimeIndex_(0),
339     startTime_(0),
340     endTime_(0),
342     stopAt_(saEndTime),
343     writeControl_(wcTimeStep),
344     writeInterval_(GREAT),
345     purgeWrite_(0),
346     subCycling_(false),
348     writeFormat_(IOstream::ASCII),
349     writeVersion_(IOstream::currentVersion),
350     writeCompression_(IOstream::UNCOMPRESSED),
351     graphFormat_("raw"),
352     runTimeModifiable_(true),
354     functionObjects_(*this, !args.optionFound("noFunctionObjects"))
356     libs_.open(controlDict_, "libs");
358     // Explicitly set read flags on objectRegistry so anything constructed
359     // from it reads as well (e.g. fvSolution).
360     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
362     setControls();
364     // Time objects not registered so do like objectRegistry::checkIn ourselves.
365     if (runTimeModifiable_)
366     {
367         monitorPtr_.reset
368         (
369             new fileMonitor
370             (
371                 regIOobject::fileModificationChecking == inotify
372              || regIOobject::fileModificationChecking == inotifyMaster
373             )
374         );
376         // File might not exist yet.
377         fileName f(controlDict_.filePath());
379         if (!f.size())
380         {
381             // We don't have this file but would like to re-read it.
382             // Possibly if in master-only reading mode. Use a non-existing
383             // file to keep fileMonitor synced.
384             f = controlDict_.objectPath();
385         }
387         controlDict_.watchIndex() = addWatch(f);
388     }
392 Foam::Time::Time
394     const dictionary& dict,
395     const fileName& rootPath,
396     const fileName& caseName,
397     const word& systemName,
398     const word& constantName,
399     const bool enableFunctionObjects
402     TimePaths
403     (
404         rootPath,
405         caseName,
406         systemName,
407         constantName
408     ),
410     objectRegistry(*this),
412     libs_(),
414     controlDict_
415     (
416         IOobject
417         (
418             controlDictName,
419             system(),
420             *this,
421             IOobject::NO_READ,
422             IOobject::NO_WRITE,
423             false
424         ),
425         dict
426     ),
428     startTimeIndex_(0),
429     startTime_(0),
430     endTime_(0),
432     stopAt_(saEndTime),
433     writeControl_(wcTimeStep),
434     writeInterval_(GREAT),
435     purgeWrite_(0),
436     subCycling_(false),
438     writeFormat_(IOstream::ASCII),
439     writeVersion_(IOstream::currentVersion),
440     writeCompression_(IOstream::UNCOMPRESSED),
441     graphFormat_("raw"),
442     runTimeModifiable_(true),
444     functionObjects_(*this, enableFunctionObjects)
446     libs_.open(controlDict_, "libs");
449     // Explicitly set read flags on objectRegistry so anything constructed
450     // from it reads as well (e.g. fvSolution).
451     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
453     // Since could not construct regIOobject with setting:
454     controlDict_.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
456     setControls();
458     // Time objects not registered so do like objectRegistry::checkIn ourselves.
459     if (runTimeModifiable_)
460     {
461         monitorPtr_.reset
462         (
463             new fileMonitor
464             (
465                 regIOobject::fileModificationChecking == inotify
466              || regIOobject::fileModificationChecking == inotifyMaster
467             )
468         );
470         // File might not exist yet.
471         fileName f(controlDict_.filePath());
473         if (!f.size())
474         {
475             // We don't have this file but would like to re-read it.
476             // Possibly if in master-only reading mode. Use a non-existing
477             // file to keep fileMonitor synced.
478             f = controlDict_.objectPath();
479         }
481         controlDict_.watchIndex() = addWatch(f);
482     }
486 Foam::Time::Time
488     const fileName& rootPath,
489     const fileName& caseName,
490     const word& systemName,
491     const word& constantName,
492     const bool enableFunctionObjects
495     TimePaths
496     (
497         rootPath,
498         caseName,
499         systemName,
500         constantName
501     ),
503     objectRegistry(*this),
505     libs_(),
507     controlDict_
508     (
509         IOobject
510         (
511             controlDictName,
512             system(),
513             *this,
514             IOobject::NO_READ,
515             IOobject::NO_WRITE,
516             false
517         )
518     ),
520     startTimeIndex_(0),
521     startTime_(0),
522     endTime_(0),
524     stopAt_(saEndTime),
525     writeControl_(wcTimeStep),
526     writeInterval_(GREAT),
527     purgeWrite_(0),
528     subCycling_(false),
530     writeFormat_(IOstream::ASCII),
531     writeVersion_(IOstream::currentVersion),
532     writeCompression_(IOstream::UNCOMPRESSED),
533     graphFormat_("raw"),
534     runTimeModifiable_(true),
536     functionObjects_(*this, enableFunctionObjects)
538     libs_.open(controlDict_, "libs");
542 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
544 Foam::Time::~Time()
546     if (controlDict_.watchIndex() != -1)
547     {
548         removeWatch(controlDict_.watchIndex());
549     }
551     // destroy function objects first
552     functionObjects_.clear();
556 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
558 Foam::label Foam::Time::addWatch(const fileName& fName) const
560     return monitorPtr_().addWatch(fName);
564 bool Foam::Time::removeWatch(const label watchIndex) const
566     return monitorPtr_().removeWatch(watchIndex);
569 const Foam::fileName& Foam::Time::getFile(const label watchIndex) const
571     return monitorPtr_().getFile(watchIndex);
575 Foam::fileMonitor::fileState Foam::Time::getState
577     const label watchFd
578 ) const
580     return monitorPtr_().getState(watchFd);
584 void Foam::Time::setUnmodified(const label watchFd) const
586     monitorPtr_().setUnmodified(watchFd);
590 Foam::word Foam::Time::timeName(const scalar t)
592     std::ostringstream buf;
593     buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
594     buf.precision(precision_);
595     buf << t;
596     return buf.str();
600 Foam::word Foam::Time::timeName() const
602     return dimensionedScalar::name();
606 // Search the construction path for times
607 Foam::instantList Foam::Time::times() const
609     return findTimes(path());
613 Foam::word Foam::Time::findInstancePath(const instant& t) const
615     instantList timeDirs = findTimes(path());
617     forAllReverse(timeDirs, timeI)
618     {
619         if (timeDirs[timeI] == t)
620         {
621             return timeDirs[timeI].name();
622         }
623     }
625     return word::null;
629 Foam::instant Foam::Time::findClosestTime(const scalar t) const
631     instantList timeDirs = findTimes(path());
633     // there is only one time (likely "constant") so return it
634     if (timeDirs.size() == 1)
635     {
636         return timeDirs[0];
637     }
639     if (t < timeDirs[1].value())
640     {
641         return timeDirs[1];
642     }
643     else if (t > timeDirs.last().value())
644     {
645         return timeDirs.last();
646     }
648     label nearestIndex = -1;
649     scalar deltaT = GREAT;
651     for (label timeI=1; timeI < timeDirs.size(); ++timeI)
652     {
653         scalar diff = mag(timeDirs[timeI].value() - t);
654         if (diff < deltaT)
655         {
656             deltaT = diff;
657             nearestIndex = timeI;
658         }
659     }
661     return timeDirs[nearestIndex];
665 // This should work too,
666 // if we don't worry about checking "constant" explicitly
668 // Foam::instant Foam::Time::findClosestTime(const scalar t) const
669 // {
670 //     instantList timeDirs = findTimes(path());
671 //     label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
672 //     return timeDirs[timeIndex];
673 // }
675 Foam::label Foam::Time::findClosestTimeIndex
677     const instantList& timeDirs,
678     const scalar t
681     label nearestIndex = -1;
682     scalar deltaT = GREAT;
684     forAll(timeDirs, timeI)
685     {
686         if (timeDirs[timeI].name() == "constant") continue;
688         scalar diff = mag(timeDirs[timeI].value() - t);
689         if (diff < deltaT)
690         {
691             deltaT = diff;
692             nearestIndex = timeI;
693         }
694     }
696     return nearestIndex;
700 Foam::label Foam::Time::startTimeIndex() const
702     return startTimeIndex_;
706 Foam::dimensionedScalar Foam::Time::startTime() const
708     return dimensionedScalar("startTime", dimTime, startTime_);
712 Foam::dimensionedScalar Foam::Time::endTime() const
714     return dimensionedScalar("endTime", dimTime, endTime_);
718 bool Foam::Time::run() const
720     bool running = value() < (endTime_ - 0.5*deltaT_);
722     if (!subCycling_)
723     {
724         // only execute when the condition is no longer true
725         // ie, when exiting the control loop
726         if (!running && timeIndex_ != startTimeIndex_)
727         {
728             // Note, end() also calls an indirect start() as required
729             functionObjects_.end();
730         }
731     }
733     if (running)
734     {
735         if (!subCycling_)
736         {
737             const_cast<Time&>(*this).readModifiedObjects();
739             if (timeIndex_ == startTimeIndex_)
740             {
741                 functionObjects_.start();
742             }
743             else
744             {
745                 functionObjects_.execute();
746             }
747         }
749         // Update the "running" status following the
750         // possible side-effects from functionObjects
751         running = value() < (endTime_ - 0.5*deltaT_);
752     }
754     return running;
758 bool Foam::Time::loop()
760     bool running = run();
762     if (running)
763     {
764         operator++();
765     }
767     return running;
771 bool Foam::Time::end() const
773     return value() > (endTime_ + 0.5*deltaT_);
777 bool Foam::Time::stopAt(const stopAtControls sa) const
779     const bool changed = (stopAt_ != sa);
780     stopAt_ = sa;
782     // adjust endTime
783     if (sa == saEndTime)
784     {
785         controlDict_.lookup("endTime") >> endTime_;
786     }
787     else
788     {
789         endTime_ = GREAT;
790     }
791     return changed;
795 void Foam::Time::setTime(const Time& t)
797     value() = t.value();
798     dimensionedScalar::name() = t.dimensionedScalar::name();
799     timeIndex_ = t.timeIndex_;
803 void Foam::Time::setTime(const instant& inst, const label newIndex)
805     value() = inst.value();
806     dimensionedScalar::name() = inst.name();
807     timeIndex_ = newIndex;
809     IOdictionary timeDict
810     (
811         IOobject
812         (
813             "time",
814             timeName(),
815             "uniform",
816             *this,
817             IOobject::READ_IF_PRESENT,
818             IOobject::NO_WRITE,
819             false
820         )
821     );
823     timeDict.readIfPresent("deltaT", deltaT_);
824     timeDict.readIfPresent("deltaT0", deltaT0_);
825     timeDict.readIfPresent("index", timeIndex_);
829 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
831     setTime(newTime.value(), newIndex);
835 void Foam::Time::setTime(const scalar newTime, const label newIndex)
837     value() = newTime;
838     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
839     timeIndex_ = newIndex;
843 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
845     setEndTime(endTime.value());
849 void Foam::Time::setEndTime(const scalar endTime)
851     endTime_ = endTime;
855 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
857     setDeltaT(deltaT.value());
861 void Foam::Time::setDeltaT(const scalar deltaT)
863     deltaT_ = deltaT;
864     deltaTchanged_ = true;
865     adjustDeltaT();
869 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
871     subCycling_ = true;
872     prevTimeState_.set(new TimeState(*this));
874     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
875     deltaT_ /= nSubCycles;
876     deltaT0_ /= nSubCycles;
877     deltaTSave_ = deltaT0_;
879     return prevTimeState();
883 void Foam::Time::endSubCycle()
885     if (subCycling_)
886     {
887         subCycling_ = false;
888         TimeState::operator=(prevTimeState());
889         prevTimeState_.clear();
890     }
894 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
896 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
898     return operator+=(deltaT.value());
902 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
904     setDeltaT(deltaT);
905     return operator++();
909 Foam::Time& Foam::Time::operator++()
911     deltaT0_ = deltaTSave_;
912     deltaTSave_ = deltaT_;
914     // Save old time name
915     const word oldTimeName = dimensionedScalar::name();
917     setTime(value() + deltaT_, timeIndex_ + 1);
919     if (!subCycling_)
920     {
921         // If the time is very close to zero reset to zero
922         if (mag(value()) < 10*SMALL*deltaT_)
923         {
924             setTime(0.0, timeIndex_);
925         }
926     }
929     // Check that new time representation differs from old one
930     if (dimensionedScalar::name() == oldTimeName)
931     {
932         int oldPrecision = precision_;
933         do
934         {
935             precision_++;
936             setTime(value(), timeIndex());
937         }
938         while (precision_ < 100 && dimensionedScalar::name() == oldTimeName);
940         WarningIn("Time::operator++()")
941             << "Increased the timePrecision from " << oldPrecision
942             << " to " << precision_
943             << " to distinguish between timeNames at time " << value()
944             << endl;
945     }
948     if (!subCycling_)
949     {
950         switch (writeControl_)
951         {
952             case wcTimeStep:
953                 outputTime_ = !(timeIndex_ % label(writeInterval_));
954             break;
956             case wcRunTime:
957             case wcAdjustableRunTime:
958             {
959                 label outputIndex = label
960                 (
961                     ((value() - startTime_) + 0.5*deltaT_)
962                   / writeInterval_
963                 );
965                 if (outputIndex > outputTimeIndex_)
966                 {
967                     outputTime_ = true;
968                     outputTimeIndex_ = outputIndex;
969                 }
970                 else
971                 {
972                     outputTime_ = false;
973                 }
974             }
975             break;
977             case wcCpuTime:
978             {
979                 label outputIndex = label
980                 (
981                     returnReduce(elapsedCpuTime(), maxOp<double>())
982                   / writeInterval_
983                 );
984                 if (outputIndex > outputTimeIndex_)
985                 {
986                     outputTime_ = true;
987                     outputTimeIndex_ = outputIndex;
988                 }
989                 else
990                 {
991                     outputTime_ = false;
992                 }
993             }
994             break;
996             case wcClockTime:
997             {
998                 label outputIndex = label
999                 (
1000                     returnReduce(label(elapsedClockTime()), maxOp<label>())
1001                   / writeInterval_
1002                 );
1003                 if (outputIndex > outputTimeIndex_)
1004                 {
1005                     outputTime_ = true;
1006                     outputTimeIndex_ = outputIndex;
1007                 }
1008                 else
1009                 {
1010                     outputTime_ = false;
1011                 }
1012             }
1013             break;
1014         }
1016         // see if endTime needs adjustment to stop at the next run()/end() check
1017         if (!end())
1018         {
1019             if (stopAt_ == saNoWriteNow)
1020             {
1021                 endTime_ = value();
1022             }
1023             else if (stopAt_ == saWriteNow)
1024             {
1025                 endTime_ = value();
1026                 outputTime_ = true;
1027             }
1028             else if (stopAt_ == saNextWrite && outputTime_ == true)
1029             {
1030                 endTime_ = value();
1031             }
1032         }
1033     }
1035     return *this;
1039 Foam::Time& Foam::Time::operator++(int)
1041     return operator++();
1045 // ************************************************************************* //