View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0002820 | OpenFOAM | Bug | public | 2018-01-27 12:51 | 2018-02-08 11:26 |
Reporter | Shorty | Assigned To | will | ||
Priority | normal | Severity | major | Reproducibility | random |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | Ubuntu | OS Version | 16.04 |
Product Version | dev | ||||
Fixed in Version | dev | ||||
Summary | 0002820: Tricky bug in »adjustableRunTime« while using functionObjects | ||||
Description | Hi all, yesterday I build two tutorials for @Will in order to check the residual control. During the evening, I was investigating into the transient case and added some more stuff. However, I realized, that after adding a functionObject, the solver behaves strange. The time step decreased (for any reason) sometimes to 1e-20 and I had to make several iterations in order to reach the next write step. First, I thought it is a problem of my set-up, or local computer. I rechecked the set-up two times. With functionObject, strange behavior, without everything is fine. The time decrease can be so extreme, that I get an floating point exception. After restarting the case, everything seems to be okay for a few more iterations till the same strange behavior occurs. Thus, I checked it on another machine too. The same situation. During the night I was debugging and found, that the function »adjustDeltaT« in the Time class makes some trouble. The problem is that we are not allowed to change the deltaT_ before calling the objectFunction.adjustTimeStep() function. If we do so, this function has a wrong deltaT information for estimating the new deltaT. Thus, if one has some » unluckily « chosen writeInterval values in the controlDict (solver) and functionObject, they end up in this strange behavior. I hope I made a clear statement. | ||||
Steps To Reproduce | Start the run script in the attachment | ||||
Additional Information | A patch is included, that changes the set-up of deltaT_ in the Time class after we have the correct estimation of the function object adjustment. | ||||
Tags | Time Control | ||||
|
|
|
Time.C (29,225 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- 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 "IOdictionary.H" #include <sstream> // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * // namespace Foam { defineTypeNameAndDebug(Time, 0); template<> const char* Foam::NamedEnum < Foam::Time::stopAtControls, 4 >::names[] = { "endTime", "noWriteNow", "writeNow", "nextWrite" }; template<> const char* Foam::NamedEnum < Foam::Time::writeControls, 5 >::names[] = { "timeStep", "runTime", "adjustableRunTime", "clockTime", "cpuTime" }; } const Foam::NamedEnum<Foam::Time::stopAtControls, 4> Foam::Time::stopAtControlNames_; const Foam::NamedEnum<Foam::Time::writeControls, 5> Foam::Time::writeControlNames_; 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"); // * * * * * * * * * * * * 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_) ); } // Estimation for controlDict settings scalar deltaT = 0.; if (adjustTime) { scalar nSteps = timeToNextWrite/deltaT_ - SMALL; // For tiny deltaT the label can overflow! if (nSteps < labelMax) { label nStepsToNextWrite = label(nSteps) + 1; 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_); } } } functionObjects_.adjustTimeStep(); // Adjust the time step based on functionObjects settings or controlDict deltaT_ = min(deltaT_, deltaT); } void Foam::Time::setControls() { // default is to resume calculation from "latestTime" const word startFrom = controlDict_.lookupOrDefault<word> ( "startFrom", "latestTime" ); if (startFrom == "startTime") { controlDict_.lookup("startTime") >> startTime_; } else { // Search directory for valid time directories instantList timeDirs = findTimes(path(), constant()); if (startFrom == "firstTime") { if (timeDirs.size()) { if (timeDirs[0].name() == constant() && timeDirs.size() >= 2) { startTime_ = timeDirs[1].value(); } else { startTime_ = timeDirs[0].value(); } } } else if (startFrom == "latestTime") { if (timeDirs.size()) { startTime_ = timeDirs.last().value(); } } else { FatalIOErrorInFunction(controlDict_) << "expected startTime, firstTime or latestTime" << " found '" << startFrom << "'" << exit(FatalIOError); } } setTime(startTime_, 0); readDict(); deltaTSave_ = deltaT_; deltaT0_ = deltaT_; // Check if time directory exists // If not increase time precision to see if it is formatted differently. if (!fileHandler().exists(timePath(), false)) { int oldPrecision = precision_; int requiredPrecision = -1; bool found = false; word oldTime(timeName()); for ( precision_ = maxPrecision_; precision_ > oldPrecision; precision_-- ) { // 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); if (found) { requiredPrecision = precision_; } } if (requiredPrecision > 0) { // Update the time precision precision_ = requiredPrecision; // Update the time formatting setTime(startTime_, 0); WarningInFunction << "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 " << 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_.lookupOrDefault<Switch>("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; } } } } // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // Foam::Time::Time ( const word& controlDictName, const fileName& rootPath, const fileName& caseName, const word& systemName, const word& constantName, const bool enableFunctionObjects ) : TimePaths ( rootPath, caseName, systemName, constantName ), objectRegistry(*this), libs_(), controlDict_ ( IOobject ( controlDictName, 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), writeOnce_(false), subCycling_(false), sigWriteNow_(true, *this), sigStopAtWriteNow_(true, *this), writeFormat_(IOstream::ASCII), writeVersion_(IOstream::currentVersion), writeCompression_(IOstream::UNCOMPRESSED), graphFormat_("raw"), runTimeModifiable_(false), functionObjects_(*this, enableFunctionObjects) { libs_.open(controlDict_, "libs"); // Explicitly set read flags on objectRegistry so anything constructed // from it reads as well (e.g. fvSolution). readOpt() = IOobject::MUST_READ_IF_MODIFIED; setControls(); // 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 controlDict_.files().clear(); } Foam::Time::Time ( const word& controlDictName, const argList& args, const word& systemName, const word& constantName ) : TimePaths ( args.parRunControl().parRun(), args.rootPath(), args.globalCaseName(), args.caseName(), systemName, constantName ), objectRegistry(*this), libs_(), controlDict_ ( IOobject ( controlDictName, 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), writeOnce_(false), subCycling_(false), sigWriteNow_(true, *this), sigStopAtWriteNow_(true, *this), writeFormat_(IOstream::ASCII), writeVersion_(IOstream::currentVersion), writeCompression_(IOstream::UNCOMPRESSED), graphFormat_("raw"), runTimeModifiable_(false), functionObjects_ ( *this, argList::validOptions.found("withFunctionObjects") ? args.optionFound("withFunctionObjects") : !args.optionFound("noFunctionObjects") ) { libs_.open(controlDict_, "libs"); // Explicitly set read flags on objectRegistry so anything constructed // from it reads as well (e.g. fvSolution). readOpt() = IOobject::MUST_READ_IF_MODIFIED; setControls(); // 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 since not needed controlDict_.files().clear(); } Foam::Time::Time ( const dictionary& dict, const fileName& rootPath, const fileName& caseName, const word& systemName, const word& constantName, const bool enableFunctionObjects ) : TimePaths ( rootPath, caseName, systemName, constantName ), objectRegistry(*this), libs_(), 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), writeOnce_(false), subCycling_(false), sigWriteNow_(true, *this), sigStopAtWriteNow_(true, *this), writeFormat_(IOstream::ASCII), writeVersion_(IOstream::currentVersion), writeCompression_(IOstream::UNCOMPRESSED), graphFormat_("raw"), runTimeModifiable_(false), functionObjects_(*this, enableFunctionObjects) { libs_.open(controlDict_, "libs"); // 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; setControls(); // 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 since not needed controlDict_.files().clear(); } Foam::Time::Time ( const fileName& rootPath, const fileName& caseName, const word& systemName, const word& constantName, const bool enableFunctionObjects ) : TimePaths ( rootPath, caseName, systemName, constantName ), objectRegistry(*this), libs_(), 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), writeOnce_(false), subCycling_(false), writeFormat_(IOstream::ASCII), writeVersion_(IOstream::currentVersion), writeCompression_(IOstream::UNCOMPRESSED), graphFormat_("raw"), runTimeModifiable_(false), functionObjects_(*this, enableFunctionObjects) { libs_.open(controlDict_, "libs"); } // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // Foam::Time::~Time() { forAllReverse(controlDict_.watchIndices(), i) { fileHandler().removeWatch(controlDict_.watchIndices()[i]); } // Destroy function objects first functionObjects_.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.precision(precision); buf << t; return buf.str(); } Foam::word Foam::Time::timeName() const { return dimensionedScalar::name(); } Foam::instantList Foam::Time::times() const { return findTimes(path(), constant()); } 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 forAllReverse(timeDirs, i) { if (t.equal(timeDirs[i].value())) { return timeDirs[i].name(); } } return word::null; } Foam::word Foam::Time::findInstancePath(const instant& t) const { return findInstancePath(path(), t); } Foam::instant Foam::Time::findClosestTime(const scalar t) const { instantList timeDirs = findTimes(path(), constant()); // There is only one time (likely "constant") so return it if (timeDirs.size() == 1) { return timeDirs[0]; } if (t < timeDirs[1].value()) { return timeDirs[1]; } else if (t > timeDirs.last().value()) { return timeDirs.last(); } label nearestIndex = -1; scalar deltaT = GREAT; for (label timei=1; timei < timeDirs.size(); ++timei) { scalar diff = mag(timeDirs[timei].value() - t); if (diff < deltaT) { deltaT = diff; nearestIndex = timei; } } return timeDirs[nearestIndex]; } Foam::label Foam::Time::findClosestTimeIndex ( const instantList& timeDirs, const scalar t, const word& constantName ) { label nearestIndex = -1; scalar deltaT = GREAT; forAll(timeDirs, timei) { if (timeDirs[timei].name() == constantName) continue; scalar diff = mag(timeDirs[timei].value() - t); if (diff < deltaT) { deltaT = diff; nearestIndex = timei; } } return nearestIndex; } 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_); } bool Foam::Time::run() const { bool running = value() < (endTime_ - 0.5*deltaT_); if (!subCycling_) { if (!running && timeIndex_ != startTimeIndex_) { functionObjects_.execute(); functionObjects_.end(); } } if (running) { if (!subCycling_) { const_cast<Time&>(*this).readModifiedObjects(); if (timeIndex_ == startTimeIndex_) { functionObjects_.start(); } else { functionObjects_.execute(); } } // Update the "running" status following the // possible side-effects from functionObjects running = value() < (endTime_ - 0.5*deltaT_); } return running; } bool Foam::Time::loop() { bool running = run(); if (running) { operator++(); } return running; } bool Foam::Time::end() const { return value() > (endTime_ + 0.5*deltaT_); } bool Foam::Time::stopAt(const stopAtControls sa) const { const bool changed = (stopAt_ != sa); stopAt_ = sa; // adjust endTime if (sa == saEndTime) { controlDict_.lookup("endTime") >> endTime_; } else { endTime_ = GREAT; } return changed; } void Foam::Time::setTime(const Time& t) { value() = t.value(); dimensionedScalar::name() = t.dimensionedScalar::name(); timeIndex_ = t.timeIndex_; fileHandler().setTime(*this); } void Foam::Time::setTime(const instant& inst, const label newIndex) { value() = inst.value(); dimensionedScalar::name() = inst.name(); timeIndex_ = newIndex; IOdictionary timeDict ( IOobject ( "time", timeName(), "uniform", *this, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false ) ); timeDict.readIfPresent("deltaT", deltaT_); timeDict.readIfPresent("deltaT0", deltaT0_); timeDict.readIfPresent("index", timeIndex_); fileHandler().setTime(*this); } void Foam::Time::setTime(const dimensionedScalar& newTime, const label newIndex) { setTime(newTime.value(), newIndex); } void Foam::Time::setTime(const scalar newTime, const label newIndex) { value() = newTime; dimensionedScalar::name() = timeName(timeToUserTime(newTime)); timeIndex_ = newIndex; fileHandler().setTime(*this); } void Foam::Time::setEndTime(const dimensionedScalar& endTime) { setEndTime(endTime.value()); } void Foam::Time::setEndTime(const scalar endTime) { endTime_ = endTime; } void Foam::Time::setDeltaT ( const dimensionedScalar& deltaT, const bool bAdjustDeltaT ) { setDeltaT(deltaT.value(), bAdjustDeltaT); } void Foam::Time::setDeltaT(const scalar deltaT, const bool bAdjustDeltaT) { deltaT_ = deltaT; deltaTchanged_ = true; if (bAdjustDeltaT) { adjustDeltaT(); } } Foam::TimeState Foam::Time::subCycle(const label nSubCycles) { subCycling_ = true; prevTimeState_.set(new TimeState(*this)); setTime(*this - deltaT(), (timeIndex() - 1)*nSubCycles); deltaT_ /= nSubCycles; deltaT0_ /= nSubCycles; deltaTSave_ = deltaT0_; return prevTimeState(); } void Foam::Time::endSubCycle() { if (subCycling_) { subCycling_ = false; TimeState::operator=(prevTimeState()); prevTimeState_.clear(); } } // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * // Foam::Time& Foam::Time::operator+=(const dimensionedScalar& deltaT) { return operator+=(deltaT.value()); } Foam::Time& Foam::Time::operator+=(const scalar deltaT) { setDeltaT(deltaT); return operator++(); } Foam::Time& Foam::Time::operator++() { deltaT0_ = deltaTSave_; deltaTSave_ = deltaT_; // Save old time value and name const scalar oldTimeValue = timeToUserTime(value()); const word oldTimeName = dimensionedScalar::name(); // Increment time setTime(value() + deltaT_, timeIndex_ + 1); if (!subCycling_) { // If the time is very close to zero reset to zero if (mag(value()) < 10*SMALL*deltaT_) { setTime(0.0, timeIndex_); } if (sigStopAtWriteNow_.active() || sigWriteNow_.active()) { // A signal might have been sent on one processor only // Reduce so all decide the same. label flag = 0; if (sigStopAtWriteNow_.active() && stopAt_ == saWriteNow) { flag += 1; } if (sigWriteNow_.active() && writeOnce_) { flag += 2; } reduce(flag, maxOp<label>()); if (flag & 1) { stopAt_ = saWriteNow; } if (flag & 2) { writeOnce_ = true; } } writeTime_ = false; switch (writeControl_) { case wcTimeStep: writeTime_ = !(timeIndex_ % label(writeInterval_)); break; case wcRunTime: case wcAdjustableRunTime: { label writeIndex = label ( ((value() - startTime_) + 0.5*deltaT_) / writeInterval_ ); if (writeIndex > writeTimeIndex_) { writeTime_ = true; writeTimeIndex_ = writeIndex; } } break; case wcCpuTime: { label writeIndex = label ( returnReduce(elapsedCpuTime(), maxOp<double>()) / writeInterval_ ); if (writeIndex > writeTimeIndex_) { writeTime_ = true; writeTimeIndex_ = writeIndex; } } break; case wcClockTime: { label writeIndex = label ( returnReduce(label(elapsedClockTime()), maxOp<label>()) / writeInterval_ ); if (writeIndex > writeTimeIndex_) { writeTime_ = true; writeTimeIndex_ = writeIndex; } } break; } // Check if endTime needs adjustment to stop at the next run()/end() if (!end()) { if (stopAt_ == saNoWriteNow) { endTime_ = value(); } else if (stopAt_ == saWriteNow) { endTime_ = value(); writeTime_ = true; } else if (stopAt_ == saNextWrite && writeTime_ == true) { endTime_ = value(); } } // Override writeTime if one-shot writing if (writeOnce_) { writeTime_ = true; writeOnce_ = false; } // Adjust the precision of the time directory name if necessary if (writeTime_) { // Tolerance used when testing time equivalence const scalar timeTol = max(min(pow(10.0, -precision_), 0.1*deltaT_), SMALL); // User-time equivalent of deltaT const scalar userDeltaT = timeToUserTime(deltaT_); // Time value obtained by reading timeName scalar timeNameValue = -VGREAT; // Check that new time representation differs from old one // reinterpretation of the word if ( readScalar(dimensionedScalar::name().c_str(), timeNameValue) && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol) ) { int oldPrecision = precision_; while ( precision_ < maxPrecision_ && readScalar(dimensionedScalar::name().c_str(), timeNameValue) && (mag(timeNameValue - oldTimeValue - userDeltaT) > timeTol) ) { precision_++; setTime(value(), timeIndex()); } if (precision_ != oldPrecision) { WarningInFunction << "Increased the timePrecision from " << oldPrecision << " to " << precision_ << " to distinguish between timeNames at time " << dimensionedScalar::name() << endl; if (precision_ == maxPrecision_) { // Reached maxPrecision limit WarningInFunction << "Current time name " << dimensionedScalar::name() << nl << " The maximum time precision has been reached" " which might result in overwriting previous" " results." << endl; } // Check if round-off error caused time-reversal scalar oldTimeNameValue = -VGREAT; if ( readScalar(oldTimeName.c_str(), oldTimeNameValue) && ( sign(timeNameValue - oldTimeNameValue) != sign(deltaT_) ) ) { WarningInFunction << "Current time name " << dimensionedScalar::name() << " is set to an instance prior to the " "previous one " << oldTimeName << nl << " This might result in temporal " "discontinuities." << endl; } } } } } return *this; } Foam::Time& Foam::Time::operator++(int) { return operator++(); } // ************************************************************************* // |
|
That was really a tricky one :) Hope you can reproduce it on your machines. |
|
Foam_Time_adjustDeltaT.diff (1,032 bytes)
diff --git a/src/OpenFOAM/db/Time/Time.C b/src/OpenFOAM/db/Time/Time.C index f002f414f..fc828e00b 100644 --- a/src/OpenFOAM/db/Time/Time.C +++ b/src/OpenFOAM/db/Time/Time.C @@ -96,6 +96,9 @@ void Foam::Time::adjustDeltaT() ); } + // Estimation for controlDict settings + scalar deltaT = 0; + if (adjustTime) { scalar nSteps = timeToNextWrite/deltaT_ - small; @@ -111,16 +114,19 @@ void Foam::Time::adjustDeltaT() // and the decrease within a factor of 5. if (newDeltaT >= deltaT_) { - deltaT_ = min(newDeltaT, 2.0*deltaT_); + deltaT = min(newDeltaT, 2.0*deltaT_); } else { - deltaT_ = max(newDeltaT, 0.2*deltaT_); + deltaT = max(newDeltaT, 0.2*deltaT_); } } } functionObjects_.adjustTimeStep(); + + // Adjust the time step based on functionObjects settings or controlDict + deltaT_ = min(deltaT_, deltaT); } |
|
I have reproduced this, and the provided fix does solve the problem for this case. However, the change breaks other function objects' modification of the time step; notably setTimeStepFunctionObject. I haven't figured a way around this yet. Also, the local deltaT variable needs initialising to the class deltaT_, rather than zero, if non-adjustable time stepping is to continue working. I've attached a patch that applies cleanly to dev (post all the great/vGreat/small/... changes), in case anyone else wants to tinker. |
|
Thanks @will for your investigation. I will check the setTimeStepFunctionObject and see if I can find any solution. |
|
I just wanted to ask something more. Do you have any defined or structured procedure in order to identify if the patch makes problems with other function objects? Or do you just try one by one? I am asking in order to do this process myself in order to save your time. |
|
The only automated check we have is ./Alltest in the tutorials directory. That's how I noticed that the deltaT initialisation was causing fixed time-step cases to fail. Sometimes you just have to run tests and/or read the code, though. In this instance I noticed the change's incompatibility with other function objects by reading through the implementations of "adjustTimeStep". There are only two that actually do anything. |
|
I had a look at it again. First, yes, deltaT should not be initialized with 0. After I changed it, the timeStep is working as expected in combination with the functionObjects. Please just change scalar deltaT = 0; to scalar deltaT = deltaT_; This fixes it. I did not figure out any other problems anymore. |
|
Your change still breaks the setTimeStepFunctionObject. Try the multiphase/reactingTwoPhaseEulerFoam/laminar/steamInjection tutorial. I've had to change how it works more fundamentally. The function object now returns the time to the next write, rather than setting the time directly. This keeps the logic within the Time class and fundamentally prevents this sort of interaction from occurring. An added benefit is that the setTimeStepFunctionObject can now be used with or without adjustable write time. Resolved in dev by commit 139ffcc1. |
Date Modified | Username | Field | Change |
---|---|---|---|
2018-01-27 12:51 | Shorty | New Issue | |
2018-01-27 12:51 | Shorty | File Added: chtResControlCheck_Transient_BugReport.tar.gz | |
2018-01-27 12:51 | Shorty | Tag Attached: Time Control | |
2018-01-27 12:51 | Shorty | File Added: Time.C | |
2018-01-27 12:52 | Shorty | Note Added: 0009225 | |
2018-02-02 15:44 | will | File Added: Foam_Time_adjustDeltaT.diff | |
2018-02-02 15:44 | will | Note Added: 0009241 | |
2018-02-02 15:55 | Shorty | Note Added: 0009244 | |
2018-02-02 15:58 | Shorty | Note Added: 0009245 | |
2018-02-02 16:39 | will | Note Added: 0009248 | |
2018-02-04 18:14 | Shorty | Note Added: 0009251 | |
2018-02-08 11:26 | will | Assigned To | => will |
2018-02-08 11:26 | will | Status | new => resolved |
2018-02-08 11:26 | will | Resolution | open => fixed |
2018-02-08 11:26 | will | Fixed in Version | => dev |
2018-02-08 11:26 | will | Note Added: 0009282 |