ISTLSolverBda.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 IRIS AS
3 Copyright 2019, 2020 Equinor ASA
4 Copyright 2020 SINTEF Digital, Mathematics and Cybernetics
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#ifndef OPM_ISTLSOLVER_WITH_BDA_HEADER_INCLUDED
23#define OPM_ISTLSOLVER_WITH_BDA_HEADER_INCLUDED
24
26
27#include <cstddef>
28#include <memory>
29#include <set>
30#include <string>
31#include <vector>
32
33namespace Opm {
34
35class Well;
36
37template<class Matrix, class Vector, int block_size> class BdaBridge;
38template<class Scalar> class WellContributions;
39namespace detail {
40
41template<class Matrix, class Vector>
43{
44 using Scalar = typename Vector::field_type;
45 using WellContribFunc = std::function<void(WellContributions<Scalar>&)>;
47
48 BdaSolverInfo(const std::string& accelerator_mode,
49 const int linear_solver_verbosity,
50 const int maxit,
51 const Scalar tolerance,
52 const int platformID,
53 const int deviceID,
54 const bool opencl_ilu_parallel,
55 const std::string& linsolver);
56
58
59 template<class Grid>
60 void prepare(const Grid& grid,
61 const Dune::CartesianIndexMapper<Grid>& cartMapper,
62 const std::vector<Well>& wellsForConn,
63 const std::unordered_map<std::string, std::set<int>>& possibleFutureConnections,
64 const std::vector<int>& cellPartition,
65 const std::size_t nonzeroes,
66 const bool useWellConn);
67
68 bool apply(Vector& rhs,
69 const bool useWellConn,
70 WellContribFunc getContribs,
71 const int rank,
72 Matrix& matrix,
73 Vector& x,
75
76 bool gpuActive();
77
79
80private:
83 template<class Grid>
84 void blockJacobiAdjacency(const Grid& grid,
85 const std::vector<int>& cell_part,
86 std::size_t nonzeroes);
87
88 void copyMatToBlockJac(const Matrix& mat, Matrix& blockJac);
89
90 std::unique_ptr<Bridge> bridge_;
91 std::string accelerator_mode_;
92 std::unique_ptr<Matrix> blockJacobiForGPUILU0_;
93 std::vector<std::set<int>> wellConnectionsGraph_;
94};
95
96}
97
102template <class TypeTag>
103class ISTLSolverBda : public ISTLSolver<TypeTag>
104{
105protected:
114 using Matrix = typename SparseMatrixAdapter::IstlMatrix;
117 using AbstractSolverType = Dune::InverseOperator<Vector, Vector>;
118 using AbstractOperatorType = Dune::AssembledLinearOperator<Matrix, Vector, Vector>;
123
124#if HAVE_MPI
125 using CommunicationType = Dune::OwnerOverlapCopyCommunication<int,int>;
126#else
127 using CommunicationType = Dune::Communication<int>;
128#endif
129
130public:
131 using AssembledLinearOperatorType = Dune::AssembledLinearOperator< Matrix, Vector, Vector >;
132
137 ISTLSolverBda(const Simulator& simulator, const FlowLinearSolverParameters& parameters)
138 : ParentType(simulator, parameters)
139 {
141 }
142
145 explicit ISTLSolverBda(const Simulator& simulator)
146 : ParentType(simulator)
147 {
149 }
150
152 {
153 OPM_TIMEBLOCK(initializeBda);
154
155 std::string accelerator_mode = Parameters::Get<Parameters::AcceleratorMode>();
156 // Force accelerator mode to none if using MPI.
157 if ((this->simulator_.vanguard().grid().comm().size() > 1) && (accelerator_mode != "none")) {
158 const bool on_io_rank = (this->simulator_.gridView().comm().rank() == 0);
159 if (on_io_rank) {
160 OpmLog::warning("Cannot use AcceleratorMode feature with MPI, setting AcceleratorMode to 'none'.");
161 }
162 accelerator_mode = "none";
163 }
164
165 if (accelerator_mode == "none") {
166 return;
167 }
168
169 // Initialize the BdaBridge
170 const int platformID = Parameters::Get<Parameters::OpenclPlatformId>();
171 const int deviceID = Parameters::Get<Parameters::BdaDeviceId>();
172 const int maxit = Parameters::Get<Parameters::LinearSolverMaxIter>();
173 const double tolerance = Parameters::Get<Parameters::LinearSolverReduction>();
174 const bool opencl_ilu_parallel = Parameters::Get<Parameters::OpenclIluParallel>();
175 const int linear_solver_verbosity = this->parameters_[0].linear_solver_verbosity_;
176 std::string linsolver = Parameters::Get<Parameters::LinearSolver>();
177 bdaBridge_ = std::make_unique<detail::BdaSolverInfo<Matrix,Vector>>(accelerator_mode,
178 linear_solver_verbosity,
179 maxit,
180 tolerance,
181 platformID,
182 deviceID,
183 opencl_ilu_parallel,
184 linsolver);
185 }
186
187 void prepare(const Matrix& M, Vector& b)
188 {
189 OPM_TIMEBLOCK(prepare);
190 [[maybe_unused]] const bool firstcall = (this->matrix_ == nullptr);
191
192 // Avoid performing the decomposition on CPU when we also do it on GPU,
193 // but we do need to initialize the pointers.
194 if (bdaBridge_) {
196 } else {
198 }
199
200#if HAVE_OPENCL || HAVE_ROCSPARSE || HAVE_CUDA
201 // update matrix entries for solvers.
202 if (firstcall && bdaBridge_) {
203 // model will not change the matrix object. Hence simply store a pointer
204 // to the original one with a deleter that does nothing.
205 // Outch! We need to be able to scale the linear system! Hence const_cast
206 // setup sparsity pattern for jacobi matrix for preconditioner (only used for openclSolver)
207 bdaBridge_->numJacobiBlocks_ = Parameters::Get<Parameters::NumJacobiBlocks>();
208 bdaBridge_->prepare(this->simulator_.vanguard().grid(),
209 this->simulator_.vanguard().cartesianIndexMapper(),
210 this->simulator_.vanguard().schedule().getWellsatEnd(),
211 this->simulator_.vanguard().schedule().getPossibleFutureConnections(),
212 this->simulator_.vanguard().cellPartition(),
213 this->getMatrix().nonzeroes(), this->useWellConn_);
214 }
215#endif
216 }
217
218
219 void setResidual(Vector& /* b */)
220 {
221 // rhs_ = &b; // Must be handled in prepare() instead.
222 }
223
224 void getResidual(Vector& b) const
225 {
226 b = *(this->rhs_);
227 }
228
229 void setMatrix(const SparseMatrixAdapter& /* M */)
230 {
231 // matrix_ = &M.istlMatrix(); // Must be handled in prepare() instead.
232 }
233
234 bool solve(Vector& x)
235 {
236 if (!bdaBridge_) {
237 return ParentType::solve(x);
238 }
239
240 OPM_TIMEBLOCK(istlSolverBdaSolve);
241 this->solveCount_ += 1;
242 // Write linear system if asked for.
243 const int verbosity = this->prm_[this->activeSolverNum_].template get<int>("verbosity", 0);
244 const bool write_matrix = verbosity > 10;
245 if (write_matrix) {
246 Helper::writeSystem(this->simulator_, //simulator is only used to get names
247 this->getMatrix(),
248 *(this->rhs_),
249 this->comm_.get());
250 }
251
252 // Solve system.
254
255 std::function<void(WellContributions<Scalar>&)> getContribs =
257 {
258 this->simulator_.problem().wellModel().getWellContributions(w);
259 };
260 if (!bdaBridge_->apply(*(this->rhs_), this->useWellConn_, getContribs,
261 this->simulator_.gridView().comm().rank(),
262 const_cast<Matrix&>(this->getMatrix()),
263 x, result))
264 {
265 if(bdaBridge_->gpuActive()){
266 // bda solve fails use istl solver setup need to be done since it is not setup in prepare
268 }
269 assert(this->flexibleSolver_[this->activeSolverNum_].solver_);
270 this->flexibleSolver_[this->activeSolverNum_].solver_->apply(x, *(this->rhs_), result);
271 }
272
273 // Check convergence, iterations etc.
274 this->checkConvergence(result);
275
276 return this->converged_;
277 }
278
279protected:
280 std::unique_ptr<detail::BdaSolverInfo<Matrix, Vector>> bdaBridge_;
281}; // end ISTLSolver
282
283} // namespace Opm
284
285#endif // OPM_ISTLSOLVER_WITH_BDA_HEADER_INCLUDED
Definition: CollectDataOnIORank.hpp:49
Interface class adding the update() method to the preconditioner interface.
Definition: PreconditionerWithUpdate.hpp:32
BdaBridge acts as interface between opm-simulators with the BdaSolvers.
Definition: BdaBridge.hpp:37
Definition: ISTLSolverBda.hpp:104
ISTLSolverBda(const Simulator &simulator, const FlowLinearSolverParameters &parameters)
Definition: ISTLSolverBda.hpp:137
ISTLSolverBda(const Simulator &simulator)
Definition: ISTLSolverBda.hpp:145
void prepare(const Matrix &M, Vector &b)
Definition: ISTLSolverBda.hpp:187
static constexpr std::size_t pressureIndex
Definition: ISTLSolverBda.hpp:122
void setResidual(Vector &)
Definition: ISTLSolverBda.hpp:219
bool solve(Vector &x)
Definition: ISTLSolverBda.hpp:234
void initializeBda()
Definition: ISTLSolverBda.hpp:151
typename SparseMatrixAdapter::IstlMatrix Matrix
Definition: ISTLSolverBda.hpp:114
GetPropType< TypeTag, Properties::GlobalEqVector > Vector
Definition: ISTLSolverBda.hpp:110
void setMatrix(const SparseMatrixAdapter &)
Definition: ISTLSolverBda.hpp:229
void getResidual(Vector &b) const
Definition: ISTLSolverBda.hpp:224
std::unique_ptr< detail::BdaSolverInfo< Matrix, Vector > > bdaBridge_
Definition: ISTLSolverBda.hpp:280
Definition: ISTLSolver.hpp:144
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: ISTLSolver.hpp:147
std::shared_ptr< CommunicationType > comm_
Definition: ISTLSolver.hpp:660
std::vector< FlowLinearSolverParameters > parameters_
Definition: ISTLSolver.hpp:656
GetPropType< TypeTag, Properties::GridView > GridView
Definition: ISTLSolver.hpp:146
Dune::InverseOperator< Vector, Vector > AbstractSolverType
Definition: ISTLSolver.hpp:156
GetPropType< TypeTag, Properties::WellModel > WellModel
Definition: ISTLSolver.hpp:151
int solveCount_
Definition: ISTLSolver.hpp:641
Matrix & getMatrix()
Definition: ISTLSolver.hpp:629
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: ISTLSolver.hpp:148
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
Definition: ISTLSolver.hpp:168
Matrix * matrix_
Definition: ISTLSolver.hpp:646
void prepareFlexibleSolver()
Definition: ISTLSolver.hpp:488
GetPropType< TypeTag, Properties::ThreadManager > ThreadManager
Definition: ISTLSolver.hpp:154
GetPropType< TypeTag, Properties::ElementMapper > ElementMapper
Definition: ISTLSolver.hpp:160
void initPrepare(const Matrix &M, Vector &b)
Definition: ISTLSolver.hpp:344
std::vector< detail::FlexibleSolverInfo< Matrix, Vector, CommunicationType > > flexibleSolver_
Definition: ISTLSolver.hpp:650
Dune::AssembledLinearOperator< Matrix, Vector, Vector > AbstractOperatorType
Definition: ISTLSolver.hpp:157
Dune::AssembledLinearOperator< Matrix, Vector, Vector > AssembledLinearOperatorType
Definition: ISTLSolver.hpp:174
void checkConvergence(const Dune::InverseOperatorResult &result) const
Definition: ISTLSolver.hpp:458
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: ISTLSolver.hpp:152
Vector * rhs_
Definition: ISTLSolver.hpp:647
int activeSolverNum_
Definition: ISTLSolver.hpp:649
GetPropType< TypeTag, Properties::Indices > Indices
Definition: ISTLSolver.hpp:150
const Simulator & simulator_
Definition: ISTLSolver.hpp:639
bool converged_
Definition: ISTLSolver.hpp:642
void prepare(const SparseMatrixAdapter &M, Vector &b)
Definition: ISTLSolver.hpp:372
bool solve(Vector &x)
Definition: ISTLSolver.hpp:410
std::vector< PropertyTree > prm_
Definition: ISTLSolver.hpp:658
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: ISTLSolver.hpp:155
Definition: WellContributions.hpp:53
Definition: WellOperators.hpp:70
void writeSystem(const SimulatorType &simulator, const MatrixType &matrix, const VectorType &rhs, const Communicator *comm)
Definition: WriteSystemMatrixHelper.hpp:34
Definition: blackoilboundaryratevector.hh:37
Dune::InverseOperatorResult InverseOperatorResult
Definition: BdaBridge.hpp:32
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition: FlowLinearSolverParameters.hpp:95
Definition: ISTLSolverBda.hpp:43
typename Vector::field_type Scalar
Definition: ISTLSolverBda.hpp:44
bool apply(Vector &rhs, const bool useWellConn, WellContribFunc getContribs, const int rank, Matrix &matrix, Vector &x, Dune::InverseOperatorResult &result)
std::function< void(WellContributions< Scalar > &)> WellContribFunc
Definition: ISTLSolverBda.hpp:45
BdaSolverInfo(const std::string &accelerator_mode, const int linear_solver_verbosity, const int maxit, const Scalar tolerance, const int platformID, const int deviceID, const bool opencl_ilu_parallel, const std::string &linsolver)
int numJacobiBlocks_
Definition: ISTLSolverBda.hpp:78
void prepare(const Grid &grid, const Dune::CartesianIndexMapper< Grid > &cartMapper, const std::vector< Well > &wellsForConn, const std::unordered_map< std::string, std::set< int > > &possibleFutureConnections, const std::vector< int > &cellPartition, const std::size_t nonzeroes, const bool useWellConn)