Provides all unmodified linear solvers from dune-istl. More...

#include <parallelistlbackend.hh>

Inheritance diagram for Opm::Linear::ParallelIstlSolverBackend< TypeTag >:
Inheritance graph

Public Member Functions

 ParallelIstlSolverBackend (const Simulator &simulator)
 
void eraseMatrix ()
 Causes the solve() method to discared the structure of the linear system of equations the next time it is called. More...
 
void prepare (const SparseMatrixAdapter &M, const Vector &)
 Set up the internal data structures required for the linear solver. More...
 
void setResidual (const Vector &b)
 Assign values to the internal data structure for the residual vector. More...
 
void getResidual (Vector &b) const
 Retrieve the synchronized internal residual vector. More...
 
void setMatrix (const SparseMatrixAdapter &M)
 Sets the values of the residual's Jacobian matrix. More...
 
bool solve (Vector &x)
 Actually solve the linear system of equations. More...
 
size_t iterations () const
 Return number of iterations used during last solve. More...
 

Static Public Member Functions

static void registerParameters ()
 Register all run-time parameters for the linear solver. More...
 

Protected Types

enum  { dimWorld = GridView::dimensionworld }
 
using Implementation = GetPropType< TypeTag, Properties::LinearSolverBackend >
 
using LinearSolverScalar = GetPropType< TypeTag, Properties::LinearSolverScalar >
 
using Vector = GetPropType< TypeTag, Properties::GlobalEqVector >
 
using BorderListCreator = GetPropType< TypeTag, Properties::BorderListCreator >
 
using GridView = GetPropType< TypeTag, Properties::GridView >
 
using Overlap = GetPropType< TypeTag, Properties::Overlap >
 
using OverlappingMatrix = GetPropType< TypeTag, Properties::OverlappingMatrix >
 
using PreconditionerWrapper = GetPropType< TypeTag, Properties::PreconditionerWrapper >
 
using SequentialPreconditioner = typename PreconditionerWrapper::SequentialPreconditioner
 

Protected Member Functions

std::shared_ptr< RawLinearSolver > prepareSolver_ (ParallelOperator &parOperator, ParallelScalarProduct &parScalarProduct, ParallelPreconditioner &parPreCond)
 
void cleanupSolver_ ()
 
std::pair< bool, int > runSolver_ (std::shared_ptr< RawLinearSolver > solver)
 
ImplementationasImp_ ()
 
const ImplementationasImp_ () const
 
void cleanup_ ()
 
std::shared_ptr< ParallelPreconditionerpreparePreconditioner_ ()
 
void cleanupPreconditioner_ ()
 
void writeOverlapToVTK_ ()
 

Protected Attributes

friend ParentType
 
LinearSolverWrapper solverWrapper_
 
const Simulatorsimulator_
 
int gridSequenceNumber_
 
size_t lastIterations_
 
OverlappingMatrixoverlappingMatrix_
 
OverlappingVectoroverlappingb_
 
OverlappingVectoroverlappingx_
 
PreconditionerWrapper precWrapper_
 

Detailed Description

template<class TypeTag>
class Opm::Linear::ParallelIstlSolverBackend< TypeTag >

Provides all unmodified linear solvers from dune-istl.

To set the linear solver, use

template<class TypeTag>
struct LinearSolverWrapper<TypeTag, TTag::YourTypeTag>
{ using type = Opm::Linear::SolverWrapper$SOLVER<TypeTag>; };

The possible choices for '$SOLVER' are:

  • Richardson: A fixpoint solver using the Richardson iteration
  • SteepestDescent: The steepest descent solver
  • ConjugatedGradients: A conjugated gradients solver
  • BiCGStab: A stabilized bi-conjugated gradients solver
  • MinRes: A solver based on the minimized residual algorithm
  • RestartedGMRes: A restarted GMRES solver

Chosing the preconditioner works in an analogous way:

template<class TypeTag>
struct PreconditionerWrapper<TypeTag, TTag::YourTypeTag>
{ using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER<TypeTag>; };
GetPropType< TypeTag, Properties::PreconditionerWrapper > PreconditionerWrapper
Definition: parallelbasebackend.hh:125

Where the choices possible for '$PRECONDITIONER' are:

  • Jacobi: A Jacobi preconditioner
  • GaussSeidel: A Gauss-Seidel preconditioner
  • SSOR: A symmetric successive overrelaxation (SSOR) preconditioner
  • SOR: A successive overrelaxation (SOR) preconditioner
  • ILUn: An ILU(n) preconditioner
  • ILU0: A specialized (and optimized) ILU(0) preconditioner

Member Typedef Documentation

◆ BorderListCreator

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::BorderListCreator = GetPropType<TypeTag, Properties::BorderListCreator>
protectedinherited

◆ GridView

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::GridView = GetPropType<TypeTag, Properties::GridView>
protectedinherited

◆ Implementation

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::Implementation = GetPropType<TypeTag, Properties::LinearSolverBackend>
protectedinherited

◆ LinearSolverScalar

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::LinearSolverScalar = GetPropType<TypeTag, Properties::LinearSolverScalar>
protectedinherited

◆ Overlap

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::Overlap = GetPropType<TypeTag, Properties::Overlap>
protectedinherited

◆ OverlappingMatrix

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::OverlappingMatrix = GetPropType<TypeTag, Properties::OverlappingMatrix>
protectedinherited

◆ PreconditionerWrapper

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::PreconditionerWrapper = GetPropType<TypeTag, Properties::PreconditionerWrapper>
protectedinherited

◆ SequentialPreconditioner

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::SequentialPreconditioner = typename PreconditionerWrapper::SequentialPreconditioner
protectedinherited

◆ Vector

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::Vector = GetPropType<TypeTag, Properties::GlobalEqVector>
protectedinherited

Member Enumeration Documentation

◆ anonymous enum

template<class TypeTag >
anonymous enum
protectedinherited
Enumerator
dimWorld 

Constructor & Destructor Documentation

◆ ParallelIstlSolverBackend()

template<class TypeTag >
Opm::Linear::ParallelIstlSolverBackend< TypeTag >::ParallelIstlSolverBackend ( const Simulator simulator)
inline

Member Function Documentation

◆ asImp_() [1/2]

◆ asImp_() [2/2]

template<class TypeTag >
const Implementation & Opm::Linear::ParallelBaseBackend< TypeTag >::asImp_ ( ) const
inlineprotectedinherited

◆ cleanup_()

◆ cleanupPreconditioner_()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::cleanupPreconditioner_ ( )
inlineprotectedinherited

◆ cleanupSolver_()

template<class TypeTag >
void Opm::Linear::ParallelIstlSolverBackend< TypeTag >::cleanupSolver_ ( )
inlineprotected

◆ eraseMatrix()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::eraseMatrix ( )
inlineinherited

Causes the solve() method to discared the structure of the linear system of equations the next time it is called.

References Opm::Linear::ParallelBaseBackend< TypeTag >::cleanup_().

◆ getResidual()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::getResidual ( Vector b) const
inlineinherited

Retrieve the synchronized internal residual vector.

This only deals with entries which are local to the current process.

References Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingb_.

◆ iterations()

template<class TypeTag >
size_t Opm::Linear::ParallelBaseBackend< TypeTag >::iterations ( ) const
inlineinherited

Return number of iterations used during last solve.

References Opm::Linear::ParallelBaseBackend< TypeTag >::lastIterations_.

◆ prepare()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::prepare ( const SparseMatrixAdapter M,
const Vector  
)
inlineinherited

Set up the internal data structures required for the linear solver.

This only specified the topology of the linear system of equations; it does does not assign the values of the residual vector and its Jacobian matrix.

References Opm::Linear::ParallelBaseBackend< TypeTag >::asImp_(), Opm::Linear::ParallelBaseBackend< TypeTag >::gridSequenceNumber_, Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingb_, Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingMatrix_, Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingx_, and Opm::Linear::ParallelBaseBackend< TypeTag >::simulator_.

◆ preparePreconditioner_()

◆ prepareSolver_()

template<class TypeTag >
std::shared_ptr< RawLinearSolver > Opm::Linear::ParallelIstlSolverBackend< TypeTag >::prepareSolver_ ( ParallelOperator parOperator,
ParallelScalarProduct parScalarProduct,
ParallelPreconditioner parPreCond 
)
inlineprotected

◆ registerParameters()

template<class TypeTag >
static void Opm::Linear::ParallelIstlSolverBackend< TypeTag >::registerParameters ( )
inlinestatic

Register all run-time parameters for the linear solver.

References Opm::Linear::ParallelBaseBackend< TypeTag >::registerParameters().

◆ runSolver_()

template<class TypeTag >
std::pair< bool, int > Opm::Linear::ParallelIstlSolverBackend< TypeTag >::runSolver_ ( std::shared_ptr< RawLinearSolver >  solver)
inlineprotected

◆ setMatrix()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::setMatrix ( const SparseMatrixAdapter M)
inlineinherited

Sets the values of the residual's Jacobian matrix.

This method also synchronizes the data structure across the processes which are involved in the simulation run.

References Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingMatrix_.

◆ setResidual()

template<class TypeTag >
void Opm::Linear::ParallelBaseBackend< TypeTag >::setResidual ( const Vector b)
inlineinherited

Assign values to the internal data structure for the residual vector.

This method also cares about synchronizing that vector with the peer processes.

References Opm::Linear::ParallelBaseBackend< TypeTag >::overlappingb_.

◆ solve()

template<class TypeTag >
bool Opm::Linear::ParallelBaseBackend< TypeTag >::solve ( Vector x)
inlineinherited

◆ writeOverlapToVTK_()

Member Data Documentation

◆ gridSequenceNumber_

template<class TypeTag >
int Opm::Linear::ParallelBaseBackend< TypeTag >::gridSequenceNumber_
protectedinherited

◆ lastIterations_

template<class TypeTag >
size_t Opm::Linear::ParallelBaseBackend< TypeTag >::lastIterations_
protectedinherited

◆ overlappingb_

◆ overlappingMatrix_

◆ overlappingx_

◆ ParentType

template<class TypeTag >
friend Opm::Linear::ParallelIstlSolverBackend< TypeTag >::ParentType
protected

◆ precWrapper_

◆ simulator_

◆ solverWrapper_


The documentation for this class was generated from the following file: