Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / db / Time / Time.C
blobe90758dafd1dc70fb12d361815c4d81a6aa56349
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2011 OpenCFD Ltd.
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_ = deltaTSave_;
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", deltaTSave_))
200     {
201         deltaT0_ = deltaTSave_;
202     }
204     if (timeDict.readIfPresent("index", startTimeIndex_))
205     {
206         timeIndex_ = startTimeIndex_;
207     }
211 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
213 Foam::Time::Time
215     const word& controlDictName,
216     const fileName& rootPath,
217     const fileName& caseName,
218     const word& systemName,
219     const word& constantName,
220     const bool enableFunctionObjects
223     TimePaths
224     (
225         rootPath,
226         caseName,
227         systemName,
228         constantName
229     ),
231     objectRegistry(*this),
233     libs_(),
235     controlDict_
236     (
237         IOobject
238         (
239             controlDictName,
240             system(),
241             *this,
242             IOobject::MUST_READ_IF_MODIFIED,
243             IOobject::NO_WRITE,
244             false
245         )
246     ),
248     startTimeIndex_(0),
249     startTime_(0),
250     endTime_(0),
252     stopAt_(saEndTime),
253     writeControl_(wcTimeStep),
254     writeInterval_(GREAT),
255     purgeWrite_(0),
256     subCycling_(false),
258     writeFormat_(IOstream::ASCII),
259     writeVersion_(IOstream::currentVersion),
260     writeCompression_(IOstream::UNCOMPRESSED),
261     graphFormat_("raw"),
262     runTimeModifiable_(true),
264     functionObjects_(*this, enableFunctionObjects)
266     libs_.open(controlDict_, "libs");
268     // Explicitly set read flags on objectRegistry so anything constructed
269     // from it reads as well (e.g. fvSolution).
270     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
272     setControls();
274     // Time objects not registered so do like objectRegistry::checkIn ourselves.
275     if (runTimeModifiable_)
276     {
277         monitorPtr_.reset
278         (
279             new fileMonitor
280             (
281                 regIOobject::fileModificationChecking == inotify
282              || regIOobject::fileModificationChecking == inotifyMaster
283             )
284         );
286         // File might not exist yet.
287         fileName f(controlDict_.filePath());
289         if (!f.size())
290         {
291             // We don't have this file but would like to re-read it.
292             // Possibly if in master-only reading mode. Use a non-existing
293             // file to keep fileMonitor synced.
294             f = controlDict_.objectPath();
295         }
297         controlDict_.watchIndex() = addWatch(f);
298     }
302 Foam::Time::Time
304     const word& controlDictName,
305     const argList& args,
306     const word& systemName,
307     const word& constantName
310     TimePaths
311     (
312         args.rootPath(),
313         args.caseName(),
314         systemName,
315         constantName
316     ),
318     objectRegistry(*this),
320     libs_(),
322     controlDict_
323     (
324         IOobject
325         (
326             controlDictName,
327             system(),
328             *this,
329             IOobject::MUST_READ_IF_MODIFIED,
330             IOobject::NO_WRITE,
331             false
332         )
333     ),
335     startTimeIndex_(0),
336     startTime_(0),
337     endTime_(0),
339     stopAt_(saEndTime),
340     writeControl_(wcTimeStep),
341     writeInterval_(GREAT),
342     purgeWrite_(0),
343     subCycling_(false),
345     writeFormat_(IOstream::ASCII),
346     writeVersion_(IOstream::currentVersion),
347     writeCompression_(IOstream::UNCOMPRESSED),
348     graphFormat_("raw"),
349     runTimeModifiable_(true),
351     functionObjects_(*this, !args.optionFound("noFunctionObjects"))
353     libs_.open(controlDict_, "libs");
355     // Explicitly set read flags on objectRegistry so anything constructed
356     // from it reads as well (e.g. fvSolution).
357     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
359     setControls();
361     // Time objects not registered so do like objectRegistry::checkIn ourselves.
362     if (runTimeModifiable_)
363     {
364         monitorPtr_.reset
365         (
366             new fileMonitor
367             (
368                 regIOobject::fileModificationChecking == inotify
369              || regIOobject::fileModificationChecking == inotifyMaster
370             )
371         );
373         // File might not exist yet.
374         fileName f(controlDict_.filePath());
376         if (!f.size())
377         {
378             // We don't have this file but would like to re-read it.
379             // Possibly if in master-only reading mode. Use a non-existing
380             // file to keep fileMonitor synced.
381             f = controlDict_.objectPath();
382         }
384         controlDict_.watchIndex() = addWatch(f);
385     }
389 Foam::Time::Time
391     const dictionary& dict,
392     const fileName& rootPath,
393     const fileName& caseName,
394     const word& systemName,
395     const word& constantName,
396     const bool enableFunctionObjects
399     TimePaths
400     (
401         rootPath,
402         caseName,
403         systemName,
404         constantName
405     ),
407     objectRegistry(*this),
409     libs_(),
411     controlDict_
412     (
413         IOobject
414         (
415             controlDictName,
416             system(),
417             *this,
418             IOobject::NO_READ,
419             IOobject::NO_WRITE,
420             false
421         ),
422         dict
423     ),
425     startTimeIndex_(0),
426     startTime_(0),
427     endTime_(0),
429     stopAt_(saEndTime),
430     writeControl_(wcTimeStep),
431     writeInterval_(GREAT),
432     purgeWrite_(0),
433     subCycling_(false),
435     writeFormat_(IOstream::ASCII),
436     writeVersion_(IOstream::currentVersion),
437     writeCompression_(IOstream::UNCOMPRESSED),
438     graphFormat_("raw"),
439     runTimeModifiable_(true),
441     functionObjects_(*this, enableFunctionObjects)
443     libs_.open(controlDict_, "libs");
446     // Explicitly set read flags on objectRegistry so anything constructed
447     // from it reads as well (e.g. fvSolution).
448     readOpt() = IOobject::MUST_READ_IF_MODIFIED;
450     // Since could not construct regIOobject with setting:
451     controlDict_.readOpt() = IOobject::MUST_READ_IF_MODIFIED;
453     setControls();
455     // Time objects not registered so do like objectRegistry::checkIn ourselves.
456     if (runTimeModifiable_)
457     {
458         monitorPtr_.reset
459         (
460             new fileMonitor
461             (
462                 regIOobject::fileModificationChecking == inotify
463              || regIOobject::fileModificationChecking == inotifyMaster
464             )
465         );
467         // File might not exist yet.
468         fileName f(controlDict_.filePath());
470         if (!f.size())
471         {
472             // We don't have this file but would like to re-read it.
473             // Possibly if in master-only reading mode. Use a non-existing
474             // file to keep fileMonitor synced.
475             f = controlDict_.objectPath();
476         }
478         controlDict_.watchIndex() = addWatch(f);
479     }
483 Foam::Time::Time
485     const fileName& rootPath,
486     const fileName& caseName,
487     const word& systemName,
488     const word& constantName,
489     const bool enableFunctionObjects
492     TimePaths
493     (
494         rootPath,
495         caseName,
496         systemName,
497         constantName
498     ),
500     objectRegistry(*this),
502     libs_(),
504     controlDict_
505     (
506         IOobject
507         (
508             controlDictName,
509             system(),
510             *this,
511             IOobject::NO_READ,
512             IOobject::NO_WRITE,
513             false
514         )
515     ),
517     startTimeIndex_(0),
518     startTime_(0),
519     endTime_(0),
521     stopAt_(saEndTime),
522     writeControl_(wcTimeStep),
523     writeInterval_(GREAT),
524     purgeWrite_(0),
525     subCycling_(false),
527     writeFormat_(IOstream::ASCII),
528     writeVersion_(IOstream::currentVersion),
529     writeCompression_(IOstream::UNCOMPRESSED),
530     graphFormat_("raw"),
531     runTimeModifiable_(true),
533     functionObjects_(*this, enableFunctionObjects)
535     libs_.open(controlDict_, "libs");
539 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
541 Foam::Time::~Time()
543     if (controlDict_.watchIndex() != -1)
544     {
545         removeWatch(controlDict_.watchIndex());
546     }
548     // destroy function objects first
549     functionObjects_.clear();
553 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
555 Foam::label Foam::Time::addWatch(const fileName& fName) const
557     return monitorPtr_().addWatch(fName);
561 bool Foam::Time::removeWatch(const label watchIndex) const
563     return monitorPtr_().removeWatch(watchIndex);
566 const Foam::fileName& Foam::Time::getFile(const label watchIndex) const
568     return monitorPtr_().getFile(watchIndex);
572 Foam::fileMonitor::fileState Foam::Time::getState
574     const label watchFd
575 ) const
577     return monitorPtr_().getState(watchFd);
581 void Foam::Time::setUnmodified(const label watchFd) const
583     monitorPtr_().setUnmodified(watchFd);
587 Foam::word Foam::Time::timeName(const scalar t)
589     std::ostringstream buf;
590     buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
591     buf.precision(precision_);
592     buf << t;
593     return buf.str();
597 Foam::word Foam::Time::timeName() const
599     return dimensionedScalar::name();
603 // Search the construction path for times
604 Foam::instantList Foam::Time::times() const
606     return findTimes(path());
610 Foam::word Foam::Time::findInstancePath(const instant& t) const
612     instantList timeDirs = findTimes(path());
614     forAllReverse(timeDirs, timeI)
615     {
616         if (timeDirs[timeI] == t)
617         {
618             return timeDirs[timeI].name();
619         }
620     }
622     return word::null;
626 Foam::instant Foam::Time::findClosestTime(const scalar t) const
628     instantList timeDirs = findTimes(path());
630     // there is only one time (likely "constant") so return it
631     if (timeDirs.size() == 1)
632     {
633         return timeDirs[0];
634     }
636     if (t < timeDirs[1].value())
637     {
638         return timeDirs[1];
639     }
640     else if (t > timeDirs.last().value())
641     {
642         return timeDirs.last();
643     }
645     label nearestIndex = -1;
646     scalar deltaT = GREAT;
648     for (label timeI=1; timeI < timeDirs.size(); ++timeI)
649     {
650         scalar diff = mag(timeDirs[timeI].value() - t);
651         if (diff < deltaT)
652         {
653             deltaT = diff;
654             nearestIndex = timeI;
655         }
656     }
658     return timeDirs[nearestIndex];
662 // This should work too,
663 // if we don't worry about checking "constant" explicitly
665 // Foam::instant Foam::Time::findClosestTime(const scalar t) const
666 // {
667 //     instantList timeDirs = findTimes(path());
668 //     label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
669 //     return timeDirs[timeIndex];
670 // }
672 Foam::label Foam::Time::findClosestTimeIndex
674     const instantList& timeDirs,
675     const scalar t
678     label nearestIndex = -1;
679     scalar deltaT = GREAT;
681     forAll(timeDirs, timeI)
682     {
683         if (timeDirs[timeI].name() == "constant") continue;
685         scalar diff = mag(timeDirs[timeI].value() - t);
686         if (diff < deltaT)
687         {
688             deltaT = diff;
689             nearestIndex = timeI;
690         }
691     }
693     return nearestIndex;
697 Foam::label Foam::Time::startTimeIndex() const
699     return startTimeIndex_;
703 Foam::dimensionedScalar Foam::Time::startTime() const
705     return dimensionedScalar("startTime", dimTime, startTime_);
709 Foam::dimensionedScalar Foam::Time::endTime() const
711     return dimensionedScalar("endTime", dimTime, endTime_);
715 bool Foam::Time::run() const
717     bool running = value() < (endTime_ - 0.5*deltaT_);
719     if (!subCycling_)
720     {
721         // only execute when the condition is no longer true
722         // ie, when exiting the control loop
723         if (!running && timeIndex_ != startTimeIndex_)
724         {
725             // Note, end() also calls an indirect start() as required
726             functionObjects_.end();
727         }
728     }
730     if (running)
731     {
732         if (!subCycling_)
733         {
734             const_cast<Time&>(*this).readModifiedObjects();
736             if (timeIndex_ == startTimeIndex_)
737             {
738                 functionObjects_.start();
739             }
740             else
741             {
742                 functionObjects_.execute();
743             }
744         }
746         // Update the "running" status following the
747         // possible side-effects from functionObjects
748         running = value() < (endTime_ - 0.5*deltaT_);
749     }
751     return running;
755 bool Foam::Time::loop()
757     bool running = run();
759     if (running)
760     {
761         operator++();
762     }
764     return running;
768 bool Foam::Time::end() const
770     return value() > (endTime_ + 0.5*deltaT_);
774 bool Foam::Time::stopAt(const stopAtControls sa) const
776     const bool changed = (stopAt_ != sa);
777     stopAt_ = sa;
779     // adjust endTime
780     if (sa == saEndTime)
781     {
782         controlDict_.lookup("endTime") >> endTime_;
783     }
784     else
785     {
786         endTime_ = GREAT;
787     }
788     return changed;
792 void Foam::Time::setTime(const Time& t)
794     value() = t.value();
795     dimensionedScalar::name() = t.dimensionedScalar::name();
796     timeIndex_ = t.timeIndex_;
800 void Foam::Time::setTime(const instant& inst, const label newIndex)
802     value() = inst.value();
803     dimensionedScalar::name() = inst.name();
804     timeIndex_ = newIndex;
806     IOdictionary timeDict
807     (
808         IOobject
809         (
810             "time",
811             timeName(),
812             "uniform",
813             *this,
814             IOobject::READ_IF_PRESENT,
815             IOobject::NO_WRITE,
816             false
817         )
818     );
820     timeDict.readIfPresent("deltaT", deltaT_);
821     timeDict.readIfPresent("deltaT0", deltaT0_);
822     timeDict.readIfPresent("index", timeIndex_);
826 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
828     setTime(newTime.value(), newIndex);
832 void Foam::Time::setTime(const scalar newTime, const label newIndex)
834     value() = newTime;
835     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
836     timeIndex_ = newIndex;
840 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
842     setEndTime(endTime.value());
846 void Foam::Time::setEndTime(const scalar endTime)
848     endTime_ = endTime;
852 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
854     setDeltaT(deltaT.value());
858 void Foam::Time::setDeltaT(const scalar deltaT)
860     deltaT_ = deltaT;
861     deltaTchanged_ = true;
862     adjustDeltaT();
866 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
868     subCycling_ = true;
869     prevTimeState_.set(new TimeState(*this));
871     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
872     deltaT_ /= nSubCycles;
873     deltaT0_ /= nSubCycles;
874     deltaTSave_ = deltaT0_;
876     return prevTimeState();
880 void Foam::Time::endSubCycle()
882     if (subCycling_)
883     {
884         subCycling_ = false;
885         TimeState::operator=(prevTimeState());
886         prevTimeState_.clear();
887     }
891 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
893 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
895     return operator+=(deltaT.value());
899 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
901     setDeltaT(deltaT);
902     return operator++();
906 Foam::Time& Foam::Time::operator++()
908     deltaT0_ = deltaTSave_;
909     deltaTSave_ = deltaT_;
911     // Save old time name
912     const word oldTimeName = dimensionedScalar::name();
914     setTime(value() + deltaT_, timeIndex_ + 1);
916     if (!subCycling_)
917     {
918         // If the time is very close to zero reset to zero
919         if (mag(value()) < 10*SMALL*deltaT_)
920         {
921             setTime(0.0, timeIndex_);
922         }
923     }
926     // Check that new time representation differs from old one
927     if (dimensionedScalar::name() == oldTimeName)
928     {
929         int oldPrecision = precision_;
930         do
931         {
932             precision_++;
933             setTime(value(), timeIndex());
934         }
935         while (precision_ < 100 && dimensionedScalar::name() == oldTimeName);
937         WarningIn("Time::operator++()")
938             << "Increased the timePrecision from " << oldPrecision
939             << " to " << precision_
940             << " to distinguish between timeNames at time " << value()
941             << endl;
942     }
945     if (!subCycling_)
946     {
947         switch (writeControl_)
948         {
949             case wcTimeStep:
950                 outputTime_ = !(timeIndex_ % label(writeInterval_));
951             break;
953             case wcRunTime:
954             case wcAdjustableRunTime:
955             {
956                 label outputIndex = label
957                 (
958                     ((value() - startTime_) + 0.5*deltaT_)
959                   / writeInterval_
960                 );
962                 if (outputIndex > outputTimeIndex_)
963                 {
964                     outputTime_ = true;
965                     outputTimeIndex_ = outputIndex;
966                 }
967                 else
968                 {
969                     outputTime_ = false;
970                 }
971             }
972             break;
974             case wcCpuTime:
975             {
976                 label outputIndex = label
977                 (
978                     returnReduce(elapsedCpuTime(), maxOp<double>())
979                   / writeInterval_
980                 );
981                 if (outputIndex > outputTimeIndex_)
982                 {
983                     outputTime_ = true;
984                     outputTimeIndex_ = outputIndex;
985                 }
986                 else
987                 {
988                     outputTime_ = false;
989                 }
990             }
991             break;
993             case wcClockTime:
994             {
995                 label outputIndex = label
996                 (
997                     returnReduce(label(elapsedClockTime()), maxOp<label>())
998                   / writeInterval_
999                 );
1000                 if (outputIndex > outputTimeIndex_)
1001                 {
1002                     outputTime_ = true;
1003                     outputTimeIndex_ = outputIndex;
1004                 }
1005                 else
1006                 {
1007                     outputTime_ = false;
1008                 }
1009             }
1010             break;
1011         }
1013         // see if endTime needs adjustment to stop at the next run()/end() check
1014         if (!end())
1015         {
1016             if (stopAt_ == saNoWriteNow)
1017             {
1018                 endTime_ = value();
1019             }
1020             else if (stopAt_ == saWriteNow)
1021             {
1022                 endTime_ = value();
1023                 outputTime_ = true;
1024             }
1025             else if (stopAt_ == saNextWrite && outputTime_ == true)
1026             {
1027                 endTime_ = value();
1028             }
1029         }
1030     }
1032     return *this;
1036 Foam::Time& Foam::Time::operator++(int)
1038     return operator++();
1042 // ************************************************************************* //