My Project
Loading...
Searching...
No Matches
Opm::Linear::ParallelBaseBackend< TypeTag > Class Template Reference

Provides the common code which is required by most linear solvers. More...

#include <parallelbasebackend.hh>

Inheritance diagram for Opm::Linear::ParallelBaseBackend< TypeTag >:
Opm::Linear::ParallelAmgBackend< TypeTag > Opm::Linear::ParallelBiCGStabSolverBackend< TypeTag > Opm::Linear::ParallelIstlSolverBackend< TypeTag >

Public Member Functions

 ParallelBaseBackend (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.
 
void prepare (const SparseMatrixAdapter &M, const Vector &)
 Set up the internal data structures required for the linear solver.
 
void setResidual (const Vector &b)
 Assign values to the internal data structure for the residual vector.
 
void getResidual (Vector &b) const
 Retrieve the synchronized internal residual vector.
 
void setMatrix (const SparseMatrixAdapter &M)
 Sets the values of the residual's Jacobian matrix.
 
bool solve (Vector &x)
 Actually solve the linear system of equations.
 
size_t iterations () const
 Return number of iterations used during last solve.
 

Static Public Member Functions

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

Protected Types

enum  { dimWorld = GridView::dimensionworld }
 
using Implementation = GetPropType<TypeTag, Properties::LinearSolverBackend>
 
using Simulator = GetPropType<TypeTag, Properties::Simulator>
 
using Scalar = GetPropType<TypeTag, Properties::Scalar>
 
using LinearSolverScalar = GetPropType<TypeTag, Properties::LinearSolverScalar>
 
using SparseMatrixAdapter = GetPropType<TypeTag, Properties::SparseMatrixAdapter>
 
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 OverlappingVector = GetPropType<TypeTag, Properties::OverlappingVector>
 
using OverlappingMatrix = GetPropType<TypeTag, Properties::OverlappingMatrix>
 
using PreconditionerWrapper = GetPropType<TypeTag, Properties::PreconditionerWrapper>
 
using SequentialPreconditioner = typename PreconditionerWrapper::SequentialPreconditioner
 
using ParallelPreconditioner = Opm::Linear::OverlappingPreconditioner<SequentialPreconditioner, Overlap>
 
using ParallelScalarProduct = Opm::Linear::OverlappingScalarProduct<OverlappingVector, Overlap>
 
using ParallelOperator
 

Protected Member Functions

Implementation & asImp_ ()
 
const Implementation & asImp_ () const
 
void cleanup_ ()
 
std::shared_ptr< ParallelPreconditionerpreparePreconditioner_ ()
 
void cleanupPreconditioner_ ()
 
void writeOverlapToVTK_ ()
 

Protected Attributes

const Simulator & simulator_
 
int gridSequenceNumber_
 
size_t lastIterations_
 
OverlappingMatrix * overlappingMatrix_
 
OverlappingVector * overlappingb_
 
OverlappingVector * overlappingx_
 
PreconditionerWrapper precWrapper_
 

Detailed Description

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

Provides the common code which is required by most linear solvers.

This class provides access to all preconditioners offered by dune-istl using the PreconditionerWrapper property:

template<class TypeTag>
struct PreconditionerWrapper<TypeTag, TTag::YourTypeTag>
{ using type = Opm::Linear::PreconditionerWrapper$PRECONDITIONER<TypeTag>; };

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: An ILU(0) preconditioner. The results of this preconditioner are the same as setting the PreconditionerOrder property to 0 and using the ILU(n) preconditioner. The reason for the existence of ILU0 is that it is computationally cheaper because it does not need to consider things which are only required for higher orders

Member Typedef Documentation

◆ ParallelOperator

template<class TypeTag >
using Opm::Linear::ParallelBaseBackend< TypeTag >::ParallelOperator
protected
Initial value:
OverlappingVector,
OverlappingVector>
An overlap aware linear operator usable by ISTL.
Definition overlappingoperator.hh:42

Member Function Documentation

◆ getResidual()

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

Retrieve the synchronized internal residual vector.

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

◆ prepare()

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

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.

◆ setMatrix()

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

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.

◆ setResidual()

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

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

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

◆ solve()

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

Actually solve the linear system of equations.

Returns
true if the residual reduction could be achieved, else false.

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