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