mirror of
https://github.com/OPM/ResInsight.git
synced 2025-02-25 18:55:39 -06:00
Rewrite of file acces helper class
RifEcliseOutputFileTools is now a class with only static functions Moved some application logic from file readers to RifEcliseOutput Added more tests for extracting metainfo related to results p4#: 20368
This commit is contained in:
@@ -39,7 +39,7 @@ RifEclipseRestartFilesetAccess::~RifEclipseRestartFilesetAccess()
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
/// Open files
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
bool RifEclipseRestartFilesetAccess::open(const QStringList& fileSet, const std::vector<size_t>& matrixActiveCellCounts, const std::vector<size_t>& fractureActiveCellCounts)
|
||||
bool RifEclipseRestartFilesetAccess::open(const QStringList& fileSet)
|
||||
{
|
||||
close();
|
||||
|
||||
@@ -47,18 +47,12 @@ bool RifEclipseRestartFilesetAccess::open(const QStringList& fileSet, const std:
|
||||
int i;
|
||||
for (i = 0; i < numFiles; i++)
|
||||
{
|
||||
cvf::ref<RifEclipseOutputFileTools> fileAccess = new RifEclipseOutputFileTools;
|
||||
if (!fileAccess->open(fileSet[i], matrixActiveCellCounts, fractureActiveCellCounts))
|
||||
{
|
||||
close();
|
||||
return false;
|
||||
}
|
||||
ecl_file_type* ecl_file = ecl_file_open(fileSet[i].toAscii().data());
|
||||
if (!ecl_file) return false;
|
||||
|
||||
m_files.push_back(fileAccess);
|
||||
m_ecl_files.push_back(ecl_file);
|
||||
}
|
||||
|
||||
m_gridCount = matrixActiveCellCounts.size();
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
@@ -67,7 +61,12 @@ bool RifEclipseRestartFilesetAccess::open(const QStringList& fileSet, const std:
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
void RifEclipseRestartFilesetAccess::close()
|
||||
{
|
||||
m_files.clear();
|
||||
for (size_t i = 0; i < m_ecl_files.size(); i++)
|
||||
{
|
||||
ecl_file_close(m_ecl_files[i]);
|
||||
}
|
||||
m_ecl_files.clear();
|
||||
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
@@ -75,7 +74,7 @@ void RifEclipseRestartFilesetAccess::close()
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
size_t RifEclipseRestartFilesetAccess::numTimeSteps()
|
||||
{
|
||||
return m_files.size();
|
||||
return m_ecl_files.size();
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
@@ -90,7 +89,8 @@ QStringList RifEclipseRestartFilesetAccess::timeStepsText()
|
||||
for (i = 0; i < numSteps; i++)
|
||||
{
|
||||
QStringList stepText;
|
||||
m_files[i]->timeStepsText(&stepText);
|
||||
RifEclipseOutputFileTools::timeStepsText(m_ecl_files[i], &stepText);
|
||||
|
||||
timeSteps.append(stepText.size() == 1 ? stepText : QStringList(QString("Step %1").arg(i+1)));
|
||||
}
|
||||
|
||||
@@ -109,7 +109,8 @@ QList<QDateTime> RifEclipseRestartFilesetAccess::timeSteps()
|
||||
for (i = 0; i < numSteps; i++)
|
||||
{
|
||||
QList<QDateTime> stepTime;
|
||||
m_files[i]->timeSteps(&stepTime);
|
||||
|
||||
RifEclipseOutputFileTools::timeSteps(m_ecl_files[i], &stepTime);
|
||||
|
||||
if (stepTime.size() == 1)
|
||||
{
|
||||
@@ -127,29 +128,31 @@ QList<QDateTime> RifEclipseRestartFilesetAccess::timeSteps()
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
/// Get list of result names
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
QStringList RifEclipseRestartFilesetAccess::resultNames(RifReaderInterface::PorosityModelResultType matrixOrFracture)
|
||||
void RifEclipseRestartFilesetAccess::resultNames(QStringList* resultNames, std::vector<size_t>* resultDataItemCounts)
|
||||
{
|
||||
CVF_ASSERT(numTimeSteps() > 0);
|
||||
|
||||
// Get the results found on the first file
|
||||
QStringList resultsList;
|
||||
m_files[0]->validKeywords(&resultsList, matrixOrFracture);
|
||||
std::vector<size_t> valueCountForOneFile;
|
||||
RifEclipseOutputFileTools::findKeywordsAndDataItemCounts(m_ecl_files[0], resultNames, &valueCountForOneFile);
|
||||
|
||||
return resultsList;
|
||||
for (size_t i = 0; i < valueCountForOneFile.size(); i++)
|
||||
{
|
||||
resultDataItemCounts->push_back(valueCountForOneFile[i] * numTimeSteps());
|
||||
}
|
||||
}
|
||||
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
/// Get result values for given time step
|
||||
//--------------------------------------------------------------------------------------------------
|
||||
bool RifEclipseRestartFilesetAccess::results(const QString& resultName, RifReaderInterface::PorosityModelResultType matrixOrFracture, size_t timeStep, std::vector<double>* values)
|
||||
bool RifEclipseRestartFilesetAccess::results(const QString& resultName, size_t timeStep, size_t gridCount, std::vector<double>* values)
|
||||
{
|
||||
size_t numOccurrences = m_files[timeStep]->numOccurrences(resultName);
|
||||
size_t numOccurrences = ecl_file_get_num_named_kw(m_ecl_files[timeStep], resultName.toAscii().data());
|
||||
|
||||
// No results for this result variable for current time step found
|
||||
if (numOccurrences == 0) return true;
|
||||
|
||||
// Result handling depends on presents of result values for all grids
|
||||
if (m_gridCount != numOccurrences)
|
||||
if (gridCount != numOccurrences)
|
||||
{
|
||||
return false;
|
||||
}
|
||||
@@ -158,7 +161,8 @@ bool RifEclipseRestartFilesetAccess::results(const QString& resultName, RifReade
|
||||
for (i = 0; i < numOccurrences; i++)
|
||||
{
|
||||
std::vector<double> partValues;
|
||||
if (!m_files[timeStep]->keywordData(resultName, i, matrixOrFracture, &partValues)) // !! don't need to append afterwards
|
||||
|
||||
if (!RifEclipseOutputFileTools::keywordData(m_ecl_files[timeStep], resultName, i, &partValues))
|
||||
{
|
||||
return false;
|
||||
}
|
||||
@@ -177,13 +181,13 @@ void RifEclipseRestartFilesetAccess::readWellData(well_info_type* well_info)
|
||||
{
|
||||
if (!well_info) return;
|
||||
|
||||
for (size_t i = 0; i < m_files.size(); i++)
|
||||
for (size_t i = 0; i < m_ecl_files.size(); i++)
|
||||
{
|
||||
const char* fileName = ecl_file_get_src_file(m_files[i]->filePointer());
|
||||
const char* fileName = ecl_file_get_src_file(m_ecl_files[i]);
|
||||
int reportNumber = ecl_util_filename_report_nr(fileName);
|
||||
if(reportNumber != -1)
|
||||
{
|
||||
well_info_add_wells(well_info, m_files[i]->filePointer(), reportNumber);
|
||||
well_info_add_wells(well_info, m_ecl_files[i], reportNumber);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user