initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / db / Time / Time.C
blobcc2c8e710da75ba158d1b8062a4ad508c50d57ee
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "Time.H"
28 #include "PstreamReduceOps.H"
30 #include <sstream>
32 // * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(Foam::Time, 0);
36 template<>
37 const char* Foam::NamedEnum<Foam::Time::stopAtControls, 4>::names[] =
39     "endTime",
40     "noWriteNow",
41     "writeNow",
42     "nextWrite"
45 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
46     Foam::Time::stopAtControlNames_;
48 template<>
49 const char* Foam::NamedEnum<Foam::Time::writeControls, 5>::names[] =
51     "timeStep",
52     "runTime",
53     "adjustableRunTime",
54     "clockTime",
55     "cpuTime"
58 const Foam::NamedEnum<Foam::Time::writeControls, 5>
59     Foam::Time::writeControlNames_;
61 Foam::Time::fmtflags Foam::Time::format_(Foam::Time::general);
62 int Foam::Time::precision_(6);
64 Foam::word Foam::Time::controlDictName("controlDict");
67 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
69 void Foam::Time::adjustDeltaT()
71     if (writeControl_ == wcAdjustableRunTime)
72     {
73         scalar timeToNextWrite = max
74         (
75             0.0,
76             (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
77         );
79         label nStepsToNextWrite = label(timeToNextWrite/deltaT_ - SMALL) + 1;
80         scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
82         // Control the increase of the time step to within a factor of 2
83         // and the decrease within a factor of 5.
84         if (newDeltaT >= deltaT_)
85         {
86             deltaT_ = min(newDeltaT, 2.0*deltaT_);
87         }
88         else
89         {
90             deltaT_ = max(newDeltaT, 0.2*deltaT_);
91         }
92     }
96 void Foam::Time::setControls()
98     // default is to resume calculation from "latestTime"
99     word startFrom = controlDict_.lookupOrDefault<word>
100     (
101         "startFrom",
102         "latestTime"
103     );
105     if (startFrom == "startTime")
106     {
107         controlDict_.lookup("startTime") >> startTime_;
108     }
109     else
110     {
111         // Search directory for valid time directories
112         instantList timeDirs = findTimes(path());
114         if (startFrom == "firstTime")
115         {
116             if (timeDirs.size())
117             {
118                 startTime_ = timeDirs[0].value();
119             }
120         }
121         else if (startFrom == "latestTime")
122         {
123             if (timeDirs.size())
124             {
125                 startTime_ = timeDirs[timeDirs.size()-1].value();
126             }
127         }
128         else
129         {
130             WarningIn("Time::setControls()")
131                 << "    expected startTime, firstTime or latestTime"
132                 << " found '" << startFrom
133                 << "' in dictionary " << controlDict_.name() << nl
134                 << "    Setting time to " << startTime_ << endl;
135         }
136     }
138     setTime(startTime_, 0);
140     readDict();
141     deltaTSave_ = deltaT_;
142     deltaT0_ = deltaTSave_;
144     if (Pstream::parRun())
145     {
146         scalar sumStartTime = startTime_;
147         reduce(sumStartTime, sumOp<scalar>());
148         if
149         (
150             mag(Pstream::nProcs()*startTime_ - sumStartTime)
151           > Pstream::nProcs()*deltaT_/10.0
152         )
153         {
154             FatalErrorIn("Time::setControls()")
155                 << "Start time is not the same for all processors" << nl
156                 << "processor " << Pstream::myProcNo() << " has startTime "
157                 << startTime_ << exit(FatalError);
158         }
159     }
161     IOdictionary timeDict
162     (
163         IOobject
164         (
165             "time",
166             timeName(),
167             "uniform",
168             *this,
169             IOobject::READ_IF_PRESENT,
170             IOobject::NO_WRITE,
171             false
172         )
173     );
175     if (timeDict.readIfPresent("deltaT", deltaTSave_))
176     {
177         deltaT0_ = deltaTSave_;
178     }
180     if (timeDict.readIfPresent("index", startTimeIndex_))
181     {
182         timeIndex_ = startTimeIndex_;
183     }
187 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
189 Foam::Time::Time
191     const word& controlDictName,
192     const fileName& rootPath,
193     const fileName& caseName,
194     const word& systemName,
195     const word& constantName
198     TimePaths
199     (
200         rootPath,
201         caseName,
202         systemName,
203         constantName
204     ),
206     objectRegistry(*this),
208     controlDict_
209     (
210         IOobject
211         (
212             controlDictName,
213             system(),
214             *this,
215             IOobject::MUST_READ,
216             IOobject::NO_WRITE,
217             false
218         )
219     ),
221     startTimeIndex_(0),
222     startTime_(0),
223     endTime_(0),
225     stopAt_(saEndTime),
226     writeControl_(wcTimeStep),
227     writeInterval_(GREAT),
228     purgeWrite_(0),
229     subCycling_(false),
231     writeFormat_(IOstream::ASCII),
232     writeVersion_(IOstream::currentVersion),
233     writeCompression_(IOstream::UNCOMPRESSED),
234     graphFormat_("raw"),
235     runTimeModifiable_(true),
237     readLibs_(controlDict_, "libs"),
238     functionObjects_(*this)
240     setControls();
244 Foam::Time::Time
246     const dictionary& dict,
247     const fileName& rootPath,
248     const fileName& caseName,
249     const word& systemName,
250     const word& constantName
253     TimePaths
254     (
255         rootPath,
256         caseName,
257         systemName,
258         constantName
259     ),
261     objectRegistry(*this),
263     controlDict_
264     (
265         IOobject
266         (
267             controlDictName,
268             system(),
269             *this,
270             IOobject::NO_READ,
271             IOobject::NO_WRITE,
272             false
273         ),
274         dict
275     ),
277     startTimeIndex_(0),
278     startTime_(0),
279     endTime_(0),
281     stopAt_(saEndTime),
282     writeControl_(wcTimeStep),
283     writeInterval_(GREAT),
284     purgeWrite_(0),
285     subCycling_(false),
287     writeFormat_(IOstream::ASCII),
288     writeVersion_(IOstream::currentVersion),
289     writeCompression_(IOstream::UNCOMPRESSED),
290     graphFormat_("raw"),
291     runTimeModifiable_(true),
293     readLibs_(controlDict_, "libs"),
294     functionObjects_(*this)
296     setControls();
300 Foam::Time::Time
302     const fileName& rootPath,
303     const fileName& caseName,
304     const word& systemName,
305     const word& constantName
308     TimePaths
309     (
310         rootPath,
311         caseName,
312         systemName,
313         constantName
314     ),
316     objectRegistry(*this),
318     controlDict_
319     (
320         IOobject
321         (
322             controlDictName,
323             system(),
324             *this,
325             IOobject::NO_READ,
326             IOobject::NO_WRITE,
327             false
328         )
329     ),
331     startTimeIndex_(0),
332     startTime_(0),
333     endTime_(0),
335     stopAt_(saEndTime),
336     writeControl_(wcTimeStep),
337     writeInterval_(GREAT),
338     purgeWrite_(0),
339     subCycling_(false),
341     writeFormat_(IOstream::ASCII),
342     writeVersion_(IOstream::currentVersion),
343     writeCompression_(IOstream::UNCOMPRESSED),
344     graphFormat_("raw"),
345     runTimeModifiable_(true),
347     readLibs_(controlDict_, "libs"),
348     functionObjects_(*this)
352 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
354 Foam::Time::~Time()
356     // destroy function objects first
357     functionObjects_.clear();
361 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
363 Foam::word Foam::Time::timeName(const scalar t)
365     std::ostringstream buf;
366     buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
367     buf.precision(precision_);
368     buf << t;
369     return buf.str();
373 Foam::word Foam::Time::timeName() const
375     return dimensionedScalar::name();
379 // Search the construction path for times
380 Foam::instantList Foam::Time::times() const
382     return findTimes(path());
386 Foam::word Foam::Time::findInstancePath(const instant& t) const
388     instantList timeDirs = findTimes(path());
390     forAllReverse(timeDirs, timeI)
391     {
392         if (timeDirs[timeI] == t)
393         {
394             return timeDirs[timeI].name();
395         }
396     }
398     return word::null;
402 Foam::instant Foam::Time::findClosestTime(const scalar t) const
404     instantList timeDirs = findTimes(path());
406     // there is only one time (likely "constant") so return it
407     if (timeDirs.size() == 1)
408     {
409         return timeDirs[0];
410     }
412     if (t < timeDirs[1].value())
413     {
414         return timeDirs[1];
415     }
416     else if (t > timeDirs[timeDirs.size()-1].value())
417     {
418         return timeDirs[timeDirs.size()-1];
419     }
421     label nearestIndex = -1;
422     scalar deltaT = GREAT;
424     for (label timeI=1; timeI < timeDirs.size(); ++timeI)
425     {
426         scalar diff = mag(timeDirs[timeI].value() - t);
427         if (diff < deltaT)
428         {
429             deltaT = diff;
430             nearestIndex = timeI;
431         }
432     }
434     return timeDirs[nearestIndex];
438 // This should work too,
439 // if we don't worry about checking "constant" explicitly
441 // Foam::instant Foam::Time::findClosestTime(const scalar t) const
442 // {
443 //     instantList timeDirs = findTimes(path());
444 //     label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
445 //     return timeDirs[timeIndex];
446 // }
448 Foam::label Foam::Time::findClosestTimeIndex
450     const instantList& timeDirs,
451     const scalar t
454     label nearestIndex = -1;
455     scalar deltaT = GREAT;
457     forAll(timeDirs, timeI)
458     {
459         if (timeDirs[timeI].name() == "constant") continue;
461         scalar diff = mag(timeDirs[timeI].value() - t);
462         if (diff < deltaT)
463         {
464             deltaT = diff;
465             nearestIndex = timeI;
466         }
467     }
469     return nearestIndex;
473 Foam::label Foam::Time::startTimeIndex() const
475     return startTimeIndex_;
479 Foam::dimensionedScalar Foam::Time::startTime() const
481     return dimensionedScalar("startTime", dimTime, startTime_);
485 Foam::dimensionedScalar Foam::Time::endTime() const
487     return dimensionedScalar("endTime", dimTime, endTime_);
491 bool Foam::Time::run() const
493     bool running = value() < (endTime_ - 0.5*deltaT_);
495     if (!subCycling_)
496     {
497         // only execute when the condition is no longer true
498         // ie, when exiting the control loop
499         if (!running && timeIndex_ != startTimeIndex_)
500         {
501             // Note, end() also calls an indirect start() as required
502             functionObjects_.end();
503         }
504     }
506     return running;
510 bool Foam::Time::loop()
512     bool running = run();
514     if (running)
515     {
516         operator++();
517     }
519     return running;
523 bool Foam::Time::end() const
525     return value() > (endTime_ + 0.5*deltaT_);
529 void Foam::Time::setTime(const Time& t)
531     value() = t.value();
532     dimensionedScalar::name() = t.dimensionedScalar::name();
533     timeIndex_ = t.timeIndex_;
537 void Foam::Time::setTime(const instant& inst, const label newIndex)
539     value() = inst.value();
540     dimensionedScalar::name() = inst.name();
541     timeIndex_ = newIndex;
543     IOdictionary timeDict
544     (
545         IOobject
546         (
547             "time",
548             timeName(),
549             "uniform",
550             *this,
551             IOobject::READ_IF_PRESENT,
552             IOobject::NO_WRITE,
553             false
554         )
555     );
557     timeDict.readIfPresent("deltaT", deltaT_);
558     timeDict.readIfPresent("deltaT0", deltaT0_);
559     timeDict.readIfPresent("index", timeIndex_);
563 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
565     setTime(newTime.value(), newIndex);
569 void Foam::Time::setTime(const scalar newTime, const label newIndex)
571     value() = newTime;
572     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
573     timeIndex_ = newIndex;
577 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
579     setEndTime(endTime.value());
583 void Foam::Time::setEndTime(const scalar endTime)
585     endTime_ = endTime;
589 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
591     setDeltaT(deltaT.value());
595 void Foam::Time::setDeltaT(const scalar deltaT)
597     deltaT_ = deltaT;
598     deltaTchanged_ = true;
599     adjustDeltaT();
603 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
605     subCycling_ = true;
607     TimeState ts = *this;
608     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
609     deltaT_ /= nSubCycles;
610     deltaT0_ /= nSubCycles;
611     deltaTSave_ = deltaT0_;
613     return ts;
617 void Foam::Time::endSubCycle(const TimeState& ts)
619     subCycling_ = false;
620     TimeState::operator=(ts);
624 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
626 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
628     return operator+=(deltaT.value());
632 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
634     setDeltaT(deltaT);
635     return operator++();
639 Foam::Time& Foam::Time::operator++()
641     readModifiedObjects();
643     if (!subCycling_)
644     {
645         if (timeIndex_ == startTimeIndex_)
646         {
647             functionObjects_.start();
648         }
649         else
650         {
651             functionObjects_.execute();
652         }
653     }
655     deltaT0_ = deltaTSave_;
656     deltaTSave_ = deltaT_;
657     setTime(value() + deltaT_, timeIndex_ + 1);
659     // If the time is very close to zero reset to zero
660     if (mag(value()) < 10*SMALL*deltaT_)
661     {
662         setTime(0.0, timeIndex_);
663     }
665     switch (writeControl_)
666     {
667         case wcTimeStep:
668             outputTime_ = !(timeIndex_ % label(writeInterval_));
669         break;
671         case wcRunTime:
672         case wcAdjustableRunTime:
673         {
674             label outputIndex =
675                 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
677             if (outputIndex > outputTimeIndex_)
678             {
679                 outputTime_ = true;
680                 outputTimeIndex_ = outputIndex;
681             }
682             else
683             {
684                 outputTime_ = false;
685             }
686         }
687         break;
689         case wcCpuTime:
690         {
691             label outputIndex = label(elapsedCpuTime()/writeInterval_);
692             if (outputIndex > outputTimeIndex_)
693             {
694                 outputTime_ = true;
695                 outputTimeIndex_ = outputIndex;
696             }
697             else
698             {
699                 outputTime_ = false;
700             }
701         }
702         break;
704         case wcClockTime:
705         {
706             label outputIndex = label(elapsedClockTime()/writeInterval_);
707             if (outputIndex > outputTimeIndex_)
708             {
709                 outputTime_ = true;
710                 outputTimeIndex_ = outputIndex;
711             }
712             else
713             {
714                 outputTime_ = false;
715             }
716         }
717         break;
718     }
720     // see if endTime needs adjustment to stop at the next run()/end() check
721     if (!end())
722     {
723         if (stopAt_ == saNoWriteNow)
724         {
725             endTime_ = value();
726         }
727         else if (stopAt_ == saWriteNow)
728         {
729             endTime_ = value();
730             outputTime_ = true;
731         }
732         else if (stopAt_ == saNextWrite && outputTime_ == true)
733         {
734             endTime_ = value();
735         }
736     }
738     return *this;
742 Foam::Time& Foam::Time::operator++(int)
744     return operator++();
748 // ************************************************************************* //