1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2008 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 * * * * * * * * * * * * * * //
36 defineTypeNameAndDebug(Time, 0);
41 const char* Foam::NamedEnum<Foam::Time::stopAtControls, 4>::names[] =
49 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
50 Foam::Time::stopAtControlNames_;
53 const char* Foam::NamedEnum<Foam::Time::writeControls, 5>::names[] =
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)
77 scalar timeToNextWrite = max
80 (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
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_)
90 deltaT_ = min(newDeltaT, 2.0*deltaT_);
94 deltaT_ = max(newDeltaT, 0.2*deltaT_);
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")
109 controlDict_.lookup("startTime") >> startTime_;
113 // Search directory for valid time directories
114 instantList Times = findTimes(path());
116 if (startFrom == "firstTime")
120 startTime_ = Times[0].value();
123 else if (startFrom == "latestTime")
127 startTime_ = Times[Times.size()-1].value();
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;
140 setTime(startTime_, 0);
143 deltaTSave_ = deltaT_;
144 deltaT0_ = deltaTSave_;
146 if (Pstream::parRun())
148 scalar sumStartTime = startTime_;
149 reduce(sumStartTime, sumOp<scalar>());
152 mag(Pstream::nProcs()*startTime_ - sumStartTime)
153 > Pstream::nProcs()*deltaT_/10.0
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);
163 IOdictionary timeDict
171 IOobject::READ_IF_PRESENT,
177 if (timeDict.readIfPresent("deltaT", deltaTSave_))
179 deltaT0_ = deltaTSave_;
182 if (timeDict.readIfPresent("index", startTimeIndex_))
184 timeIndex_ = startTimeIndex_;
189 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193 const word& controlDictName,
194 const fileName& rootPath,
195 const fileName& caseName,
196 const word& systemName,
197 const word& constantName
208 objectRegistry(*this),
228 writeControl_(wcTimeStep),
229 writeInterval_(GREAT),
233 writeFormat_(IOstream::ASCII),
234 writeVersion_(IOstream::currentVersion),
235 writeCompression_(IOstream::UNCOMPRESSED),
237 runTimeModifiable_(true),
239 readLibs_(controlDict_, "libs"),
240 functionObjects_(*this)
248 const dictionary& dict,
249 const fileName& rootPath,
250 const fileName& caseName,
251 const word& systemName,
252 const word& constantName
263 objectRegistry(*this),
284 writeControl_(wcTimeStep),
285 writeInterval_(GREAT),
289 writeFormat_(IOstream::ASCII),
290 writeVersion_(IOstream::currentVersion),
291 writeCompression_(IOstream::UNCOMPRESSED),
293 runTimeModifiable_(true),
295 readLibs_(controlDict_, "libs"),
296 functionObjects_(*this)
304 const fileName& rootPath,
305 const fileName& caseName,
306 const word& systemName,
307 const word& constantName
318 objectRegistry(*this),
338 writeControl_(wcTimeStep),
339 writeInterval_(GREAT),
343 writeFormat_(IOstream::ASCII),
344 writeVersion_(IOstream::currentVersion),
345 writeCompression_(IOstream::UNCOMPRESSED),
347 runTimeModifiable_(true),
349 readLibs_(controlDict_, "libs"),
350 functionObjects_(*this)
354 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
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_);
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)
394 return times[i].name();
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)
412 if (t < times[1].value())
416 else if (t > times[times.size() - 1].value())
418 return times[times.size() - 1];
421 label nearestIndex = -1;
422 scalar deltaT = GREAT;
424 for (label i=1; i<times.size(); i++)
426 scalar diff = mag(times[i].value() - t);
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
443 // instantList times = Time::findTimes(path());
444 // label timeIndex = min(findClosestTimeIndex(times, t), 0);
445 // return times[timeIndex];
448 Foam::label Foam::Time::findClosestTimeIndex
450 const instantList& times,
454 label nearestIndex = -1;
455 scalar deltaT = GREAT;
459 if (times[i].name() == "constant") continue;
461 scalar diff = fabs(times[i].value() - t);
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_)
497 const_cast<functionObjectList&>(functionObjects_).execute();
504 bool Foam::Time::end() const
506 return (value() > (endTime_ + 0.5*deltaT_));
510 void Foam::Time::setTime(const Time& t)
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
532 IOobject::READ_IF_PRESENT,
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)
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)
570 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
572 setDeltaT(deltaT.value());
576 void Foam::Time::setDeltaT(const scalar deltaT)
579 deltaTchanged_ = true;
584 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
588 TimeState ts = *this;
589 setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
590 deltaT_ /= nSubCycles;
591 deltaT0_ /= nSubCycles;
592 deltaTSave_ = deltaT0_;
598 void Foam::Time::endSubCycle(const TimeState& ts)
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();
619 if (timeIndex_ == startTimeIndex_)
621 functionObjects_.start();
625 functionObjects_.execute();
636 Foam::Time& Foam::Time::operator++()
638 readModifiedObjects();
642 if (timeIndex_ == startTimeIndex_)
644 functionObjects_.start();
648 functionObjects_.execute();
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_)
659 setTime(0.0, timeIndex_);
662 switch(writeControl_)
665 outputTime_ = !(timeIndex_%label(writeInterval_));
669 case wcAdjustableRunTime:
671 label outputTimeIndex =
672 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
674 if (outputTimeIndex > outputTimeIndex_)
677 outputTimeIndex_ = outputTimeIndex;
688 label outputTimeIndex =
689 label(elapsedCpuTime()/writeInterval_);
691 if (outputTimeIndex > outputTimeIndex_)
694 outputTimeIndex_ = outputTimeIndex;
705 label outputTimeIndex = label(elapsedClockTime()/writeInterval_);
706 if (outputTimeIndex > outputTimeIndex_)
709 outputTimeIndex_ = outputTimeIndex;
721 if (stopAt_ == saNoWriteNow)
725 else if (stopAt_ == saWriteNow)
730 else if (stopAt_ == saNextWrite && outputTime_ == true)
740 Foam::Time& Foam::Time::operator++(int)
746 // ************************************************************************* //