initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / OpenFOAM / db / Time / Time.C
blobef4fb9c21fa23b78a973f2bfe90028534bf36a53
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 namespace Foam
36     defineTypeNameAndDebug(Time, 0);
40 template<>
41 const char* Foam::NamedEnum<Foam::Time::stopAtControls, 4>::names[] =
43     "endTime",
44     "noWriteNow",
45     "writeNow",
46     "nextWrite"
49 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
50     Foam::Time::stopAtControlNames_;
52 template<>
53 const char* Foam::NamedEnum<Foam::Time::writeControls, 5>::names[] =
55     "timeStep",
56     "runTime",
57     "adjustableRunTime",
58     "clockTime",
59     "cpuTime"
62 const Foam::NamedEnum<Foam::Time::writeControls, 5>
63     Foam::Time::writeControlNames_;
65 Foam::Time::fmtflags Foam::Time::format_(Foam::Time::general);
66 int Foam::Time::precision_(6);
68 Foam::word Foam::Time::controlDictName("controlDict");
71 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
73 void Foam::Time::adjustDeltaT()
75     if (writeControl_ == wcAdjustableRunTime)
76     {
77         scalar timeToNextWrite = max
78         (
79             0.0,
80             (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
81         );
83         label nStepsToNextWrite = label(timeToNextWrite/deltaT_ - SMALL) + 1;
84         scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
86         // Control the increase of the time step to within a factor of 2
87         // and the decrease within a factor of 5.
88         if (newDeltaT >= deltaT_)
89         {
90             deltaT_ = min(newDeltaT, 2.0*deltaT_);
91         }
92         else
93         {
94             deltaT_ = max(newDeltaT, 0.2*deltaT_);
95         }
96     }
100 void Foam::Time::setControls()
102     // default is to resume calculation from "latestTime"
103     word startFrom("latestTime");
105     controlDict_.readIfPresent("startFrom", startFrom);
107     if (startFrom == "startTime")
108     {
109         controlDict_.lookup("startTime") >> startTime_;
110     }
111     else
112     {
113         // Search directory for valid time directories
114         instantList Times = findTimes(path());
116         if (startFrom == "firstTime")
117         {
118             if (Times.size())
119             {
120                 startTime_ = Times[0].value();
121             }
122         }
123         else if (startFrom == "latestTime")
124         {
125             if (Times.size())
126             {
127                 startTime_ = Times[Times.size()-1].value();
128             }
129         }
130         else
131         {
132             WarningIn("Time::setControls()")
133                 << "    expected startTime, firstTime or latestTime"
134                 << " found '" << startFrom
135                 << "' in dictionary " << controlDict_.name() << nl
136                 << "    Setting time to " << startTime_ << endl;
137         }
138     }
140     setTime(startTime_, 0);
142     readDict();
143     deltaTSave_ = deltaT_;
144     deltaT0_ = deltaTSave_;
146     if (Pstream::parRun())
147     {
148         scalar sumStartTime = startTime_;
149         reduce(sumStartTime, sumOp<scalar>());
150         if
151         (
152             mag(Pstream::nProcs()*startTime_ - sumStartTime)
153           > Pstream::nProcs()*deltaT_/10.0
154         )
155         {
156             FatalErrorIn("Time::setControls()")
157                 << "Start time is not the same for all processors" << nl
158                 << "processor " << Pstream::myProcNo() << " has startTime "
159                 << startTime_ << exit(FatalError);
160         }
161     }
163     IOdictionary timeDict
164     (
165         IOobject
166         (
167             "time",
168             timeName(),
169             "uniform",
170             *this,
171             IOobject::READ_IF_PRESENT,
172             IOobject::NO_WRITE,
173             false
174         )
175     );
177     if (timeDict.readIfPresent("deltaT", deltaTSave_))
178     {
179         deltaT0_ = deltaTSave_;
180     }
182     if (timeDict.readIfPresent("index", startTimeIndex_))
183     {
184         timeIndex_ = startTimeIndex_;
185     }
189 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
191 Foam::Time::Time
193     const word& controlDictName,
194     const fileName& rootPath,
195     const fileName& caseName,
196     const word& systemName,
197     const word& constantName
200     TimePaths
201     (
202         rootPath,
203         caseName,
204         systemName,
205         constantName
206     ),
208     objectRegistry(*this),
210     controlDict_
211     (
212         IOobject
213         (
214             controlDictName,
215             system(),
216             *this,
217             IOobject::MUST_READ,
218             IOobject::NO_WRITE,
219             false
220         )
221     ),
223     startTimeIndex_(0),
224     startTime_(0),
225     endTime_(0),
227     stopAt_(saEndTime),
228     writeControl_(wcTimeStep),
229     writeInterval_(GREAT),
230     purgeWrite_(0),
231     subCycling_(false),
233     writeFormat_(IOstream::ASCII),
234     writeVersion_(IOstream::currentVersion),
235     writeCompression_(IOstream::UNCOMPRESSED),
236     graphFormat_("raw"),
237     runTimeModifiable_(true),
239     readLibs_(controlDict_, "libs"),
240     functionObjects_(*this)
242     setControls();
246 Foam::Time::Time
248     const dictionary& dict,
249     const fileName& rootPath,
250     const fileName& caseName,
251     const word& systemName,
252     const word& constantName
255     TimePaths
256     (
257         rootPath,
258         caseName,
259         systemName,
260         constantName
261     ),
263     objectRegistry(*this),
265     controlDict_
266     (
267         IOobject
268         (
269             controlDictName,
270             system(),
271             *this,
272             IOobject::NO_READ,
273             IOobject::NO_WRITE,
274             false
275         ),
276         dict
277     ),
279     startTimeIndex_(0),
280     startTime_(0),
281     endTime_(0),
283     stopAt_(saEndTime),
284     writeControl_(wcTimeStep),
285     writeInterval_(GREAT),
286     purgeWrite_(0),
287     subCycling_(false),
289     writeFormat_(IOstream::ASCII),
290     writeVersion_(IOstream::currentVersion),
291     writeCompression_(IOstream::UNCOMPRESSED),
292     graphFormat_("raw"),
293     runTimeModifiable_(true),
295     readLibs_(controlDict_, "libs"),
296     functionObjects_(*this)
298     setControls();
302 Foam::Time::Time
304     const fileName& rootPath,
305     const fileName& caseName,
306     const word& systemName,
307     const word& constantName
310     TimePaths
311     (
312         rootPath,
313         caseName,
314         systemName,
315         constantName
316     ),
318     objectRegistry(*this),
320     controlDict_
321     (
322         IOobject
323         (
324             controlDictName,
325             system(),
326             *this,
327             IOobject::NO_READ,
328             IOobject::NO_WRITE,
329             false
330         )
331     ),
333     startTimeIndex_(0),
334     startTime_(0),
335     endTime_(0),
337     stopAt_(saEndTime),
338     writeControl_(wcTimeStep),
339     writeInterval_(GREAT),
340     purgeWrite_(0),
341     subCycling_(false),
343     writeFormat_(IOstream::ASCII),
344     writeVersion_(IOstream::currentVersion),
345     writeCompression_(IOstream::UNCOMPRESSED),
346     graphFormat_("raw"),
347     runTimeModifiable_(true),
349     readLibs_(controlDict_, "libs"),
350     functionObjects_(*this)
354 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
356 Foam::Time::~Time()
360 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
362 Foam::word Foam::Time::timeName(const scalar t)
364     std::ostringstream osBuffer;
365     osBuffer.setf(ios_base::fmtflags(format_), ios_base::floatfield);
366     osBuffer.precision(precision_);
367     osBuffer << t;
368     return osBuffer.str();
372 Foam::word Foam::Time::timeName() const
374     return dimensionedScalar::name();
375     //return timeName(timeOutputValue());
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 times = Time::findTimes(path());
390     forAllReverse(times, i)
391     {
392         if (times[i] == t)
393         {
394             return times[i].name();
395         }
396     }
398     return word::null;
402 Foam::instant Foam::Time::findClosestTime(const scalar t) const
404     instantList times = Time::findTimes(path());
406     // If there is only one time it is "constant" so return it
407     if (times.size() == 1)
408     {
409         return times[0];
410     }
412     if (t < times[1].value())
413     {
414         return times[1];
415     }
416     else if (t > times[times.size() - 1].value())
417     {
418         return times[times.size() - 1];
419     }
421     label nearestIndex = -1;
422     scalar deltaT = GREAT;
424     for (label i=1; i<times.size(); i++)
425     {
426         scalar diff = mag(times[i].value() - t);
427         if (diff < deltaT)
428         {
429             deltaT = diff;
430             nearestIndex = i;
431         }
432     }
434     return times[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 times = Time::findTimes(path());
444 //     label timeIndex = min(findClosestTimeIndex(times, t), 0);
445 //     return times[timeIndex];
446 // }
448 Foam::label Foam::Time::findClosestTimeIndex
450     const instantList& times,
451     const scalar t
454     label nearestIndex = -1;
455     scalar deltaT = GREAT;
457     forAll (times, i)
458     {
459         if (times[i].name() == "constant") continue;
461         scalar diff = fabs(times[i].value() - t);
462         if (diff < deltaT)
463         {
464             deltaT = diff;
465             nearestIndex = i;
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_ && !running && timeIndex_ != startTimeIndex_)
496     {
497         const_cast<functionObjectList&>(functionObjects_).execute();
498     }
500     return running;
504 bool Foam::Time::end() const
506     return (value() > (endTime_ + 0.5*deltaT_));
510 void Foam::Time::setTime(const Time& t)
512     value() = t.value();
513     dimensionedScalar::name() = t.dimensionedScalar::name();
514     timeIndex_ = t.timeIndex_;
518 void Foam::Time::setTime(const instant& inst, const label newIndex)
520     value() = inst.value();
521     dimensionedScalar::name() = inst.name();
522     timeIndex_ = newIndex;
524     IOdictionary timeDict
525     (
526         IOobject
527         (
528             "time",
529             timeName(),
530             "uniform",
531             *this,
532             IOobject::READ_IF_PRESENT,
533             IOobject::NO_WRITE,
534             false
535         )
536     );
538     timeDict.readIfPresent("deltaT", deltaT_);
539     timeDict.readIfPresent("deltaT0", deltaT0_);
540     timeDict.readIfPresent("index", timeIndex_);
544 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
546     setTime(newTime.value(), newIndex);
550 void Foam::Time::setTime(const scalar newTime, const label newIndex)
552     value() = newTime;
553     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
554     timeIndex_ = newIndex;
558 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
560     setEndTime(endTime.value());
564 void Foam::Time::setEndTime(const scalar endTime)
566     endTime_ = endTime;
570 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
572     setDeltaT(deltaT.value());
576 void Foam::Time::setDeltaT(const scalar deltaT)
578     deltaT_ = deltaT;
579     deltaTchanged_ = true;
580     adjustDeltaT();
584 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
586     subCycling_ = true;
588     TimeState ts = *this;
589     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
590     deltaT_ /= nSubCycles;
591     deltaT0_ /= nSubCycles;
592     deltaTSave_ = deltaT0_;
594     return ts;
598 void Foam::Time::endSubCycle(const TimeState& ts)
600     subCycling_ = false;
601     TimeState::operator=(ts);
605 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
607 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
609     return operator+=(deltaT.value());
613 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
615     readModifiedObjects();
617     if (!subCycling_)
618     {
619         if (timeIndex_ == startTimeIndex_)
620         {
621             functionObjects_.start();
622         }
623         else
624         {
625             functionObjects_.execute();
626         }
627     }
629     setDeltaT(deltaT);
630     operator++();
632     return *this;
636 Foam::Time& Foam::Time::operator++()
638     readModifiedObjects();
640     if (!subCycling_)
641     {
642         if (timeIndex_ == startTimeIndex_)
643         {
644             functionObjects_.start();
645         }
646         else
647         {
648             functionObjects_.execute();
649         }
650     }
652     deltaT0_ = deltaTSave_;
653     deltaTSave_ = deltaT_;
654     setTime(value() + deltaT_, timeIndex_ + 1);
656     // If the time is very close to zero reset to zero
657     if (mag(value()) < 10*SMALL*deltaT_)
658     {
659         setTime(0.0, timeIndex_);
660     }
662     switch(writeControl_)
663     {
664         case wcTimeStep:
665             outputTime_ = !(timeIndex_%label(writeInterval_));
666         break;
668         case wcRunTime:
669         case wcAdjustableRunTime:
670         {
671             label outputTimeIndex =
672                 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
674             if (outputTimeIndex > outputTimeIndex_)
675             {
676                 outputTime_ = true;
677                 outputTimeIndex_ = outputTimeIndex;
678             }
679             else
680             {
681                 outputTime_ = false;
682             }
683         }
684         break;
686         case wcCpuTime:
687         {
688             label outputTimeIndex =
689                 label(elapsedCpuTime()/writeInterval_);
691             if (outputTimeIndex > outputTimeIndex_)
692             {
693                 outputTime_ = true;
694                 outputTimeIndex_ = outputTimeIndex;
695             }
696             else
697             {
698                 outputTime_ = false;
699             }
700         }
701         break;
703         case wcClockTime:
704         {
705             label outputTimeIndex = label(elapsedClockTime()/writeInterval_);
706             if (outputTimeIndex > outputTimeIndex_)
707             {
708                 outputTime_ = true;
709                 outputTimeIndex_ = outputTimeIndex;
710             }
711             else
712             {
713                 outputTime_ = false;
714             }
715         }
716         break;
717     };
719     if (!end())
720     {
721         if (stopAt_ == saNoWriteNow)
722         {
723             endTime_ = value();
724         }
725         else if (stopAt_ == saWriteNow)
726         {
727             endTime_ = value();
728             outputTime_ = true;
729         }
730         else if (stopAt_ == saNextWrite && outputTime_ == true)
731         {
732             endTime_ = value();
733         }
734     }
736     return *this;
740 Foam::Time& Foam::Time::operator++(int)
742     return operator++();
746 // ************************************************************************* //