ENH: Automatically increase the time precision to make timeNames differ.
[OpenFOAM-1.6.x.git] / src / OpenFOAM / db / Time / Time.C
blob427b5c454d764004484f82c62ccdd9775bab3dda
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             FatalIOErrorIn("Time::setControls()", controlDict_)
131                 << "expected startTime, firstTime or latestTime"
132                 << " found '" << startFrom << "'"
133                 << exit(FatalIOError);
134         }
135     }
137     setTime(startTime_, 0);
139     readDict();
140     deltaTSave_ = deltaT_;
141     deltaT0_ = deltaTSave_;
143     if (Pstream::parRun())
144     {
145         scalar sumStartTime = startTime_;
146         reduce(sumStartTime, sumOp<scalar>());
147         if
148         (
149             mag(Pstream::nProcs()*startTime_ - sumStartTime)
150           > Pstream::nProcs()*deltaT_/10.0
151         )
152         {
153             FatalIOErrorIn("Time::setControls()", controlDict_)
154                 << "Start time is not the same for all processors" << nl
155                 << "processor " << Pstream::myProcNo() << " has startTime "
156                 << startTime_ << exit(FatalIOError);
157         }
158     }
160     IOdictionary timeDict
161     (
162         IOobject
163         (
164             "time",
165             timeName(),
166             "uniform",
167             *this,
168             IOobject::READ_IF_PRESENT,
169             IOobject::NO_WRITE,
170             false
171         )
172     );
174     if (timeDict.readIfPresent("deltaT", deltaTSave_))
175     {
176         deltaT0_ = deltaTSave_;
177     }
179     if (timeDict.readIfPresent("index", startTimeIndex_))
180     {
181         timeIndex_ = startTimeIndex_;
182     }
186 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
188 Foam::Time::Time
190     const word& controlDictName,
191     const fileName& rootPath,
192     const fileName& caseName,
193     const word& systemName,
194     const word& constantName
197     TimePaths
198     (
199         rootPath,
200         caseName,
201         systemName,
202         constantName
203     ),
205     objectRegistry(*this),
207     controlDict_
208     (
209         IOobject
210         (
211             controlDictName,
212             system(),
213             *this,
214             IOobject::MUST_READ,
215             IOobject::NO_WRITE,
216             false
217         )
218     ),
220     startTimeIndex_(0),
221     startTime_(0),
222     endTime_(0),
224     stopAt_(saEndTime),
225     writeControl_(wcTimeStep),
226     writeInterval_(GREAT),
227     purgeWrite_(0),
228     subCycling_(false),
230     writeFormat_(IOstream::ASCII),
231     writeVersion_(IOstream::currentVersion),
232     writeCompression_(IOstream::UNCOMPRESSED),
233     graphFormat_("raw"),
234     runTimeModifiable_(true),
236     readLibs_(controlDict_, "libs"),
237     functionObjects_(*this)
239     setControls();
243 Foam::Time::Time
245     const dictionary& dict,
246     const fileName& rootPath,
247     const fileName& caseName,
248     const word& systemName,
249     const word& constantName
252     TimePaths
253     (
254         rootPath,
255         caseName,
256         systemName,
257         constantName
258     ),
260     objectRegistry(*this),
262     controlDict_
263     (
264         IOobject
265         (
266             controlDictName,
267             system(),
268             *this,
269             IOobject::NO_READ,
270             IOobject::NO_WRITE,
271             false
272         ),
273         dict
274     ),
276     startTimeIndex_(0),
277     startTime_(0),
278     endTime_(0),
280     stopAt_(saEndTime),
281     writeControl_(wcTimeStep),
282     writeInterval_(GREAT),
283     purgeWrite_(0),
284     subCycling_(false),
286     writeFormat_(IOstream::ASCII),
287     writeVersion_(IOstream::currentVersion),
288     writeCompression_(IOstream::UNCOMPRESSED),
289     graphFormat_("raw"),
290     runTimeModifiable_(true),
292     readLibs_(controlDict_, "libs"),
293     functionObjects_(*this)
295     setControls();
299 Foam::Time::Time
301     const fileName& rootPath,
302     const fileName& caseName,
303     const word& systemName,
304     const word& constantName
307     TimePaths
308     (
309         rootPath,
310         caseName,
311         systemName,
312         constantName
313     ),
315     objectRegistry(*this),
317     controlDict_
318     (
319         IOobject
320         (
321             controlDictName,
322             system(),
323             *this,
324             IOobject::NO_READ,
325             IOobject::NO_WRITE,
326             false
327         )
328     ),
330     startTimeIndex_(0),
331     startTime_(0),
332     endTime_(0),
334     stopAt_(saEndTime),
335     writeControl_(wcTimeStep),
336     writeInterval_(GREAT),
337     purgeWrite_(0),
338     subCycling_(false),
340     writeFormat_(IOstream::ASCII),
341     writeVersion_(IOstream::currentVersion),
342     writeCompression_(IOstream::UNCOMPRESSED),
343     graphFormat_("raw"),
344     runTimeModifiable_(true),
346     readLibs_(controlDict_, "libs"),
347     functionObjects_(*this)
351 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
353 Foam::Time::~Time()
355     // destroy function objects first
356     functionObjects_.clear();
360 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
362 Foam::word Foam::Time::timeName(const scalar t)
364     std::ostringstream buf;
365     buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
366     buf.precision(precision_);
367     buf << t;
368     return buf.str();
372 Foam::word Foam::Time::timeName() const
374     return dimensionedScalar::name();
378 // Search the construction path for times
379 Foam::instantList Foam::Time::times() const
381     return findTimes(path());
385 Foam::word Foam::Time::findInstancePath(const instant& t) const
387     instantList timeDirs = findTimes(path());
389     forAllReverse(timeDirs, timeI)
390     {
391         if (timeDirs[timeI] == t)
392         {
393             return timeDirs[timeI].name();
394         }
395     }
397     return word::null;
401 Foam::instant Foam::Time::findClosestTime(const scalar t) const
403     instantList timeDirs = findTimes(path());
405     // there is only one time (likely "constant") so return it
406     if (timeDirs.size() == 1)
407     {
408         return timeDirs[0];
409     }
411     if (t < timeDirs[1].value())
412     {
413         return timeDirs[1];
414     }
415     else if (t > timeDirs[timeDirs.size()-1].value())
416     {
417         return timeDirs[timeDirs.size()-1];
418     }
420     label nearestIndex = -1;
421     scalar deltaT = GREAT;
423     for (label timeI=1; timeI < timeDirs.size(); ++timeI)
424     {
425         scalar diff = mag(timeDirs[timeI].value() - t);
426         if (diff < deltaT)
427         {
428             deltaT = diff;
429             nearestIndex = timeI;
430         }
431     }
433     return timeDirs[nearestIndex];
437 // This should work too,
438 // if we don't worry about checking "constant" explicitly
440 // Foam::instant Foam::Time::findClosestTime(const scalar t) const
441 // {
442 //     instantList timeDirs = findTimes(path());
443 //     label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
444 //     return timeDirs[timeIndex];
445 // }
447 Foam::label Foam::Time::findClosestTimeIndex
449     const instantList& timeDirs,
450     const scalar t
453     label nearestIndex = -1;
454     scalar deltaT = GREAT;
456     forAll(timeDirs, timeI)
457     {
458         if (timeDirs[timeI].name() == "constant") continue;
460         scalar diff = mag(timeDirs[timeI].value() - t);
461         if (diff < deltaT)
462         {
463             deltaT = diff;
464             nearestIndex = timeI;
465         }
466     }
468     return nearestIndex;
472 Foam::label Foam::Time::startTimeIndex() const
474     return startTimeIndex_;
478 Foam::dimensionedScalar Foam::Time::startTime() const
480     return dimensionedScalar("startTime", dimTime, startTime_);
484 Foam::dimensionedScalar Foam::Time::endTime() const
486     return dimensionedScalar("endTime", dimTime, endTime_);
490 bool Foam::Time::run() const
492     bool running = value() < (endTime_ - 0.5*deltaT_);
494     if (!subCycling_)
495     {
496         // only execute when the condition is no longer true
497         // ie, when exiting the control loop
498         if (!running && timeIndex_ != startTimeIndex_)
499         {
500             // Note, end() also calls an indirect start() as required
501             functionObjects_.end();
502         }
503     }
505     return running;
509 bool Foam::Time::loop()
511     bool running = run();
513     if (running)
514     {
515         operator++();
516     }
518     return running;
522 bool Foam::Time::end() const
524     return value() > (endTime_ + 0.5*deltaT_);
528 void Foam::Time::setTime(const Time& t)
530     value() = t.value();
531     dimensionedScalar::name() = t.dimensionedScalar::name();
532     timeIndex_ = t.timeIndex_;
536 void Foam::Time::setTime(const instant& inst, const label newIndex)
538     value() = inst.value();
539     dimensionedScalar::name() = inst.name();
540     timeIndex_ = newIndex;
542     IOdictionary timeDict
543     (
544         IOobject
545         (
546             "time",
547             timeName(),
548             "uniform",
549             *this,
550             IOobject::READ_IF_PRESENT,
551             IOobject::NO_WRITE,
552             false
553         )
554     );
556     timeDict.readIfPresent("deltaT", deltaT_);
557     timeDict.readIfPresent("deltaT0", deltaT0_);
558     timeDict.readIfPresent("index", timeIndex_);
562 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
564     setTime(newTime.value(), newIndex);
568 void Foam::Time::setTime(const scalar newTime, const label newIndex)
570     value() = newTime;
571     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
572     timeIndex_ = newIndex;
576 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
578     setEndTime(endTime.value());
582 void Foam::Time::setEndTime(const scalar endTime)
584     endTime_ = endTime;
588 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
590     setDeltaT(deltaT.value());
594 void Foam::Time::setDeltaT(const scalar deltaT)
596     deltaT_ = deltaT;
597     deltaTchanged_ = true;
598     adjustDeltaT();
602 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
604     subCycling_ = true;
606     TimeState ts = *this;
607     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
608     deltaT_ /= nSubCycles;
609     deltaT0_ /= nSubCycles;
610     deltaTSave_ = deltaT0_;
612     return ts;
616 void Foam::Time::endSubCycle(const TimeState& ts)
618     subCycling_ = false;
619     TimeState::operator=(ts);
623 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
625 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
627     return operator+=(deltaT.value());
631 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
633     setDeltaT(deltaT);
634     return operator++();
638 Foam::Time& Foam::Time::operator++()
640     readModifiedObjects();
642     if (!subCycling_)
643     {
644         if (timeIndex_ == startTimeIndex_)
645         {
646             functionObjects_.start();
647         }
648         else
649         {
650             functionObjects_.execute();
651         }
652     }
654     deltaT0_ = deltaTSave_;
655     deltaTSave_ = deltaT_;
657     const word oldTimeName = dimensionedScalar::name();
659     setTime(value() + deltaT_, timeIndex_ + 1);
661     // If the time is very close to zero reset to zero
662     if (mag(value()) < 10*SMALL*deltaT_)
663     {
664         setTime(0.0, timeIndex_);
665     }
667     // Check that new time representation differs from old one
668     if (dimensionedScalar::name() == oldTimeName)
669     {
670         int oldPrecision = precision_;
671         do
672         {
673             precision_++;
674             setTime(value(), timeIndex());
675         }
676         while (precision_ < 100 && dimensionedScalar::name() == oldTimeName);
678         WarningIn("Time::operator++()")
679             << "Increased the timePrecision from " << oldPrecision
680             << " to " << precision_
681             << " to distinguish between timeNames at time " << value()
682             << endl;
683     }
685     switch (writeControl_)
686     {
687         case wcTimeStep:
688             outputTime_ = !(timeIndex_ % label(writeInterval_));
689         break;
691         case wcRunTime:
692         case wcAdjustableRunTime:
693         {
694             label outputIndex =
695                 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
697             if (outputIndex > outputTimeIndex_)
698             {
699                 outputTime_ = true;
700                 outputTimeIndex_ = outputIndex;
701             }
702             else
703             {
704                 outputTime_ = false;
705             }
706         }
707         break;
709         case wcCpuTime:
710         {
711             label outputIndex = label(elapsedCpuTime()/writeInterval_);
712             if (outputIndex > outputTimeIndex_)
713             {
714                 outputTime_ = true;
715                 outputTimeIndex_ = outputIndex;
716             }
717             else
718             {
719                 outputTime_ = false;
720             }
721         }
722         break;
724         case wcClockTime:
725         {
726             label outputIndex = label(elapsedClockTime()/writeInterval_);
727             if (outputIndex > outputTimeIndex_)
728             {
729                 outputTime_ = true;
730                 outputTimeIndex_ = outputIndex;
731             }
732             else
733             {
734                 outputTime_ = false;
735             }
736         }
737         break;
738     }
740     // see if endTime needs adjustment to stop at the next run()/end() check
741     if (!end())
742     {
743         if (stopAt_ == saNoWriteNow)
744         {
745             endTime_ = value();
746         }
747         else if (stopAt_ == saWriteNow)
748         {
749             endTime_ = value();
750             outputTime_ = true;
751         }
752         else if (stopAt_ == saNextWrite && outputTime_ == true)
753         {
754             endTime_ = value();
755         }
756     }
758     return *this;
762 Foam::Time& Foam::Time::operator++(int)
764     return operator++();
768 // ************************************************************************* //