PressureTransferPolicy.hpp
Go to the documentation of this file.
1/*
2 Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
3 Copyright 2019 Dr. Blatt - HPC-Simulation-Software & Services.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
22#define OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
23
24
27#include <opm/simulators/linalg/matrixblock.hh>
28
29#include <cstddef>
30
31namespace Opm
32{
33
34namespace Details
35{
36 using PressureMatrixType = Dune::BCRSMatrix<Opm::MatrixBlock<double, 1, 1>>;
37 using PressureVectorType = Dune::BlockVector<Dune::FieldVector<double, 1>>;
38 using SeqCoarseOperatorType = Dune::MatrixAdapter<PressureMatrixType, PressureVectorType, PressureVectorType>;
39 template <class Comm>
41 = Dune::OverlappingSchwarzOperator<PressureMatrixType, PressureVectorType, PressureVectorType, Comm>;
42 template <class Comm>
43 using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
45 ParCoarseOperatorType<Comm>>;
46} // namespace Details
47
48
49
50template <class FineOperator, class Communication, bool transpose = false>
52 : public Dune::Amg::LevelTransferPolicyCpr<FineOperator, Details::CoarseOperatorType<Communication>>
53{
54public:
58 typedef typename FineOperator::domain_type FineVectorType;
59
60public:
62 const FineVectorType& weights,
63 const Opm::PropertyTree& /*prm*/,
64 int pressure_var_index)
65 : communication_(&const_cast<Communication&>(comm))
66 , weights_(weights)
67 , pressure_var_index_(pressure_var_index)
68 {
69 }
70
71 virtual void createCoarseLevelSystem(const FineOperator& fineOperator) override
72 {
73 using CoarseMatrix = typename CoarseOperator::matrix_type;
74 const auto& fineLevelMatrix = fineOperator.getmat();
75 coarseLevelMatrix_.reset(new CoarseMatrix(fineLevelMatrix.N(), fineLevelMatrix.M(), CoarseMatrix::row_wise));
76 auto createIter = coarseLevelMatrix_->createbegin();
77
78 for (const auto& row : fineLevelMatrix) {
79 for (auto col = row.begin(), cend = row.end(); col != cend; ++col) {
80 createIter.insert(col.index());
81 }
82 ++createIter;
83 }
84
85 calculateCoarseEntries(fineOperator);
86 coarseLevelCommunication_.reset(communication_, [](Communication*) {});
87
88 this->lhs_.resize(this->coarseLevelMatrix_->M());
89 this->rhs_.resize(this->coarseLevelMatrix_->N());
90 using OperatorArgs = typename Dune::Amg::ConstructionTraits<CoarseOperator>::Arguments;
91 OperatorArgs oargs(coarseLevelMatrix_, *coarseLevelCommunication_);
92 this->operator_ = Dune::Amg::ConstructionTraits<CoarseOperator>::construct(oargs);
93 }
94
95 virtual void calculateCoarseEntries(const FineOperator& fineOperator) override
96 {
97 const auto& fineMatrix = fineOperator.getmat();
98 *coarseLevelMatrix_ = 0;
99 auto rowCoarse = coarseLevelMatrix_->begin();
100 for (auto row = fineMatrix.begin(), rowEnd = fineMatrix.end(); row != rowEnd; ++row, ++rowCoarse) {
101 assert(row.index() == rowCoarse.index());
102 auto entryCoarse = rowCoarse->begin();
103 for (auto entry = row->begin(), entryEnd = row->end(); entry != entryEnd; ++entry, ++entryCoarse) {
104 assert(entry.index() == entryCoarse.index());
105 double matrix_el = 0;
106 if (transpose) {
107 const auto& bw = weights_[entry.index()];
108 for (std::size_t i = 0; i < bw.size(); ++i) {
109 matrix_el += (*entry)[pressure_var_index_][i] * bw[i];
110 }
111 } else {
112 const auto& bw = weights_[row.index()];
113 for (std::size_t i = 0; i < bw.size(); ++i) {
114 matrix_el += (*entry)[i][pressure_var_index_] * bw[i];
115 }
116 }
117 (*entryCoarse) = matrix_el;
118 }
119 }
120 assert(rowCoarse == coarseLevelMatrix_->end());
121 }
122
123 virtual void moveToCoarseLevel(const typename ParentType::FineRangeType& fine) override
124 {
125 // Set coarse vector to zero
126 this->rhs_ = 0;
127
128 auto end = fine.end(), begin = fine.begin();
129
130 for (auto block = begin; block != end; ++block) {
131 const auto& bw = weights_[block.index()];
132 double rhs_el = 0.0;
133 if (transpose) {
134 rhs_el = (*block)[pressure_var_index_];
135 } else {
136 for (std::size_t i = 0; i < block->size(); ++i) {
137 rhs_el += (*block)[i] * bw[i];
138 }
139 }
140 this->rhs_[block - begin] = rhs_el;
141 }
142
143 this->lhs_ = 0;
144 }
145
146 virtual void moveToFineLevel(typename ParentType::FineDomainType& fine) override
147 {
148 auto end = fine.end(), begin = fine.begin();
149
150 for (auto block = begin; block != end; ++block) {
151 if (transpose) {
152 const auto& bw = weights_[block.index()];
153 for (std::size_t i = 0; i < block->size(); ++i) {
154 (*block)[i] = this->lhs_[block - begin] * bw[i];
155 }
156 } else {
157 (*block)[pressure_var_index_] = this->lhs_[block - begin];
158 }
159 }
160 }
161
162 virtual PressureTransferPolicy* clone() const override
163 {
164 return new PressureTransferPolicy(*this);
165 }
166
168 {
169 return *coarseLevelCommunication_;
170 }
171
172 std::size_t getPressureIndex() const
173 {
174 return pressure_var_index_;
175 }
176
177private:
178 Communication* communication_;
179 const FineVectorType& weights_;
180 const std::size_t pressure_var_index_;
181 std::shared_ptr<Communication> coarseLevelCommunication_;
182 std::shared_ptr<typename CoarseOperator::matrix_type> coarseLevelMatrix_;
183};
184
185} // namespace Opm
186
187#endif // OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
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: PressureTransferPolicy.hpp:53
virtual PressureTransferPolicy * clone() const override
Clone the current object.
Definition: PressureTransferPolicy.hpp:162
virtual void createCoarseLevelSystem(const FineOperator &fineOperator) override
Algebraically creates the coarse level system.
Definition: PressureTransferPolicy.hpp:71
virtual void calculateCoarseEntries(const FineOperator &fineOperator) override
???.
Definition: PressureTransferPolicy.hpp:95
PressureTransferPolicy(const Communication &comm, const FineVectorType &weights, const Opm::PropertyTree &, int pressure_var_index)
Definition: PressureTransferPolicy.hpp:61
Dune::Amg::LevelTransferPolicyCpr< FineOperator, CoarseOperator > ParentType
Definition: PressureTransferPolicy.hpp:56
FineOperator::domain_type FineVectorType
Definition: PressureTransferPolicy.hpp:58
Communication ParallelInformation
Definition: PressureTransferPolicy.hpp:57
const Communication & getCoarseLevelCommunication() const
Definition: PressureTransferPolicy.hpp:167
virtual void moveToCoarseLevel(const typename ParentType::FineRangeType &fine) override
Definition: PressureTransferPolicy.hpp:123
Details::CoarseOperatorType< Communication > CoarseOperator
Definition: PressureTransferPolicy.hpp:55
std::size_t getPressureIndex() const
Definition: PressureTransferPolicy.hpp:172
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:146
Definition: PropertyTree.hpp:37
std::conditional_t< std::is_same< Comm, Dune::Amg::SequentialInformation >::value, SeqCoarseOperatorType, ParCoarseOperatorType< Comm > > CoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:88
Dune::BlockVector< Dune::FieldVector< double, 1 > > PressureVectorType
Definition: PressureBhpTransferPolicy.hpp:80
Dune::MatrixAdapter< PressureMatrixType, PressureVectorType, PressureVectorType > SeqCoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:81
Dune::OverlappingSchwarzOperator< PressureMatrixType, PressureVectorType, PressureVectorType, Comm > ParCoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:84
Dune::BCRSMatrix< Opm::MatrixBlock< double, 1, 1 > > PressureMatrixType
Definition: PressureBhpTransferPolicy.hpp:79
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27
Algebraic twolevel methods.