Template Numerical Library version\ main:8b8c8226
Loading...
Searching...
No Matches
Namespaces | Classes | Functions
TNL::Matrices Namespace Reference

Namespace for different matrix types. More...

Namespaces

namespace  Sandbox
 Namespace for sandbox matrices.
 

Classes

class  DenseMatrix
 Implementation of dense matrix, i.e. matrix storing explicitly all of its elements including zeros. More...
 
class  DenseMatrixElement
 Accessor for dense matrix elements. More...
 
class  DenseMatrixRowView
 RowView is a simple structure for accessing rows of dense matrix. More...
 
class  DenseMatrixView
 Implementation of dense matrix view. More...
 
class  DistributedMatrix
 
struct  GeneralMatrix
 General non-symmetric matrix type. More...
 
class  GinkgoOperator
 Wraps a general TNL matrix as a Ginkgo LinOp. More...
 
class  HypreCSRMatrix
 Wrapper for Hypre's sequential CSR matrix. More...
 
class  HypreParCSRMatrix
 Wrapper for Hypre's sequential CSR matrix. More...
 
class  LambdaMatrix
 "Matrix-free matrix" based on lambda functions. More...
 
class  LambdaMatrixElement
 Accessor for elements of lambda matrix. More...
 
struct  LambdaMatrixFactory
 Helper class for creating instances of LambdaMatrix. More...
 
class  LambdaMatrixRowView
 RowView is a simple structure for accessing rows of Lambda matrix. More...
 
class  LambdaMatrixRowViewIterator
 
class  Matrix
 Base class for other matrix types. More...
 
struct  MatrixInfo
 
class  MatrixOperations
 
class  MatrixOperations< Devices::Cuda >
 
class  MatrixReader
 Helper class for importing of matrices from different input formats. More...
 
class  MatrixRowViewIterator
 
class  MatrixSetter
 
class  MatrixSetterTraverserUserData
 
struct  MatrixType
 Structure for specifying type of sparse matrix. More...
 
class  MatrixView
 Base class for other matrix types views. More...
 
class  MatrixWriter
 Helper class for exporting of matrices to different output formats. More...
 
class  MultidiagonalMatrix
 Implementation of sparse multidiagonal matrix. More...
 
class  MultidiagonalMatrixElement
 Accessor for multidiagonal matrix elements. More...
 
class  MultidiagonalMatrixRowView
 RowView is a simple structure for accessing rows of multidiagonal matrix. More...
 
class  MultidiagonalMatrixView
 Implementation of sparse multidiagonal matrix. More...
 
class  SparseMatrix
 Implementation of sparse matrix, i.e. matrix storing only non-zero elements. More...
 
class  SparseMatrixElement
 Accessor for sparse matrix elements. More...
 
class  SparseMatrixRowView
 RowView is a simple structure for accessing rows of sparse matrix. More...
 
class  SparseMatrixView
 Implementation of sparse matrix view. More...
 
class  StaticMatrix
 
struct  SymmetricMatrix
 Symmetric matrix type. More...
 
class  TridiagonalMatrix
 Implementation of sparse tridiagonal matrix. More...
 
class  TridiagonalMatrixRowView
 RowView is a simple structure for accessing rows of tridiagonal matrix. More...
 
class  TridiagonalMatrixView
 Implementation of sparse tridiagonal matrix. More...
 

Functions

template<typename Matrix , typename InVector , typename OutVector >
__global__ void ColumnMajorDenseMatrixViewVectorMultiplicationKernel (const Matrix matrix, const InVector inVector, OutVector outVector, const int begin, const int end, int gridIdx)
 
template<typename Matrix , typename AdjacencyMatrix >
void copyAdjacencyStructure (const Matrix &A, AdjacencyMatrix &B, bool has_symmetric_pattern=false, bool ignore_diagonal=true)
 
template<typename Matrix1 , typename Matrix2 >
void copySparseMatrix (Matrix1 &A, const Matrix2 &B)
 
template<typename Matrix1 , typename Matrix2 >
std::enable_if< std::is_same< typenameMatrix1::DeviceType, typenameMatrix2::DeviceType >::value >::type copySparseMatrix_impl (Matrix1 &A, const Matrix2 &B)
 
template<typename Matrix1 , typename Matrix2 >
std::enable_if<!std::is_same< typenameMatrix1::DeviceType, typenameMatrix2::DeviceType >::value &&std::is_same< typenameMatrix2::DeviceType, Devices::Host >::value >::type copySparseMatrix_impl (Matrix1 &A, const Matrix2 &B)
 
template<typename Matrix1 , typename Matrix2 >
std::enable_if<!std::is_same< typenameMatrix1::DeviceType, typenameMatrix2::DeviceType >::value &&std::is_same< typenameMatrix2::DeviceType, Devices::Cuda >::value >::type copySparseMatrix_impl (Matrix1 &A, const Matrix2 &B)
 
template<int tileDim, int tileRowBlockSize, typename ResultMatrix , typename Matrix1 , typename Matrix2 >
__global__ void DenseMatrixProductKernel (ResultMatrix resultMatrix, const Matrix1 matrixA, const Matrix2 matrixB, const typename ResultMatrix::RealType matrixMultiplicator, const typename ResultMatrix::IndexType gridIdx_x, const typename ResultMatrix::IndexType gridIdx_y)
 
template<int tileDim, int tileRowBlockSize, typename OutputMatrix , typename InputMatrix , typename Real , typename Index >
__global__ void DenseTranspositionAlignedKernel (OutputMatrix resultMatrix, const InputMatrix inputMatrix, const Real matrixMultiplicator, const Index gridIdx_x, const Index gridIdx_y)
 
template<int tileDim, int tileRowBlockSize, typename OutputMatrix , typename InputMatrix , typename Real , typename Index >
__global__ void DenseTranspositionNonAlignedKernel (OutputMatrix resultMatrix, const InputMatrix inputMatrix, const Real matrixMultiplicator, const Index gridIdx_x, const Index gridIdx_y)
 
template<typename RealType , typename IndexType >
__global__ void GeamCudaKernel (const IndexType m, const IndexType n, const RealType alpha, const RealType *A, const IndexType lda, const RealType beta, const RealType *B, const IndexType ldb, RealType *C, const IndexType ldc)
 
template<typename RealType , typename IndexType >
__global__ void GemvCudaKernel (const IndexType m, const IndexType n, const RealType alpha, const RealType *A, const IndexType lda, const RealType *x, const RealType beta, RealType *y)
 
template<typename Matrix >
auto getGinkgoMatrixCsr (std::shared_ptr< const gko::Executor > exec, Matrix &matrix) -> std::unique_ptr< gko::matrix::Csr< typename Matrix::RealType, typename Matrix::IndexType > >
 Converts any TNL sparse matrix to a Ginkgo Csr matrix. More...
 
template<typename Matrix >
auto getGinkgoMatrixCsrView (std::shared_ptr< const gko::Executor > exec, Matrix &matrix) -> std::unique_ptr< gko::matrix::Csr< typename Matrix::RealType, typename Matrix::IndexType > >
 Creates a Ginkgo Csr matrix view from a TNL CSR matrix. More...
 
template<typename InMatrixView , typename OutMatrixView , typename Real , typename Index >
__global__ void MultidiagonalMatrixTranspositionCudaKernel (const InMatrixView inMatrix, OutMatrixView outMatrix, const Real matrixMultiplicator, const Index gridIdx)
 
template<typename Real , typename Device , typename Index , typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator >
bool operator!= (const DenseMatrixView< Real, Device, Index, Organization > &leftMatrix, const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator > &rightMatrix)
 Comparison operator with another dense matrix view. More...
 
template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
std::ostreamoperator<< (std::ostream &str, const DenseMatrix< Real, Device, Index, Organization, RealAllocator > &matrix)
 Insertion operator for dense matrix and output stream. More...
 
template<typename MatrixElementsLambda , typename CompressedRowLengthsLambda , typename Real , typename Device , typename Index >
std::ostreamoperator<< (std::ostream &str, const LambdaMatrix< MatrixElementsLambda, CompressedRowLengthsLambda, Real, Device, Index > &matrix)
 Insertion operator for lambda matrix and output stream. More...
 
template<typename MatrixElementsLambda , typename CompressedRowLengthsLambda , typename Real , typename Index >
std::ostreamoperator<< (std::ostream &str, const LambdaMatrixRowView< MatrixElementsLambda, CompressedRowLengthsLambda, Real, Index > &row)
 Insertion operator for a Lambda matrix row. More...
 
template<typename Real , typename Device , typename Index >
std::ostreamoperator<< (std::ostream &str, const Matrix< Real, Device, Index > &matrix)
 Overloaded insertion operator for printing a matrix to output stream. More...
 
template<typename Real , typename Device , typename Index >
std::ostreamoperator<< (std::ostream &str, const MatrixView< Real, Device, Index > &matrix)
 Overloaded insertion operator for printing a matrix to output stream. More...
 
template<typename SegmentView , typename ValuesView , typename ColumnsIndexesView >
std::ostreamoperator<< (std::ostream &str, const SparseMatrixRowView< SegmentView, ValuesView, ColumnsIndexesView > &row)
 Insertion operator for a sparse matrix row. More...
 
template<typename Real , typename Device , typename Index , typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator >
bool operator== (const DenseMatrixView< Real, Device, Index, Organization > &leftMatrix, const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator > &rightMatrix)
 Comparison operator with another dense matrix view. More...
 
template<typename Matrix , typename PermutationArray >
void permuteMatrixColumns (Matrix &matrix, const PermutationArray &iperm)
 
template<typename Matrix , typename PermutationArray >
void permuteMatrixRows (Matrix &matrix, const PermutationArray &perm)
 
template<typename Array1 , typename Array2 , typename PermutationArray >
void reorderArray (const Array1 &src, Array2 &dest, const PermutationArray &perm)
 
template<typename Matrix1 , typename Matrix2 , typename PermutationArray >
void reorderSparseMatrix (const Matrix1 &matrix1, Matrix2 &matrix2, const PermutationArray &perm, const PermutationArray &iperm)
 
template<typename Matrix , typename InVector , typename OutVector >
__global__ void RowMajorDenseMatrixViewVectorMultiplicationKernel (const Matrix matrix, const InVector inVector, OutVector outVector, const int first, const int last, int gridIdx)
 
template<typename Matrix1 , typename Matrix2 >
__global__ void SparseMatrixCopyKernel (Matrix1 *A, const Matrix2 *B, const typename Matrix2::IndexType *rowLengths, typename Matrix2::IndexType rows)
 
template<typename Vector , typename Matrix >
__global__ void SparseMatrixSetRowLengthsVectorKernel (Vector *rowLengths, const Matrix *matrix, typename Matrix::IndexType rows, typename Matrix::IndexType cols)
 
template<typename InMatrixView , typename OutMatrixView , typename Real , typename Index >
__global__ void TridiagonalMatrixTranspositionCudaKernel (const InMatrixView inMatrix, OutMatrixView outMatrix, Real matrixMultiplicator, Index gridIdx)
 
template<int BlockSize, int ThreadsPerRow, typename Matrix , typename InVector , typename OutVector >
__global__ void VectorColumnMajorDenseMatrixViewVectorMultiplicationKernel (const Matrix matrix, const InVector inVector, OutVector outVector, const int begin, const int end, int gridIdx)
 
template<typename Device , typename Real , typename Index >
SparseMatrixView< Real, Device, Index, GeneralMatrix, Algorithms::Segments::CSRViewDefaultwrapCSRMatrix (const Index &rows, const Index &columns, Index *rowPointers, Real *values, Index *columnIndexes)
 Function for wrapping of arrays defining CSR format into a sparse matrix view. More...
 
template<typename Device , typename Real , typename Index , ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization()>
DenseMatrixView< Real, Device, Index, Organization > wrapDenseMatrix (const Index &rows, const Index &columns, Real *values)
 Function for wrapping an array of values into a dense matrix view. More...
 
template<typename Device , ElementsOrganization Organization, typename Real , typename Index , int Alignment = 1>
auto wrapEllpackMatrix (const Index rows, const Index columns, const Index nonzerosPerRow, Real *values, Index *columnIndexes) -> decltype(EllpackMatrixWrapper< Device, Organization, Real, Index, Alignment >::wrap(rows, columns, nonzerosPerRow, values, columnIndexes))
 Function for wrapping of arrays defining Ellpack format into a sparse matrix view. More...
 

Detailed Description

Namespace for different matrix types.

Namespace for different matrix formats.

Namespace for matrix formats.

Function Documentation

◆ operator!=()

template<typename Real , typename Device , typename Index , typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator >
bool TNL::Matrices::operator!= ( const DenseMatrixView< Real, Device, Index, Organization > &  leftMatrix,
const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator > &  rightMatrix 
)

Comparison operator with another dense matrix view.

Parameters
leftMatrixis the left-hand side matrix view.
rightMatrixis the right-hand side matrix.
Returns
false if the both matrices are is equal, true otherwise.

◆ operator<<() [1/6]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const DenseMatrix< Real, Device, Index, Organization, RealAllocator > &  matrix 
)

Insertion operator for dense matrix and output stream.

Parameters
stris the output stream.
matrixis the dense matrix.
Returns
reference to the stream.

◆ operator<<() [2/6]

template<typename MatrixElementsLambda , typename CompressedRowLengthsLambda , typename Real , typename Device , typename Index >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const LambdaMatrix< MatrixElementsLambda, CompressedRowLengthsLambda, Real, Device, Index > &  matrix 
)

Insertion operator for lambda matrix and output stream.

Parameters
stris the output stream.
matrixis the lambda matrix.
Returns
reference to the stream.

◆ operator<<() [3/6]

template<typename MatrixElementsLambda , typename CompressedRowLengthsLambda , typename Real , typename Index >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const LambdaMatrixRowView< MatrixElementsLambda, CompressedRowLengthsLambda, Real, Index > &  row 
)

Insertion operator for a Lambda matrix row.

Parameters
stris an output stream.
rowis an input Lambda matrix row.
Returns
reference to the output stream.

◆ operator<<() [4/6]

template<typename Real , typename Device , typename Index >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const Matrix< Real, Device, Index > &  matrix 
)

Overloaded insertion operator for printing a matrix to output stream.

Template Parameters
Realis a type of the matrix elements.
Deviceis a device where the matrix is allocated.
Indexis a type used for the indexing of the matrix elements.
Parameters
stris a output stream.
matrixis the matrix to be printed.
Returns
a reference to the output stream std::ostream.

◆ operator<<() [5/6]

template<typename Real , typename Device , typename Index >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const MatrixView< Real, Device, Index > &  matrix 
)

Overloaded insertion operator for printing a matrix to output stream.

Template Parameters
Realis a type of the matrix elements.
Deviceis a device where the matrix is allocated.
Indexis a type used for the indexing of the matrix elements.
Parameters
stris a output stream.
matrixis the matrix to be printed.
Returns
a reference to the output stream std::ostream.

◆ operator<<() [6/6]

template<typename SegmentView , typename ValuesView , typename ColumnsIndexesView >
std::ostream & TNL::Matrices::operator<< ( std::ostream str,
const SparseMatrixRowView< SegmentView, ValuesView, ColumnsIndexesView > &  row 
)

Insertion operator for a sparse matrix row.

Parameters
stris an output stream.
rowis an input sparse matrix row.
Returns
reference to the output stream.

◆ operator==()

template<typename Real , typename Device , typename Index , typename Real_ , typename Device_ , typename Index_ , ElementsOrganization Organization, typename RealAllocator >
bool TNL::Matrices::operator== ( const DenseMatrixView< Real, Device, Index, Organization > &  leftMatrix,
const DenseMatrix< Real_, Device_, Index_, Organization, RealAllocator > &  rightMatrix 
)

Comparison operator with another dense matrix view.

Parameters
leftMatrixis the left-hand side matrix view.
rightMatrixis the right-hand side matrix.
Returns
true if the both matrices are is equal, false otherwise.

◆ VectorColumnMajorDenseMatrixViewVectorMultiplicationKernel()

template<int BlockSize, int ThreadsPerRow, typename Matrix , typename InVector , typename OutVector >
__global__ void TNL::Matrices::VectorColumnMajorDenseMatrixViewVectorMultiplicationKernel ( const Matrix  matrix,
const InVector  inVector,
OutVector  outVector,
const int  begin,
const int  end,
int  gridIdx 
)

The following kernel is an attempt to map more CUDA threads to one matrix row.

◆ wrapCSRMatrix()

template<typename Device , typename Real , typename Index >
SparseMatrixView< Real, Device, Index, GeneralMatrix, Algorithms::Segments::CSRViewDefault > TNL::Matrices::wrapCSRMatrix ( const Index &  rows,
const Index &  columns,
Index *  rowPointers,
Real values,
Index *  columnIndexes 
)

Function for wrapping of arrays defining CSR format into a sparse matrix view.

Template Parameters
Deviceis a device on which the arrays are allocated.
Realis a type of matrix elements values.
Indexis a type for matrix elements indexing.
Parameters
rowsis a number of matrix rows.
columnsis a number of matrix columns.
rowPointersis an array holding row pointers of the CSR format ( ROW_INDEX here)
valuesis an array with values of matrix elements ( V here)
columnIndexesis an array with column indexes of matrix elements ( COL_INDEX here)
Returns
instance of SparseMatrixView with CSR format.

The size of array rowPointers must be equal to number of rows + 1. The last element of the array equals to the number of all nonzero matrix elements. The sizes of arrays values and columnIndexes must be equal to this number.

Example
#include <iostream>
#include <TNL/Matrices/DenseMatrix.h>
#include <TNL/Matrices/MatrixWrapping.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void wrapMatrixView()
{
/***
* Encode the following matrix to CSR format...
*
* / 1 2 0 0 \.
* | 0 6 0 0 |
* | 9 0 0 0 |
* \ 0 0 15 16 /
*/
const int rows( 4 ), columns( 4 );
TNL::Containers::Vector< double, Device > valuesVector { 1, 2, 6, 9, 15, 16 };
TNL::Containers::Vector< int, Device > columnIndexesVector { 0, 1, 1, 0, 2, 3 };
TNL::Containers::Vector< int, Device > rowPointersVector { 0, 2, 3, 4, 6 };
double* values = valuesVector.getData();
int* columnIndexes = columnIndexesVector.getData();
int* rowPointers = rowPointersVector.getData();
/***
* Wrap the arrays `rowPointers, `values` and `columnIndexes` to sparse matrix view
*/
auto matrix = TNL::Matrices::wrapCSRMatrix< Device >( rows, columns, rowPointers, values, columnIndexes );
std::cout << "Matrix reads as: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Wraping matrix view on host: " << std::endl;
wrapMatrixView< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Wraping matrix view on CUDA device: " << std::endl;
wrapMatrixView< TNL::Devices::Cuda >();
#endif
}
__cuda_callable__ const Value * getData() const
Returns a const-qualified raw pointer to the data.
Definition: Array.hpp:324
Vector extends Array with algebraic operations.
Definition: Vector.h:40
T endl(T... args)
Output
Wraping matrix view on host:
Matrix reads as:
Row: 0 -> 0:1 1:2
Row: 1 -> 1:6
Row: 2 -> 0:9
Row: 3 -> 2:15 3:16
Wraping matrix view on CUDA device:
Matrix reads as:
Row: 0 -> 0:1 1:2
Row: 1 -> 1:6
Row: 2 -> 0:9
Row: 3 -> 2:15 3:16

◆ wrapDenseMatrix()

template<typename Device , typename Real , typename Index , ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization()>
DenseMatrixView< Real, Device, Index, Organization > TNL::Matrices::wrapDenseMatrix ( const Index &  rows,
const Index &  columns,
Real values 
)

Function for wrapping an array of values into a dense matrix view.

Template Parameters
Deviceis a device on which the array is allocated.
Realis a type of array elements.
Indexis a type for indexing of matrix elements.
Organizationis matrix elements organization - see TNL::Algorithms::Segments::ElementsOrganization.
Parameters
rowsis a number of matrix rows.
columnsis a number of matrix columns.
valuesis the array with matrix elements values.
Returns
instance of DenseMatrixView wrapping the array.

The array size must be equal to product of rows and columns. The dense matrix view does not deallocate the input array at the end of its lifespan.

Example
#include <iostream>
#include <TNL/Matrices/DenseMatrix.h>
#include <TNL/Matrices/MatrixWrapping.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void wrapMatrixView()
{
const int rows( 3 ), columns( 4 );
1, 2, 3, 4,
5, 6, 7, 8,
9, 10, 11, 12 };
double* values = valuesVector.getData();
/***
* Wrap the array `values` to dense matrix view
*/
auto matrix = TNL::Matrices::wrapDenseMatrix< Device >( rows, columns, values );
std::cout << "Matrix reads as: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Wraping matrix view on host: " << std::endl;
wrapMatrixView< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Wraping matrix view on CUDA device: " << std::endl;
wrapMatrixView< TNL::Devices::Cuda >();
#endif
}
Output
Wraping matrix view on host:
Matrix reads as:
Row: 0 -> 0:1 1:2 2:3 3:4
Row: 1 -> 0:5 1:6 2:7 3:8
Row: 2 -> 0:9 1:10 2:11 3:12
Wraping matrix view on CUDA device:
Matrix reads as:
Row: 0 -> 0:1 1:4 2:7 3:10
Row: 1 -> 0:2 1:5 2:8 3:11
Row: 2 -> 0:3 1:6 2:9 3:12

◆ wrapEllpackMatrix()

template<typename Device , ElementsOrganization Organization, typename Real , typename Index , int Alignment = 1>
auto TNL::Matrices::wrapEllpackMatrix ( const Index  rows,
const Index  columns,
const Index  nonzerosPerRow,
Real values,
Index *  columnIndexes 
) -> decltype( EllpackMatrixWrapper< Device, Organization, Real, Index, Alignment >::wrap( rows, columns, nonzerosPerRow, values, columnIndexes ) )

Function for wrapping of arrays defining Ellpack format into a sparse matrix view.

This is to prevent from appearing in Doxygen documentation.

Template Parameters
Deviceis a device on which the arrays are allocated.
Realis a type of matrix elements values.
Indexis a type for matrix elements indexing.
Alignmentdefines alignment of data. The number of matrix rows is rounded to a multiple of this number. It it usefull mainly for GPUs.
Parameters
rowsis a number of matrix rows.
columnsis a number of matrix columns.
nonzerosPerRowis number of nonzero matrix elements in each row.
valuesis an array with values of matrix elements.
columnIndexesis an array with column indexes of matrix elements.
Returns
instance of SparseMatrixView with CSR format.

The sizes of arrays values and columnIndexes must be equal to rows * nonzerosPerRow. Use -1 as a column index for padding zeros.

Example
#include <iostream>
#include <TNL/Matrices/DenseMatrix.h>
#include <TNL/Matrices/MatrixWrapping.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void wrapMatrixView()
{
/***
* Encode the following matrix to Ellpack format...
*
* / 1 2 0 0 \.
* | 0 6 0 0 |
* | 9 0 0 0 |
* \ 0 0 15 16 /
*/
const int rows( 4 ), columns( 4 );
TNL::Containers::Vector< double, Device > valuesVector { 1, 2, 6, 0, 9, 0, 15, 16 };
TNL::Containers::Vector< int, Device > columnIndexesVector { 0, 1, 1, -1, 0, -1, 2, 3 };
double* values = valuesVector.getData();
int* columnIndexes = columnIndexesVector.getData();
/***
* Wrap the arrays `values` and `columnIndexes` to sparse matrix view
*/
auto matrix = TNL::Matrices::wrapEllpackMatrix< Device, TNL::Algorithms::Segments::RowMajorOrder >( rows, columns, 2, values, columnIndexes );
std::cout << "Matrix reads as: " << std::endl << matrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Wraping matrix view on host: " << std::endl;
wrapMatrixView< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Wraping matrix view on CUDA device: " << std::endl;
wrapMatrixView< TNL::Devices::Cuda >();
#endif
}
Output
Wraping matrix view on host:
Matrix reads as:
Row: 0 -> 0:1 1:2
Row: 1 -> 1:6
Row: 2 -> 0:9
Row: 3 -> 2:15 3:16
Wraping matrix view on CUDA device:
Matrix reads as:
Row: 0 -> 0:1 1:2
Row: 1 -> 1:6
Row: 2 -> 0:9
Row: 3 -> 2:15 3:16