Added new option (report steps only) for compareECL for comparing summary results

This commit is contained in:
Torbjørn Skille 2019-06-30 17:08:20 +02:00
parent 597bbde151
commit 51adb23728
10 changed files with 91 additions and 7 deletions

View File

@ -317,6 +317,7 @@ if(ENABLE_ECL_OUTPUT)
tests/SUMMARY_EFF_FAC.DATA
tests/SPE1CASE1.DATA
tests/SPE1CASE1.SMSPEC
tests/SPE1CASE1A.SMSPEC
tests/SPE9_CP_PACKED.DATA
tests/SOFR_TEST.DATA
)
@ -335,6 +336,7 @@ if(ENABLE_ECL_INPUT)
tests/SPE1_TESTCASE.F0025
tests/SPE1_TESTCASE.X0025
tests/SPE1CASE1.UNSMRY
tests/SPE1CASE1A.UNSMRY
tests/SPE1CASE1_RST60.SMSPEC
tests/SPE1CASE1_RST60.UNSMRY
)

View File

@ -34,6 +34,9 @@ public:
bool hasKey(const std::string& key) const;
const std::vector<float>& get(const std::string& name) const;
std::vector<float> get_at_rstep(const std::string& name) const;
const std::vector<std::string>& keywordList() const { return keyword; }
private:

View File

@ -26,6 +26,7 @@
#include <algorithm>
#include <unistd.h>
#include <limits>
#include <limits.h>
#include <set>
#include <opm/io/eclipse/EclFile.hpp>
@ -57,8 +58,11 @@ ESmry::ESmry(const std::string &filename, bool loadBaseRunData)
std::string rootN;
bool formatted=false;
char buff[1000];
getcwd( buff, 1000 );
char buff[PATH_MAX];
if (getcwd( buff, PATH_MAX )==NULL){
throw std::invalid_argument("failed when trying to get current working directory");
}
std::string currentWorkingDir(buff);
@ -424,4 +428,18 @@ const std::vector<float>& ESmry::get(const std::string& name) const
return param[ind];
}
std::vector<float> ESmry::get_at_rstep(const std::string& name) const
{
std::vector<float> full_vector= this->get(name);
std::vector<float> rstep_vector;
rstep_vector.reserve(seqIndex.size());
for (auto ind : seqIndex){
rstep_vector.push_back(full_vector[ind]);
}
return rstep_vector;
}
}} // namespace Opm::ecl

View File

@ -824,7 +824,11 @@ void ECLRegressionTest::results_smry()
std::string reference = "Summary file";
std::cout << "\nComparing summary files " << std::endl;
if (reportStepOnly){
std::cout << " -- Values at report steps will be compared. Time steps in between reports are ignored " << std::endl;
}
std::vector<std::string> keywords1 = smry1.keywordList();
std::vector<std::string> keywords2 = smry2.keywordList();
@ -870,9 +874,18 @@ void ECLRegressionTest::results_smry()
std::cout << "\nChecking " << keywords1.size() << " vectors ... ";
for (size_t i = 0; i < keywords1.size(); i++) {
std::vector<float> vect1 = smry1.get( keywords1[i]);
std::vector<float> vect2 = smry2.get( keywords1[i]);
std::vector<float> vect1;
std::vector<float> vect2;
if (reportStepOnly){
vect1 = smry1.get_at_rstep( keywords1[i]);
vect2 = smry2.get_at_rstep( keywords1[i]);
} else {
vect1 = smry1.get( keywords1[i]);
vect2 = smry2.get( keywords1[i]);
}
if (vect1.size() != vect2.size()) {
OPM_THROW(std::runtime_error, "\nKeyword " << keywords1[i] << " summary vector of different length");
}

View File

@ -63,6 +63,10 @@ public:
// will ignore extra keywords which are only present
// in the new simulation.
void setReportStepOnly(bool reportStepOnlyArg) {
this->reportStepOnly = reportStepOnlyArg;
}
void setAcceptExtraKeywords(bool acceptExtraKeywordsArg) {
this->acceptExtraKeywords = acceptExtraKeywordsArg;
}
@ -159,6 +163,8 @@ private:
"TRANX", "TRANY", "TRANZ", "TRANNNC", "SGRP", "SCON", "DOUBHEAD"
};
bool reportStepOnly = false;
// Only compare last occurrence
bool onlyLastSequence = false;

View File

@ -34,6 +34,7 @@ static void printHelp() {
<< "In addition, the program takes these options (which must be given before the arguments):\n\n"
<< "-a Run a full analysis of errors.\n"
<< "-h Print help and exit.\n"
<< "-d Use report steps only when comparing results from summary files.\n"
<< "-i Execute integration test (regression test is default).\n"
<< " The integration test compares SGAS, SWAT and PRESSURE in unified restart files, and WOPR, WGPR, WWPR and WBHP (all wells) in summary file. \n"
<< "-k Specify specific keyword to compare (capitalized), for examples -k PRESSURE or -k WOPR:A-1H \n"
@ -66,6 +67,7 @@ static void printHelp() {
int main(int argc, char** argv) {
bool integrationTest = false;
bool onlyLastSequence = false;
bool reportStepOnly = false;
bool printKeywords = false;
bool specificKeyword = false;
bool specificReportStepNumber = false;
@ -79,7 +81,7 @@ int main(int argc, char** argv) {
int reportStepNumber = -1;
std::string fileTypeString;
while ((c = getopt(argc, argv, "hik:alnpt:Rr:x")) != -1) {
while ((c = getopt(argc, argv, "hik:alnpt:Rr:x:d")) != -1) {
switch (c) {
case 'a':
analysis = true;
@ -87,6 +89,9 @@ int main(int argc, char** argv) {
case 'h':
printHelp();
return 0;
case 'd':
reportStepOnly = true;
break;
case 'i':
integrationTest = true;
break;
@ -168,7 +173,11 @@ int main(int argc, char** argv) {
if (onlyLastSequence) {
comparator.setOnlyLastReportNumber(true);
}
if (reportStepOnly) {
comparator.setReportStepOnly(true);
}
if (specificKeyword) {
comparator.compareSpesificKeyword(keyword);
}

BIN
tests/SPE1CASE1A.SMSPEC Normal file

Binary file not shown.

BIN
tests/SPE1CASE1A.UNSMRY Normal file

Binary file not shown.

View File

@ -144,6 +144,7 @@ std::vector<float> getFrom(const std::vector<float> &ref_vect,int from){
return vect;
}
BOOST_AUTO_TEST_CASE(TestESmry_1) {
std::vector <float> time_ref, wgpr_prod_ref, wbhp_prod_ref, wbhp_inj_ref, fgor_ref, bpr_111_ref, bpr_10103_ref;
@ -337,4 +338,19 @@ BOOST_AUTO_TEST_CASE(TestESmry_3) {
}
}
BOOST_AUTO_TEST_CASE(TestESmry_4) {
std::vector<float> time_ref = {31, 59, 90, 120, 151, 181, 212, 243, 273, 304, 334, 365, 396, 424, 455, 485, 516, 546, 577, 608, 638, 669, 699, 730, 761, 789, 820, 850, 881, 911, 942, 973, 1003, 1034, 1064, 1095, 1126, 1154, 1185, 1215, 1246, 1276, 1307, 1338, 1368, 1399, 1429, 1460, 1491, 1519, 1550, 1580, 1611, 1641, 1672, 1703, 1733, 1764, 1794, 1825, 1856, 1884, 1915, 1945, 1976, 2006, 2037, 2068, 2098, 2129, 2159, 2190, 2221, 2249, 2280, 2310, 2341, 2371, 2402, 2433, 2463, 2494, 2524, 2555, 2586, 2614, 2645, 2675, 2706, 2736, 2767, 2798, 2828, 2859, 2889, 2920, 2951, 2979, 3010, 3040, 3071, 3101, 3132, 3163, 3193, 3224, 3254, 3285, 3316, 3344, 3375, 3405, 3436, 3466, 3497, 3528, 3558, 3589, 3619, 3650};
ESmry smry1("SPE1CASE1.SMSPEC");
std::vector<float> smryVect = smry1.get("TIME");
std::vector<float> smryVect_rstep = smry1.get_at_rstep("TIME");
BOOST_CHECK_EQUAL(smryVect==time_ref, false);
BOOST_CHECK_EQUAL(smryVect_rstep==time_ref, true);
}

View File

@ -241,6 +241,7 @@ void makeRftFile(const std::string &fileName,
}
BOOST_AUTO_TEST_CASE(gridCompare) {
std::vector<float> coord = {2000,2000,2000,1999.9127,1999.8691,2009.9951,2099.9849,2000,2001.7452,2099.8975,1999.8691,2011.7404,2199.9695,2000,2003.4905,2199.8823,
@ -1105,6 +1106,22 @@ BOOST_AUTO_TEST_CASE(results_unsmry_2) {
}
BOOST_AUTO_TEST_CASE(results_unsmry_3) {
ECLRegressionTest test1("SPE1CASE1", "SPE1CASE1A", 1e-3, 1e-3);
BOOST_CHECK_THROW(test1.results_smry(),std::runtime_error);
test1.setReportStepOnly(true);
BOOST_CHECK_NO_THROW(test1.results_smry());
}
BOOST_AUTO_TEST_CASE(results_rft_1) {
using Date = std::tuple<int, int, int>;