From 6ce3f7a385d0154ad661443f2904ca5540da887f Mon Sep 17 00:00:00 2001 From: Lisa Julia Nebel Date: Wed, 16 Oct 2024 12:27:29 +0200 Subject: [PATCH] Add pw_info_ to the MultisegmentWellEquations and use it to fill the matrices B and C correctly --- opm/simulators/wells/MultisegmentWellEquations.cpp | 13 ++++++++++--- opm/simulators/wells/MultisegmentWellEquations.hpp | 5 ++++- opm/simulators/wells/MultisegmentWellEval.cpp | 2 +- 3 files changed, 15 insertions(+), 5 deletions(-) diff --git a/opm/simulators/wells/MultisegmentWellEquations.cpp b/opm/simulators/wells/MultisegmentWellEquations.cpp index ac30d7958..a8ca05a43 100644 --- a/opm/simulators/wells/MultisegmentWellEquations.cpp +++ b/opm/simulators/wells/MultisegmentWellEquations.cpp @@ -48,8 +48,9 @@ namespace Opm { template MultisegmentWellEquations:: -MultisegmentWellEquations(const MultisegmentWellGeneric& well) +MultisegmentWellEquations(const MultisegmentWellGeneric& well, const ParallelWellInfo& pw_info) : well_(well) + , pw_info_(pw_info) { } @@ -107,7 +108,10 @@ init(const int numPerfs, end = duneC_.createend(); row != end; ++row) { // the number of the row corresponds to the segment number now. for (const int& perf : perforations[row.index()]) { - row.insert(perf); + const int local_perf_index = pw_info_.globalToLocal(perf); + if (local_perf_index < 0) // then the perforation is not on this process + continue; + row.insert(local_perf_index); } } @@ -116,7 +120,10 @@ init(const int numPerfs, end = duneB_.createend(); row != end; ++row) { // the number of the row corresponds to the segment number now. for (const int& perf : perforations[row.index()]) { - row.insert(perf); + const int local_perf_index = pw_info_.globalToLocal(perf); + if (local_perf_index < 0) // then the perforation is not on this process + continue; + row.insert(local_perf_index); } } diff --git a/opm/simulators/wells/MultisegmentWellEquations.hpp b/opm/simulators/wells/MultisegmentWellEquations.hpp index fceaf8549..fbc50469d 100644 --- a/opm/simulators/wells/MultisegmentWellEquations.hpp +++ b/opm/simulators/wells/MultisegmentWellEquations.hpp @@ -22,6 +22,7 @@ #ifndef OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED #define OPM_MULTISEGMENTWELL_EQUATIONS_HEADER_INCLUDED +#include #include #include #include @@ -67,7 +68,7 @@ public: using OffDiagMatrixBlockWellType = Dune::FieldMatrix; using OffDiagMatWell = Dune::BCRSMatrix; - MultisegmentWellEquations(const MultisegmentWellGeneric& well); + MultisegmentWellEquations(const MultisegmentWellGeneric& well, const ParallelWellInfo& pw_info); //! \brief Setup sparsity pattern for the matrices. //! \param numPerfs Number of perforations @@ -146,6 +147,8 @@ public: // Store the global index of well perforated cells std::vector cells_; + + const ParallelWellInfo& pw_info_; }; } diff --git a/opm/simulators/wells/MultisegmentWellEval.cpp b/opm/simulators/wells/MultisegmentWellEval.cpp index 20439500b..aadf14810 100644 --- a/opm/simulators/wells/MultisegmentWellEval.cpp +++ b/opm/simulators/wells/MultisegmentWellEval.cpp @@ -57,7 +57,7 @@ MultisegmentWellEval(WellInterfaceIndices& baseif, const Pa : MultisegmentWellGeneric(baseif) , pw_info_(pw_info) , baseif_(baseif) - , linSys_(*this) + , linSys_(*this, pw_info) , primary_variables_(baseif) , segments_(this->numberOfSegments(), pw_info.communication().sum(baseif.numPerfs()), baseif) , cell_perforation_depth_diffs_(baseif_.numPerfs(), 0.0)