Add small class WellContainer to manage well data in WellState

This commit is contained in:
Joakim Hove 2021-05-05 10:27:51 +02:00
parent 66932936b3
commit 36cc9e8567
3 changed files with 187 additions and 0 deletions

View File

@ -188,6 +188,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/core/props/phaseUsageFromDeck.hpp
opm/core/props/satfunc/RelpermDiagnostics.hpp
opm/simulators/timestepping/SimulatorReport.hpp
opm/simulators/wells/WellContainer.hpp
opm/simulators/wells/WellState.hpp
opm/simulators/aquifers/AquiferInterface.hpp
opm/simulators/aquifers/AquiferCarterTracy.hpp

View File

@ -0,0 +1,125 @@
/*
Copyright 2021 Equinor ASA
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef OPM_WELL_CONTAINER_HEADER_INCLUDED
#define OPM_WELL_CONTAINER_HEADER_INCLUDED
#include <stdexcept>
#include <string>
#include <unordered_map>
#include <vector>
namespace Opm {
template <class T>
class WellContainer {
public:
std::size_t size() const {
return this->data.size();
}
void add(const std::string& name, T&& value) {
if (index_map.count(name) != 0)
throw std::logic_error("An object with name: " + name + " already exists in container");
this->index_map.emplace(name, this->data.size());
this->data.push_back(std::forward<T>(value));
}
bool has(const std::string& name) const {
return (index_map.count(name) != 0);
}
void update(const std::string& name, T&& value) {
auto index = this->index_map.at(name);
this->data[index] = std::forward<T>(value);
}
void update(std::size_t index, T&& value) {
this->data.at(index) = std::forward<T>(value);
}
void copy_welldata(const WellContainer<T>& other) {
if (this->index_map == other.index_map)
this->data = other.data;
else {
for (const auto& [name, index] : this->index_map)
this->update_if(index, name, other);
}
}
void copy_welldata(const WellContainer<T>& other, const std::string& name) {
auto this_index = this->index_map.at(name);
auto other_index = other.index_map.at(name);
this->data[this_index] = other.data[other_index];
}
T& operator[](std::size_t index) {
return this->data.at(index);
}
const T& operator[](std::size_t index) const {
return this->data.at(index);
}
T& operator[](const std::string& name) {
auto index = this->index_map.at(name);
return this->data[index];
}
const T& operator[](const std::string& name) const {
auto index = this->index_map.at(name);
return this->data[index];
}
void clear() {
this->data.clear();
this->index_map.clear();
}
typename std::vector<T>::const_iterator begin() const {
return this->data.begin();
}
typename std::vector<T>::const_iterator end() const {
return this->data.end();
}
private:
void update_if(std::size_t index, const std::string& name, const WellContainer<T>& other) {
auto other_iter = other.index_map.find(name);
if (other_iter == other.index_map.end())
return;
auto other_index = other_iter->second;
this->data[index] = other.data[other_index];
}
std::vector<T> data;
std::unordered_map<std::string, std::size_t> index_map;
};
}
#endif

View File

@ -24,6 +24,7 @@
#include "MpiFixture.hpp"
#include <opm/simulators/wells/GlobalWellInfo.hpp>
#include <opm/simulators/wells/WellStateFullyImplicitBlackoil.hpp>
#include <opm/simulators/wells/WellContainer.hpp>
#include <opm/parser/eclipse/Python/Python.hpp>
#include <boost/test/unit_test.hpp>
@ -422,7 +423,67 @@ BOOST_AUTO_TEST_CASE(GlobalWellInfo_TEST) {
BOOST_CHECK(!gwi.in_producing_group("PROD01"));
}
BOOST_AUTO_TEST_CASE(TESTWellContainer) {
Opm::WellContainer<int> wc;
BOOST_CHECK_EQUAL(wc.size(), 0);
wc.add("W1", 1);
wc.add("W2", 2);
BOOST_CHECK_EQUAL(wc.size(), 2);
BOOST_CHECK_THROW(wc.add("W1", 1), std::exception);
BOOST_CHECK_THROW(wc[10], std::exception);
BOOST_CHECK_EQUAL(wc[0], 1);
BOOST_CHECK_EQUAL(wc[1], 2);
BOOST_CHECK_THROW(wc["INVALID_WELL"], std::exception);
BOOST_CHECK_EQUAL(wc["W1"], 1);
BOOST_CHECK_EQUAL(wc["W2"], 2);
Opm::WellContainer<int> wc2;
wc2.copy_welldata(wc);
BOOST_CHECK_EQUAL(wc2.size() , 0);
wc2.add("W1", 100);
BOOST_CHECK_EQUAL(wc2["W1"], 100);
wc2.copy_welldata(wc);
BOOST_CHECK_EQUAL(wc2["W1"], 1);
Opm::WellContainer<int> wc3;
wc3.add("W2", 100);
wc3.copy_welldata(wc);
BOOST_CHECK_EQUAL(wc3["W2"], 2);
BOOST_CHECK_EQUAL(wc3[0], 2);
wc3["W2"] = 200;
wc3.add("W3", 300);
wc3.copy_welldata(wc);
BOOST_CHECK_EQUAL(wc3["W2"], 2);
BOOST_CHECK_EQUAL(wc3[0], 2);
BOOST_CHECK_EQUAL(wc3["W3"], 300);
BOOST_CHECK_EQUAL(wc3[1], 300);
BOOST_CHECK_THROW(wc3.copy_welldata(wc, "W1"), std::exception);
BOOST_CHECK_THROW(wc3.copy_welldata(wc, "W3"), std::exception);
wc.clear();
BOOST_CHECK_EQUAL(wc.size(), 0);
BOOST_CHECK_THROW(wc3.update("NO_SUCH_WELL", -1), std::exception);
BOOST_CHECK_THROW(wc3.update(100, -1), std::exception);
BOOST_CHECK(wc3.has("W2"));
BOOST_CHECK(!wc3.has("NO_SUCH_WELL"));
std::vector<int> vec_copy(wc3.begin(), wc3.end());
BOOST_CHECK_EQUAL(vec_copy.size(), wc3.size());
for (std::size_t i = 0; i < wc3.size(); i++)
BOOST_CHECK_EQUAL(vec_copy[i], wc3[i]);
}
BOOST_AUTO_TEST_SUITE_END()