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
29
30#include <cstddef>
31
32namespace Opm { namespace Details {
33 template<class Scalar> using PressureMatrixType = Dune::BCRSMatrix<MatrixBlock<Scalar, 1, 1>>;
34 template<class Scalar> using PressureVectorType = Dune::BlockVector<Dune::FieldVector<Scalar, 1>>;
35 template<class Scalar> using SeqCoarseOperatorType = Dune::MatrixAdapter<PressureMatrixType<Scalar>,
36 PressureVectorType<Scalar>,
37 PressureVectorType<Scalar>>;
38 template<class Scalar, class Comm>
41 PressureVectorType<Scalar>,
42 PressureVectorType<Scalar>,
43 Comm>;
44 template<class Scalar, class Comm>
45 using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
46 SeqCoarseOperatorType<Scalar>,
47 ParCoarseOperatorType<Scalar,Comm>>;
48} // namespace Details
49
50
51
52template <class FineOperator, class Communication, class Scalar, bool transpose = false>
54 : public Dune::Amg::LevelTransferPolicyCpr<FineOperator, Details::CoarseOperatorType<Scalar,Communication>>
55{
56public:
60 using FineVectorType = typename FineOperator::domain_type;
61
62public:
64 const FineVectorType& weights,
65 const PropertyTree& /*prm*/,
66 int pressure_var_index)
67 : communication_(&const_cast<Communication&>(comm))
68 , weights_(weights)
69 , pressure_var_index_(pressure_var_index)
70 {
71 }
72
73 void createCoarseLevelSystem(const FineOperator& fineOperator) override
74 {
75 using CoarseMatrix = typename CoarseOperator::matrix_type;
76 const auto& fineLevelMatrix = fineOperator.getmat();
77 coarseLevelMatrix_.reset(new CoarseMatrix(fineLevelMatrix.N(), fineLevelMatrix.M(), CoarseMatrix::row_wise));
78 auto createIter = coarseLevelMatrix_->createbegin();
79
80 for (const auto& row : fineLevelMatrix) {
81 for (auto col = row.begin(), cend = row.end(); col != cend; ++col) {
82 createIter.insert(col.index());
83 }
84 ++createIter;
85 }
86
87 calculateCoarseEntries(fineOperator);
88 coarseLevelCommunication_.reset(communication_, [](Communication*) {});
89
90 this->lhs_.resize(this->coarseLevelMatrix_->M());
91 this->rhs_.resize(this->coarseLevelMatrix_->N());
92 using OperatorArgs = typename Dune::Amg::ConstructionTraits<CoarseOperator>::Arguments;
93 OperatorArgs oargs(coarseLevelMatrix_, *coarseLevelCommunication_);
94 this->operator_ = Dune::Amg::ConstructionTraits<CoarseOperator>::construct(oargs);
95 }
96
97 void calculateCoarseEntries(const FineOperator& fineOperator) override
98 {
99 const auto& fineMatrix = fineOperator.getmat();
100 *coarseLevelMatrix_ = 0;
101 auto rowCoarse = coarseLevelMatrix_->begin();
102 for (auto row = fineMatrix.begin(), rowEnd = fineMatrix.end(); row != rowEnd; ++row, ++rowCoarse) {
103 assert(row.index() == rowCoarse.index());
104 auto entryCoarse = rowCoarse->begin();
105 for (auto entry = row->begin(), entryEnd = row->end(); entry != entryEnd; ++entry, ++entryCoarse) {
106 assert(entry.index() == entryCoarse.index());
107 Scalar matrix_el = 0;
108 if (transpose) {
109 const auto& bw = weights_[entry.index()];
110 for (std::size_t i = 0; i < bw.size(); ++i) {
111 matrix_el += (*entry)[pressure_var_index_][i] * bw[i];
112 }
113 } else {
114 const auto& bw = weights_[row.index()];
115 for (std::size_t i = 0; i < bw.size(); ++i) {
116 matrix_el += (*entry)[i][pressure_var_index_] * bw[i];
117 }
118 }
119 (*entryCoarse) = matrix_el;
120 }
121 }
122 assert(rowCoarse == coarseLevelMatrix_->end());
123 }
124
125 void moveToCoarseLevel(const typename ParentType::FineRangeType& fine) override
126 {
127 // Set coarse vector to zero
128 this->rhs_ = 0;
129
130 auto end = fine.end(), begin = fine.begin();
131
132 for (auto block = begin; block != end; ++block) {
133 const auto& bw = weights_[block.index()];
134 Scalar rhs_el = 0.0;
135 if (transpose) {
136 rhs_el = (*block)[pressure_var_index_];
137 } else {
138 for (std::size_t i = 0; i < block->size(); ++i) {
139 rhs_el += (*block)[i] * bw[i];
140 }
141 }
142 this->rhs_[block - begin] = rhs_el;
143 }
144
145 this->lhs_ = 0;
146 }
147
148 void moveToFineLevel(typename ParentType::FineDomainType& fine) override
149 {
150 auto end = fine.end(), begin = fine.begin();
151
152 for (auto block = begin; block != end; ++block) {
153 if (transpose) {
154 const auto& bw = weights_[block.index()];
155 for (std::size_t i = 0; i < block->size(); ++i) {
156 (*block)[i] = this->lhs_[block - begin] * bw[i];
157 }
158 } else {
159 (*block)[pressure_var_index_] = this->lhs_[block - begin];
160 }
161 }
162 }
163
165 {
166 return new PressureTransferPolicy(*this);
167 }
168
170 {
171 return *coarseLevelCommunication_;
172 }
173
174 std::size_t getPressureIndex() const
175 {
176 return pressure_var_index_;
177 }
178
179private:
180 Communication* communication_;
181 const FineVectorType& weights_;
182 const std::size_t pressure_var_index_;
183 std::shared_ptr<Communication> coarseLevelCommunication_;
184 std::shared_ptr<typename CoarseOperator::matrix_type> coarseLevelMatrix_;
185};
186
187} // namespace Opm
188
189#endif // OPM_PRESSURE_TRANSFER_POLICY_HEADER_INCLUDED
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:285
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
Dune linear operator that assumes ghost rows are ordered after interior rows. Avoids some computation...
Definition: WellOperators.hpp:340
Definition: PressureTransferPolicy.hpp:55
void moveToCoarseLevel(const typename ParentType::FineRangeType &fine) override
Definition: PressureTransferPolicy.hpp:125
void moveToFineLevel(typename ParentType::FineDomainType &fine) override
Updates the fine level linear system after the correction of the coarse levels system.
Definition: PressureTransferPolicy.hpp:148
const Communication & getCoarseLevelCommunication() const
Definition: PressureTransferPolicy.hpp:169
typename FineOperator::domain_type FineVectorType
Definition: PressureTransferPolicy.hpp:60
PressureTransferPolicy(const Communication &comm, const FineVectorType &weights, const PropertyTree &, int pressure_var_index)
Definition: PressureTransferPolicy.hpp:63
PressureTransferPolicy * clone() const override
Clone the current object.
Definition: PressureTransferPolicy.hpp:164
typename Details::CoarseOperatorType< Scalar, Communication > CoarseOperator
Definition: PressureTransferPolicy.hpp:57
void calculateCoarseEntries(const FineOperator &fineOperator) override
???.
Definition: PressureTransferPolicy.hpp:97
Communication ParallelInformation
Definition: PressureTransferPolicy.hpp:59
void createCoarseLevelSystem(const FineOperator &fineOperator) override
Algebraically creates the coarse level system.
Definition: PressureTransferPolicy.hpp:73
std::size_t getPressureIndex() const
Definition: PressureTransferPolicy.hpp:174
Definition: PropertyTree.hpp:37
Dune::BCRSMatrix< MatrixBlock< Scalar, 1, 1 > > PressureMatrixType
Definition: PressureBhpTransferPolicy.hpp:80
Opm::GhostLastMatrixAdapter< PressureMatrixType< Scalar >, PressureVectorType< Scalar >, PressureVectorType< Scalar >, Comm > ParCoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:90
std::conditional_t< std::is_same< Comm, Dune::Amg::SequentialInformation >::value, SeqCoarseOperatorType< Scalar >, ParCoarseOperatorType< Scalar, Comm > > CoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:94
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > PressureVectorType
Definition: PressureBhpTransferPolicy.hpp:81
Dune::MatrixAdapter< PressureMatrixType< Scalar >, PressureVectorType< Scalar >, PressureVectorType< Scalar > > SeqCoarseOperatorType
Definition: PressureBhpTransferPolicy.hpp:84
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37
Algebraic twolevel methods.