Internalize EDITNNCR and make it available with the NNC information.

If there is an EDITNNCR entry and a NNC entry for the same cell pair
then the EDITNNCR entry is not be represented inernally but we simply
overwrite the transmissibility of the corresponding NNC entry.

If there is an EDITNNC entry and an EDITNNCR entry for the same cell
pair then the EDITNNC entry is removed while internalizing EDITNNCR.

Order matters for EDITNNCR entries in the sense that later specified
values will overwrite previous entries when we internalize.

Note that similar to EDITNNC only the first 7 options are represented
and the rest is ignored by OPM flow.

Note that EDITNNCR entries for neighboring cells are ignored (like for
EDITNNC).
This commit is contained in:
Markus Blatt 2023-03-31 00:03:57 +02:00
parent fb3c2ab5a1
commit f2ffebc814
3 changed files with 218 additions and 0 deletions

View File

@ -111,8 +111,10 @@ public:
bool addNNC(const size_t cell1, const size_t cell2, const double trans);
const std::vector<NNCdata>& input() const { return m_input; }
const std::vector<NNCdata>& edit() const { return m_edit; }
const std::vector<NNCdata>& editr() const { return m_editr; }
KeywordLocation input_location(const NNCdata& nnc) const;
KeywordLocation edit_location(const NNCdata& nnc) const;
KeywordLocation editr_location(const NNCdata& nnc) const;
bool operator==(const NNC& data) const;
@ -122,20 +124,25 @@ public:
{
serializer(m_input);
serializer(m_edit);
serializer(m_editr);
serializer(m_nnc_location);
serializer(m_edit_location);
serializer(m_editr_location);
}
private:
void load_input(const EclipseGrid& grid, const Deck& deck);
void load_edit(const EclipseGrid& grid, const Deck& deck);
void load_editr(const EclipseGrid& grid, const Deck& deck);
void add_edit(const NNCdata& edit_node);
std::vector<NNCdata> m_input;
std::vector<NNCdata> m_edit;
std::vector<NNCdata> m_editr;
std::optional<KeywordLocation> m_nnc_location;
std::optional<KeywordLocation> m_edit_location;
std::optional<KeywordLocation> m_editr_location;
};

View File

@ -90,6 +90,7 @@ bool is_neighbor(const EclipseGrid& grid, std::size_t g1, std::size_t g2) {
NNC::NNC(const EclipseGrid& grid, const Deck& deck) {
this->load_input(grid, deck);
this->load_edit(grid, deck);
this->load_editr(grid, deck);
}
@ -138,6 +139,9 @@ bool is_neighbor(const EclipseGrid& grid, std::size_t g1, std::size_t g2) {
std::sort(nnc_edit.begin(), nnc_edit.end());
// If we have a corresponding NNC already, then we apply
// the multiplier from EDITNNC to it. Otherwise we internalize
// it into m_edit
auto current_input = this->m_input.begin();
for (const auto& current_edit : nnc_edit) {
if (current_input == this->m_input.end()) {
@ -171,14 +175,127 @@ bool is_neighbor(const EclipseGrid& grid, std::size_t g1, std::size_t g2) {
}
}
void NNC::load_editr(const EclipseGrid& grid, const Deck& deck) {
using NNCinsert = std::tuple<std::array<std::size_t,3>,double>;
std::vector<NNCinsert> nnc_editr;
std::size_t insertion_index = 0;
const auto& keyword_list = deck.getKeywordList<ParserKeywords::EDITNNCR>();
if (keyword_list.empty()) {
return;
}
for (const auto& keyword_ptr : deck.getKeywordList<ParserKeywords::EDITNNCR>()) {
const auto& records = *keyword_ptr;
if (records.empty()) {
continue;
}
for (const auto& record : records) {
auto index_pair = make_index_pair(grid, record);
if (!index_pair)
continue;
auto [g1, g2] = index_pair.value();
if (is_neighbor(grid, g1, g2))
continue;
double trans = record.getItem(6).getSIDouble(0);
nnc_editr.push_back( {{g1, g2, insertion_index }, trans});
}
if (!this->m_editr_location)
this->m_editr_location = keyword_ptr->location();
}
if (nnc_editr.empty()) {
return;
}
// Make the entries unique (i.e. only one entry for each index pair)
// the entry surviving should be last one (i.e. later enrtries overwrite previous entries)
std::sort(nnc_editr.begin(), nnc_editr.end(),
[](const NNCinsert& d1, NNCinsert& d2) {
return std::get<0>(d1) < std::get<0>(d2);
}
);
std::size_t g1 = -1, g2 = -1;
auto target = nnc_editr.begin() - 1;
for(auto current = nnc_editr.begin(); current != nnc_editr.end(); ++current)
{
const auto& indices = std::get<0>(*current);
if ( indices[0] == g1 && indices[1] == g2){
// overwrites previous entry
*target = std::move(*current);
}
else{
++target;
if (target != current) {
*target = std::move(*current);
}
g1 = indices[0];
g2 = indices[1];
}
}
nnc_editr.resize((++target) - nnc_editr.begin());
// Remove corresponding EDITNNC entries in m_edit as EDIT_NNCR
// will overwrite transmissibilities anyway
std::vector<NNCdata> slim_edit;
slim_edit.reserve(m_edit.size());
struct Comp
{
bool operator()(const NNCinsert& d1, NNCdata& d2){
return std::get<0>(d1)[0] < d2.cell1 && std::get<0>(d1)[1] < d2.cell2 ;
}
bool operator()(const NNCdata& d1, NNCinsert& d2){
return d1.cell1 < std::get<0>(d2)[0] && d1.cell2 < std::get<0>(d2)[1];
}
};
std::set_difference(m_edit.begin(), m_edit.end(),
nnc_editr.begin(), nnc_editr.end(),
std::back_inserter(slim_edit),
Comp());
m_edit = std::move(slim_edit);
// If we have a corresponding NNC already, then we overwrite
// its transmissibility from EDITNNCR. Otherwise we internalize it
// in m_editr
auto current_input = this->m_input.begin();
for (const auto& current_editr : nnc_editr) {
if (current_input == this->m_input.end()) {
m_editr.push_back({std::get<0>(current_editr)[0], std::get<0>(current_editr)[1], std::get<double>(current_editr)});
continue;
}
if (current_input->cell1 != std::get<0>(current_editr)[0] || current_input->cell2 != std::get<0>(current_editr)[1]) {
current_input = std::lower_bound(current_input,
this->m_input.end(),
NNCdata(std::get<0>(current_editr)[0], std::get<0>(current_editr)[1], 0));
}
if (current_input != this->m_input.end()
&& current_input->cell1 == std::get<0>(current_editr)[0]
&& current_input->cell2 == std::get<0>(current_editr)[1]) {
current_input->trans = std::get<double>(current_editr);
}
else {
m_editr.push_back({std::get<0>(current_editr)[0], std::get<0>(current_editr)[1], std::get<double>(current_editr)});
}
}
}
NNC NNC::serializationTestObject()
{
NNC result;
result.m_input= {{1,2,1.0},{2,3,2.0}};
result.m_edit= {{1,2,1.0},{2,3,2.0}};
result.m_editr= {{1,2,1.0},{2,3,2.0}};
result.m_nnc_location = {"NNC?", "File", 123};
result.m_edit_location = {"EDITNNC?", "File", 123};
result.m_edit_location = {"EDITNNCR?", "File", 123};
return result;
}
@ -237,5 +354,11 @@ bool is_neighbor(const EclipseGrid& grid, std::size_t g1, std::size_t g2) {
return {};
}
KeywordLocation NNC::editr_location([[maybe_unused]] const NNCdata& nnc) const {
if (this->m_editr_location)
return *this->m_editr_location;
else
return {};
}
} // namespace Opm

View File

@ -220,6 +220,65 @@ EDITNNC
/
)";
const std::string editnncr_input = R"(
RUNSPEC
OIL
GAS
WATER
DIMENS
10 10 1 /
GRID
NNC
7 1 1 3 1 1 0.1 / -- Transmissibility will be overwritten with 10 by EDITNNCR
3 1 1 5 1 1 0.1 /
/
DXV
10*1000.0
/
DYV
10*1000.0
/
DZ
100*20.0
/
TOPS
100*10
/
PORO
100*0.15 /
EDIT
EDITNNC
5 1 1 3 1 1 2.0 / -- Will be removed after internalization as there is an NNC
3 1 1 1 1 1 0.1 / -- stays internalized
2 1 1 2 3 1 0.3 / -- Will be removed after internalization as there is an EDTNNCR
/
EDITNNCR
8 1 1 3 1 1 2.0 /
1 1 1 1 2 1 0.01 / -- This is ignored because the two cells are ijk neighbours
2 1 1 2 3 1 0.2 / -- Overwritten with 0.3 by next entry
7 1 1 3 1 1 10 / -- Will overwrite entry in NNC
/
EDITNNCR
1 1 1 2 1 1 0.1 / -- Ignored
2 1 1 2 3 1 0.3 / -- Overwrites previous value
/
)";
std::optional<NNCdata> find_nnc(const std::vector<NNCdata>& v, std::size_t c1, std::size_t c2) {
if (c1 > c2)
@ -270,6 +329,7 @@ void check_order(const std::vector<NNCdata>& d) {
void check_order(const NNC& nnc) {
check_order(nnc.input());
check_order(nnc.edit());
check_order(nnc.editr());
}
@ -314,6 +374,7 @@ BOOST_AUTO_TEST_CASE(noNNC_EDIT)
EclipseState eclipseState(deck);
const auto& editnnc = eclipseState.getInputNNC();
BOOST_CHECK(editnnc.edit().empty());
BOOST_CHECK(editnnc.editr().empty());
}
@ -339,6 +400,33 @@ BOOST_AUTO_TEST_CASE(readDeck_EDIT)
check_nnc(input, grid.getGlobalIndex(2,0,0), grid.getGlobalIndex(6,0,0), 0.10);
}
BOOST_AUTO_TEST_CASE(readDeck_EDITR)
{
Parser parser;
auto deck = parser.parseString(editnncr_input);
EclipseGrid grid(10,10,10);
NNC nnc(grid, deck);
const std::vector<NNCdata>& input = nnc.input();
const std::vector<NNCdata>& edit = nnc.edit();
const std::vector<NNCdata>& editr = nnc.editr();
check_order(nnc);
BOOST_CHECK_EQUAL(input.size(), 2);
check_nnc(input, grid.getGlobalIndex(6,0,0), grid.getGlobalIndex(2,0,0), 10);
check_nnc(input, grid.getGlobalIndex(2,0,0), grid.getGlobalIndex(4,0,0), 0.1*2.0);
BOOST_CHECK_EQUAL(edit.size(), 1);
BOOST_CHECK(!has_nnc(edit, grid.getGlobalIndex(4,0,0), grid.getGlobalIndex(2,0,0)));
check_edit_nnc(edit, grid.getGlobalIndex(2,0,0), grid.getGlobalIndex(0,0,0), 0.1);
BOOST_CHECK_EQUAL(editr.size(), 2);
BOOST_CHECK(!has_nnc(editr, grid.getGlobalIndex(0,0,0), grid.getGlobalIndex(0,1,0)));
BOOST_CHECK(!has_nnc(editr, grid.getGlobalIndex(6,0,0), grid.getGlobalIndex(2,0,0)));
BOOST_CHECK(!has_nnc(editr, grid.getGlobalIndex(0,0,0), grid.getGlobalIndex(1,0,0)));
check_nnc(editr, grid.getGlobalIndex(1,0,0), grid.getGlobalIndex(1,2,0), 0.3);
check_nnc(editr, grid.getGlobalIndex(2,0,0), grid.getGlobalIndex(7,0,0), 2.0);
}
BOOST_AUTO_TEST_CASE(ACTNUM)
{