using AquiferConfig to handle AQUFLUX aquifers

This commit is contained in:
Kai Bao
2023-02-13 13:38:46 +01:00
parent d412158cf7
commit f081ee5826
5 changed files with 42 additions and 14 deletions

View File

@@ -23,6 +23,7 @@
#include <opm/input/eclipse/EclipseState/Aquifer/Aquancon.hpp>
#include <opm/input/eclipse/EclipseState/Aquifer/Aquifetp.hpp>
#include <opm/input/eclipse/EclipseState/Aquifer/AquiferCT.hpp>
#include <opm/input/eclipse/EclipseState/Aquifer/AquiferFlux.hpp>
#include <opm/input/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquifers.hpp>
#include <cstddef>
@@ -43,10 +44,12 @@ namespace Opm {
class AquiferConfig {
public:
using AquFluxs = std::unordered_map<int, AquiferFlux>;
AquiferConfig() = default;
AquiferConfig(const TableManager& tables, const EclipseGrid& grid,
const Deck& deck, const FieldPropsManager& field_props);
AquiferConfig(const Aquifetp& fetp, const AquiferCT& ct, const Aquancon& conn);
AquiferConfig(const Aquifetp& fetp, const AquiferCT& ct, const AquFluxs& aqufluxs, const Aquancon& conn);
void load_connections(const Deck& deck, const EclipseGrid& grid);
void pruneDeactivatedAquiferConnections(const std::vector<std::size_t>& deactivated_cells);
@@ -58,6 +61,7 @@ public:
bool active() const;
const AquiferCT& ct() const;
const Aquifetp& fetp() const;
const AquFluxs& aquflux() const;
const Aquancon& connections() const;
bool operator==(const AquiferConfig& other) const;
bool hasAquifer(const int aquID) const;
@@ -74,12 +78,15 @@ public:
serializer(aquifetp);
serializer(aquiferct);
serializer(aqconn);
// TODO:
// serializer(aquiferflux);
serializer(numerical_aquifers);
}
private:
Aquifetp aquifetp{};
AquiferCT aquiferct{};
AquFluxs aquiferflux{};
mutable NumericalAquifers numerical_aquifers{};
Aquancon aqconn{};
};

View File

@@ -40,6 +40,8 @@ namespace Opm {
// to work with ScheduleState::map_member
int name() const;
bool operator==(const AquiferFlux& other) const;
static std::unordered_map<int, AquiferFlux> aqufluxFromKeywords(const std::vector<const DeckKeyword*>& keywords);
};
} // end of namespace Opm

View File

@@ -18,6 +18,7 @@
*/
#include <opm/input/eclipse/Deck/DeckSection.hpp>
#include <opm/input/eclipse/EclipseState/Aquifer/AquiferConfig.hpp>
#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
#include <opm/input/eclipse/EclipseState/Grid/EclipseGrid.hpp>
@@ -33,14 +34,17 @@ AquiferConfig::AquiferConfig(const TableManager& tables,
const FieldPropsManager& field_props)
: aquifetp(tables, deck)
, aquiferct(tables, deck)
, aquiferflux(AquiferFlux::aqufluxFromKeywords(SOLUTIONSection(deck).getKeywordList("AQUFLUX")) )
, numerical_aquifers(deck, grid, field_props)
{}
AquiferConfig::AquiferConfig(const Aquifetp& fetp,
const AquiferCT& ct,
const AquFluxs& aqufluxs,
const Aquancon& conn)
: aquifetp(fetp)
, aquiferct(ct)
, aquiferflux(aqufluxs)
, aqconn(conn)
{}
@@ -62,6 +66,7 @@ void AquiferConfig::loadFromRestart(const RestartIO::RstAquifer& aquifers,
this->aquifetp.loadFromRestart(aquifers, tables);
this->aquiferct.loadFromRestart(aquifers, tables);
this->aqconn.loadFromRestart(aquifers);
// TODO: loadFromrestart for AQUFLUX;
}
AquiferConfig AquiferConfig::serializationTestObject()
@@ -71,6 +76,7 @@ AquiferConfig AquiferConfig::serializationTestObject()
result.aquiferct = AquiferCT::serializationTestObject();
result.aqconn = Aquancon::serializationTestObject();
result.numerical_aquifers = NumericalAquifers::serializationTestObject();
// TODO: serializationTestObject for AQUFLUX
return result;
}
@@ -83,6 +89,7 @@ bool AquiferConfig::active() const {
bool AquiferConfig::operator==(const AquiferConfig& other) const {
return this->aquifetp == other.aquifetp &&
this->aquiferct == other.aquiferct &&
this->aquiferflux == other.aquiferflux &&
this->aqconn == other.aqconn &&
this->numerical_aquifers == other.numerical_aquifers;
}
@@ -99,6 +106,10 @@ const Aquancon& AquiferConfig::connections() const {
return this->aqconn;
}
const AquiferConfig::AquFluxs& AquiferConfig::aquflux() const {
return this->aquiferflux;
}
bool AquiferConfig::hasAquifer(const int aquID) const {
return this->hasAnalyticalAquifer(aquID) ||
numerical_aquifers.hasAquifer(aquID);
@@ -106,7 +117,8 @@ bool AquiferConfig::hasAquifer(const int aquID) const {
bool AquiferConfig::hasAnalyticalAquifer(const int aquID) const {
return aquifetp.hasAquifer(aquID) ||
aquiferct.hasAquifer(aquID);
aquiferct.hasAquifer(aquID) ||
aquiferflux.count(aquID) > 0;
}
bool AquiferConfig::hasNumericalAquifer() const {
@@ -123,7 +135,8 @@ NumericalAquifers& AquiferConfig::mutableNumericalAquifers() const {
bool AquiferConfig::hasAnalyticalAquifer() const {
return (this->aquiferct.size() > std::size_t{0})
|| (this->aquifetp.size() > std::size_t{0});
|| (this->aquifetp.size() > std::size_t{0})
|| !this->aquiferflux.empty();
}
}
@@ -141,6 +154,9 @@ std::vector<int> Opm::analyticAquiferIDs(const AquiferConfig& cfg)
for (const auto& aquifer : cfg.fetp())
aquiferIDs.push_back(aquifer.aquiferID);
for ([[maybe_unused]] const auto& [id, aqu] : cfg.aquflux())
aquiferIDs.push_back(id);
std::sort(aquiferIDs.begin(), aquiferIDs.end());
return aquiferIDs;

View File

@@ -38,6 +38,14 @@ namespace Opm {
}
}
bool AquiferFlux::operator==(const AquiferFlux& other) const {
return this->id == other.id &&
this->flux == other.flux &&
this->salt_concentration == other.salt_concentration &&
this->temperature == other.temperature &&
this->datum_pressure == other.datum_pressure;
}
int AquiferFlux::name() const
{
return this->id;

View File

@@ -883,7 +883,8 @@ Opm::RestartIO::getSimulationTimePoint(const std::time_t start,
namespace {
int getNumberOfAnalyticAquifers(const Opm::AquiferConfig& cfg)
{
const auto numAnalyticAquifers = cfg.ct().size() + cfg.fetp().size();
// TODO: this can be function for AquiferConfig
const auto numAnalyticAquifers = cfg.ct().size() + cfg.fetp().size() + cfg.aquflux().size();
return static_cast<int>(numAnalyticAquifers);
}
@@ -912,16 +913,10 @@ namespace {
int getMaximumAnalyticAquiferID(const Opm::AquiferConfig& cfg)
{
auto maxAquID = [](const auto& aquiferCollection) -> int
{
return std::accumulate(aquiferCollection.begin(), aquiferCollection.end(), 0,
[](const int maxID, const auto& aquiferData)
{
return std::max(maxID, aquiferData.aquiferID);
});
};
return std::max(maxAquID(cfg.ct()), maxAquID(cfg.fetp()));
// TODO: this can be based on AQUCT and AQUFETP and AQUFLUX
// but there are pros and cons
const auto& aquifer_ids = analyticAquiferIDs(cfg);
return *max_element(std::begin(aquifer_ids), std::end(aquifer_ids));
}
}