Add Method for Identifying Numerical Aquifer Cell IDs
This commit adds a new, focused, member function NumericalAquifers::allAquiferCellIds() which returns a vector of those Cartesian/global cells that have been marked as defining the model's numerical aquifers through the AQUNUM keyword. We intend to use this to identify those NNCs that go to numerical aquifers--or between numerical aquifer cells--as those may need special treatment when processing the MULTREGT keyword.
This commit is contained in:
parent
3ab1f65485
commit
b84c5f46bf
@ -20,13 +20,15 @@
|
||||
#ifndef OPM_NUMERICALAQUIFERS_HPP
|
||||
#define OPM_NUMERICALAQUIFERS_HPP
|
||||
|
||||
#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/SingleNumericalAquifer.hpp>
|
||||
|
||||
#include <cstddef>
|
||||
#include <map>
|
||||
#include <unordered_map>
|
||||
#include <vector>
|
||||
|
||||
#include <stddef.h>
|
||||
|
||||
#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/SingleNumericalAquifer.hpp>
|
||||
|
||||
namespace Opm {
|
||||
class Deck;
|
||||
class EclipseGrid;
|
||||
@ -46,6 +48,7 @@ namespace Opm {
|
||||
bool operator==(const NumericalAquifers& other) const;
|
||||
|
||||
std::unordered_map<size_t, const NumericalAquiferCell*> allAquiferCells() const;
|
||||
std::vector<std::size_t> allAquiferCellIds() const;
|
||||
|
||||
std::unordered_map<size_t, double> aquiferCellVolumes() const;
|
||||
|
||||
|
@ -35,8 +35,12 @@
|
||||
|
||||
#include <fmt/format.h>
|
||||
|
||||
#include <algorithm>
|
||||
#include <cstddef>
|
||||
#include <unordered_set>
|
||||
#include <utility>
|
||||
#include <stdexcept>
|
||||
#include <vector>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
@ -155,6 +159,22 @@ namespace Opm {
|
||||
return cells;
|
||||
}
|
||||
|
||||
std::vector<std::size_t> NumericalAquifers::allAquiferCellIds() const {
|
||||
auto cellIds = std::vector<std::size_t>{};
|
||||
|
||||
for (const auto& aquifer_pair : this->m_aquifers) {
|
||||
const auto nc = aquifer_pair.second.numCells();
|
||||
for (auto aquCell = 0*nc; aquCell < nc; ++aquCell) {
|
||||
cellIds.push_back(aquifer_pair.second.getCellPrt(aquCell)->global_index);
|
||||
}
|
||||
}
|
||||
|
||||
std::sort(cellIds.begin(), cellIds.end());
|
||||
cellIds.erase(std::unique(cellIds.begin(), cellIds.end()), cellIds.end());
|
||||
|
||||
return cellIds;
|
||||
}
|
||||
|
||||
const std::map<size_t, SingleNumericalAquifer>& NumericalAquifers::aquifers() const {
|
||||
return this->m_aquifers;
|
||||
}
|
||||
|
@ -795,6 +795,14 @@ BOOST_AUTO_TEST_CASE(NumericalAquiferTest)
|
||||
BOOST_CHECK_EQUAL(c3->global_index, 3);
|
||||
}
|
||||
|
||||
{
|
||||
const auto numAquCells = num_aqu.allAquiferCellIds();
|
||||
const auto expect = std::vector { 0, 2, 3, };
|
||||
|
||||
BOOST_CHECK_EQUAL_COLLECTIONS(numAquCells.begin(), numAquCells.end(),
|
||||
expect .begin(), expect .end());
|
||||
}
|
||||
|
||||
// using processed actnum for numerical aquifer connection generation
|
||||
std::vector<int> new_actnum(360, 1);
|
||||
new_actnum[0] = 0;
|
||||
|
Loading…
Reference in New Issue
Block a user