Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l > Class Template Reference

Converts the field type (eg. double to float) to benchmark single precision preconditioners. More...

#include <PreconditionerConvertFieldTypeAdapter.hpp>

Inheritance diagram for Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >:
Inheritance graph

Public Types

using matrix_type = typename std::remove_const< M >::type
 The matrix type the preconditioner is for. More...
 
using domain_type = X
 The domain type of the preconditioner. More...
 
using range_type = Y
 The range type of the preconditioner. More...
 
using field_type = typename X::field_type
 The field type of the preconditioner. More...
 
using domain_type_to = typename CudaPreconditionerType::domain_type
 
using range_type_to = typename CudaPreconditionerType::range_type
 The range type of the preconditioner. More...
 
using field_type_to = typename domain_type_to::field_type
 The field type of the preconditioner. More...
 
using block_type = typename domain_type::block_type
 
using XTo = Dune::BlockVector< Dune::FieldVector< field_type_to, block_type::dimension > >
 
using YTo = Dune::BlockVector< Dune::FieldVector< field_type_to, block_type::dimension > >
 
using matrix_type_to = typename Dune::BCRSMatrix< Dune::FieldMatrix< field_type_to, block_type::dimension, block_type::dimension > >
 

Public Member Functions

 PreconditionerConvertFieldTypeAdapter (const M &matrix)
 Constructor. More...
 
virtual void pre (X &x, Y &b) override
 Not used at the moment. More...
 
virtual void apply (X &v, const Y &d) override
 Apply the preconditoner. More...
 
virtual void post (X &x) override
 Not used at the moment. More...
 
virtual Dune::SolverCategory::Category category () const override
 Category of the preconditioner (see SolverCategory::Category) More...
 
virtual void update () override
 
const matrix_type_togetConvertedMatrix () const
 
void setUnderlyingPreconditioner (const std::shared_ptr< CudaPreconditionerType > &conditioner)
 

Detailed Description

template<class CudaPreconditionerType, class M, class X, class Y, int l = 1>
class Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >

Converts the field type (eg. double to float) to benchmark single precision preconditioners.

Note
This is not a fast conversion, it is simply meant to benchmark the potential of some preconditioners on consumer grade GPUs where the double precision performance is often artificially limited.
In theory this can do any field_type conversion that is meaningful, but it is only tested on double to float conversion.
Remember to set the underlying preconditioner with setUnderlyingPreconditioner (should use the matrix from getConvertedMatrix())
One could simply change the constructor design by accepting a creator function for the underlying preconditioner. For the current use cases this is however not needed.

To use this, use something like the following code:

using XDouble = Dune::BlockVector<Dune::FieldVector<double, 2>>;
using MDouble = Dune::FieldMatrix<double, 2, 2>;
using SpMatrixDouble = Dune::BCRSMatrix<MDouble>;
using XFloat = Dune::BlockVector<Dune::FieldVector<float, 2>>;
using MFloat = Dune::FieldMatrix<float, 2, 2>;
using SpMatrixFloat = Dune::BCRSMatrix<MFloat>;
template<class ParallelInfo>
void applyILU0AsFloat(const MDouble& matrix, const XDouble& x, XDouble& y) {
using DoubleToFloatConverter = typename Opm::cuistl::PreconditionerConvertFieldTypeAdapter<FloatILU0, MDouble,
XDouble, XDouble>;
// Note that we do not need to make a new instance for every invocation, this
// is just done for this example
auto doubleToFloatConverter = DoubleToFloatConverter(matrix);
const auto& convertedMatrix = doubleToFloatConverter.getConvertedMatrix();
auto floatILU0 = std::make_shared<FloatILU0>(convertedMatrix, 0, 1.0, 0);
doubleToFloatConverter.setUnderlyingPreconditioner(floatILU0);
// This will convert x and y to float, then call floatILU0.apply on the converted arguments
doubleToFloatConverter.apply(x, y);
}
A two-step version of an overlapping Schwarz preconditioner using one step ILU0 as.
Definition: ParallelOverlappingILU0.hpp:131
Converts the field type (eg. double to float) to benchmark single precision preconditioners.
Definition: PreconditionerConvertFieldTypeAdapter.hpp:86

Member Typedef Documentation

◆ block_type

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::block_type = typename domain_type::block_type

◆ domain_type

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::domain_type = X

The domain type of the preconditioner.

◆ domain_type_to

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::domain_type_to = typename CudaPreconditionerType::domain_type

◆ field_type

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::field_type = typename X::field_type

The field type of the preconditioner.

◆ field_type_to

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::field_type_to = typename domain_type_to::field_type

The field type of the preconditioner.

◆ matrix_type

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::matrix_type = typename std::remove_const<M>::type

The matrix type the preconditioner is for.

◆ matrix_type_to

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<field_type_to, block_type::dimension, block_type::dimension> >

◆ range_type

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::range_type = Y

The range type of the preconditioner.

◆ range_type_to

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::range_type_to = typename CudaPreconditionerType::range_type

The range type of the preconditioner.

◆ XTo

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::XTo = Dune::BlockVector<Dune::FieldVector<field_type_to, block_type::dimension> >

◆ YTo

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
using Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::YTo = Dune::BlockVector<Dune::FieldVector<field_type_to, block_type::dimension> >

Constructor & Destructor Documentation

◆ PreconditionerConvertFieldTypeAdapter()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::PreconditionerConvertFieldTypeAdapter ( const M &  matrix)
inlineexplicit

Constructor.

Parameters
AThe matrix to operate on.
Note
After the PreconditionerConvertFieldTypeAdapter you can get the converted matrix by calling getConvertedMatrix(), which in turn can be used to create the underlying preconditioner. Once the underlying preconditioner has been called, this must be supplied to setUnderlyingPreconditioner.

Member Function Documentation

◆ apply()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
virtual void Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::apply ( X &  v,
const Y &  d 
)
inlineoverridevirtual

Apply the preconditoner.

◆ category()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
virtual Dune::SolverCategory::Category Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::category ( ) const
inlineoverridevirtual

Category of the preconditioner (see SolverCategory::Category)

◆ getConvertedMatrix()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
const matrix_type_to & Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::getConvertedMatrix ( ) const
inline

◆ post()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
virtual void Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::post ( X &  x)
inlineoverridevirtual

Not used at the moment.

◆ pre()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
virtual void Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::pre ( X &  x,
Y &  b 
)
inlineoverridevirtual

Not used at the moment.

◆ setUnderlyingPreconditioner()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
void Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::setUnderlyingPreconditioner ( const std::shared_ptr< CudaPreconditionerType > &  conditioner)
inline

◆ update()

template<class CudaPreconditionerType , class M , class X , class Y , int l = 1>
virtual void Opm::cuistl::PreconditionerConvertFieldTypeAdapter< CudaPreconditionerType, M, X, Y, l >::update ( )
inlineoverridevirtual

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