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 WarningIn("Time::setControls()")
131 << " expected startTime, firstTime or latestTime"
132 << " found '" << startFrom
133 << "' in dictionary " << controlDict_.name() << nl
134 << " Setting time to " << startTime_ << endl;
138 setTime(startTime_, 0);
141 deltaTSave_ = deltaT_;
142 deltaT0_ = deltaTSave_;
144 if (Pstream::parRun())
146 scalar sumStartTime = startTime_;
147 reduce(sumStartTime, sumOp<scalar>());
150 mag(Pstream::nProcs()*startTime_ - sumStartTime)
151 > Pstream::nProcs()*deltaT_/10.0
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);
161 IOdictionary timeDict
169 IOobject::READ_IF_PRESENT,
175 if (timeDict.readIfPresent("deltaT", deltaTSave_))
177 deltaT0_ = deltaTSave_;
180 if (timeDict.readIfPresent("index", startTimeIndex_))
182 timeIndex_ = startTimeIndex_;
187 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191 const word& controlDictName,
192 const fileName& rootPath,
193 const fileName& caseName,
194 const word& systemName,
195 const word& constantName
206 objectRegistry(*this),
226 writeControl_(wcTimeStep),
227 writeInterval_(GREAT),
231 writeFormat_(IOstream::ASCII),
232 writeVersion_(IOstream::currentVersion),
233 writeCompression_(IOstream::UNCOMPRESSED),
235 runTimeModifiable_(true),
237 readLibs_(controlDict_, "libs"),
238 functionObjects_(*this)
246 const dictionary& dict,
247 const fileName& rootPath,
248 const fileName& caseName,
249 const word& systemName,
250 const word& constantName
261 objectRegistry(*this),
282 writeControl_(wcTimeStep),
283 writeInterval_(GREAT),
287 writeFormat_(IOstream::ASCII),
288 writeVersion_(IOstream::currentVersion),
289 writeCompression_(IOstream::UNCOMPRESSED),
291 runTimeModifiable_(true),
293 readLibs_(controlDict_, "libs"),
294 functionObjects_(*this)
302 const fileName& rootPath,
303 const fileName& caseName,
304 const word& systemName,
305 const word& constantName
316 objectRegistry(*this),
336 writeControl_(wcTimeStep),
337 writeInterval_(GREAT),
341 writeFormat_(IOstream::ASCII),
342 writeVersion_(IOstream::currentVersion),
343 writeCompression_(IOstream::UNCOMPRESSED),
345 runTimeModifiable_(true),
347 readLibs_(controlDict_, "libs"),
348 functionObjects_(*this)
352 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
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_);
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)
392 if (timeDirs[timeI] == t)
394 return timeDirs[timeI].name();
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)
412 if (t < timeDirs[1].value())
416 else if (t > timeDirs[timeDirs.size()-1].value())
418 return timeDirs[timeDirs.size()-1];
421 label nearestIndex = -1;
422 scalar deltaT = GREAT;
424 for (label timeI=1; timeI < timeDirs.size(); ++timeI)
426 scalar diff = mag(timeDirs[timeI].value() - t);
430 nearestIndex = timeI;
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
443 // instantList timeDirs = findTimes(path());
444 // label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
445 // return timeDirs[timeIndex];
448 Foam::label Foam::Time::findClosestTimeIndex
450 const instantList& timeDirs,
454 label nearestIndex = -1;
455 scalar deltaT = GREAT;
457 forAll(timeDirs, timeI)
459 if (timeDirs[timeI].name() == "constant") continue;
461 scalar diff = mag(timeDirs[timeI].value() - t);
465 nearestIndex = timeI;
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_);
497 // only execute when the condition is no longer true
498 // ie, when exiting the control loop
499 if (!running && timeIndex_ != startTimeIndex_)
501 // Note, end() also calls an indirect start() as required
502 functionObjects_.end();
510 bool Foam::Time::loop()
512 bool running = run();
523 bool Foam::Time::end() const
525 return value() > (endTime_ + 0.5*deltaT_);
529 void Foam::Time::setTime(const Time& t)
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
551 IOobject::READ_IF_PRESENT,
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)
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)
589 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
591 setDeltaT(deltaT.value());
595 void Foam::Time::setDeltaT(const scalar deltaT)
598 deltaTchanged_ = true;
603 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
607 TimeState ts = *this;
608 setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
609 deltaT_ /= nSubCycles;
610 deltaT0_ /= nSubCycles;
611 deltaTSave_ = deltaT0_;
617 void Foam::Time::endSubCycle(const TimeState& ts)
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)
639 Foam::Time& Foam::Time::operator++()
641 readModifiedObjects();
645 if (timeIndex_ == startTimeIndex_)
647 functionObjects_.start();
651 functionObjects_.execute();
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_)
662 setTime(0.0, timeIndex_);
665 switch (writeControl_)
668 outputTime_ = !(timeIndex_ % label(writeInterval_));
672 case wcAdjustableRunTime:
675 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
677 if (outputIndex > outputTimeIndex_)
680 outputTimeIndex_ = outputIndex;
691 label outputIndex = label(elapsedCpuTime()/writeInterval_);
692 if (outputIndex > outputTimeIndex_)
695 outputTimeIndex_ = outputIndex;
706 label outputIndex = label(elapsedClockTime()/writeInterval_);
707 if (outputIndex > outputTimeIndex_)
710 outputTimeIndex_ = outputIndex;
720 // see if endTime needs adjustment to stop at the next run()/end() check
723 if (stopAt_ == saNoWriteNow)
727 else if (stopAt_ == saWriteNow)
732 else if (stopAt_ == saNextWrite && outputTime_ == true)
742 Foam::Time& Foam::Time::operator++(int)
748 // ************************************************************************* //