FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

Time.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2011 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "Time.H"
00027 #include <OpenFOAM/PstreamReduceOps.H>
00028 
00029 #include <sstream>
00030 
00031 // * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * * * //
00032 
00033 defineTypeNameAndDebug(Foam::Time, 0);
00034 
00035 template<>
00036 const char* Foam::NamedEnum<Foam::Time::stopAtControls, 4>::names[] =
00037 {
00038     "endTime",
00039     "noWriteNow",
00040     "writeNow",
00041     "nextWrite"
00042 };
00043 
00044 const Foam::NamedEnum<Foam::Time::stopAtControls, 4>
00045     Foam::Time::stopAtControlNames_;
00046 
00047 template<>
00048 const char* Foam::NamedEnum<Foam::Time::writeControls, 5>::names[] =
00049 {
00050     "timeStep",
00051     "runTime",
00052     "adjustableRunTime",
00053     "clockTime",
00054     "cpuTime"
00055 };
00056 
00057 const Foam::NamedEnum<Foam::Time::writeControls, 5>
00058     Foam::Time::writeControlNames_;
00059 
00060 Foam::Time::fmtflags Foam::Time::format_(Foam::Time::general);
00061 int Foam::Time::precision_(6);
00062 
00063 Foam::word Foam::Time::controlDictName("controlDict");
00064 
00065 
00066 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00067 
00068 void Foam::Time::adjustDeltaT()
00069 {
00070     if (writeControl_ == wcAdjustableRunTime)
00071     {
00072         scalar timeToNextWrite = max
00073         (
00074             0.0,
00075             (outputTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
00076         );
00077 
00078         scalar nSteps = timeToNextWrite/deltaT_ - SMALL;
00079 
00080         // For tiny deltaT the label can overflow!
00081         if (nSteps < labelMax)
00082         {
00083             label nStepsToNextWrite = label(nSteps) + 1;
00084 
00085             scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;
00086 
00087             // Control the increase of the time step to within a factor of 2
00088             // and the decrease within a factor of 5.
00089             if (newDeltaT >= deltaT_)
00090             {
00091                 deltaT_ = min(newDeltaT, 2.0*deltaT_);
00092             }
00093             else
00094             {
00095                 deltaT_ = max(newDeltaT, 0.2*deltaT_);
00096             }
00097         }
00098     }
00099 }
00100 
00101 
00102 void Foam::Time::setControls()
00103 {
00104     // default is to resume calculation from "latestTime"
00105     word startFrom = controlDict_.lookupOrDefault<word>
00106     (
00107         "startFrom",
00108         "latestTime"
00109     );
00110 
00111     if (startFrom == "startTime")
00112     {
00113         controlDict_.lookup("startTime") >> startTime_;
00114     }
00115     else
00116     {
00117         // Search directory for valid time directories
00118         instantList timeDirs = findTimes(path());
00119 
00120         if (startFrom == "firstTime")
00121         {
00122             if (timeDirs.size())
00123             {
00124                 startTime_ = timeDirs[0].value();
00125             }
00126         }
00127         else if (startFrom == "latestTime")
00128         {
00129             if (timeDirs.size())
00130             {
00131                 startTime_ = timeDirs[timeDirs.size()-1].value();
00132             }
00133         }
00134         else
00135         {
00136             FatalIOErrorIn("Time::setControls()", controlDict_)
00137                 << "expected startTime, firstTime or latestTime"
00138                 << " found '" << startFrom << "'"
00139                 << exit(FatalIOError);
00140         }
00141     }
00142 
00143     setTime(startTime_, 0);
00144 
00145     readDict();
00146     deltaTSave_ = deltaT_;
00147     deltaT0_ = deltaTSave_;
00148 
00149     if (Pstream::parRun())
00150     {
00151         scalar sumStartTime = startTime_;
00152         reduce(sumStartTime, sumOp<scalar>());
00153         if
00154         (
00155             mag(Pstream::nProcs()*startTime_ - sumStartTime)
00156           > Pstream::nProcs()*deltaT_/10.0
00157         )
00158         {
00159             FatalIOErrorIn("Time::setControls()", controlDict_)
00160                 << "Start time is not the same for all processors" << nl
00161                 << "processor " << Pstream::myProcNo() << " has startTime "
00162                 << startTime_ << exit(FatalIOError);
00163         }
00164     }
00165 
00166     IOdictionary timeDict
00167     (
00168         IOobject
00169         (
00170             "time",
00171             timeName(),
00172             "uniform",
00173             *this,
00174             IOobject::READ_IF_PRESENT,
00175             IOobject::NO_WRITE,
00176             false
00177         )
00178     );
00179 
00180     if (timeDict.readIfPresent("deltaT", deltaTSave_))
00181     {
00182         deltaT0_ = deltaTSave_;
00183     }
00184 
00185     if (timeDict.readIfPresent("index", startTimeIndex_))
00186     {
00187         timeIndex_ = startTimeIndex_;
00188     }
00189 }
00190 
00191 
00192 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00193 
00194 Foam::Time::Time
00195 (
00196     const word& controlDictName,
00197     const fileName& rootPath,
00198     const fileName& caseName,
00199     const word& systemName,
00200     const word& constantName
00201 )
00202 :
00203     TimePaths
00204     (
00205         rootPath,
00206         caseName,
00207         systemName,
00208         constantName
00209     ),
00210 
00211     objectRegistry(*this),
00212 
00213     controlDict_
00214     (
00215         IOobject
00216         (
00217             controlDictName,
00218             system(),
00219             *this,
00220             IOobject::MUST_READ,
00221             IOobject::NO_WRITE,
00222             false
00223         )
00224     ),
00225 
00226     startTimeIndex_(0),
00227     startTime_(0),
00228     endTime_(0),
00229 
00230     stopAt_(saEndTime),
00231     writeControl_(wcTimeStep),
00232     writeInterval_(GREAT),
00233     purgeWrite_(0),
00234     subCycling_(false),
00235 
00236     writeFormat_(IOstream::ASCII),
00237     writeVersion_(IOstream::currentVersion),
00238     writeCompression_(IOstream::UNCOMPRESSED),
00239     graphFormat_("raw"),
00240     runTimeModifiable_(true),
00241 
00242     readLibs_(controlDict_, "libs"),
00243     functionObjects_(*this)
00244 {
00245     setControls();
00246 }
00247 
00248 
00249 Foam::Time::Time
00250 (
00251     const dictionary& dict,
00252     const fileName& rootPath,
00253     const fileName& caseName,
00254     const word& systemName,
00255     const word& constantName
00256 )
00257 :
00258     TimePaths
00259     (
00260         rootPath,
00261         caseName,
00262         systemName,
00263         constantName
00264     ),
00265 
00266     objectRegistry(*this),
00267 
00268     controlDict_
00269     (
00270         IOobject
00271         (
00272             controlDictName,
00273             system(),
00274             *this,
00275             IOobject::NO_READ,
00276             IOobject::NO_WRITE,
00277             false
00278         ),
00279         dict
00280     ),
00281 
00282     startTimeIndex_(0),
00283     startTime_(0),
00284     endTime_(0),
00285 
00286     stopAt_(saEndTime),
00287     writeControl_(wcTimeStep),
00288     writeInterval_(GREAT),
00289     purgeWrite_(0),
00290     subCycling_(false),
00291 
00292     writeFormat_(IOstream::ASCII),
00293     writeVersion_(IOstream::currentVersion),
00294     writeCompression_(IOstream::UNCOMPRESSED),
00295     graphFormat_("raw"),
00296     runTimeModifiable_(true),
00297 
00298     readLibs_(controlDict_, "libs"),
00299     functionObjects_(*this)
00300 {
00301     setControls();
00302 }
00303 
00304 
00305 Foam::Time::Time
00306 (
00307     const fileName& rootPath,
00308     const fileName& caseName,
00309     const word& systemName,
00310     const word& constantName
00311 )
00312 :
00313     TimePaths
00314     (
00315         rootPath,
00316         caseName,
00317         systemName,
00318         constantName
00319     ),
00320 
00321     objectRegistry(*this),
00322 
00323     controlDict_
00324     (
00325         IOobject
00326         (
00327             controlDictName,
00328             system(),
00329             *this,
00330             IOobject::NO_READ,
00331             IOobject::NO_WRITE,
00332             false
00333         )
00334     ),
00335 
00336     startTimeIndex_(0),
00337     startTime_(0),
00338     endTime_(0),
00339 
00340     stopAt_(saEndTime),
00341     writeControl_(wcTimeStep),
00342     writeInterval_(GREAT),
00343     purgeWrite_(0),
00344     subCycling_(false),
00345 
00346     writeFormat_(IOstream::ASCII),
00347     writeVersion_(IOstream::currentVersion),
00348     writeCompression_(IOstream::UNCOMPRESSED),
00349     graphFormat_("raw"),
00350     runTimeModifiable_(true),
00351 
00352     readLibs_(controlDict_, "libs"),
00353     functionObjects_(*this)
00354 {}
00355 
00356 
00357 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00358 
00359 Foam::Time::~Time()
00360 {
00361     // destroy function objects first
00362     functionObjects_.clear();
00363 }
00364 
00365 
00366 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00367 
00368 Foam::word Foam::Time::timeName(const scalar t)
00369 {
00370     std::ostringstream buf;
00371     buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
00372     buf.precision(precision_);
00373     buf << t;
00374     return buf.str();
00375 }
00376 
00377 
00378 Foam::word Foam::Time::timeName() const
00379 {
00380     return dimensionedScalar::name();
00381 }
00382 
00383 
00384 // Search the construction path for times
00385 Foam::instantList Foam::Time::times() const
00386 {
00387     return findTimes(path());
00388 }
00389 
00390 
00391 Foam::word Foam::Time::findInstancePath(const instant& t) const
00392 {
00393     instantList timeDirs = findTimes(path());
00394 
00395     forAllReverse(timeDirs, timeI)
00396     {
00397         if (timeDirs[timeI] == t)
00398         {
00399             return timeDirs[timeI].name();
00400         }
00401     }
00402 
00403     return word::null;
00404 }
00405 
00406 
00407 Foam::instant Foam::Time::findClosestTime(const scalar t) const
00408 {
00409     instantList timeDirs = findTimes(path());
00410 
00411     // there is only one time (likely "constant") so return it
00412     if (timeDirs.size() == 1)
00413     {
00414         return timeDirs[0];
00415     }
00416 
00417     if (t < timeDirs[1].value())
00418     {
00419         return timeDirs[1];
00420     }
00421     else if (t > timeDirs[timeDirs.size()-1].value())
00422     {
00423         return timeDirs[timeDirs.size()-1];
00424     }
00425 
00426     label nearestIndex = -1;
00427     scalar deltaT = GREAT;
00428 
00429     for (label timeI=1; timeI < timeDirs.size(); ++timeI)
00430     {
00431         scalar diff = mag(timeDirs[timeI].value() - t);
00432         if (diff < deltaT)
00433         {
00434             deltaT = diff;
00435             nearestIndex = timeI;
00436         }
00437     }
00438 
00439     return timeDirs[nearestIndex];
00440 }
00441 
00442 
00443 // This should work too,
00444 // if we don't worry about checking "constant" explicitly
00445 //
00446 // Foam::instant Foam::Time::findClosestTime(const scalar t) const
00447 // {
00448 //     instantList timeDirs = findTimes(path());
00449 //     label timeIndex = min(findClosestTimeIndex(timeDirs, t), 0);
00450 //     return timeDirs[timeIndex];
00451 // }
00452 
00453 Foam::label Foam::Time::findClosestTimeIndex
00454 (
00455     const instantList& timeDirs,
00456     const scalar t
00457 )
00458 {
00459     label nearestIndex = -1;
00460     scalar deltaT = GREAT;
00461 
00462     forAll(timeDirs, timeI)
00463     {
00464         if (timeDirs[timeI].name() == "constant") continue;
00465 
00466         scalar diff = mag(timeDirs[timeI].value() - t);
00467         if (diff < deltaT)
00468         {
00469             deltaT = diff;
00470             nearestIndex = timeI;
00471         }
00472     }
00473 
00474     return nearestIndex;
00475 }
00476 
00477 
00478 Foam::label Foam::Time::startTimeIndex() const
00479 {
00480     return startTimeIndex_;
00481 }
00482 
00483 
00484 Foam::dimensionedScalar Foam::Time::startTime() const
00485 {
00486     return dimensionedScalar("startTime", dimTime, startTime_);
00487 }
00488 
00489 
00490 Foam::dimensionedScalar Foam::Time::endTime() const
00491 {
00492     return dimensionedScalar("endTime", dimTime, endTime_);
00493 }
00494 
00495 
00496 bool Foam::Time::run() const
00497 {
00498     bool running = value() < (endTime_ - 0.5*deltaT_);
00499 
00500     if (!subCycling_)
00501     {
00502         // only execute when the condition is no longer true
00503         // ie, when exiting the control loop
00504         if (!running && timeIndex_ != startTimeIndex_)
00505         {
00506             // Note, end() also calls an indirect start() as required
00507             functionObjects_.end();
00508         }
00509     }
00510 
00511     return running;
00512 }
00513 
00514 
00515 bool Foam::Time::loop()
00516 {
00517     bool running = run();
00518 
00519     if (running)
00520     {
00521         operator++();
00522     }
00523 
00524     return running;
00525 }
00526 
00527 
00528 bool Foam::Time::end() const
00529 {
00530     return value() > (endTime_ + 0.5*deltaT_);
00531 }
00532 
00533 
00534 void Foam::Time::setTime(const Time& t)
00535 {
00536     value() = t.value();
00537     dimensionedScalar::name() = t.dimensionedScalar::name();
00538     timeIndex_ = t.timeIndex_;
00539 }
00540 
00541 
00542 void Foam::Time::setTime(const instant& inst, const label newIndex)
00543 {
00544     value() = inst.value();
00545     dimensionedScalar::name() = inst.name();
00546     timeIndex_ = newIndex;
00547 
00548     IOdictionary timeDict
00549     (
00550         IOobject
00551         (
00552             "time",
00553             timeName(),
00554             "uniform",
00555             *this,
00556             IOobject::READ_IF_PRESENT,
00557             IOobject::NO_WRITE,
00558             false
00559         )
00560     );
00561 
00562     timeDict.readIfPresent("deltaT", deltaT_);
00563     timeDict.readIfPresent("deltaT0", deltaT0_);
00564     timeDict.readIfPresent("index", timeIndex_);
00565 }
00566 
00567 
00568 void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex)
00569 {
00570     setTime(newTime.value(), newIndex);
00571 }
00572 
00573 
00574 void Foam::Time::setTime(const scalar newTime, const label newIndex)
00575 {
00576     value() = newTime;
00577     dimensionedScalar::name() = timeName(timeToUserTime(newTime));
00578     timeIndex_ = newIndex;
00579 }
00580 
00581 
00582 void Foam::Time::setEndTime(const dimensionedScalar& endTime)
00583 {
00584     setEndTime(endTime.value());
00585 }
00586 
00587 
00588 void Foam::Time::setEndTime(const scalar endTime)
00589 {
00590     endTime_ = endTime;
00591 }
00592 
00593 
00594 void Foam::Time::setDeltaT(const dimensionedScalar& deltaT)
00595 {
00596     setDeltaT(deltaT.value());
00597 }
00598 
00599 
00600 void Foam::Time::setDeltaT(const scalar deltaT)
00601 {
00602     deltaT_ = deltaT;
00603     deltaTchanged_ = true;
00604     adjustDeltaT();
00605 }
00606 
00607 
00608 Foam::TimeState Foam::Time::subCycle(const label nSubCycles)
00609 {
00610     subCycling_ = true;
00611     prevTimeState_.set(new TimeState(*this));
00612 
00613     setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles);
00614     deltaT_ /= nSubCycles;
00615     deltaT0_ /= nSubCycles;
00616     deltaTSave_ = deltaT0_;
00617 
00618     return prevTimeState();
00619 }
00620 
00621 
00622 void Foam::Time::endSubCycle()
00623 {
00624     if (subCycling_)
00625     {
00626         subCycling_ = false;
00627         TimeState::operator=(prevTimeState());
00628         prevTimeState_.clear();
00629     }
00630 }
00631 
00632 
00633 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00634 
00635 Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT)
00636 {
00637     return operator+=(deltaT.value());
00638 }
00639 
00640 
00641 Foam::Time& Foam::Time::operator+=(const scalar deltaT)
00642 {
00643     setDeltaT(deltaT);
00644     return operator++();
00645 }
00646 
00647 
00648 Foam::Time& Foam::Time::operator++()
00649 {
00650     readModifiedObjects();
00651 
00652     if (!subCycling_)
00653     {
00654         if (timeIndex_ == startTimeIndex_)
00655         {
00656             functionObjects_.start();
00657         }
00658         else
00659         {
00660             functionObjects_.execute();
00661         }
00662     }
00663 
00664     deltaT0_ = deltaTSave_;
00665     deltaTSave_ = deltaT_;
00666 
00667     const word oldTimeName = dimensionedScalar::name();
00668 
00669     setTime(value() + deltaT_, timeIndex_ + 1);
00670 
00671     // If the time is very close to zero reset to zero
00672     if (mag(value()) < 10*SMALL*deltaT_)
00673     {
00674         setTime(0.0, timeIndex_);
00675     }
00676 
00677     // Check that new time representation differs from old one
00678     if (dimensionedScalar::name() == oldTimeName)
00679     {
00680         int oldPrecision = precision_;
00681         do
00682         {
00683             precision_++;
00684             setTime(value(), timeIndex());
00685         }
00686         while (precision_ < 100 && dimensionedScalar::name() == oldTimeName);
00687 
00688         WarningIn("Time::operator++()")
00689             << "Increased the timePrecision from " << oldPrecision
00690             << " to " << precision_
00691             << " to distinguish between timeNames at time " << value()
00692             << endl;
00693     }
00694 
00695     switch (writeControl_)
00696     {
00697         case wcTimeStep:
00698             outputTime_ = !(timeIndex_ % label(writeInterval_));
00699         break;
00700 
00701         case wcRunTime:
00702         case wcAdjustableRunTime:
00703         {
00704             label outputIndex =
00705                 label(((value() - startTime_) + 0.5*deltaT_)/writeInterval_);
00706 
00707             if (outputIndex > outputTimeIndex_)
00708             {
00709                 outputTime_ = true;
00710                 outputTimeIndex_ = outputIndex;
00711             }
00712             else
00713             {
00714                 outputTime_ = false;
00715             }
00716         }
00717         break;
00718 
00719         case wcCpuTime:
00720         {
00721             label outputIndex = label
00722             (
00723                 returnReduce(elapsedCpuTime(), maxOp<double>())
00724               / writeInterval_
00725             );
00726             if (outputIndex > outputTimeIndex_)
00727             {
00728                 outputTime_ = true;
00729                 outputTimeIndex_ = outputIndex;
00730             }
00731             else
00732             {
00733                 outputTime_ = false;
00734             }
00735         }
00736         break;
00737 
00738         case wcClockTime:
00739         {
00740             label outputIndex = label
00741             (
00742                 returnReduce(label(elapsedClockTime()), maxOp<label>())
00743               / writeInterval_
00744             );
00745             if (outputIndex > outputTimeIndex_)
00746             {
00747                 outputTime_ = true;
00748                 outputTimeIndex_ = outputIndex;
00749             }
00750             else
00751             {
00752                 outputTime_ = false;
00753             }
00754         }
00755         break;
00756     }
00757 
00758     // see if endTime needs adjustment to stop at the next run()/end() check
00759     if (!end())
00760     {
00761         if (stopAt_ == saNoWriteNow)
00762         {
00763             endTime_ = value();
00764         }
00765         else if (stopAt_ == saWriteNow)
00766         {
00767             endTime_ = value();
00768             outputTime_ = true;
00769         }
00770         else if (stopAt_ == saNextWrite && outputTime_ == true)
00771         {
00772             endTime_ = value();
00773         }
00774     }
00775 
00776     return *this;
00777 }
00778 
00779 
00780 Foam::Time& Foam::Time::operator++(int)
00781 {
00782     return operator++();
00783 }
00784 
00785 
00786 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines