21#ifndef OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
22#define OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
26#include <opm/simulators/linalg/PropertyTree.hpp>
27#include <opm/simulators/linalg/matrixblock.hh>
35 using PressureMatrixType = Dune::BCRSMatrix<Opm::MatrixBlock<double, 1, 1>>;
36 using PressureVectorType = Dune::BlockVector<Dune::FieldVector<double, 1>>;
37 using SeqCoarseOperatorType = Dune::MatrixAdapter<PressureMatrixType, PressureVectorType, PressureVectorType>;
39 using ParCoarseOperatorType
40 = Dune::OverlappingSchwarzOperator<PressureMatrixType, PressureVectorType, PressureVectorType, Comm>;
42 using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
43 SeqCoarseOperatorType,
44 ParCoarseOperatorType<Comm>>;
49template <
class FineOperator,
class Communication,
bool transpose = false>
54 typedef typename Details::CoarseOperatorType<Communication> CoarseOperator;
56 typedef Communication ParallelInformation;
57 typedef typename FineOperator::domain_type FineVectorType;
61 const FineVectorType& weights,
64 : communication_(&
const_cast<Communication&
>(comm))
72 using CoarseMatrix =
typename CoarseOperator::matrix_type;
75 auto createIter = coarseLevelMatrix_->createbegin();
85 coarseLevelCommunication_.reset(communication_, [](Communication*) {});
87 this->
lhs_.resize(this->coarseLevelMatrix_->M());
88 this->
rhs_.resize(this->coarseLevelMatrix_->N());
89 using OperatorArgs =
typename Dune::Amg::ConstructionTraits<CoarseOperator>::Arguments;
91 this->
operator_ = Dune::Amg::ConstructionTraits<CoarseOperator>::construct(
oargs);
97 *coarseLevelMatrix_ = 0;
98 auto rowCoarse = coarseLevelMatrix_->begin();
106 const auto&
bw = weights_[
entry.index()];
107 for (
size_t i = 0;
i <
bw.size(); ++
i) {
111 const auto&
bw = weights_[
row.index()];
112 for (
size_t i = 0;
i <
bw.size(); ++
i) {
127 auto end =
fine.end(), begin =
fine.begin();
130 const auto&
bw = weights_[
block.index()];
133 rhs_el = (*block)[pressure_var_index_];
135 for (
size_t i = 0; i < block->size(); ++i) {
136 rhs_el += (*block)[i] * bw[i];
139 this->
rhs_[block - begin] = rhs_el;
147 auto end =
fine.end(), begin =
fine.begin();
151 const auto&
bw = weights_[
block.index()];
152 for (
size_t i = 0;
i <
block->size(); ++
i) {
156 (*block)[pressure_var_index_] = this->
lhs_[
block - begin];
166 const Communication& getCoarseLevelCommunication()
const
168 return *coarseLevelCommunication_;
171 std::size_t getPressureIndex()
const
173 return pressure_var_index_;
176 Communication* communication_;
177 const FineVectorType& weights_;
178 const std::size_t pressure_var_index_;
179 std::shared_ptr<Communication> coarseLevelCommunication_;
180 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 PressureTransferPolicy.hpp:52
virtual PressureTransferPolicy * clone() const override
Clone the current object.
Definition PressureTransferPolicy.hpp:161
virtual void createCoarseLevelSystem(const FineOperator &fineOperator) override
Algebraically creates the coarse level system.
Definition PressureTransferPolicy.hpp:70
virtual void calculateCoarseEntries(const FineOperator &fineOperator) override
???.
Definition PressureTransferPolicy.hpp:94
virtual void moveToFineLevel(typename ParentType::FineDomainType &fine) override
Updates the fine level linear system after the correction of the coarse levels system.
Definition PressureTransferPolicy.hpp:145
Definition PropertyTree.hpp:37
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition BlackoilPhases.hpp:27
Algebraic twolevel methods.