ISTLSolver.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_HEADER_INCLUDED
23#define OPM_ISTLSOLVER_HEADER_INCLUDED
24
25#include <dune/istl/owneroverlapcopy.hh>
26#include <dune/istl/solver.hh>
27
28#include <opm/common/CriticalError.hpp>
29#include <opm/common/ErrorMacros.hpp>
30#include <opm/common/Exceptions.hpp>
31#include <opm/common/TimingMacros.hpp>
32
33#include <opm/grid/utility/ElementChunks.hpp>
34
54
55#include <fmt/format.h>
56
57#include <any>
58#include <cstddef>
59#include <functional>
60#include <memory>
61#include <set>
62#include <sstream>
63#include <string>
64#include <tuple>
65#include <vector>
66
67namespace Opm::Properties {
68
69namespace TTag {
71 using InheritsFrom = std::tuple<FlowIstlSolverParams>;
72};
73}
74
75template <class TypeTag, class MyTypeTag>
76struct WellModel;
77
80template<class TypeTag>
81struct SparseMatrixAdapter<TypeTag, TTag::FlowIstlSolver>
82{
83private:
85 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
87
88public:
90};
91
92} // namespace Opm::Properties
93
94namespace Opm
95{
96
97
98namespace detail
99{
100
101template<class Matrix, class Vector, class Comm>
103{
104 using AbstractSolverType = Dune::InverseOperator<Vector, Vector>;
105 using AbstractOperatorType = Dune::AssembledLinearOperator<Matrix, Vector, Vector>;
107
108 void create(const Matrix& matrix,
109 bool parallel,
110 const PropertyTree& prm,
111 std::size_t pressureIndex,
112 std::function<Vector()> weightCalculator,
113 const bool forceSerial,
114 Comm* comm);
115
116 std::unique_ptr<AbstractSolverType> solver_;
117 std::unique_ptr<AbstractOperatorType> op_;
118 std::unique_ptr<LinearOperatorExtra<Vector,Vector>> wellOperator_;
120 std::size_t interiorCellNum_ = 0;
121};
122
123
124#ifdef HAVE_MPI
126void copyParValues(std::any& parallelInformation, std::size_t size,
127 Dune::OwnerOverlapCopyCommunication<int,int>& comm);
128#endif
129
132template<class Matrix>
133void makeOverlapRowsInvalid(Matrix& matrix,
134 const std::vector<int>& overlapRows);
135
138template<class Matrix, class Grid>
139std::unique_ptr<Matrix> blockJacobiAdjacency(const Grid& grid,
140 const std::vector<int>& cell_part,
141 std::size_t nonzeroes,
142 const std::vector<std::set<int>>& wellConnectionsGraph);
143}
144
149 template <class TypeTag>
150 class ISTLSolver : public AbstractISTLSolver<GetPropType<TypeTag, Properties::SparseMatrixAdapter>,
151 GetPropType<TypeTag, Properties::GlobalEqVector>>
152 {
153 protected:
161 using Matrix = typename SparseMatrixAdapter::IstlMatrix;
164 using AbstractSolverType = Dune::InverseOperator<Vector, Vector>;
165 using AbstractOperatorType = Dune::AssembledLinearOperator<Matrix, Vector, Vector>;
169 using ElementChunksType = ElementChunks<GridView, Dune::Partitions::All>;
170
172
173 enum { enablePolymerMolarWeight = getPropValue<TypeTag, Properties::EnablePolymerMW>() };
175
176#if HAVE_MPI
177 using CommunicationType = Dune::OwnerOverlapCopyCommunication<int,int>;
178#else
179 using CommunicationType = Dune::Communication<int>;
180#endif
181
182 public:
183 using AssembledLinearOperatorType = Dune::AssembledLinearOperator< Matrix, Vector, Vector >;
184
185 static void registerParameters()
186 {
188 }
189
197 ISTLSolver(const Simulator& simulator,
198 const FlowLinearSolverParameters& parameters,
199 bool forceSerial = false)
200 : simulator_(simulator),
201 iterations_( 0 ),
202 matrix_(nullptr),
203 parameters_{parameters},
204 forceSerial_(forceSerial)
205 {
206 initialize();
207 }
208
211 explicit ISTLSolver(const Simulator& simulator)
212 : simulator_(simulator),
213 iterations_( 0 ),
214 solveCount_(0),
215 matrix_(nullptr)
216 {
217 parameters_.resize(1);
218 parameters_[0].init(simulator_.vanguard().eclState().getSimulationConfig().useCPR());
219 initialize();
220 }
221
223 {
224 OPM_TIMEBLOCK(IstlSolver);
225
227 // Polymer injectivity is incompatible with the CPRW linear solver.
228 // Instead of aborting the run, emit a warning and fall back to ILU0
229 // so that polymer injectivity cases work with default parameters.
230 if (parameters_[0].linsolver_ == "cprw" || parameters_[0].linsolver_ == "hybrid") {
231 const bool on_io_rank = (simulator_.gridView().comm().rank() == 0);
232 const std::string msg =
233 fmt::format("The polymer injectivity model is incompatible with the '{}' "
234 "linear solver. Falling back to --linear-solver=ilu0.",
235 parameters_[0].linsolver_);
236 if (on_io_rank) {
237 OpmLog::warning(msg);
238 }
239 parameters_[0].linsolver_ = "ilu0";
240 }
241 }
242
243 if (parameters_[0].linsolver_ == "hybrid") {
244 // Experimental hybrid configuration.
245 // When chosen, will set up two solvers, one with CPRW
246 // and the other with ILU0 preconditioner. More general
247 // options may be added later.
248 prm_.clear();
249 parameters_.clear();
250 {
252 para.init(false);
253 para.linsolver_ = "cprw";
254 parameters_.push_back(para);
256 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
257 Parameters::IsSet<Parameters::LinearSolverReduction>()));
258 }
259 {
261 para.init(false);
262 para.linsolver_ = "ilu0";
263 parameters_.push_back(para);
265 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
266 Parameters::IsSet<Parameters::LinearSolverReduction>()));
267 }
268 // ------------
269 } else {
270 assert(parameters_.size() == 1);
271 assert(prm_.empty());
272
273 // Do a normal linear solver setup.
274 if (parameters_[0].is_nldd_local_solver_) {
276 Parameters::IsSet<Parameters::NlddLocalLinearSolverMaxIter>(),
277 Parameters::IsSet<Parameters::NlddLocalLinearSolverReduction>()));
278 }
279 else {
281 Parameters::IsSet<Parameters::LinearSolverMaxIter>(),
282 Parameters::IsSet<Parameters::LinearSolverReduction>()));
283 }
284 }
285 flexibleSolver_.resize(prm_.size());
286
287 const bool on_io_rank = (simulator_.gridView().comm().rank() == 0);
288#if HAVE_MPI
289 comm_.reset( new CommunicationType( simulator_.vanguard().grid().comm() ) );
290#endif
292
293 // For some reason simulator_.model().elementMapper() is not initialized at this stage
294 //const auto& elemMapper = simulator_.model().elementMapper(); //does not work.
295 // Set it up manually
296 ElementMapper elemMapper(simulator_.vanguard().gridView(), Dune::mcmgElementLayout());
298 useWellConn_ = Parameters::Get<Parameters::MatrixAddWellContributions>();
299 const bool ownersFirst = Parameters::Get<Parameters::OwnerCellsFirst>();
300 if (!ownersFirst) {
301 const std::string msg = "The linear solver no longer supports --owner-cells-first=false.";
302 if (on_io_rank) {
303 OpmLog::error(msg);
304 }
305 OPM_THROW_NOLOG(std::runtime_error, msg);
306 }
307
308 const int interiorCellNum_ = detail::numMatrixRowsToUseInSolver(simulator_.vanguard().grid(), true);
309 for (auto& f : flexibleSolver_) {
310 f.interiorCellNum_ = interiorCellNum_;
311 }
312
313#if HAVE_MPI
314 if (isParallel()) {
315 const std::size_t size = simulator_.vanguard().grid().leafGridView().size(0);
317 }
318#endif
319
320 // Print parameters to PRT/DBG logs.
322
323 element_chunks_ = std::make_unique<ElementChunksType>(simulator_.vanguard().gridView(), Dune::Partitions::all, ThreadManager::maxThreads());
324 }
325
326 // nothing to clean here
327 void eraseMatrix() override
328 {
329 }
330
331 void setActiveSolver(const int num) override
332 {
333 if (num > static_cast<int>(prm_.size()) - 1) {
334 OPM_THROW(std::logic_error, "Solver number " + std::to_string(num) + " not available.");
335 }
336 activeSolverNum_ = num;
337 if (simulator_.gridView().comm().rank() == 0) {
338 OpmLog::debug("Active solver = " + std::to_string(activeSolverNum_)
339 + " (" + parameters_[activeSolverNum_].linsolver_ + ")");
340 }
341 }
342
343 int numAvailableSolvers() const override
344 {
345 return flexibleSolver_.size();
346 }
347
348 void initPrepare(const Matrix& M, Vector& b)
349 {
350 const bool firstcall = (matrix_ == nullptr);
351
352 // update matrix entries for solvers.
353 if (firstcall) {
354 // model will not change the matrix object. Hence simply store a pointer
355 // to the original one with a deleter that does nothing.
356 // Outch! We need to be able to scale the linear system! Hence const_cast
357 matrix_ = const_cast<Matrix*>(&M);
358
359 useWellConn_ = Parameters::Get<Parameters::MatrixAddWellContributions>();
360 // setup sparsity pattern for jacobi matrix for preconditioner (only used for openclSolver)
361 } else {
362 // Pointers should not change
363 if ( &M != matrix_ ) {
364 OPM_THROW(std::logic_error,
365 "Matrix objects are expected to be reused when reassembling!");
366 }
367 }
368 rhs_ = &b;
369
370 // TODO: check all solvers, not just one.
371 // We use lower case as the internal canonical representation of solver names
372 std::string type = prm_[activeSolverNum_].template get<std::string>("preconditioner.type", "paroverilu0");
373 std::ranges::transform(type, type.begin(), ::tolower);
374 if (isParallel() && type != "paroverilu0") {
376 }
377 }
378
379 void prepare(const SparseMatrixAdapter& M, Vector& b) override
380 {
381 prepare(M.istlMatrix(), b);
382 }
383
384 void prepare(const Matrix& M, Vector& b) override
385 {
386 OPM_TIMEBLOCK(istlSolverPrepare);
387 try {
388 initPrepare(M,b);
389
391 } OPM_CATCH_AND_RETHROW_AS_CRITICAL_ERROR("This is likely due to a faulty linear solver JSON specification. Check for errors related to missing nodes.");
392 }
393
394
395 void setResidual(Vector& /* b */) override
396 {
397 // rhs_ = &b; // Must be handled in prepare() instead.
398 }
399
400 void getResidual(Vector& b) const override
401 {
402 b = *rhs_;
403 }
404
405 void setMatrix(const SparseMatrixAdapter& /* M */) override
406 {
407 // matrix_ = &M.istlMatrix(); // Must be handled in prepare() instead.
408 }
409
410 int getSolveCount() const override {
411 return solveCount_;
412 }
413
415 solveCount_ = 0;
416 }
417
418 bool solve(Vector& x) override
419 {
420 OPM_TIMEBLOCK(istlSolverSolve);
421 ++solveCount_;
422 // Write linear system if asked for.
423 const int verbosity = prm_[activeSolverNum_].get("verbosity", 0);
424 const bool write_matrix = verbosity > 10;
425 if (write_matrix) {
426 Helper::writeSystem(simulator_, //simulator is only used to get names
427 getMatrix(),
428 *rhs_,
429 comm_.get());
430 }
431
432 // Solve system.
434 {
435 OPM_TIMEBLOCK(flexibleSolverApply);
436 assert(flexibleSolver_[activeSolverNum_].solver_);
437 flexibleSolver_[activeSolverNum_].solver_->apply(x, *rhs_, result);
438 }
439
440 iterations_ = result.iterations;
441
442 // Check convergence, iterations etc.
443 return checkConvergence(result);
444 }
445
446
452
454 int iterations () const override { return iterations_; }
455
457 const std::any& parallelInformation() const { return parallelInformation_; }
458
459 const CommunicationType* comm() const override { return comm_.get(); }
460
461 void setDomainIndex(const int index)
462 {
463 domainIndex_ = index;
464 }
465
466 bool isNlddLocalSolver() const
467 {
468 return parameters_[activeSolverNum_].is_nldd_local_solver_;
469 }
470
471 protected:
472#if HAVE_MPI
473 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
474#endif
475
477 {
479 }
480
481 bool isParallel() const {
482#if HAVE_MPI
483 return !forceSerial_ && comm_->communicator().size() > 1;
484#else
485 return false;
486#endif
487 }
488
490 {
491 OPM_TIMEBLOCK(flexibleSolverPrepare);
492 if (shouldCreateSolver()) {
493 if (!useWellConn_) {
494 if (isNlddLocalSolver()) {
495 auto wellOp = std::make_unique<DomainWellModelAsLinearOperator<WellModel, Vector, Vector>>(simulator_.problem().wellModel());
496 wellOp->setDomainIndex(domainIndex_);
497 flexibleSolver_[activeSolverNum_].wellOperator_ = std::move(wellOp);
498 }
499 else {
500 auto wellOp = std::make_unique<WellModelOperator>(simulator_.problem().wellModel());
501 flexibleSolver_[activeSolverNum_].wellOperator_ = std::move(wellOp);
502 }
503 }
504 std::function<Vector()> weightCalculator = this->getWeightsCalculator(prm_[activeSolverNum_], getMatrix(), pressureIndex);
505 OPM_TIMEBLOCK(flexibleSolverCreate);
507 isParallel(),
510 weightCalculator,
512 comm_.get());
513 }
514 else
515 {
516 OPM_TIMEBLOCK(flexibleSolverUpdate);
517 flexibleSolver_[activeSolverNum_].pre_->update();
518 }
519 }
520
521
525 {
526 // Decide if we should recreate the solver or just do
527 // a minimal preconditioner update.
528 if (flexibleSolver_.empty()) {
529 return true;
530 }
531 if (!flexibleSolver_[activeSolverNum_].solver_) {
532 return true;
533 }
534
535 if (flexibleSolver_[activeSolverNum_].pre_->hasPerfectUpdate()) {
536 return false;
537 }
538
539 // For AMG based preconditioners, the hierarchy depends on the matrix values
540 // so it is recreated at certain intervals
541 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 0) {
542 // Always recreate solver.
543 return true;
544 }
545 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 1) {
546 // Recreate solver on the first iteration of every timestep.
547 return this->simulator_.problem().iterationContext().isFirstGlobalIteration();
548 }
549 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 2) {
550 // Recreate solver if the last solve used more than 10 iterations.
551 return this->iterations() > 10;
552 }
553 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 3) {
554 // Never recreate the solver
555 return false;
556 }
557 if (this->parameters_[activeSolverNum_].cpr_reuse_setup_ == 4) {
558 // Recreate solver every 'step' solve calls.
559 const int step = this->parameters_[activeSolverNum_].cpr_reuse_interval_;
560 const bool create = ((solveCount_ % step) == 0);
561 return create;
562 }
563 // If here, we have an invalid parameter.
564 const bool on_io_rank = (simulator_.gridView().comm().rank() == 0);
565 std::string msg = "Invalid value: " + std::to_string(this->parameters_[activeSolverNum_].cpr_reuse_setup_)
566 + " for --cpr-reuse-setup parameter, run with --help to see allowed values.";
567 if (on_io_rank) {
568 OpmLog::error(msg);
569 }
570 throw std::runtime_error(msg);
571
572 return false;
573 }
574
575
576 // Weights to make approximate pressure equations.
577 // Calculated from the storage terms (only) of the
578 // conservation equations, ignoring all other terms.
579 std::function<Vector()> getWeightsCalculator(const PropertyTree& prm,
580 const Matrix& matrix,
581 std::size_t pressIndex) const
582 {
583 std::function<Vector()> weightsCalculator;
584
585 using namespace std::string_literals;
586
587 auto preconditionerType = prm.get("preconditioner.type"s, "cpr"s);
588 // We use lower case as the internal canonical representation of solver names
589 std::ranges::transform(preconditionerType, preconditionerType.begin(), ::tolower);
590 if (preconditionerType == "cpr" || preconditionerType == "cprt"
591 || preconditionerType == "cprw" || preconditionerType == "cprwt") {
592 const bool transpose = preconditionerType == "cprt" || preconditionerType == "cprwt";
593 const bool enableThreadParallel = this->parameters_[0].cpr_weights_thread_parallel_;
594 const auto weightsType = prm.get("preconditioner.weight_type"s, "quasiimpes"s);
595 if (weightsType == "quasiimpes") {
596 // weights will be created as default in the solver
597 // assignment p = pressureIndex prevent compiler warning about
598 // capturing variable with non-automatic storage duration
599 weightsCalculator = [matrix, transpose, pressIndex, enableThreadParallel]() {
600 return Amg::getQuasiImpesWeights<Matrix, Vector>(matrix,
601 pressIndex,
602 transpose,
603 enableThreadParallel);
604 };
605 } else if ( weightsType == "trueimpes" ) {
606 weightsCalculator =
607 [this, pressIndex, enableThreadParallel]
608 {
609 Vector weights(rhs_->size());
610 ElementContext elemCtx(simulator_);
611 Amg::getTrueImpesWeights(pressIndex,
612 weights,
613 elemCtx,
614 simulator_.model(),
616 enableThreadParallel
617 );
618 return weights;
619 };
620 } else if (weightsType == "trueimpesanalytic" ) {
621 weightsCalculator =
622 [this, pressIndex, enableThreadParallel]
623 {
624 Vector weights(rhs_->size());
625 ElementContext elemCtx(simulator_);
627 weights,
628 elemCtx,
629 simulator_.model(),
631 enableThreadParallel
632 );
633 return weights;
634 };
635 } else {
636 OPM_THROW(std::invalid_argument,
637 "Weights type " + weightsType +
638 "not implemented for cpr."
639 " Please use quasiimpes, trueimpes or trueimpesanalytic.");
640 }
641 }
642 return weightsCalculator;
643 }
644
645
647 {
648 return *matrix_;
649 }
650
651 const Matrix& getMatrix() const
652 {
653 return *matrix_;
654 }
655
657 mutable int iterations_;
658 mutable int solveCount_;
660
661 // non-const to be able to scale the linear system
664
666 std::vector<detail::FlexibleSolverInfo<Matrix,Vector,CommunicationType>> flexibleSolver_;
667 std::vector<int> overlapRows_;
668 std::vector<int> interiorRows_;
669
670 int domainIndex_ = -1;
671
673
674 std::vector<FlowLinearSolverParameters> parameters_;
675 bool forceSerial_ = false;
676 std::vector<PropertyTree> prm_;
677
678 std::shared_ptr< CommunicationType > comm_;
679 std::unique_ptr<ElementChunksType> element_chunks_;
680 }; // end ISTLSolver
681
682} // namespace Opm
683
684#endif // OPM_ISTLSOLVER_HEADER_INCLUDED
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:325
Interface class adding the update() method to the preconditioner interface.
Definition: PreconditionerWithUpdate.hpp:34
Abstract interface for ISTL solvers.
Definition: AbstractISTLSolver.hpp:45
static bool checkConvergence(const Dune::InverseOperatorResult &result, const FlowLinearSolverParameters &parameters)
Check the convergence of the linear solver.
Definition: AbstractISTLSolver.hpp:192
Definition: ISTLSolver.hpp:152
void getResidual(Vector &b) const override
Definition: ISTLSolver.hpp:400
void initialize()
Definition: ISTLSolver.hpp:222
const Matrix & getMatrix() const
Definition: ISTLSolver.hpp:651
ElementChunks< GridView, Dune::Partitions::All > ElementChunksType
Definition: ISTLSolver.hpp:169
ISTLSolver(const Simulator &simulator, const FlowLinearSolverParameters &parameters, bool forceSerial=false)
Definition: ISTLSolver.hpp:197
void setDomainIndex(const int index)
Definition: ISTLSolver.hpp:461
int iterations() const override
Definition: ISTLSolver.hpp:454
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: ISTLSolver.hpp:155
std::shared_ptr< CommunicationType > comm_
Definition: ISTLSolver.hpp:678
static constexpr bool isIncompatibleWithCprw
Definition: ISTLSolver.hpp:174
std::vector< FlowLinearSolverParameters > parameters_
Definition: ISTLSolver.hpp:674
void setActiveSolver(const int num) override
Set the active solver by its index.
Definition: ISTLSolver.hpp:331
GetPropType< TypeTag, Properties::GridView > GridView
Definition: ISTLSolver.hpp:154
Dune::InverseOperator< Vector, Vector > AbstractSolverType
Definition: ISTLSolver.hpp:164
void setMatrix(const SparseMatrixAdapter &) override
Definition: ISTLSolver.hpp:405
void prepare(const Matrix &M, Vector &b) override
Definition: ISTLSolver.hpp:384
typename SparseMatrixAdapter::IstlMatrix Matrix
Definition: ISTLSolver.hpp:161
GetPropType< TypeTag, Properties::WellModel > WellModel
Definition: ISTLSolver.hpp:159
int solveCount_
Definition: ISTLSolver.hpp:658
bool isNlddLocalSolver() const
Definition: ISTLSolver.hpp:466
Matrix & getMatrix()
Definition: ISTLSolver.hpp:646
int numAvailableSolvers() const override
Get the number of available solvers.
Definition: ISTLSolver.hpp:343
@ enablePolymerMolarWeight
Definition: ISTLSolver.hpp:173
GetPropType< TypeTag, Properties::SparseMatrixAdapter > SparseMatrixAdapter
Definition: ISTLSolver.hpp:156
void prepare(const SparseMatrixAdapter &M, Vector &b) override
Definition: ISTLSolver.hpp:379
Dune::OwnerOverlapCopyCommunication< int, int > CommunicationType
Definition: ISTLSolver.hpp:177
int getSolveCount() const override
Get the count of how many times the solver has been called.
Definition: ISTLSolver.hpp:410
Matrix * matrix_
Definition: ISTLSolver.hpp:662
bool useWellConn_
Definition: ISTLSolver.hpp:672
bool shouldCreateSolver() const
Definition: ISTLSolver.hpp:524
bool checkConvergence(const Dune::InverseOperatorResult &result) const
Definition: ISTLSolver.hpp:476
static constexpr std::size_t pressureIndex
Definition: ISTLSolver.hpp:171
void prepareFlexibleSolver()
Definition: ISTLSolver.hpp:489
GetPropType< TypeTag, Properties::ThreadManager > ThreadManager
Definition: ISTLSolver.hpp:162
GetPropType< TypeTag, Properties::ElementMapper > ElementMapper
Definition: ISTLSolver.hpp:168
void resetSolveCount()
Definition: ISTLSolver.hpp:414
GetPropType< TypeTag, Properties::GlobalEqVector > Vector
Definition: ISTLSolver.hpp:157
const std::any & parallelInformation() const
Definition: ISTLSolver.hpp:457
void initPrepare(const Matrix &M, Vector &b)
Definition: ISTLSolver.hpp:348
std::vector< detail::FlexibleSolverInfo< Matrix, Vector, CommunicationType > > flexibleSolver_
Definition: ISTLSolver.hpp:666
void eraseMatrix() override
Signals that the memory for the matrix internally in the solver could be erased.
Definition: ISTLSolver.hpp:327
Dune::AssembledLinearOperator< Matrix, Vector, Vector > AbstractOperatorType
Definition: ISTLSolver.hpp:165
std::function< Vector()> getWeightsCalculator(const PropertyTree &prm, const Matrix &matrix, std::size_t pressIndex) const
Definition: ISTLSolver.hpp:579
Dune::AssembledLinearOperator< Matrix, Vector, Vector > AssembledLinearOperatorType
Definition: ISTLSolver.hpp:183
ISTLSolver(const Simulator &simulator)
Definition: ISTLSolver.hpp:211
int iterations_
Definition: ISTLSolver.hpp:657
int domainIndex_
Definition: ISTLSolver.hpp:670
std::any parallelInformation_
Definition: ISTLSolver.hpp:659
GetPropType< TypeTag, Properties::Simulator > Simulator
Definition: ISTLSolver.hpp:160
Vector * rhs_
Definition: ISTLSolver.hpp:663
int activeSolverNum_
Definition: ISTLSolver.hpp:665
std::vector< int > overlapRows_
Definition: ISTLSolver.hpp:667
std::unique_ptr< ElementChunksType > element_chunks_
Definition: ISTLSolver.hpp:679
GetPropType< TypeTag, Properties::Indices > Indices
Definition: ISTLSolver.hpp:158
const Simulator & simulator_
Definition: ISTLSolver.hpp:656
std::vector< int > interiorRows_
Definition: ISTLSolver.hpp:668
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: ISTLSolver.hpp:473
bool forceSerial_
Definition: ISTLSolver.hpp:675
bool solve(Vector &x) override
Definition: ISTLSolver.hpp:418
static void registerParameters()
Definition: ISTLSolver.hpp:185
std::vector< PropertyTree > prm_
Definition: ISTLSolver.hpp:676
GetPropType< TypeTag, Properties::ElementContext > ElementContext
Definition: ISTLSolver.hpp:163
void setResidual(Vector &) override
Definition: ISTLSolver.hpp:395
const CommunicationType * comm() const override
Get the communication object used by the solver.
Definition: ISTLSolver.hpp:459
bool isParallel() const
Definition: ISTLSolver.hpp:481
A sparse matrix interface backend for BCRSMatrix from dune-istl.
Definition: istlsparsematrixadapter.hh:43
Definition: matrixblock.hh:229
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
T get(const std::string &key) const
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hpp:66
Definition: WellOperators.hpp:70
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
void getTrueImpesWeights(int pressureVarIndex, Vector &weights, const ElementContext &elemCtx, const Model &model, const ElementChunksType &element_chunks, bool enable_thread_parallel)
Definition: getQuasiImpesWeights.hpp:161
void getTrueImpesWeightsAnalytic(int, Vector &weights, const ElementContext &elemCtx, const Model &model, const ElementChunksType &element_chunks, bool enable_thread_parallel)
Definition: getQuasiImpesWeights.hpp:228
void writeSystem(const SimulatorType &simulator, const MatrixType &matrix, const VectorType &rhs, const std::string &sysName, const Communicator *comm)
Definition: WriteSystemMatrixHelper.hpp:197
Definition: blackoilmodel.hh:80
std::unique_ptr< Matrix > blockJacobiAdjacency(const Grid &grid, const std::vector< int > &cell_part, std::size_t nonzeroes, const std::vector< std::set< int > > &wellConnectionsGraph)
void copyParValues(std::any &parallelInformation, std::size_t size, Dune::OwnerOverlapCopyCommunication< int, int > &comm)
Copy values in parallel.
std::size_t numMatrixRowsToUseInSolver(const Grid &grid, bool ownerFirst)
If ownerFirst=true, returns the number of interior cells in grid, else just numCells().
Definition: findOverlapRowsAndColumns.hpp:122
void makeOverlapRowsInvalid(Matrix &matrix, const std::vector< int > &overlapRows)
void findOverlapAndInterior(const Grid &grid, const Mapper &mapper, std::vector< int > &overlapRows, std::vector< int > &interiorRows)
Find the rows corresponding to overlap cells.
Definition: findOverlapRowsAndColumns.hpp:92
void printLinearSolverParameters(const FlowLinearSolverParameters &parameters, const VectorOrSingle &prm, const Comm &comm)
Print the linear solver parameters to the log if requested.
Definition: printlinearsolverparameter.hpp:61
Definition: blackoilbioeffectsmodules.hh:45
Dune::InverseOperatorResult InverseOperatorResult
Definition: GpuBridge.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:233
void extractParallelGridInformationToISTL(const Dune::CpGrid &grid, std::any &anyComm)
Extracts the information about the data decomposition from the grid for dune-istl.
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
PropertyTree setupPropertyTree(FlowLinearSolverParameters p, bool linearSolverMaxIterSet, bool linearSolverReductionSet, bool tpsaSetup=false)
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
This class carries all parameters for the NewtonIterationBlackoilInterleaved class.
Definition: FlowLinearSolverParameters.hpp:98
void init(bool cprRequestedInDataFile)
std::string linsolver_
Definition: FlowLinearSolverParameters.hpp:113
typename Linear::IstlSparseMatrixAdapter< Block > type
Definition: ISTLSolver.hpp:89
The class that allows to manipulate sparse matrices.
Definition: linalgproperties.hh:50
Definition: ISTLSolver.hpp:70
std::tuple< FlowIstlSolverParams > InheritsFrom
Definition: ISTLSolver.hpp:71
Definition: FlowBaseProblemProperties.hpp:98
Definition: ISTLSolver.hpp:103
std::unique_ptr< AbstractSolverType > solver_
Definition: ISTLSolver.hpp:116
std::size_t interiorCellNum_
Definition: ISTLSolver.hpp:120
Dune::InverseOperator< Vector, Vector > AbstractSolverType
Definition: ISTLSolver.hpp:104
AbstractPreconditionerType * pre_
Definition: ISTLSolver.hpp:119
Dune::AssembledLinearOperator< Matrix, Vector, Vector > AbstractOperatorType
Definition: ISTLSolver.hpp:105
void create(const Matrix &matrix, bool parallel, const PropertyTree &prm, std::size_t pressureIndex, std::function< Vector()> weightCalculator, const bool forceSerial, Comm *comm)
std::unique_ptr< LinearOperatorExtra< Vector, Vector > > wellOperator_
Definition: ISTLSolver.hpp:118
std::unique_ptr< AbstractOperatorType > op_
Definition: ISTLSolver.hpp:117