23#include <opm/common/TimingMacros.hpp>
25#include <opm/simulators/linalg/matrixblock.hh>
26#include <opm/simulators/linalg/PropertyTree.hpp>
29#include <dune/istl/paamg/pinfo.hh>
33 template <
class Communication>
34 void extendCommunicatorWithWells(
const Communication& comm,
35 std::shared_ptr<Communication>& commRW,
38 OPM_TIMEBLOCK(extendCommunicatorWithWells);
40 using IndexSet =
typename Communication::ParallelIndexSet;
41 using LocalIndex =
typename IndexSet::LocalIndex;
42 const IndexSet& indset = comm.indexSet();
43 IndexSet& indset_rw = commRW->indexSet();
44 const int max_nw = comm.communicator().max(nw);
45 const int rank = comm.communicator().rank();
48 indset_rw.beginResize();
49 for (
auto ind = indset.begin(), indend = indset.end(); ind != indend; ++ind) {
50 indset_rw.add(ind->global(), LocalIndex(ind->local(), ind->local().attribute(),
true));
51 const int glo = ind->global();
52 const size_t loc = ind->local().local();
53 glo_max = std::max(glo_max, glo);
54 loc_max = std::max(loc_max, loc);
56 const int global_max = comm.communicator().max(glo_max);
58 assert(loc_max + 1 == indset.size());
59 size_t local_ind = loc_max + 1;
60 for (
int i = 0; i < nw; ++i) {
62 const size_t v = global_max + max_nw * rank + i + 1;
64 indset_rw.add(v, LocalIndex(local_ind, Dune::OwnerOverlapCopyAttributeSet::owner,
true));
67 indset_rw.endResize();
68 assert(indset_rw.size() == indset.size() + nw);
70 commRW->remoteIndices().setNeighbours(comm.remoteIndices().getNeighbours());
71 commRW->remoteIndices().template rebuild<true>();
77 using PressureMatrixType = Dune::BCRSMatrix<Opm::MatrixBlock<double, 1, 1>>;
78 using PressureVectorType = Dune::BlockVector<Dune::FieldVector<double, 1>>;
79 using SeqCoarseOperatorType = Dune::MatrixAdapter<PressureMatrixType, PressureVectorType, PressureVectorType>;
81 using ParCoarseOperatorType
82 = Dune::OverlappingSchwarzOperator<PressureMatrixType, PressureVectorType, PressureVectorType, Comm>;
84 using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
85 SeqCoarseOperatorType,
86 ParCoarseOperatorType<Comm>>;
89 template <
class FineOperator,
class Communication,
bool transpose = false>
93 typedef typename Details::CoarseOperatorType<Communication> CoarseOperator;
95 typedef Communication ParallelInformation;
96 typedef typename FineOperator::domain_type FineVectorType;
100 const FineVectorType& weights,
102 const std::size_t pressureIndex)
103 : communication_(&
const_cast<Communication&
>(comm))
106 , pressure_var_index_(pressureIndex)
113 using CoarseMatrix =
typename CoarseOperator::matrix_type;
116 if (prm_.get<
bool>(
"add_wells")) {
124 CoarseMatrix::implicit));
128 coarseLevelMatrix_->entry(
rownum,
col.index()) = 0.0;
133 coarseLevelMatrix_.reset(
135 auto createIter = coarseLevelMatrix_->createbegin();
143 if constexpr (std::is_same_v<Communication, Dune::Amg::SequentialInformation>) {
144 coarseLevelCommunication_ = std::make_shared<Communication>();
146 coarseLevelCommunication_ = std::make_shared<Communication>(
147 communication_->communicator(), communication_->category(),
false);
149 if (prm_.get<
bool>(
"add_wells")) {
150 fineOperator.addWellPressureEquationsStruct(*coarseLevelMatrix_);
151 coarseLevelMatrix_->compress();
152 if constexpr (!std::is_same_v<Communication, Dune::Amg::SequentialInformation>) {
153 extendCommunicatorWithWells(*communication_, coarseLevelCommunication_,
nw);
158 this->
lhs_.resize(this->coarseLevelMatrix_->M());
159 this->
rhs_.resize(this->coarseLevelMatrix_->N());
160 using OperatorArgs =
typename Dune::Amg::ConstructionTraits<CoarseOperator>::Arguments;
162 this->
operator_ = Dune::Amg::ConstructionTraits<CoarseOperator>::construct(
oargs);
169 *coarseLevelMatrix_ = 0;
170 auto rowCoarse = coarseLevelMatrix_->begin();
178 const auto&
bw = weights_[
entry.index()];
179 for (
size_t i = 0;
i <
bw.size(); ++
i) {
183 const auto&
bw = weights_[
row.index()];
184 for (
size_t i = 0;
i <
bw.size(); ++
i) {
191 if (prm_.get<
bool>(
"add_wells")) {
211 auto end =
fine.end(), begin =
fine.begin();
214 const auto&
bw = weights_[
block.index()];
217 rhs_el = (*block)[pressure_var_index_];
219 for (
size_t i = 0; i < block->size(); ++i) {
220 rhs_el += (*block)[i] * bw[i];
223 this->
rhs_[block - begin] = rhs_el;
233 auto end =
fine.end(), begin =
fine.begin();
237 const auto&
bw = weights_[
block.index()];
238 for (
size_t i = 0;
i <
block->size(); ++
i) {
242 (*block)[pressure_var_index_] = this->
lhs_[
block - begin];
252 const Communication& getCoarseLevelCommunication()
const
254 return *coarseLevelCommunication_;
258 Communication* communication_;
259 const FineVectorType& weights_;
261 const int pressure_var_index_;
262 std::shared_ptr<Communication> coarseLevelCommunication_;
263 std::shared_ptr<typename CoarseOperator::matrix_type> coarseLevelMatrix_;
Abstract base class for transfer between levels and creation of the coarse level system.
Definition twolevelmethodcpr.hh:44
FineOperatorType::range_type FineRangeType
The type of the range of the fine level operator.
Definition twolevelmethodcpr.hh:54
CoarseDomainType lhs_
The coarse level lhs.
Definition twolevelmethodcpr.hh:141
CoarseRangeType rhs_
The coarse level rhs.
Definition twolevelmethodcpr.hh:139
FineOperatorType::domain_type FineDomainType
The type of the domain of the fine level operator.
Definition twolevelmethodcpr.hh:58
std::shared_ptr< CoarseOperatorType > operator_
the coarse level linear operator.
Definition twolevelmethodcpr.hh:143
Definition AquiferInterface.hpp:35
Definition PressureBhpTransferPolicy.hpp:91
virtual void moveToFineLevel(typename ParentType::FineDomainType &fine) override
Updates the fine level linear system after the correction of the coarse levels system.
Definition PressureBhpTransferPolicy.hpp:229
virtual void calculateCoarseEntries(const FineOperator &fineOperator) override
???.
Definition PressureBhpTransferPolicy.hpp:165
virtual PressureBhpTransferPolicy * clone() const override
Clone the current object.
Definition PressureBhpTransferPolicy.hpp:247
virtual void createCoarseLevelSystem(const FineOperator &fineOperator) override
Algebraically creates the coarse level system.
Definition PressureBhpTransferPolicy.hpp:110
Definition PropertyTree.hpp:37
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Algebraic twolevel methods.