Template Numerical Library version\ main:1437bf49
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator > Class Template Reference

Implementation of dense matrix, i.e. matrix storing explicitly all of its elements including zeros. More...

#include <TNL/Matrices/DenseMatrix.h>

Inheritance diagram for TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >:
Inheritance graph
[legend]
Collaboration diagram for TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >:
Collaboration graph
[legend]

Public Types

using ConstViewType = typename DenseMatrixView< Real, Device, Index, Organization >::ConstViewType
 Matrix view type for constant instances.
 
using RealAllocatorType = RealAllocator
 The allocator for matrix elements.
 
template<typename _Real = Real, typename _Device = Device, typename _Index = Index, ElementsOrganization _Organization = Algorithms::Segments::DefaultElementsOrganization< _Device >::getOrganization(), typename _RealAllocator = typename Allocators::Default< _Device >::template Allocator< _Real >>
using Self = DenseMatrix< _Real, _Device, _Index, _Organization, _RealAllocator >
 Helper type for getting self type or its modifications.
 
using ValuesVectorType = Containers::Vector< Real, Device, Index, RealAllocator >
 Type of vector holding values of matrix elements.
 
using ViewType = DenseMatrixView< Real, Device, Index, Organization >
 Type of related matrix view.
 
- Public Types inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
using ConstRowView = typename RowView::ConstRowView
 Type for accessing immutable matrix row.
 
using DeviceType = Device
 The device where the matrix is allocated.
 
using IndexType = Index
 The type used for matrix elements indexing.
 
using RealType = Real
 The type of matrix elements.
 
using RowView = DenseMatrixRowView< SegmentViewType, typename Base::ValuesViewType >
 Type for accessing matrix row.
 
- Public Types inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
using ConstValuesViewType
 Type of constant vector view holding values of matrix elements.
 
using DeviceType
 The device where the matrix is allocated.
 
using IndexType
 The type used for matrix elements indexing.
 
using RealType
 The type of matrix elements.
 
using RowCapacitiesType
 
using ValuesViewType
 Type of vector view holding values of matrix elements.
 

Public Member Functions

 DenseMatrix (const DenseMatrix &matrix)
 Copy constructor.
 
 DenseMatrix (const RealAllocatorType &allocator=RealAllocatorType())
 Constructor only with values allocator.
 
 DenseMatrix (DenseMatrix &&matrix) noexcept=default
 Move constructor.
 
 DenseMatrix (Index rows, Index columns, const RealAllocatorType &allocator=RealAllocatorType())
 Constructor with matrix dimensions.
 
template<typename Value >
 DenseMatrix (std::initializer_list< std::initializer_list< Value > > data, const RealAllocatorType &allocator=RealAllocatorType())
 Constructor with 2D initializer list.
 
ConstViewType getConstView () const
 Returns a non-modifiable view of the dense matrix.
 
template<typename Matrix1 , typename Matrix2 , int tileDim = 32>
void getMatrixProduct (const Matrix1 &matrix1, const Matrix2 &matrix2, Real matrixMultiplicator=1.0)
 
template<typename Matrix , int tileDim = 32>
void getTransposition (const Matrix &matrix, Real matrixMultiplicator=1.0)
 
ViewType getView ()
 Returns a modifiable view of the dense matrix.
 
void load (const String &fileName)
 Method for loading the matrix from the file with given filename.
 
void load (File &file) override
 Method for loading the matrix from a file.
 
DenseMatrixoperator= (const DenseMatrix &matrix)
 Copy-assignment of exactly the same matrix type.
 
template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const DenseMatrix< Real, Device, Index, Organization, RealAllocator > &matrix)
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , typename RHSRealAllocator >
DenseMatrixoperator= (const DenseMatrix< RHSReal, RHSDevice, RHSIndex, Organization, RHSRealAllocator > &matrix)
 Assignment operator with the same organization.
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , typename RHSRealAllocator >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const DenseMatrix< RHSReal, RHSDevice, RHSIndex, Organization, RHSRealAllocator > &matrix)
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization, typename RHSRealAllocator >
DenseMatrixoperator= (const DenseMatrix< RHSReal, RHSDevice, RHSIndex, RHSOrganization, RHSRealAllocator > &matrix)
 Assignment operator with other dense matrices.
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization, typename RHSRealAllocator >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const DenseMatrix< RHSReal, RHSDevice, RHSIndex, RHSOrganization, RHSRealAllocator > &matrix)
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex >
DenseMatrixoperator= (const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, Organization > &matrix)
 Assignment operator with matrix view having the same elements organization.
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, Organization > &matrix)
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization>
DenseMatrixoperator= (const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, RHSOrganization > &matrix)
 Assignment operator with other dense matrices.
 
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization>
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, RHSOrganization > &matrix)
 
template<typename RHSMatrix >
DenseMatrixoperator= (const RHSMatrix &matrix)
 Assignment operator with other (sparse) types of matrices.
 
template<typename RHSMatrix >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (const RHSMatrix &matrix)
 
DenseMatrixoperator= (DenseMatrix &&matrix) noexcept(false)
 Move-assignment of exactly the same matrix type.
 
template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
DenseMatrix< Real, Device, Index, Organization, RealAllocator > & operator= (DenseMatrix< Real, Device, Index, Organization, RealAllocator > &&matrix) noexcept(false)
 
void reset ()
 Resets the matrix to zero dimensions.
 
void save (const String &fileName) const
 Method for saving the matrix to the file with given filename.
 
void save (File &file) const override
 Method for saving the matrix to a file.
 
void setDimensions (Index rows, Index columns)
 Set number of rows and columns of this matrix.
 
template<typename MapIndex , typename MapValue >
void setElements (const std::map< std::pair< MapIndex, MapIndex >, MapValue > &map)
 This method sets the dense matrix elements from std::map.
 
template<typename Value >
void setElements (std::initializer_list< std::initializer_list< Value > > data)
 This method recreates the dense matrix from 2D initializer list.
 
template<typename Matrix >
void setLike (const Matrix &matrix)
 Set the number of matrix rows and columns by the given matrix.
 
template<typename Matrix_ >
void setLike (const Matrix_ &matrix)
 
template<typename RowCapacitiesVector >
void setRowCapacities (const RowCapacitiesVector &rowCapacities)
 This method is only for the compatibility with the sparse matrices.
 
- Public Member Functions inherited from TNL::Object
virtual ~Object ()=default
 Destructor.
 
virtual std::string getSerializationTypeVirtual () const
 
void load (const String &fileName)
 Method for restoring the object from a file.
 
void save (const String &fileName) const
 Method for saving the object to a file as a binary data.
 
- Public Member Functions inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
__cuda_callable__ DenseMatrixBase ()=default
 Constructor without parameters.
 
__cuda_callable__ DenseMatrixBase (const DenseMatrixBase &matrix)=default
 Copy constructor.
 
__cuda_callable__ DenseMatrixBase (DenseMatrixBase &&matrix) noexcept=default
 Move constructor.
 
__cuda_callable__ DenseMatrixBase (IndexType rows, IndexType columns, typename Base::ValuesViewType values)
 Constructor with matrix dimensions and values.
 
__cuda_callable__ void addElement (IndexType row, IndexType column, const RealType &value, const RealType &thisElementMultiplicator=1.0)
 Add element at given row and column to given value.
 
template<typename Matrix >
void addMatrix (const Matrix &matrix, const RealType &matrixMultiplicator=1.0, const RealType &thisMatrixMultiplicator=1.0)
 
template<typename Function >
void forAllElements (Function &&function)
 This method calls forElements for all matrix rows.
 
template<typename Function >
void forAllElements (Function &&function) const
 This method calls forElements for all matrix rows.
 
template<typename Function >
void forAllRows (Function &&function)
 Method for parallel iteration over all matrix rows.
 
template<typename Function >
void forAllRows (Function &&function) const
 Method for parallel iteration over all matrix rows for constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &&function)
 Method for iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void forElements (IndexType begin, IndexType end, Function &&function) const
 Method for iteration over all matrix rows for constant instances.
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function)
 Method for parallel iteration over matrix rows from interval [begin, end).
 
template<typename Function >
void forRows (IndexType begin, IndexType end, Function &&function) const
 Method for parallel iteration over matrix rows from interval [begin, end) for constant instances.
 
template<typename Vector >
void getCompressedRowLengths (Vector &rowLengths) const
 Computes number of non-zeros in each row.
 
__cuda_callable__ Real getElement (IndexType row, IndexType column) const
 Returns value of matrix element at position given by its row and column index.
 
__cuda_callable__ RowView getRow (IndexType rowIdx)
 Non-constant getter of simple structure for accessing given matrix row.
 
__cuda_callable__ ConstRowView getRow (IndexType rowIdx) const
 Constant getter of simple structure for accessing given matrix row.
 
template<typename Vector >
void getRowCapacities (Vector &rowCapacities) const
 Compute capacities of all rows.
 
template<typename Real_ , typename Device_ , typename Index_ >
bool operator!= (const DenseMatrixBase< Real_, Device_, Index_, Organization > &matrix) const
 Comparison operator with another dense matrix view.
 
__cuda_callable__ Realoperator() (IndexType row, IndexType column)
 Returns non-constant reference to element at row row and column column.
 
__cuda_callable__ const Realoperator() (IndexType row, IndexType column) const
 Returns constant reference to element at row row and column column.
 
DenseMatrixBaseoperator= (const DenseMatrixBase &)=delete
 Copy-assignment operator.
 
DenseMatrixBaseoperator= (DenseMatrixBase &&)=delete
 Move-assignment operator.
 
template<typename Real_ , typename Device_ , typename Index_ >
bool operator== (const DenseMatrixBase< Real_, Device_, Index_, Organization > &matrix) const
 Comparison operator with another dense matrix view.
 
void print (std::ostream &str) const
 Method for printing the matrix to output stream.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceAllRows (Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on ALL matrix rows for constant instances.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchReal >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchReal &identity) const
 Method for performing general reduction on matrix rows for constant instances.
 
template<typename Fetch , typename Reduce , typename Keep , typename FetchValue >
void reduceRows (IndexType begin, IndexType end, Fetch &fetch, const Reduce &reduce, Keep &keep, const FetchValue &identity) const
 
template<typename Function >
void sequentialForAllRows (Function &&function)
 This method calls sequentialForRows for all matrix rows.
 
template<typename Function >
void sequentialForAllRows (Function &&function) const
 This method calls sequentialForRows for all matrix rows (for constant instances).
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &&function)
 Method for sequential iteration over all matrix rows for non-constant instances.
 
template<typename Function >
void sequentialForRows (IndexType begin, IndexType end, Function &&function) const
 Method for sequential iteration over all matrix rows for constant instances.
 
__cuda_callable__ void setElement (IndexType row, IndexType column, const RealType &value)
 Sets element at given row and column to given value.
 
void setValue (const RealType &v)
 Sets all matrix elements to value v.
 
template<typename InVector , typename OutVector >
void vectorProduct (const InVector &inVector, OutVector &outVector, const RealType &matrixMultiplicator=1.0, const RealType &outVectorMultiplicator=0.0, IndexType begin=0, IndexType end=0) const
 Computes product of matrix and vector.
 
- Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__cuda_callable__ MatrixBase ()=default
 Basic constructor with no parameters.
 
__cuda_callable__ MatrixBase (const MatrixBase &view)=default
 Shallow copy constructor.
 
__cuda_callable__ MatrixBase (IndexType rows, IndexType columns, ValuesViewType values)
 Constructor with matrix dimensions and matrix elements values.
 
__cuda_callable__ MatrixBase (MatrixBase &&view) noexcept=default
 Move constructor.
 
IndexType getAllocatedElementsCount () const
 Tells the number of allocated matrix elements.
 
__cuda_callable__ IndexType getColumns () const
 Returns number of matrix columns.
 
virtual IndexType getNonzeroElementsCount () const
 Computes a current number of nonzero matrix elements.
 
__cuda_callable__ IndexType getRows () const
 Returns number of matrix rows.
 
__cuda_callable__ ValuesViewTypegetValues ()
 Returns a reference to a vector with the matrix elements values.
 
__cuda_callable__ const ValuesViewTypegetValues () const
 Returns a constant reference to a vector with the matrix elements values.
 
bool operator!= (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator!= (const MatrixT &matrix) const
 
__cuda_callable__ MatrixBaseoperator= (const MatrixBase &)=delete
 Copy-assignment operator.
 
__cuda_callable__ MatrixBaseoperator= (MatrixBase &&)=delete
 Move-assignment operator.
 
bool operator== (const Matrix &matrix) const
 Comparison operator with another arbitrary matrix view type.
 
bool operator== (const MatrixT &matrix) const
 

Static Public Member Functions

static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Object
static std::string getSerializationType ()
 Static serialization type getter.
 
- Static Public Member Functions inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
static std::string getSerializationType ()
 Returns string with serialization type.
 
- Static Public Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
static constexpr ElementsOrganization getOrganization ()
 Matrix elements organization getter.
 
static constexpr bool isBinary ()
 Test of binary matrix type.
 
static constexpr bool isMatrix ()
 Test of matrix type.
 
static constexpr bool isSymmetric ()
 Test of symmetric matrix type.
 

Protected Attributes

Base::SegmentsType segments
 Instance of the segments used for indexing in the dense matrix.
 
ValuesVectorType values
 Vector containing the allocated matrix elements.
 
- Protected Attributes inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
SegmentsViewType segments
 
- Protected Attributes inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
IndexType columns
 
IndexType rows
 
ValuesViewType values
 

Additional Inherited Members

- Protected Types inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
using Base = MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
 
using SegmentsReductionKernel = Algorithms::SegmentsReductionKernels::EllpackKernel< Index, Device >
 
using SegmentsType
 
using SegmentsViewType = typename SegmentsType::ViewType
 
using SegmentViewType = typename SegmentsType::SegmentViewType
 
- Protected Member Functions inherited from TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >
__cuda_callable__ void bind (IndexType rows, IndexType columns, typename Base::ValuesViewType values, SegmentsViewType segments)
 Re-initializes the internal attributes of the base class.
 
__cuda_callable__ IndexType getElementIndex (IndexType row, IndexType column) const
 
- Protected Member Functions inherited from TNL::Matrices::MatrixBase< Real, Device, Index, GeneralMatrix, Organization >
__cuda_callable__ void bind (IndexType rows, IndexType columns, ValuesViewType values)
 Re-initializes the internal attributes of the base class.
 

Detailed Description

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
class TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >

Implementation of dense matrix, i.e. matrix storing explicitly all of its elements including zeros.

Template Parameters
Realis a type of matrix elements.
Deviceis a device where the matrix is allocated.
Indexis a type for indexing of the matrix elements.
Organizationtells the ordering of matrix elements. It is either RowMajorOrder or ColumnMajorOrder.
RealAllocatoris allocator for the matrix elements.

Member Typedef Documentation

◆ ConstViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
using TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::ConstViewType = typename DenseMatrixView< Real, Device, Index, Organization >::ConstViewType

Matrix view type for constant instances.

See DenseMatrixView.

◆ ViewType

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
using TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::ViewType = DenseMatrixView< Real, Device, Index, Organization >

Type of related matrix view.

See DenseMatrixView.

Constructor & Destructor Documentation

◆ DenseMatrix() [1/5]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::DenseMatrix ( const RealAllocatorType & allocator = RealAllocatorType())

Constructor only with values allocator.

Parameters
allocatoris used for allocation of matrix elements values.

◆ DenseMatrix() [2/5]

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

Copy constructor.

Parameters
matrixis the source matrix

◆ DenseMatrix() [3/5]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::DenseMatrix ( DenseMatrix< Real, Device, Index, Organization, RealAllocator > && matrix)
defaultnoexcept

Move constructor.

Parameters
matrixis the source matrix

◆ DenseMatrix() [4/5]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::DenseMatrix ( Index rows,
Index columns,
const RealAllocatorType & allocator = RealAllocatorType() )

Constructor with matrix dimensions.

Parameters
rowsis number of matrix rows.
columnsis number of matrix columns.
allocatoris used for allocation of matrix elements values.

◆ DenseMatrix() [5/5]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
template<typename Value >
TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::DenseMatrix ( std::initializer_list< std::initializer_list< Value > > data,
const RealAllocatorType & allocator = RealAllocatorType() )

Constructor with 2D initializer list.

The number of matrix rows is set to the outer list size and the number of matrix columns is set to maximum size of inner lists. Missing elements are filled in with zeros.

Parameters
datais a initializer list of initializer lists representing list of matrix rows.
allocatoris used for allocation of matrix elements values.
Example
#include <iostream>
#include <TNL/Matrices/DenseMatrix.h>
#include <TNL/Devices/Host.h>
template< typename Device >
void initializerListExample()
{
{ 1, 2, 3, 4, 5, 6 },
{ 7, 8, 9, 10, 11, 12 },
{ 13, 14, 15, 16, 17, 18 }
};
std::cout << "General dense matrix: " << std::endl << matrix << std::endl;
{ 1 },
{ 2, 3 },
{ 4, 5, 6 },
{ 7, 8, 9, 10 },
{ 11, 12, 13, 14, 15 }
};
std::cout << "Triangular dense matrix: " << std::endl << triangularMatrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Creating matrices on CPU ... " << std::endl;
initializerListExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Creating matrices on CUDA GPU ... " << std::endl;
initializerListExample< TNL::Devices::Cuda >();
#endif
}
Implementation of dense matrix, i.e. matrix storing explicitly all of its elements including zeros.
Definition DenseMatrix.h:31
T endl(T... args)
Output
Creating matrices on CPU ...
General dense matrix:
Row: 0 -> 0:1 1:2 2:3 3:4 4:5 5:6
Row: 1 -> 0:7 1:8 2:9 3:10 4:11 5:12
Row: 2 -> 0:13 1:14 2:15 3:16 4:17 5:18
Triangular dense matrix:
Row: 0 -> 0:1 1:0 2:0 3:0 4:0
Row: 1 -> 0:2 1:3 2:0 3:0 4:0
Row: 2 -> 0:4 1:5 2:6 3:0 4:0
Row: 3 -> 0:7 1:8 2:9 3:10 4:0
Row: 4 -> 0:11 1:12 2:13 3:14 4:15
Creating matrices on CUDA GPU ...
General dense matrix:
Row: 0 -> 0:1 1:2 2:3 3:4 4:5 5:6
Row: 1 -> 0:7 1:8 2:9 3:10 4:11 5:12
Row: 2 -> 0:13 1:14 2:15 3:16 4:17 5:18
Triangular dense matrix:
Row: 0 -> 0:1 1:0 2:0 3:0 4:0
Row: 1 -> 0:2 1:3 2:0 3:0 4:0
Row: 2 -> 0:4 1:5 2:6 3:0 4:0
Row: 3 -> 0:7 1:8 2:9 3:10 4:0
Row: 4 -> 0:11 1:12 2:13 3:14 4:15

Member Function Documentation

◆ getConstView()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
auto TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::getConstView ( ) const

Returns a non-modifiable view of the dense matrix.

See DenseMatrixView.

Returns
dense matrix view.

◆ getSerializationType()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
std::string TNL::Matrices::DenseMatrixBase< Real, Device, Index, Organization >::getSerializationType ( )
static

Returns string with serialization type.

The string has a form `MatricesDenseMatrix< RealType, [any_device], IndexType, [any_allocator], true/false >`.

Returns
String with the serialization type.

◆ getView()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
auto TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::getView ( )

Returns a modifiable view of the dense matrix.

See DenseMatrixView.

Returns
dense matrix view.

◆ load() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::load ( const String & fileName)

Method for loading the matrix from the file with given filename.

Parameters
fileNameis name of the file.

◆ load() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::load ( File & file)
overridevirtual

Method for loading the matrix from a file.

Parameters
fileis the file from which the matrix will be loaded.

Reimplemented from TNL::Object.

◆ operator=() [1/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const DenseMatrix< Real, Device, Index, Organization, RealAllocator > & matrix)

Copy-assignment of exactly the same matrix type.

Parameters
matrixis input matrix for the assignment.
Returns
reference to this matrix.

◆ operator=() [2/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename RHSReal , typename RHSDevice , typename RHSIndex , typename RHSRealAllocator >
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const DenseMatrix< RHSReal, RHSDevice, RHSIndex, Organization, RHSRealAllocator > & matrix)

Assignment operator with the same organization.

Parameters
matrixis the right-hand side matrix.
Returns
reference to this matrix.

◆ operator=() [3/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization, typename RHSRealAllocator >
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const DenseMatrix< RHSReal, RHSDevice, RHSIndex, RHSOrganization, RHSRealAllocator > & matrix)

Assignment operator with other dense matrices.

Parameters
matrixis the right-hand side matrix.
Returns
reference to this matrix.

◆ operator=() [4/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename RHSReal , typename RHSDevice , typename RHSIndex >
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, Organization > & matrix)

Assignment operator with matrix view having the same elements organization.

Parameters
matrixis the right-hand side matrix.
Returns
reference to this matrix.

◆ operator=() [5/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename RHSReal , typename RHSDevice , typename RHSIndex , ElementsOrganization RHSOrganization>
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const DenseMatrixView< RHSReal, RHSDevice, RHSIndex, RHSOrganization > & matrix)

Assignment operator with other dense matrices.

Parameters
matrixis the right-hand side matrix.
Returns
reference to this matrix.

◆ operator=() [6/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename RHSMatrix >
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( const RHSMatrix & matrix)

Assignment operator with other (sparse) types of matrices.

Parameters
matrixis the right-hand side matrix.
Returns
reference to this matrix.

◆ operator=() [7/7]

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
DenseMatrix & TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::operator= ( DenseMatrix< Real, Device, Index, Organization, RealAllocator > && matrix)

Move-assignment of exactly the same matrix type.

Parameters
matrixis input matrix for the assignment.
Returns
reference to this matrix.

◆ save() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::save ( const String & fileName) const

Method for saving the matrix to the file with given filename.

Parameters
fileNameis name of the file.

◆ save() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::save ( File & file) const
overridevirtual

Method for saving the matrix to a file.

Parameters
fileis the file where the matrix will be saved.

Reimplemented from TNL::Object.

◆ setDimensions()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::setDimensions ( Index rows,
Index columns )

Set number of rows and columns of this matrix.

Parameters
rowsis the number of matrix rows.
columnsis the number of matrix columns.

◆ setElements() [1/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
template<typename MapIndex , typename MapValue >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::setElements ( const std::map< std::pair< MapIndex, MapIndex >, MapValue > & map)

This method sets the dense matrix elements from std::map.

This is intended for compatibility with SparseMatrix, the method is used e.g. in MatrixReader.

The matrix elements values are given as a map data where keys are std::pair of matrix coordinates ( {row, column} ) and value is the matrix element value.

Template Parameters
MapIndexis a type for indexing rows and columns.
MapValueis a type for matrix elements values in the map.
Parameters
mapis std::map containing matrix elements.

◆ setElements() [2/2]

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
template<typename Value >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::setElements ( std::initializer_list< std::initializer_list< Value > > data)

This method recreates the dense matrix from 2D initializer list.

The number of matrix rows is set to the outer list size and the number of matrix columns is set to maximum size of inner lists. Missing elements are filled in with zeros.

Parameters
datais a initializer list of initializer lists representing list of matrix rows.
Example
#include <iostream>
#include <TNL/Matrices/DenseMatrix.h>
#include <TNL/Devices/Host.h>
#include <TNL/Devices/Cuda.h>
template< typename Device >
void setElementsExample()
{
matrix.setElements( {
{ 1, 2, 3, 4, 5, 6 },
{ 7, 8, 9, 10, 11, 12 },
{ 13, 14, 15, 16, 17, 18 }
} );
std::cout << matrix << std::endl;
triangularMatrix.setElements( {
{ 1 },
{ 2, 3 },
{ 4, 5, 6 },
{ 7, 8, 9, 10 },
{ 11, 12, 13, 14, 15 }
} );
std::cout << triangularMatrix << std::endl;
}
int main( int argc, char* argv[] )
{
std::cout << "Setting matrix elements on host: " << std::endl;
setElementsExample< TNL::Devices::Host >();
#ifdef __CUDACC__
std::cout << "Setting matrix elements on CUDA device: " << std::endl;
setElementsExample< TNL::Devices::Cuda >();
#endif
}
void setElements(std::initializer_list< std::initializer_list< Value > > data)
This method recreates the dense matrix from 2D initializer list.
Definition DenseMatrix.hpp:48
Output
Setting matrix elements on host:
Row: 0 -> 0:1 1:2 2:3 3:4 4:5 5:6
Row: 1 -> 0:7 1:8 2:9 3:10 4:11 5:12
Row: 2 -> 0:13 1:14 2:15 3:16 4:17 5:18
Row: 0 -> 0:1 1:0 2:0 3:0 4:0
Row: 1 -> 0:2 1:3 2:0 3:0 4:0
Row: 2 -> 0:4 1:5 2:6 3:0 4:0
Row: 3 -> 0:7 1:8 2:9 3:10 4:0
Row: 4 -> 0:11 1:12 2:13 3:14 4:15
Setting matrix elements on CUDA device:
Row: 0 -> 0:1 1:2 2:3 3:4 4:5 5:6
Row: 1 -> 0:7 1:8 2:9 3:10 4:11 5:12
Row: 2 -> 0:13 1:14 2:15 3:16 4:17 5:18
Row: 0 -> 0:1 1:0 2:0 3:0 4:0
Row: 1 -> 0:2 1:3 2:0 3:0 4:0
Row: 2 -> 0:4 1:5 2:6 3:0 4:0
Row: 3 -> 0:7 1:8 2:9 3:10 4:0
Row: 4 -> 0:11 1:12 2:13 3:14 4:15

◆ setLike()

template<typename Real = double, typename Device = Devices::Host, typename Index = int, ElementsOrganization Organization = Algorithms::Segments::DefaultElementsOrganization< Device >::getOrganization(), typename RealAllocator = typename Allocators::Default< Device >::template Allocator< Real >>
template<typename Matrix >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::setLike ( const Matrix & matrix)

Set the number of matrix rows and columns by the given matrix.

Template Parameters
Matrixis matrix type. This can be any matrix having methods getRows and getColumns.
Parameters
matrixin the input matrix dimensions of which are to be adopted.

◆ setRowCapacities()

template<typename Real , typename Device , typename Index , ElementsOrganization Organization, typename RealAllocator >
template<typename RowCapacitiesVector >
void TNL::Matrices::DenseMatrix< Real, Device, Index, Organization, RealAllocator >::setRowCapacities ( const RowCapacitiesVector & rowCapacities)

This method is only for the compatibility with the sparse matrices.

This method does nothing. In debug mode it contains assertions checking that given rowCapacities are compatible with the current matrix dimensions.


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