Gyselalib++
 
Loading...
Searching...
No Matches
MatrixBatchTridiag< ExecSpace > Class Template Reference

A structure for solving a set of independant tridiagonal systems using a direct method. More...

Inheritance diagram for MatrixBatchTridiag< ExecSpace >:
MatrixBatch< ExecSpace >

Public Member Functions

 MatrixBatchTridiag (const int batch_size, const int mat_size, DKokkosView2D const aa, DKokkosView2D const bb, DKokkosView2D const cc)
 Creates an instance of the MatrixBatchTridiag class.
 
bool check_stability () const
 Check if the matrices are in the stability area of the solver.
 
void setup_solver () final
 Perform a pre-process operation on the solver.
 
void solve (BatchedRHS const b) const final
 Solve the batched linear problem Ax=b.
 
- Public Member Functions inherited from MatrixBatch< ExecSpace >
virtual ~MatrixBatch ()=default
 Destruct.
 
std::size_t size () const
 Get the size of the square matrix corresponding to a single batch in one of its dimensions.
 
std::size_t batch_size () const
 Get the batch size of the linear problem.
 

Additional Inherited Members

- Public Types inherited from MatrixBatch< ExecSpace >
using BatchedRHS = Kokkos::View< double **, Kokkos::LayoutRight, ExecSpace >
 The type of a Kokkos::View storing batched right-hand sides. Second dimenion is batch dimension.
 
- Protected Member Functions inherited from MatrixBatch< ExecSpace >
 MatrixBatch (const std::size_t batch_size, const std::size_t mat_size)
 The constructor for MatrixBatch class.
 

Detailed Description

template<class ExecSpace>
class MatrixBatchTridiag< ExecSpace >

A structure for solving a set of independant tridiagonal systems using a direct method.

The parallelism operates on the whole collection by dispatching to threads. Each problem is treated sequentially, by the tridiagonal matrix algorithm (TDMA). This solver is stable for tridiagonal matrices which satisfy one of the following conditions:

  • Diagonally Dominant.
  • Symmetric positive-definite. Diagonally Dominant property is fully checked. Only symmetry property is checked, positivity-definiteness is not.
    Template Parameters
    ExecSpaceThe execution space related to Kokkos.

Constructor & Destructor Documentation

◆ MatrixBatchTridiag()

template<class ExecSpace >
MatrixBatchTridiag< ExecSpace >::MatrixBatchTridiag ( const int  batch_size,
const int  mat_size,
DKokkosView2D const  aa,
DKokkosView2D const  bb,
DKokkosView2D const  cc 
)
inlineexplicit

Creates an instance of the MatrixBatchTridiag class.

First dimension is the batch, second one refers to matrix entries indexed by line. The entries aa,bb,cc are 2D Kokkos views and have the same dimensions. LayoutRight: means that the "last" dimension is the contiguous one. aa(batch_idx,0) and cc(batch_idx,mat_size-1) are not used for any values of batch_idx.

Parameters
[in]batch_sizeThe size of the set of linear problems.
[in]mat_sizeThe common size of each individual matrix .
[in]aa2d Kokkos View which stores subdiagonal components for all matrices.
[in]bb2d Kokkos View which stores diagonal components for all matrices.
[in]cc2d Kokkos View which stores upper diagonal components for all matrices.

Member Function Documentation

◆ check_stability()

template<class ExecSpace >
bool MatrixBatchTridiag< ExecSpace >::check_stability ( ) const
inline

Check if the matrices are in the stability area of the solver.

It checks if each matrix of the batch has one of the following structures:

  • Diagonally Dominant.
  • Symmetric. If assertion fails, there is at least one of the matrices which does not verify these conditions. Positivity-definiteness is not checked, it is needed to fully achieve TDMA algorithm requirement: symmetric positive definite.
Returns
A boolean which indicates if the stability condition may be verified.

◆ setup_solver()

template<class ExecSpace >
void MatrixBatchTridiag< ExecSpace >::setup_solver ( )
inlinefinalvirtual

Perform a pre-process operation on the solver.

Must be called after filling the matrix.

It calls check_stability function to verify if the matrices data is in range of validity of the solver.

The stopping criterion is a reduction factor ||Ax-b||/||b||<tol with max_iter maximum iterations.

Implements MatrixBatch< ExecSpace >.

◆ solve()

template<class ExecSpace >
void MatrixBatchTridiag< ExecSpace >::solve ( BatchedRHS const  b) const
inlinefinalvirtual

Solve the batched linear problem Ax=b.

Parameters
[in,out]bA 2D Kokkos::View storing the batched right-hand sides of the problem and receiving the corresponding solutions.

Implements MatrixBatch< ExecSpace >.


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