Skip to content
Snippets Groups Projects
Time.C 33 KiB
Newer Older
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
OpenFOAM bot's avatar
OpenFOAM bot committed
    \\  /    A nd           | www.openfoam.com
     \\/     M anipulation  |
-------------------------------------------------------------------------------
OpenFOAM bot's avatar
OpenFOAM bot committed
    Copyright (C) 2011-2017 OpenFOAM Foundation
    Copyright (C) 2015-2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

\*---------------------------------------------------------------------------*/

#include "Time.H"
#include "PstreamReduceOps.H"
#include "argList.H"
#include "HashSet.H"
#include "profiling.H"
#include <sstream>

// * * * * * * * * * * * * * Static Member Data  * * * * * * * * * * * * * * //

    defineTypeNameAndDebug(Time, 0);
const Foam::Enum
<
    Foam::Time::stopAtControls
>
Mark OLESEN's avatar
Mark OLESEN committed
Foam::Time::stopAtControlNames
    { stopAtControls::saEndTime, "endTime" },
    { stopAtControls::saNoWriteNow, "noWriteNow" },
    { stopAtControls::saWriteNow, "writeNow" },
    { stopAtControls::saNextWrite, "nextWrite" },
    // Leave saUnknown untabulated - fallback to flag unknown settings
Mark OLESEN's avatar
Mark OLESEN committed
});


const Foam::Enum
<
    Foam::Time::writeControls
>
Mark OLESEN's avatar
Mark OLESEN committed
Foam::Time::writeControlNames
    { writeControls::wcNone, "none" },
    { writeControls::wcTimeStep, "timeStep" },
    { writeControls::wcRunTime, "runTime" },
    { writeControls::wcAdjustableRunTime, "adjustable" },
    { writeControls::wcAdjustableRunTime, "adjustableRunTime" },
    { writeControls::wcClockTime, "clockTime" },
    { writeControls::wcCpuTime, "cpuTime" },
    // Leave wcUnknown untabulated - fallback to flag unknown settings
Mark OLESEN's avatar
Mark OLESEN committed
});


Foam::Time::fmtflags Foam::Time::format_(Foam::Time::general);
int Foam::Time::precision_(6);

const int Foam::Time::maxPrecision_(3 - log10(SMALL));

Foam::word Foam::Time::controlDictName("controlDict");

int Foam::Time::printExecutionFormat_
(
    Foam::debug::infoSwitch("printExecutionFormat", 0)
);

registerInfoSwitch
(
    "printExecutionFormat",
    int,
    Foam::Time::printExecutionFormat_
);

// * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * * //

void Foam::Time::adjustDeltaT()
{
    bool adjustTime = false;
    scalar timeToNextWrite = VGREAT;

    if (writeControl_ == wcAdjustableRunTime)
    {
        adjustTime = true;
        timeToNextWrite = max
        (
            0.0,
            (writeTimeIndex_ + 1)*writeInterval_ - (value() - startTime_)
        scalar nSteps = timeToNextWrite/deltaT_;
        // For tiny deltaT the label can overflow!
        if (nSteps < scalar(labelMax))
            // nSteps can be < 1 so make sure at least 1
            label nStepsToNextWrite = max(1, round(nSteps));

            scalar newDeltaT = timeToNextWrite/nStepsToNextWrite;

            // Control the increase of the time step to within a factor of 2
            // and the decrease within a factor of 5.
            if (newDeltaT >= deltaT_)
            {
                deltaT_ = min(newDeltaT, 2.0*deltaT_);
            }
            else
            {
                deltaT_ = max(newDeltaT, 0.2*deltaT_);
            }
}


void Foam::Time::setControls()
{
    // default is to resume calculation from "latestTime"
    const word startFrom = controlDict_.getOrDefault<word>
Mark Olesen's avatar
Mark Olesen committed
    (
        "startFrom",
        "latestTime"
    );

    if (startFrom == "startTime")
    {
        controlDict_.readEntry("startTime", startTime_);
    }
    else
    {
        // Search directory for valid time directories
        instantList timeDirs = findTimes(path(), constant());
        const label nTimes = timeDirs.size();

        if (startFrom == "firstTime")
        {
            if (nTimes > 1 && timeDirs.first().name() == constant())
            {
                startTime_ = timeDirs[1].value();
            }
            else if (nTimes)
                startTime_ = timeDirs.first().value();
            }
        }
        else if (startFrom == "latestTime")
        {
                startTime_ = timeDirs.last().value();
            FatalIOErrorInFunction(controlDict_)
mattijs's avatar
mattijs committed
                << "expected startTime, firstTime or latestTime"
                << " found '" << startFrom << "'"
                << exit(FatalIOError);
        }
    }

    setTime(startTime_, 0);

    readDict();
    deltaTSave_ = deltaT_;
    // Check if time directory exists
    // If not increase time precision to see if it is formatted differently.
    // Note: do not use raw fileHandler().exists(...) since does not check
    //       alternative processorsDDD directories naming
    if (fileHandler().filePath(timePath()).empty())
    {
        int oldPrecision = precision_;
        int requiredPrecision = -1;
        bool found = false;
        for
        (
            precision_ = maxPrecision_;
            precision_ > oldPrecision;
        )
        {
            // Update the time formatting
            setTime(startTime_, 0);

            word newTime(timeName());
            if (newTime == oldTime)
            {
                break;
            }
            oldTime = newTime;

            // Check the existence of the time directory with the new format
            //found = fileHandler().exists(timePath(), false);
            const fileName dirName(fileHandler().filePath(timePath()));
            found = !dirName.empty();

            if (found)
            {
                requiredPrecision = precision_;
            }
        }

        if (requiredPrecision > 0)
        {
            // Update the time precision
            precision_ = requiredPrecision;

            // Update the time formatting
            setTime(startTime_, 0);

                << "Increasing the timePrecision from " << oldPrecision
                << " to " << precision_
                << " to support the formatting of the current time directory "
                << timeName() << nl << endl;
        }
        else
        {
            // Could not find time directory so assume it is not present
            precision_ = oldPrecision;

            // Revert the time formatting
            setTime(startTime_, 0);
        }
    }

    if (Pstream::parRun())
    {
        scalar sumStartTime = startTime_;
        reduce(sumStartTime, sumOp<scalar>());
        if
        (
            mag(Pstream::nProcs()*startTime_ - sumStartTime)
          > Pstream::nProcs()*deltaT_/10.0
        )
        {
            FatalIOErrorInFunction(controlDict_)
                << "Start time is not the same for all processors" << nl
                << "processor " << Pstream::myProcNo() << " has startTime "
mattijs's avatar
mattijs committed
                << startTime_ << exit(FatalIOError);
        }
    }

    IOdictionary timeDict
    (
        IOobject
        (
            "time",
            timeName(),
            "uniform",
            *this,
            IOobject::READ_IF_PRESENT,
            IOobject::NO_WRITE,
            false
        )
    );

    // Read and set the deltaT only if time-step adjustment is active
    // otherwise use the deltaT from the controlDict
    if (controlDict_.getOrDefault("adjustTimeStep", false))
        if (timeDict.readIfPresent("deltaT", deltaT_))
        {
            deltaTSave_ = deltaT_;
            deltaT0_ = deltaT_;
        }
    timeDict.readIfPresent("deltaT0", deltaT0_);

    if (timeDict.readIfPresent("index", startTimeIndex_))
    {
        timeIndex_ = startTimeIndex_;
    }

    // Check if values stored in time dictionary are consistent

    // 1. Based on time name
    bool checkValue = true;

    string storedTimeName;
    if (timeDict.readIfPresent("name", storedTimeName))
        if (storedTimeName == timeName())
        {
            // Same time. No need to check stored value
            checkValue = false;
        }
    }
    // 2. Based on time value
    //    (consistent up to the current time writing precision so it won't
    //     trigger if we just change the write precision)
    if (checkValue)
    {
        scalar storedTimeValue;
        if (timeDict.readIfPresent("value", storedTimeValue))
            word storedTimeName(timeName(storedTimeValue));

            if (storedTimeName != timeName())
            {
                IOWarningInFunction(timeDict)
                    << "Time read from time dictionary " << storedTimeName
                    << " differs from actual time " << timeName() << '.' << nl
                    << "    This may cause unexpected database behaviour."
                    << " If you are not interested" << nl
                    << "    in preserving time state delete"
                    << " the time dictionary."
                    << endl;
            }
void Foam::Time::setMonitoring(const bool forceProfiling)
    const dictionary* profilingDict = controlDict_.findDict("profiling");
    if (!profilingDict)
    {
        // ... or from etc/controlDict
        profilingDict = debug::controlDict().findDict("profiling");
    // initialize profiling on request
    // otherwise rely on profiling entry within controlDict
    // and skip if 'active' keyword is explicitly set to false
    if (forceProfiling)
    {
        profiling::initialize
        (
            IOobject
            (
                "profiling",
                timeName(),
                "uniform",
                *this,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE
            ),
            *this
        );
    }
    else if
     && profilingDict->getOrDefault("active", true)
        profiling::initialize
            IOobject
            (
                "profiling",
                timeName(),
                "uniform",
                *this,
                IOobject::NO_READ,
                IOobject::AUTO_WRITE
            ),
            *this
        );
    }

    // Time objects not registered so do like objectRegistry::checkIn ourselves.
    if (runTimeModifiable_)
    {
        // Monitor all files that controlDict depends on
        fileHandler().addWatches(controlDict_, controlDict_.files());
    }

    // Clear dependent files - not needed now
    controlDict_.files().clear();
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::Time::Time
(
    const word& ctrlDictName,
    const fileName& rootPath,
    const fileName& caseName,
    const word& systemName,
    const word& constantName,
    const bool enableFunctionObjects,
    const bool enableLibs
)
:
    TimePaths
    (
        rootPath,
        caseName,
        systemName,
        constantName
    ),

    objectRegistry(*this),

    loopProfiling_(nullptr),
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE,
            false
        )
    ),

    startTimeIndex_(0),
    startTime_(0),
    endTime_(0),

    stopAt_(saEndTime),
    writeControl_(wcTimeStep),
    writeInterval_(GREAT),
    purgeWrite_(0),
Mark OLESEN's avatar
Mark OLESEN committed
    subCycling_(0),
    writeOnce_(false),
    sigWriteNow_(*this, true),
    sigStopAtWriteNow_(*this, true),
    writeStreamOption_(IOstream::ASCII),
    graphFormat_("raw"),
    runTimeModifiable_(false),
    functionObjects_(*this, false)
    if (enableFunctionObjects)
    {
        functionObjects_.on();
    }

    if (enableLibs)
    {
        libs_.open("libs", controlDict_);

    // Explicitly set read flags on objectRegistry so anything constructed
    // from it reads as well (e.g. fvSolution).
    readOpt(IOobject::MUST_READ_IF_MODIFIED);
    setMonitoring();
    const word& ctrlDictName,
    const argList& args,
    const word& systemName,
    const word& constantName,
    const bool enableFunctionObjects,
    const bool enableLibs
    TimePaths(args, systemName, constantName),
    loopProfiling_(nullptr),
            system(),
            *this,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE,
            false
        )
    ),

    startTimeIndex_(0),
    startTime_(0),
    endTime_(0),

    stopAt_(saEndTime),
    writeControl_(wcTimeStep),
    writeInterval_(GREAT),
    purgeWrite_(0),
Mark OLESEN's avatar
Mark OLESEN committed
    subCycling_(0),
    writeOnce_(false),
    sigWriteNow_(*this, true),
    sigStopAtWriteNow_(*this, true),
    writeStreamOption_(IOstream::ASCII),
    graphFormat_("raw"),
    runTimeModifiable_(false),
    functionObjects_(*this, false)
{
    // * '-withFunctionObjects' exists and used = enable
    // * '-noFunctionObjects' exists and used = disable
    // * default: no functions if there is no way to enable/disable them
    (
        argList::validOptions.found("withFunctionObjects")
      ? args.found("withFunctionObjects")
      : argList::validOptions.found("noFunctionObjects")
      ? !args.found("noFunctionObjects")
        if (enableFunctionObjects)
        {
            functionObjects_.on();
        }
    // * enable by default unless '-no-libs' option was used
    if (enableLibs && !args.found("no-libs"))
        libs_.open("libs", controlDict_);
    // Explicitly set read flags on objectRegistry so anything constructed
    // from it reads as well (e.g. fvSolution).
    readOpt(IOobject::MUST_READ_IF_MODIFIED);

    // '-profiling' = force profiling, ignore controlDict entry
    setMonitoring(args.found("profiling"));
}


Foam::Time::Time
(
    const dictionary& dict,
    const fileName& rootPath,
    const fileName& caseName,
    const word& systemName,
    const word& constantName,
    const bool enableFunctionObjects,
    const bool enableLibs
)
:
    TimePaths
    (
        rootPath,
        caseName,
        systemName,
        constantName
    ),

    objectRegistry(*this),

    loopProfiling_(nullptr),
    controlDict_
    (
        IOobject
        (
            controlDictName,
            system(),
            *this,
            IOobject::NO_READ,
            IOobject::NO_WRITE,
            false
        ),
        dict
    ),

    startTimeIndex_(0),
    startTime_(0),
    endTime_(0),

    stopAt_(saEndTime),
    writeControl_(wcTimeStep),
    writeInterval_(GREAT),
    purgeWrite_(0),
Mark OLESEN's avatar
Mark OLESEN committed
    subCycling_(0),
    writeOnce_(false),
    sigWriteNow_(*this, true),
    sigStopAtWriteNow_(*this, true),
    writeStreamOption_(IOstream::ASCII),
    graphFormat_("raw"),
    runTimeModifiable_(false),
    functionObjects_(*this, false)
    if (enableFunctionObjects)
    {
        functionObjects_.on();
    }

    if (enableLibs)
    {
        libs_.open("libs", controlDict_);
    // Explicitly set read flags on objectRegistry so anything constructed
    // from it reads as well (e.g. fvSolution).
    readOpt(IOobject::MUST_READ_IF_MODIFIED);
    // Since could not construct regIOobject with setting:
    controlDict_.readOpt(IOobject::MUST_READ_IF_MODIFIED);
    setMonitoring();
}


Foam::Time::Time
(
    const fileName& rootPath,
    const fileName& caseName,
    const word& systemName,
    const word& constantName,
    const bool enableFunctionObjects,
    const bool enableLibs
)
:
    TimePaths
    (
        rootPath,
        caseName,
        systemName,
        constantName
    ),

    objectRegistry(*this),

    loopProfiling_(nullptr),
    controlDict_
    (
        IOobject
        (
            controlDictName,
            system(),
            *this,
            IOobject::NO_READ,
            IOobject::NO_WRITE,
            false
        )
    ),

    startTimeIndex_(0),
    startTime_(0),
    endTime_(0),

    stopAt_(saEndTime),
    writeControl_(wcTimeStep),
    writeInterval_(GREAT),
    purgeWrite_(0),
Mark OLESEN's avatar
Mark OLESEN committed
    subCycling_(0),
    writeOnce_(false),
    writeStreamOption_(IOstream::ASCII),
    graphFormat_("raw"),
    runTimeModifiable_(false),
    functionObjects_(*this, false)
    if (enableFunctionObjects)
    {
        functionObjects_.on();
    }

    if (enableLibs)
    {
        libs_.open("libs", controlDict_);
    setMonitoring(); // for profiling etc
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //

Foam::autoPtr<Foam::Time> Foam::Time::New()
{
    return
        autoPtr<Time>::New
        (
            fileName("."),  // root-path
            fileName("."),  // case-name
            false,          // No enableFunctionObjects
            false           // No enableLibs
        );
}


Foam::autoPtr<Foam::Time> Foam::Time::New(const fileName& caseDir)
{
    return
        autoPtr<Time>::New
        (
            caseDir.path(), // root-path
            caseDir.name(), // case-name
            false,          // No enableFunctionObjects
            false           // No enableLibs
        );
}


Foam::autoPtr<Foam::Time> Foam::Time::New(const argList& args)
{
    return
        autoPtr<Time>::New
        (
            Time::controlDictName,
            args,
            false,          // No enableFunctionObjects
            false           // No enableLibs
        );
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::Time::~Time()
    loopProfiling_.reset(nullptr);
    forAllReverse(controlDict_.watchIndices(), i)
        fileHandler().removeWatch(controlDict_.watchIndices()[i]);
    functionObjects_.clear();
    // Clean up profiling
    profiling::stop(*this);

    // Ensure all owned objects are also cleaned up now
    objectRegistry::clear();


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

Foam::word Foam::Time::timeName(const scalar t, const int precision)
    std::ostringstream buf;
    buf.setf(ios_base::fmtflags(format_), ios_base::floatfield);
    buf << t;
    return buf.str();
}


Foam::word Foam::Time::timeName() const
{
    return dimensionedScalar::name();
}


Foam::word Foam::Time::findInstance
(
    const fileName& dir,
    const word& name,
    const IOobject::readOption rOpt,
    const word& stopInstance
) const
{
    IOobject startIO
    (
        name,           // name might be empty!
        timeName(),
        dir,
        *this,
        rOpt
    );

    IOobject io
    (
        fileHandler().findInstance
        (
            startIO,
            timeOutputValue(),
            stopInstance
        )
    );
    return io.instance();
}


Foam::word Foam::Time::findInstancePath
(
    const fileName& directory,
    const instant& t
) const
    // Simplified version: use findTimes (readDir + sort). The expensive
    // bit is the readDir, not the sorting. Tbd: avoid calling findInstancePath
    // from filePath.
    instantList timeDirs = findTimes(path(), constant());
    // Note:
    // - times will include constant (with value 0) as first element.
    //   For backwards compatibility make sure to find 0 in preference
    //   to constant.
    // - list is sorted so could use binary search
        if (t.equal(timeDirs[i].value()))
Foam::word Foam::Time::findInstancePath(const instant& t) const
{
    return findInstancePath(path(), t);
}


Foam::label Foam::Time::startTimeIndex() const
{
    return startTimeIndex_;
}


Foam::dimensionedScalar Foam::Time::startTime() const
{
    return dimensionedScalar("startTime", dimTime, startTime_);
}


Foam::dimensionedScalar Foam::Time::endTime() const
{
    return dimensionedScalar("endTime", dimTime, endTime_);
}


Mark OLESEN's avatar
Mark OLESEN committed
Foam::Time::stopAtControls Foam::Time::stopAt() const
{
    return stopAt_;
}


bool Foam::Time::run() const
{
    loopProfiling_.reset(nullptr);
    bool isRunning = value() < (endTime_ - 0.5*deltaT_);
        // Only execute when the condition is no longer true
        // ie, when exiting the control loop
        if (!isRunning && timeIndex_ != startTimeIndex_)
            // Ensure functionObjects execute on last time step
            // (and hence write uptodate functionObjectProperties)
            {
                addProfiling(fo, "functionObjects.execute()");
                functionObjects_.execute();
            }
            {
                addProfiling(fo, "functionObjects.end()");
                functionObjects_.end();
            }
            const_cast<Time&>(*this).readModifiedObjects();

            if (timeIndex_ == startTimeIndex_)
            {
                addProfiling(functionObjects, "functionObjects.start()");
                functionObjects_.start();
            }
            else
            {
                addProfiling(functionObjects, "functionObjects.execute()");
                functionObjects_.execute();
            }

            // Check if the execution of functionObjects require re-reading
            // any files. This moves effect of e.g. 'timeActivatedFileUpdate'
            // one time step forward. Note that we cannot call
            // readModifiedObjects from within timeActivatedFileUpdate since
            // it might re-read the functionObjects themselves (and delete
            // the timeActivatedFileUpdate one)
            if (functionObjects_.filesModified())
            {
                const_cast<Time&>(*this).readModifiedObjects();
            }
        // Update the "is-running" status following the
        // possible side-effects from functionObjects
        isRunning = value() < (endTime_ - 0.5*deltaT_);

        // (re)trigger profiling
        if (profiling::active())
        {
            loopProfiling_.reset
            (
                new profilingTrigger("time.run() " + objectRegistry::name())
            );
    const bool isRunning = run();
bool Foam::Time::end() const
{
    return value() > (endTime_ + 0.5*deltaT_);
Mark OLESEN's avatar
Mark OLESEN committed
bool Foam::Time::stopAt(const stopAtControls stopCtrl) const
Mark OLESEN's avatar
Mark OLESEN committed
    if (stopCtrl == stopAtControls::saUnknown)
Mark OLESEN's avatar
Mark OLESEN committed
        return false;
Mark OLESEN's avatar
Mark OLESEN committed

    const bool changed = (stopAt_ != stopCtrl);
    stopAt_ = stopCtrl;
    endTime_ = GREAT;

    // Adjust endTime
    if (stopCtrl == stopAtControls::saEndTime)
        controlDict_.readEntry("endTime", endTime_);
bool Foam::Time::isAdjustTimeStep() const
{
    return controlDict_.getOrDefault("adjustTimeStep", false);