BlackOilAquiferModel: use aquifers through base interface

This commit is contained in:
Arne Morten Kvarving 2022-08-10 14:28:00 +02:00
parent 271d5bee3a
commit 4917c860fe
2 changed files with 31 additions and 135 deletions

View File

@ -113,28 +113,17 @@ protected:
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>; using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>; using Scalar = GetPropType<TypeTag, Properties::Scalar>;
typedef AquiferCarterTracy<TypeTag> AquiferCarterTracy_object;
typedef AquiferFetkovich<TypeTag> AquiferFetkovich_object;
Simulator& simulator_; Simulator& simulator_;
// TODO: probably we can use one variable to store both types of aquifers, because std::vector<std::unique_ptr<AquiferInterface<TypeTag>>> aquifers;
// they share the base class
mutable std::vector<AquiferCarterTracy_object> aquifers_CarterTracy;
mutable std::vector<AquiferFetkovich_object> aquifers_Fetkovich;
std::vector<AquiferNumerical<TypeTag>> aquifers_numerical;
// This initialization function is used to connect the parser objects with the ones needed by AquiferCarterTracy // This initialization function is used to connect the parser objects with the ones needed by AquiferCarterTracy
void init(); void init();
bool aquiferActive() const;
bool aquiferCarterTracyActive() const;
bool aquiferFetkovichActive() const;
bool aquiferNumericalActive() const;
}; };
} // namespace Opm } // namespace Opm
#include "BlackoilAquiferModel_impl.hpp" #include "BlackoilAquiferModel_impl.hpp"
#endif #endif

View File

@ -19,7 +19,6 @@
*/ */
#include <opm/grid/utility/cartesianToCompressed.hpp> #include <opm/grid/utility/cartesianToCompressed.hpp>
#include "BlackoilAquiferModel.hpp"
namespace Opm namespace Opm
{ {
@ -39,45 +38,16 @@ template <typename TypeTag>
void void
BlackoilAquiferModel<TypeTag>::initialSolutionApplied() BlackoilAquiferModel<TypeTag>::initialSolutionApplied()
{ {
if (aquiferCarterTracyActive()) { for (auto& aquifer : aquifers)
for (auto& aquifer : aquifers_CarterTracy) { aquifer->initialSolutionApplied();
aquifer.initialSolutionApplied();
}
}
if (aquiferFetkovichActive()) {
for (auto& aquifer : aquifers_Fetkovich) {
aquifer.initialSolutionApplied();
}
}
if (this->aquiferNumericalActive()) {
for (auto& aquifer : this->aquifers_numerical) {
aquifer.initialSolutionApplied();
}
}
} }
template <typename TypeTag> template <typename TypeTag>
void void
BlackoilAquiferModel<TypeTag>::initFromRestart(const data::Aquifers& aquiferSoln) BlackoilAquiferModel<TypeTag>::initFromRestart(const data::Aquifers& aquiferSoln)
{ {
if (this->aquiferCarterTracyActive()) { for (auto& aquifer : this->aquifers)
for (auto& aquifer : this->aquifers_CarterTracy) { aquifer->initFromRestart(aquiferSoln);
aquifer.initFromRestart(aquiferSoln);
}
}
if (this->aquiferFetkovichActive()) {
for (auto& aquifer : this->aquifers_Fetkovich) {
aquifer.initFromRestart(aquiferSoln);
}
}
if (this->aquiferNumericalActive()) {
for (auto& aquifer : this->aquifers_numerical) {
aquifer.initFromRestart(aquiferSoln);
}
}
} }
template <typename TypeTag> template <typename TypeTag>
@ -94,16 +64,8 @@ template <typename TypeTag>
void void
BlackoilAquiferModel<TypeTag>::beginTimeStep() BlackoilAquiferModel<TypeTag>::beginTimeStep()
{ {
if (aquiferCarterTracyActive()) { for (auto& aquifer : aquifers)
for (auto& aquifer : aquifers_CarterTracy) { aquifer->beginTimeStep();
aquifer.beginTimeStep();
}
}
if (aquiferFetkovichActive()) {
for (auto& aquifer : aquifers_Fetkovich) {
aquifer.beginTimeStep();
}
}
} }
template <typename TypeTag> template <typename TypeTag>
@ -114,16 +76,8 @@ BlackoilAquiferModel<TypeTag>::addToSource(RateVector& rates,
unsigned spaceIdx, unsigned spaceIdx,
unsigned timeIdx) const unsigned timeIdx) const
{ {
if (aquiferCarterTracyActive()) { for (auto& aquifer : aquifers)
for (auto& aquifer : aquifers_CarterTracy) { aquifer->addToSource(rates, context, spaceIdx, timeIdx);
aquifer.addToSource(rates, context, spaceIdx, timeIdx);
}
}
if (aquiferFetkovichActive()) {
for (auto& aquifer : aquifers_Fetkovich) {
aquifer.addToSource(rates, context, spaceIdx, timeIdx);
}
}
} }
template <typename TypeTag> template <typename TypeTag>
@ -132,16 +86,8 @@ BlackoilAquiferModel<TypeTag>::addToSource(RateVector& rates,
unsigned globalSpaceIdx, unsigned globalSpaceIdx,
unsigned timeIdx) const unsigned timeIdx) const
{ {
if (aquiferCarterTracyActive()) { for (auto& aquifer : aquifers)
for (auto& aquifer : aquifers_CarterTracy) { aquifer->addToSource(rates, globalSpaceIdx, timeIdx);
aquifer.addToSource(rates, globalSpaceIdx, timeIdx);
}
}
if (aquiferFetkovichActive()) {
for (auto& aquifer : aquifers_Fetkovich) {
aquifer.addToSource(rates, globalSpaceIdx, timeIdx);
}
}
} }
template <typename TypeTag> template <typename TypeTag>
@ -153,21 +99,12 @@ template <typename TypeTag>
void void
BlackoilAquiferModel<TypeTag>::endTimeStep() BlackoilAquiferModel<TypeTag>::endTimeStep()
{ {
if (aquiferCarterTracyActive()) { for (auto& aquifer : aquifers) {
for (auto& aquifer : aquifers_CarterTracy) { aquifer->endTimeStep();
aquifer.endTimeStep(); using NumAq = AquiferNumerical<TypeTag>;
} NumAq* num = dynamic_cast<NumAq*>(aquifer.get());
} if (num)
if (aquiferFetkovichActive()) {
for (auto& aquifer : aquifers_Fetkovich) {
aquifer.endTimeStep();
}
}
if (aquiferNumericalActive()) {
for (auto& aquifer : this->aquifers_numerical) {
aquifer.endTimeStep();
this->simulator_.vanguard().grid().comm().barrier(); this->simulator_.vanguard().grid().comm().barrier();
}
} }
} }
@ -207,13 +144,15 @@ BlackoilAquiferModel<TypeTag>::init()
const auto& connections = aquifer.connections(); const auto& connections = aquifer.connections();
for (const auto& aq : aquifer.ct()) { for (const auto& aq : aquifer.ct()) {
aquifers_CarterTracy.emplace_back(connections[aq.aquiferID], auto aqf = std::make_unique<AquiferCarterTracy<TypeTag>>(connections[aq.aquiferID],
this->simulator_, aq); this->simulator_, aq);
aquifers.push_back(std::move(aqf));
} }
for (const auto& aq : aquifer.fetp()) { for (const auto& aq : aquifer.fetp()) {
aquifers_Fetkovich.emplace_back(connections[aq.aquiferID], auto aqf = std::make_unique<AquiferFetkovich<TypeTag>>(connections[aq.aquiferID],
this->simulator_, aq); this->simulator_, aq);
aquifers.push_back(std::move(aqf));
} }
if (aquifer.hasNumericalAquifer()) { if (aquifer.hasNumericalAquifer()) {
@ -223,55 +162,23 @@ BlackoilAquiferModel<TypeTag>::init()
const std::unordered_map<int, int> cartesian_to_compressed = cartesianToCompressed(number_of_cells, const std::unordered_map<int, int> cartesian_to_compressed = cartesianToCompressed(number_of_cells,
global_cell); global_cell);
for ([[maybe_unused]]const auto& [id, aqu] : num_aquifers) { for ([[maybe_unused]]const auto& [id, aqu] : num_aquifers) {
this->aquifers_numerical.emplace_back(aqu, auto aqf = std::make_unique<AquiferNumerical<TypeTag>>(aqu,
cartesian_to_compressed, this->simulator_, global_cell); cartesian_to_compressed,
this->simulator_,
global_cell);
aquifers.push_back(std::move(aqf));
} }
} }
} }
template <typename TypeTag>
bool
BlackoilAquiferModel<TypeTag>::aquiferCarterTracyActive() const
{
return !aquifers_CarterTracy.empty();
}
template <typename TypeTag>
bool
BlackoilAquiferModel<TypeTag>::aquiferFetkovichActive() const
{
return !aquifers_Fetkovich.empty();
}
template<typename TypeTag>
bool
BlackoilAquiferModel<TypeTag>::aquiferNumericalActive() const
{
return !(this->aquifers_numerical.empty());
}
template<typename TypeTag> template<typename TypeTag>
data::Aquifers BlackoilAquiferModel<TypeTag>::aquiferData() const data::Aquifers BlackoilAquiferModel<TypeTag>::aquiferData() const
{ {
data::Aquifers data; data::Aquifers data;
if (this->aquiferCarterTracyActive()) { for (const auto& aqu : this->aquifers)
for (const auto& aqu : this->aquifers_CarterTracy) { data.insert_or_assign(aqu->aquiferID(), aqu->aquiferData());
data.insert_or_assign(aqu.aquiferID(), aqu.aquiferData());
}
}
if (this->aquiferFetkovichActive()) {
for (const auto& aqu : this->aquifers_Fetkovich) {
data.insert_or_assign(aqu.aquiferID(), aqu.aquiferData());
}
}
if (this->aquiferNumericalActive()) {
for (const auto& aqu : this->aquifers_numerical) {
data.insert_or_assign(aqu.aquiferID(), aqu.aquiferData());
}
}
return data; return data;
} }
} // namespace Opm } // namespace Opm