1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2011 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
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
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 \*---------------------------------------------------------------------------*/
27 #include "PstreamReduceOps.H"
32 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(Foam::Time, 0);
39 const char* Foam::NamedEnum
41 Foam::Time::stopAtControls,
52 const char* Foam::NamedEnum
54 Foam::Time::writeControls,
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)
84 scalar timeToNextWrite = max
87 (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
90 scalar nSteps = timeToNextWrite/deltaT_ - SMALL;
92 // For tiny deltaT the label can overflow!
93 if (nSteps < labelMax)
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_)
103 deltaT_ = min(newDeltaT, 2.0*deltaT_);
107 deltaT_ = max(newDeltaT, 0.2*deltaT_);
114 void Foam::Time::setControls()
116 // default is to resume calculation from "latestTime"
117 const word startFrom = controlDict_.lookupOrDefault<word>
123 if (startFrom == "startTime")
125 controlDict_.lookup("startTime") >> startTime_;
129 // Search directory for valid time directories
130 instantList timeDirs = findTimes(path());
132 if (startFrom == "firstTime")
136 if (timeDirs[0].name() == constant() && timeDirs.size() >= 2)
138 startTime_ = timeDirs[1].value();
142 startTime_ = timeDirs[0].value();
146 else if (startFrom == "latestTime")
150 startTime_ = timeDirs.last().value();
155 FatalIOErrorIn("Time::setControls()", controlDict_)
156 << "expected startTime, firstTime or latestTime"
157 << " found '" << startFrom << "'"
158 << exit(FatalIOError);
162 setTime(startTime_, 0);
165 deltaTSave_ = deltaT_;
166 deltaT0_ = deltaTSave_;
168 if (Pstream::parRun())
170 scalar sumStartTime = startTime_;
171 reduce(sumStartTime, sumOp<scalar>());
174 mag(Pstream::nProcs()*startTime_ - sumStartTime)
175 > Pstream::nProcs()*deltaT_/10.0
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);
185 IOdictionary timeDict
193 IOobject::READ_IF_PRESENT,
199 if (timeDict.readIfPresent("deltaT", deltaTSave_))
201 deltaT0_ = deltaTSave_;
204 if (timeDict.readIfPresent("index", startTimeIndex_))
206 timeIndex_ = startTimeIndex_;
211 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
215 const word& controlDictName,
216 const fileName& rootPath,
217 const fileName& caseName,
218 const word& systemName,
219 const word& constantName,
220 const bool enableFunctionObjects
231 objectRegistry(*this),
242 IOobject::MUST_READ_IF_MODIFIED,
253 writeControl_(wcTimeStep),
254 writeInterval_(GREAT),
258 writeFormat_(IOstream::ASCII),
259 writeVersion_(IOstream::currentVersion),
260 writeCompression_(IOstream::UNCOMPRESSED),
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;
274 // Time objects not registered so do like objectRegistry::checkIn ourselves.
275 if (runTimeModifiable_)
281 regIOobject::fileModificationChecking == inotify
282 || regIOobject::fileModificationChecking == inotifyMaster
286 // File might not exist yet.
287 fileName f(controlDict_.filePath());
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();
297 controlDict_.watchIndex() = addWatch(f);
304 const word& controlDictName,
306 const word& systemName,
307 const word& constantName
318 objectRegistry(*this),
329 IOobject::MUST_READ_IF_MODIFIED,
340 writeControl_(wcTimeStep),
341 writeInterval_(GREAT),
345 writeFormat_(IOstream::ASCII),
346 writeVersion_(IOstream::currentVersion),
347 writeCompression_(IOstream::UNCOMPRESSED),
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;
361 // Time objects not registered so do like objectRegistry::checkIn ourselves.
362 if (runTimeModifiable_)
368 regIOobject::fileModificationChecking == inotify
369 || regIOobject::fileModificationChecking == inotifyMaster
373 // File might not exist yet.
374 fileName f(controlDict_.filePath());
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();
384 controlDict_.watchIndex() = addWatch(f);
391 const dictionary& dict,
392 const fileName& rootPath,
393 const fileName& caseName,
394 const word& systemName,
395 const word& constantName,
396 const bool enableFunctionObjects
407 objectRegistry(*this),
430 writeControl_(wcTimeStep),
431 writeInterval_(GREAT),
435 writeFormat_(IOstream::ASCII),
436 writeVersion_(IOstream::currentVersion),
437 writeCompression_(IOstream::UNCOMPRESSED),
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;
455 // Time objects not registered so do like objectRegistry::checkIn ourselves.
456 if (runTimeModifiable_)
462 regIOobject::fileModificationChecking == inotify
463 || regIOobject::fileModificationChecking == inotifyMaster
467 // File might not exist yet.
468 fileName f(controlDict_.filePath());
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();
478 controlDict_.watchIndex() = addWatch(f);
485 const fileName& rootPath,
486 const fileName& caseName,
487 const word& systemName,
488 const word& constantName,
489 const bool enableFunctionObjects
500 objectRegistry(*this),
522 writeControl_(wcTimeStep),
523 writeInterval_(GREAT),
527 writeFormat_(IOstream::ASCII),
528 writeVersion_(IOstream::currentVersion),
529 writeCompression_(IOstream::UNCOMPRESSED),
531 runTimeModifiable_(true),
533 functionObjects_(*this, enableFunctionObjects)
535 libs_.open(controlDict_, "libs");
539 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
543 if (controlDict_.watchIndex() != -1)
545 removeWatch(controlDict_.watchIndex());
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
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_);
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)
616 if (timeDirs[timeI] == t)
618 return timeDirs[timeI].name();
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)
636 if (t < timeDirs[1].value())
640 else if (t > timeDirs.last().value())
642 return timeDirs.last();
645 label nearestIndex = -1;
646 scalar deltaT = GREAT;
648 for (label timeI=1; timeI < timeDirs.size(); ++timeI)
650 scalar diff = mag(timeDirs[timeI].value() - t);
654 nearestIndex = timeI;
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
667 // instantList timeDirs = findTimes(path());
668 // label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
669 // return timeDirs[timeIndex];
672 Foam::label Foam::Time::findClosestTimeIndex
674 const instantList& timeDirs,
678 label nearestIndex = -1;
679 scalar deltaT = GREAT;
681 forAll(timeDirs, timeI)
683 if (timeDirs[timeI].name() == "constant") continue;
685 scalar diff = mag(timeDirs[timeI].value() - t);
689 nearestIndex = timeI;
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_);
721 // only execute when the condition is no longer true
722 // ie, when exiting the control loop
723 if (!running && timeIndex_ != startTimeIndex_)
725 // Note, end() also calls an indirect start() as required
726 functionObjects_.end();
734 const_cast<Time&>(*this).readModifiedObjects();
736 if (timeIndex_ == startTimeIndex_)
738 functionObjects_.start();
742 functionObjects_.execute();
746 // Update the "running" status following the
747 // possible side-effects from functionObjects
748 running = value() < (endTime_ - 0.5*deltaT_);
755 bool Foam::Time::loop()
757 bool running = run();
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);
782 controlDict_.lookup("endTime") >> endTime_;
792 void Foam::Time::setTime(const Time& t)
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
814 IOobject::READ_IF_PRESENT,
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)
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)
852 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
854 setDeltaT(deltaT.value());
858 void Foam::Time::setDeltaT(const scalar deltaT)
861 deltaTchanged_ = true;
866 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
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()
885 TimeState::operator=(prevTimeState());
886 prevTimeState_.clear();
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)
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);
918 // If the time is very close to zero reset to zero
919 if (mag(value()) < 10*SMALL*deltaT_)
921 setTime(0.0, timeIndex_);
926 // Check that new time representation differs from old one
927 if (dimensionedScalar::name() == oldTimeName)
929 int oldPrecision = precision_;
933 setTime(value(), timeIndex());
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()
947 switch (writeControl_)
950 outputTime_ = !(timeIndex_ % label(writeInterval_));
954 case wcAdjustableRunTime:
956 label outputIndex = label
958 ((value() - startTime_) + 0.5*deltaT_)
962 if (outputIndex > outputTimeIndex_)
965 outputTimeIndex_ = outputIndex;
976 label outputIndex = label
978 returnReduce(elapsedCpuTime(), maxOp<double>())
981 if (outputIndex > outputTimeIndex_)
984 outputTimeIndex_ = outputIndex;
995 label outputIndex = label
997 returnReduce(label(elapsedClockTime()), maxOp<label>())
1000 if (outputIndex > outputTimeIndex_)
1003 outputTimeIndex_ = outputIndex;
1007 outputTime_ = false;
1013 // see if endTime needs adjustment to stop at the next run()/end() check
1016 if (stopAt_ == saNoWriteNow)
1020 else if (stopAt_ == saWriteNow)
1025 else if (stopAt_ == saNextWrite && outputTime_ == true)
1036 Foam::Time& Foam::Time::operator++(int)
1038 return operator++();
1042 // ************************************************************************* //