#11092 Fix import of time since simulation start

This commit is contained in:
Magne Sjaastad 2024-01-22 11:02:28 +01:00
parent e2e861ed95
commit 62609953a0
2 changed files with 18 additions and 13 deletions

View File

@ -40,7 +40,6 @@
#include "opm/io/eclipse/ERst.hpp"
#include "opm/io/eclipse/RestartFileView.hpp"
#include "opm/io/eclipse/rst/state.hpp"
#include "opm/output/eclipse/VectorItems/doubhead.hpp"
#include "opm/output/eclipse/VectorItems/group.hpp"
#include "opm/output/eclipse/VectorItems/intehead.hpp"
#include "opm/output/eclipse/VectorItems/well.hpp"
@ -273,7 +272,7 @@ void RifReaderOpmCommon::buildMetaData( RigEclipseCaseData* eclipseCase )
QDate date( timeStep.year, timeStep.month, timeStep.day );
QDateTime dateTime = RiaQDateTimeTools::createDateTime( date );
dateTimes.push_back( dateTime );
timeStepInfos.emplace_back( dateTime, timeStep.sequenceNumber, 0.0 );
timeStepInfos.emplace_back( dateTime, timeStep.sequenceNumber, timeStep.simulationTimeFromStart );
}
m_timeSteps = dateTimes;
@ -464,8 +463,14 @@ std::vector<RifReaderOpmCommon::TimeDataFile> RifReaderOpmCommon::readTimeSteps(
if ( restartFile->hasArray( doubheadString, seqNum ) )
{
auto doubhead = restartFile->getRestartData<double>( doubheadString, seqNum );
daySinceSimStart = doubhead[VI::doubhead::TsInit];
auto doubhead = restartFile->getRestartData<double>( doubheadString, seqNum );
if ( !doubhead.empty() )
{
// Read out the simulation time from start from DOUBHEAD. There is no enum defined to access this value.
// https://github.com/OPM/ResInsight/issues/11092
daySinceSimStart = doubhead[0];
}
}
reportTimeData.emplace_back(

View File

@ -35,15 +35,6 @@ class ERst;
class RifReaderOpmCommon : public RifReaderInterface
{
public:
struct TimeDataFile
{
int sequenceNumber;
int year;
int month;
int day;
double simulationTimeFromStart;
};
RifReaderOpmCommon();
~RifReaderOpmCommon() override;
@ -55,6 +46,15 @@ public:
private:
void buildMetaData( RigEclipseCaseData* eclipseCase );
struct TimeDataFile
{
int sequenceNumber;
int year;
int month;
int day;
double simulationTimeFromStart;
};
static std::vector<TimeDataFile> readTimeSteps( std::shared_ptr<Opm::EclIO::ERst> restartFile );
static void readWellCells( std::shared_ptr<Opm::EclIO::ERst> restartFile,
RigEclipseCaseData* eclipseCase,