1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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 \*---------------------------------------------------------------------------*/
28 #include "PstreamReduceOps.H"
32 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(Foam::Time, 0);
37 const char* Foam::NamedEnum<Foam::Time::stopAtControls, 4>::names[] =
45 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
46 Foam::Time::stopAtControlNames_;
49 const char* Foam::NamedEnum<Foam::Time::writeControls, 5>::names[] =
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)
73 scalar timeToNextWrite = max
76 (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
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_)
86 deltaT_ = min(newDeltaT, 2.0*deltaT_);
90 deltaT_ = max(newDeltaT, 0.2*deltaT_);
96 void Foam::Time::setControls()
98 // default is to resume calculation from "latestTime"
99 word startFrom = controlDict_.lookupOrDefault<word>
105 if (startFrom == "startTime")
107 controlDict_.lookup("startTime") >> startTime_;
111 // Search directory for valid time directories
112 instantList timeDirs = findTimes(path());
114 if (startFrom == "firstTime")
118 startTime_ = timeDirs[0].value();
121 else if (startFrom == "latestTime")
125 startTime_ = timeDirs[timeDirs.size()-1].value();
130 FatalIOErrorIn("Time::setControls()", controlDict_)
131 << "expected startTime, firstTime or latestTime"
132 << " found '" << startFrom << "'"
133 << exit(FatalIOError);
137 setTime(startTime_, 0);
140 deltaTSave_ = deltaT_;
141 deltaT0_ = deltaTSave_;
143 if (Pstream::parRun())
145 scalar sumStartTime = startTime_;
146 reduce(sumStartTime, sumOp<scalar>());
149 mag(Pstream::nProcs()*startTime_ - sumStartTime)
150 > Pstream::nProcs()*deltaT_/10.0
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);
160 IOdictionary timeDict
168 IOobject::READ_IF_PRESENT,
174 if (timeDict.readIfPresent("deltaT", deltaTSave_))
176 deltaT0_ = deltaTSave_;
179 if (timeDict.readIfPresent("index", startTimeIndex_))
181 timeIndex_ = startTimeIndex_;
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
190 const word& controlDictName,
191 const fileName& rootPath,
192 const fileName& caseName,
193 const word& systemName,
194 const word& constantName
205 objectRegistry(*this),
225 writeControl_(wcTimeStep),
226 writeInterval_(GREAT),
230 writeFormat_(IOstream::ASCII),
231 writeVersion_(IOstream::currentVersion),
232 writeCompression_(IOstream::UNCOMPRESSED),
234 runTimeModifiable_(true),
236 readLibs_(controlDict_, "libs"),
237 functionObjects_(*this)
245 const dictionary& dict,
246 const fileName& rootPath,
247 const fileName& caseName,
248 const word& systemName,
249 const word& constantName
260 objectRegistry(*this),
281 writeControl_(wcTimeStep),
282 writeInterval_(GREAT),
286 writeFormat_(IOstream::ASCII),
287 writeVersion_(IOstream::currentVersion),
288 writeCompression_(IOstream::UNCOMPRESSED),
290 runTimeModifiable_(true),
292 readLibs_(controlDict_, "libs"),
293 functionObjects_(*this)
301 const fileName& rootPath,
302 const fileName& caseName,
303 const word& systemName,
304 const word& constantName
315 objectRegistry(*this),
335 writeControl_(wcTimeStep),
336 writeInterval_(GREAT),
340 writeFormat_(IOstream::ASCII),
341 writeVersion_(IOstream::currentVersion),
342 writeCompression_(IOstream::UNCOMPRESSED),
344 runTimeModifiable_(true),
346 readLibs_(controlDict_, "libs"),
347 functionObjects_(*this)
351 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
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_);
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)
391 if (timeDirs[timeI] == t)
393 return timeDirs[timeI].name();
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)
411 if (t < timeDirs[1].value())
415 else if (t > timeDirs[timeDirs.size()-1].value())
417 return timeDirs[timeDirs.size()-1];
420 label nearestIndex = -1;
421 scalar deltaT = GREAT;
423 for (label timeI=1; timeI < timeDirs.size(); ++timeI)
425 scalar diff = mag(timeDirs[timeI].value() - t);
429 nearestIndex = timeI;
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
442 // instantList timeDirs = findTimes(path());
443 // label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
444 // return timeDirs[timeIndex];
447 Foam::label Foam::Time::findClosestTimeIndex
449 const instantList& timeDirs,
453 label nearestIndex = -1;
454 scalar deltaT = GREAT;
456 forAll(timeDirs, timeI)
458 if (timeDirs[timeI].name() == "constant") continue;
460 scalar diff = mag(timeDirs[timeI].value() - t);
464 nearestIndex = timeI;
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_);
496 // only execute when the condition is no longer true
497 // ie, when exiting the control loop
498 if (!running && timeIndex_ != startTimeIndex_)
500 // Note, end() also calls an indirect start() as required
501 functionObjects_.end();
509 bool Foam::Time::loop()
511 bool running = run();
522 bool Foam::Time::end() const
524 return value() > (endTime_ + 0.5*deltaT_);
528 void Foam::Time::setTime(const Time& t)
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
550 IOobject::READ_IF_PRESENT,
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)
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)
588 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
590 setDeltaT(deltaT.value());
594 void Foam::Time::setDeltaT(const scalar deltaT)
597 deltaTchanged_ = true;
602 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
606 TimeState ts = *this;
607 setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
608 deltaT_ /= nSubCycles;
609 deltaT0_ /= nSubCycles;
610 deltaTSave_ = deltaT0_;
616 void Foam::Time::endSubCycle(const TimeState& ts)
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)
638 Foam::Time& Foam::Time::operator++()
640 readModifiedObjects();
644 if (timeIndex_ == startTimeIndex_)
646 functionObjects_.start();
650 functionObjects_.execute();
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_)
664 setTime(0.0, timeIndex_);
667 // Check that new time representation differs from old one
668 if (dimensionedScalar::name() == oldTimeName)
670 int oldPrecision = precision_;
674 setTime(value(), timeIndex());
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()
685 switch (writeControl_)
688 outputTime_ = !(timeIndex_ % label(writeInterval_));
692 case wcAdjustableRunTime:
695 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
697 if (outputIndex > outputTimeIndex_)
700 outputTimeIndex_ = outputIndex;
711 label outputIndex = label(elapsedCpuTime()/writeInterval_);
712 if (outputIndex > outputTimeIndex_)
715 outputTimeIndex_ = outputIndex;
726 label outputIndex = label(elapsedClockTime()/writeInterval_);
727 if (outputIndex > outputTimeIndex_)
730 outputTimeIndex_ = outputIndex;
740 // see if endTime needs adjustment to stop at the next run()/end() check
743 if (stopAt_ == saNoWriteNow)
747 else if (stopAt_ == saWriteNow)
752 else if (stopAt_ == saNextWrite && outputTime_ == true)
762 Foam::Time& Foam::Time::operator++(int)
768 // ************************************************************************* //