Clone wiki

blaze / Matrix Operations

<<toc>>


Constructors

Matrices are just as easy and intuitive to create as vectors. Still, there are a few rules to be aware of:

  • In case the last template parameter (the storage order) is omitted, the matrix is per default stored in row-major order.
  • The elements of a StaticMatrix or HybridMatrix are default initialized (i.e. built-in data types are initialized to 0, class types are initialized via the default constructor).
  • Newly allocated elements of a DynamicMatrix or CompressedMatrix remain uninitialized if they are of built-in type and are default constructed if they are of class type.

Default Construction

using blaze::StaticMatrix;
using blaze::DynamicMatrix;
using blaze::CompressedMatrix;

// All matrices can be default constructed. Whereas the size of a StaticMatrix is fixed via the second and
// third template parameter, the initial size of a constructed DynamicMatrix or CompressedMatrix is 0.
StaticMatrix<int,2UL,2UL> M1;             // Instantiation of a 2x2 integer row-major
                                          // matrix. All elements are initialized to 0.
DynamicMatrix<float> M2;                  // Instantiation of a single precision dynamic
                                          // row-major matrix with 0 rows and 0 columns.
DynamicMatrix<double,columnMajor> M3;     // Instantiation of a double precision dynamic
                                          // column-major matrix with 0 rows and 0 columns.
CompressedMatrix<int> M4;                 // Instantiation of a compressed integer
                                          // row-major matrix of size 0x0.
CompressedMatrix<double,columnMajor> M5;  // Instantiation of a compressed double precision
                                          // column-major matrix of size 0x0.

Construction with Specific Size

The DynamicMatrix, HybridMatrix and CompressedMatrix classes offer a constructor that allows to immediately give the matrices a specific number of rows and columns:

DynamicMatrix<int> M6( 5UL, 4UL );                   // Instantiation of a 5x4 dynamic row-major
                                                     // matrix. The elements are not initialized.
HybridMatrix<double,5UL,9UL> M7( 3UL, 7UL );         // Instantiation of a 3x7 hybrid row-major
                                                     // matrix. The elements are not initialized.
CompressedMatrix<float,columnMajor> M8( 8UL, 6UL );  // Instantiation of an empty 8x6 compressed
                                                     // column-major matrix.

Note that dense matrices (in this case DynamicMatrix and HybridMatrix) immediately allocate enough capacity for all matrix elements. Sparse matrices on the other hand (in this example CompressedMatrix) merely acquire the size, but don't necessarily allocate memory.

Initialization Constructors

All dense matrix classes offer a constructor for a direct, homogeneous initialization of all matrix elements. In contrast, for sparse matrices the predicted number of non-zero elements can be specified.

StaticMatrix<int,4UL,3UL,columnMajor> M9( 7 );  // Instantiation of a 4x3 integer column-major
                                                // matrix. All elements are initialized to 7.
DynamicMatrix<float> M10( 2UL, 5UL, 2.0F );     // Instantiation of a 2x5 single precision row-major
                                                // matrix. All elements are initialized to 2.0F.
CompressedMatrix<int> M11( 3UL, 4UL, 4 );       // Instantiation of a 3x4 integer row-major
                                                // matrix with capacity for 4 non-zero elements.

Array Construction

Alternatively, all dense matrix classes offer a constructor for an initialization with a dynamic or static array. If the matrix is initialized from a dynamic array, the constructor expects the dimensions of values provided by the array as first and second argument, the array as third argument. In case of a static array, the fixed size of the array is used:

const std::unique_ptr<double[]> array1( new double[6] );
// ... Initialization of the dynamic array
blaze::StaticMatrix<double,2UL,3UL> M12( 2UL, 3UL, array1.get() );

int array2[2][2] = { { 4, -5 }, { -6, 7 } };
blaze::StaticMatrix<int,2UL,2UL,rowMajor> M13( array2 );

Initializer List Construction

In addition, all dense and sparse matrix classes can be directly initialized by means of an initializer list:

blaze::DynamicMatrix<float,columnMajor> M14{ {  3.1F,  6.4F },
                                                { -0.9F, -1.2F },
                                                {  4.8F,  0.6F } };
blaze::CompressedMatrix<int,rowMajor> M15{ { 3 },
                                           { 1 },
                                           { 0, 2 } };

Copy Construction

All dense and sparse matrices can be created as a copy of another dense or sparse matrix.

StaticMatrix<int,5UL,4UL,rowMajor> M16( M6 );    // Instantiation of the dense row-major matrix M16
                                                 // as copy of the dense row-major matrix M6.
DynamicMatrix<float,columnMajor> M17( M8 );      // Instantiation of the dense column-major matrix M17
                                                 // as copy of the sparse column-major matrix M8.
CompressedMatrix<double,columnMajor> M18( M7 );  // Instantiation of the compressed column-major matrix
                                                 // M18 as copy of the dense row-major matrix M7.
CompressedMatrix<float,rowMajor> M19( M8 );      // Instantiation of the compressed row-major matrix
                                                 // M19 as copy of the compressed column-major matrix M8.

Note that it is not possible to create a StaticMatrix as a copy of a matrix with a different number of rows and/or columns:

StaticMatrix<int,4UL,5UL,rowMajor> M20( M6 );     // Runtime error: Number of rows and columns does not match!
StaticMatrix<int,4UL,4UL,columnMajor> M21( M9 );  // Compile time error: Number of columns does not match!

Assignment

There are several types of assignment to dense and sparse matrices: Homogeneous Assignment, Array Assignment, Copy Assignment, and Compound Assignment.

Homogeneous Assignment

It is possible to assign the same value to all elements of a dense matrix. All dense matrix classes provide an according assignment operator:

blaze::StaticMatrix<int,3UL,2UL> M1;
blaze::DynamicMatrix<double> M2;

// Setting all integer elements of the StaticMatrix to 4
M1 = 4;

// Setting all double precision elements of the DynamicMatrix to 3.5
M2 = 3.5

Array Assignment

Dense matrices can also be assigned a static array:

blaze::StaticMatrix<int,2UL,2UL,rowMajor> M1;
blaze::StaticMatrix<int,2UL,2UL,columnMajor> M2;
blaze::DynamicMatrix<double> M3;

int array1[2][2] = { { 1, 2 }, { 3, 4 } };
double array2[3][2] = { { 3.1, 6.4 }, { -0.9, -1.2 }, { 4.8, 0.6 } };

M1 = array1;
M2 = array1;
M3 = array2;

Note that the dimensions of the static array have to match the size of a StaticMatrix, whereas a DynamicMatrix is resized according to the array dimensions:

     (  3.1  6.4 )
M3 = ( -0.9 -1.2 )
     (  4.8  0.6 )

Initializer List Assignment

Alternatively, it is possible to directly assign an initializer list to a dense or sparse matrix:

blaze::DynamicMatrix<double> M1;
blaze::CompressedMatrix<int> M2;

M1 = { { 3.1, 6.4 }, { -0.9, -1.2 }, { 4.8, 0.6 } };
M2 = { { 1, 0 }, {}, { 0, 1 }, { 2 } };

In case of sparse matrices, only the non-zero elements are considered. Missing values are considered to be default values.

Copy Assignment

All kinds of matrices can be assigned to each other. The only restriction is that since a StaticMatrix cannot change its size, the assigned matrix must match both in the number of rows and in the number of columns.

blaze::StaticMatrix<int,3UL,2UL,rowMajor> M1;
blaze::DynamicMatrix<int,rowMajor> M2( 3UL, 2UL );
blaze::DynamicMatrix<float,rowMajor> M3( 5UL, 2UL );
blaze::CompressedMatrix<int,rowMajor> M4( 3UL, 2UL );
blaze::CompressedMatrix<float,columnMajor> M5( 3UL, 2UL );

// ... Initialization of the matrices

M1 = M2;  // OK: Assignment of a 3x2 dense row-major matrix to another 3x2 dense row-major matrix
M1 = M4;  // OK: Assignment of a 3x2 sparse row-major matrix to a 3x2 dense row-major matrix
M1 = M3;  // Runtime error: Cannot assign a 5x2 matrix to a 3x2 static matrix
M1 = M5;  // OK: Assignment of a 3x2 sparse column-major matrix to a 3x2 dense row-major matrix

Compound Assignment

Compound assignment is also available for matrices: addition assignment, subtraction assignment, and multiplication assignment. In contrast to plain assignment, however, the number of rows and columns of the two operands have to match according to the arithmetic operation.

blaze::StaticMatrix<int,2UL,3UL,rowMajor> M1;
blaze::DynamicMatrix<int,rowMajor> M2( 2UL, 3UL );
blaze::CompressedMatrix<float,columnMajor> M3( 2UL, 3UL );
blaze::CompressedMatrix<float,rowMajor> M4( 2UL, 4UL );
blaze::StaticMatrix<float,2UL,4UL,rowMajor> M5;
blaze::CompressedMatrix<float,rowMajor> M6( 3UL, 2UL );

// ... Initialization of the matrices

M1 += M2 ; // OK: Addition assignment between two row-major matrices of the same dimensions
M1 -= M3;  // OK: Subtraction assignment between between a row-major and a column-major matrix
M1 += M4 ; // Runtime error: No compound assignment between matrices of different size
M1 -= M5;  // Compilation error: No compound assignment between matrices of different size
M2 *= M6;  // OK: Multiplication assignment between two row-major matrices

Note that the multiplication assignment potentially changes the number of columns of the target matrix:

( 2 0 1 )     ( 4 0 )     ( 8 3 )
( 0 3 2 )  *  ( 1 0 )  =  ( 3 6 )
              ( 0 3 )

Since a StaticMatrix cannot change its size, only a square StaticMatrix can be used in a multiplication assignment with other square matrices of the same dimensions.


Element Access

Function Call Operator

The easiest way to access a specific dense or sparse matrix element is via the function call operator. The indices to access a matrix are zero-based:

blaze::DynamicMatrix<int> M1( 4UL, 6UL );
M1(0,0) = 1;
M1(0,1) = 3;
// ...

blaze::CompressedMatrix<double> M2( 5UL, 3UL );
M2(0,2) =  4.1;
M2(1,1) = -6.3;

Since dense matrices allocate enough memory for all contained elements, using the function call operator on a dense matrix directly returns a reference to the accessed value. In case of a sparse matrix, if the accessed value is currently not contained in the matrix, the value is inserted into the matrix prior to returning a reference to the value, which can be much more expensive than the direct access to a dense matrix. Consider the following example:

blaze::CompressedMatrix<int> M1( 4UL, 4UL );

for( size_t i=0UL; i<M1.rows(); ++i ) {
   for( size_t j=0UL; j<M1.columns(); ++j ) {
      ... = M1(i,j);
   }
}

Although the compressed matrix is only used for read access within the for loop, using the function call operator temporarily inserts 16 non-zero elements into the matrix. Therefore the preferred way to traverse the non-zero elements of a sparse matrix is to use iterators.

Iterators

All matrices (sparse as well as dense) offer an alternate way via the begin(), cbegin(), end() and cend() functions to traverse all contained elements by iterator. Note that it is not possible to traverse all elements of the matrix, but that it is only possible to traverse elements in a row/column-wise fashion. In case of a non-const matrix, begin() and end() return an Iterator, which allows a manipulation of the non-zero value, in case of a constant matrix or in case cbegin() or cend() are used a ConstIterator is returned:

using blaze::CompressedMatrix;

CompressedMatrix<int,rowMajor> M1( 4UL, 6UL );

// Traversing the matrix by Iterator
for( size_t i=0UL; i<A.rows(); ++i ) {
   for( CompressedMatrix<int,rowMajor>::Iterator it=A.begin(i); it!=A.end(i); ++it ) {
      it->value() = ...;  // OK: Write access to the value of the non-zero element.
      ... = it->value();  // OK: Read access to the value of the non-zero element.
      it->index() = ...;  // Compilation error: The index of a non-zero element cannot be changed.
      ... = it->index();  // OK: Read access to the index of the non-zero element.
   }
}

// Traversing the matrix by ConstIterator
for( size_t i=0UL; i<A.rows(); ++i ) {
   for( CompressedMatrix<int,rowMajor>::ConstIterator it=A.cbegin(i); it!=A.cend(i); ++it ) {
      it->value() = ...;  // Compilation error: Assignment to the value via a ConstIterator is invalid.
      ... = it->value();  // OK: Read access to the value of the non-zero element.
      it->index() = ...;  // Compilation error: The index of a non-zero element cannot be changed.
      ... = it->index();  // OK: Read access to the index of the non-zero element.
   }
}

Note that begin(), cbegin(), end(), and cend() are also available as free functions:

for( size_t i=0UL; i<A.rows(); ++i ) {
   for( CompressedMatrix<int,rowMajor>::Iterator it=begin( A, i ); it!=end( A, i ); ++it ) {
      // ...
   }
}

for( size_t i=0UL; i<A.rows(); ++i ) {
   for( CompressedMatrix<int,rowMajor>::ConstIterator it=cbegin( A, i ); it!=cend( A, i ); ++it ) {
      // ...
   }
}

Element Insertion

Whereas a dense matrix always provides enough capacity to store all matrix elements, a sparse matrix only stores the non-zero elements. Therefore it is necessary to explicitly add elements to the matrix.

Function Call Operator

The first possibility to add elements to a sparse matrix is the function call operator:

using blaze::CompressedMatrix;

CompressedMatrix<int> M1( 3UL, 4UL );
M1(1,2) = 9;

In case the element at the given position is not yet contained in the sparse matrix, it is automatically inserted. Otherwise the old value is replaced by the new value 2. The operator returns a reference to the sparse vector element.

.set()

An alternative to the function call operator is the set() function: In case the element is not yet contained in the matrix the element is inserted, else the element's value is modified:

// Insert or modify the value at position (2,0)
M1.set( 2, 0, 1 );

.insert()

The insertion of elements can be better controlled via the insert() function. In contrast to the function call operator and the set() function it emits an exception in case the element is already contained in the matrix. In order to check for this case, the find() function can be used:

// In case the element at position (2,3) is not yet contained in the matrix it is inserted with a value of 4.
if( M1.find( 2, 3 ) == M1.end( 2 ) )
   M1.insert( 2, 3, 4 );

.append()

Although the insert() function is very flexible, due to performance reasons it is not suited for the setup of large sparse matrices. A very efficient, yet also very low-level way to fill a sparse matrix is the append() function. It requires the sparse matrix to provide enough capacity to insert a new element in the specified row/column. Additionally, the index of the new element must be larger than the index of the previous element in the same row/column. Violating these conditions results in undefined behavior!

M1.reserve( 0, 3 );     // Reserving space for three non-zero elements in row 0
M1.append( 0, 1,  2 );  // Appending the element 2 in row 0 at column index 1
M1.append( 0, 2, -4 );  // Appending the element -4 in row 0 at column index 2
// ...

The most efficient way to fill a sparse matrix with elements, however, is a combination of reserve(), append(), and the finalize() function:

// Setup of the compressed row-major matrix
//
//       ( 0 1 0 2 0 )
//   A = ( 0 0 0 0 0 )
//       ( 3 0 0 0 0 )
//
blaze::CompressedMatrix<int> M1( 3UL, 5UL );
M1.reserve( 3 );       // Reserving enough space for 3 non-zero elements
M1.append( 0, 1, 1 );  // Appending the value 1 in row 0 with column index 1
M1.append( 0, 3, 2 );  // Appending the value 2 in row 0 with column index 3
M1.finalize( 0 );      // Finalizing row 0
M1.finalize( 1 );      // Finalizing the empty row 1 to prepare row 2
M1.append( 2, 0, 3 );  // Appending the value 3 in row 2 with column index 0
M1.finalize( 2 );      // Finalizing row 2

Please note that the finalize() function has to be explicitly called for each row or column, even for empty ones! Also note that although append() does not allocate new memory, it still invalidates all iterators returned by the end() functions!


Element Removal

.erase()

The erase() member functions can be used to remove elements from a sparse matrix. The following example gives an impression of the five different flavors of erase():

using blaze::CompressedMatrix;

CompressedMatrix<int,rowMajor> A( 42, 53 );
// ... Initialization of the matrix

// Erasing the element at position (21,23)
A.erase( 21, 23 );

// Erasing a single element in row 17 via iterator
A.erase( 17, A.find( 4 ) );

// Erasing all non-zero elements in the range [7..24] of row 33
A.erase( 33, A.lowerBound( 33, 7 ), A.upperBound( 33, 24 ) );

// Erasing all non-zero elements with a value larger than 9 by passing a unary predicate
A.erase( []( int i ){ return i > 9; } );

// Erasing all non-zero elements in the range [30..40] of row 37 with a value larger than 5
CompressedMatrix<int,rowMajor>::Iterator pos1( A.lowerBound( 37, 30 ) );
CompressedMatrix<int,rowMajor>::Iterator pos2( A.upperBound( 37, 40 ) );
A.erase( 37, pos1, pos2, []( int i ){ return i > 5; } );

Element Lookup

A sparse matrix only stores the non-zero elements contained in the matrix. Therefore, whenever accessing a matrix element at a specific position a lookup operation is required. Whereas the function call operator is performing this lookup automatically, it is also possible to use the find(), lowerBound(), and upperBound() member functions for a manual lookup.

.find()

The find() function can be used to check whether a specific element is contained in the sparse matrix. It specifically searches for the element at the specified position. In case the element is found, the function returns an iterator to the element. Otherwise an iterator just past the last non-zero element of the according row or column (the end() iterator) is returned. Note that the returned iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!

using blaze::CompressedMatrix;

CompressedMatrix<int,rowMajor> A( 42, 53 );
// ... Initialization of the matrix

// Searching the element at position (7,17). In case the element is not
// contained in the vector, the end() iterator of row 7 is returned.
CompressedMatrix<int,rowMajor>::Iterator pos( A.find( 7, 17 ) );

if( pos != A.end( 7 ) ) {
   // ...
}

.lowerBound()

In case of a row-major matrix, this function returns a row iterator to the first element with an index not less then the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index not less then the given row index. In combination with the upperBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!

using blaze::CompressedMatrix;

CompressedMatrix<int,rowMajor> A( 42, 53 );
// ... Initialization of the matrix

// Searching the lower bound of column index 17 in row 7.
CompressedMatrix<int,rowMajor>::Iterator pos1( A.lowerBound( 7, 17 ) );

// Searching the upper bound of column index 28 in row 7
CompressedMatrix<int,rowMajor>::Iterator pos2( A.upperBound( 7, 28 ) );

// Erasing all elements in the specified range
A.erase( 7, pos1, pos2 );

.upperBound()

In case of a row-major matrix, this function returns a row iterator to the first element with an index greater then the given column index. In case of a column-major matrix, the function returns a column iterator to the first element with an index greater then the given row index. In combination with the lowerBound() function this function can be used to create a pair of iterators specifying a range of indices. Note that the returned iterator is subject to invalidation due to inserting operations via the function call operator, the set() function or the insert() function!

using blaze::CompressedMatrix;

CompressedMatrix<int,columnMajor> A( 42, 53 );
// ... Initialization of the matrix

// Searching the lower bound of row index 17 in column 9.
CompressedMatrix<int,columnMajor>::Iterator pos1( A.lowerBound( 17, 9 ) );

// Searching the upper bound of row index 28 in column 9
CompressedMatrix<int,columnMajor>::Iterator pos2( A.upperBound( 28, 9 ) );

// Erasing all elements in the specified range
A.erase( 9, pos1, pos2 );

Non-Modifying Operations

.rows()

The current number of rows of a matrix can be acquired via the rows() member function:

// Instantiating a dynamic matrix with 10 rows and 8 columns
blaze::DynamicMatrix<int> M1( 10UL, 8UL );
M1.rows();  // Returns 10

// Instantiating a compressed matrix with 8 rows and 12 columns
blaze::CompressedMatrix<double> M2( 8UL, 12UL );
M2.rows();  // Returns 8

Alternatively, the free functions rows() can be used to query the current number of rows of a matrix. In contrast to the member function, the free function can also be used to query the number of rows of a matrix expression:

rows( M1 );  // Returns 10, i.e. has the same effect as the member function
rows( M2 );  // Returns 8, i.e. has the same effect as the member function

rows( M1 * M2 );  // Returns 10, i.e. the number of rows of the resulting matrix

.columns()

The current number of columns of a matrix can be acquired via the columns() member function:

// Instantiating a dynamic matrix with 6 rows and 8 columns
blaze::DynamicMatrix<int> M1( 6UL, 8UL );
M1.columns();   // Returns 8

// Instantiating a compressed matrix with 8 rows and 7 columns
blaze::CompressedMatrix<double> M2( 8UL, 7UL );
M2.columns();   // Returns 7

There is also a free function columns() available, which can also be used to query the number of columns of a matrix expression:

columns( M1 );  // Returns 8, i.e. has the same effect as the member function
columns( M2 );  // Returns 7, i.e. has the same effect as the member function

columns( M1 * M2 );  // Returns 7, i.e. the number of columns of the resulting matrix

.capacity()

The capacity() member function returns the internal capacity of a dense or sparse matrix. Note that the capacity of a matrix doesn't have to be equal to the size of a matrix. In case of a dense matrix the capacity will always be greater or equal than the total number of elements of the matrix. In case of a sparse matrix, the capacity will usually be much less than the total number of elements.

blaze::DynamicMatrix<float> M1( 5UL, 7UL );
blaze::StaticMatrix<float,7UL,4UL> M2;
M1.capacity();  // Returns at least 35
M2.capacity();  // Returns at least 28

There is also a free function capacity() available to query the capacity. However, please note that this function cannot be used to query the capacity of a matrix expression:

capacity( M1 );  // Returns at least 35, i.e. has the same effect as the member function
capacity( M2 );  // Returns at least 28, i.e. has the same effect as the member function

capacity( M1 * M2 );  // Compilation error!

.nonZeros()

For both dense and sparse matrices the current number of non-zero elements can be queried via the nonZeros() member function. In case of matrices there are two flavors of the nonZeros() function: One returns the total number of non-zero elements in the matrix, the second returns the number of non-zero elements in a specific row (in case of a row-major matrix) or column (in case of a column-major matrix). Sparse matrices directly return their number of non-zero elements, dense matrices traverse their elements and count the number of non-zero elements.

blaze::DynamicMatrix<int,rowMajor> M1( 3UL, 5UL );

// ... Initializing the dense matrix

M1.nonZeros();     // Returns the total number of non-zero elements in the dense matrix
M1.nonZeros( 2 );  // Returns the number of non-zero elements in row 2
blaze::CompressedMatrix<double,columnMajor> M2( 4UL, 7UL );

// ... Initializing the sparse matrix

M2.nonZeros();     // Returns the total number of non-zero elements in the sparse matrix
M2.nonZeros( 3 );  // Returns the number of non-zero elements in column 3

The free nonZeros() function can also be used to query the number of non-zero elements in a matrix expression. However, the result is not the exact number of non-zero elements, but may be a rough estimation:

nonZeros( M1 );     // Has the same effect as the member function
nonZeros( M1, 2 );  // Has the same effect as the member function

nonZeros( M2 );     // Has the same effect as the member function
nonZeros( M2, 3 );  // Has the same effect as the member function

nonZeros( M1 * M2 );  // Estimates the number of non-zero elements in the matrix expression

isnan()

The isnan() function provides the means to check a dense or sparse matrix for non-a-number elements:

blaze::DynamicMatrix<double> A( 3UL, 4UL );
// ... Initialization
if( isnan( A ) ) { ... }
blaze::CompressedMatrix<double> A( 3UL, 4UL );
// ... Initialization
if( isnan( A ) ) { ... }

If at least one element of the matrix is not-a-number, the function returns true, otherwise it returns false. Please note that this function only works for matrices with floating point elements. The attempt to use it for a matrix with a non-floating point element type results in a compile time error.

isDefault()

The isDefault() function returns whether the given dense or sparse matrix is in default state:

blaze::HybridMatrix<int,5UL,4UL> A;
// ... Resizing and initialization
if( isDefault( A ) ) { ... }

A matrix is in default state if it appears to just have been default constructed. A resizable matrix (HybridMatrix, DynamicMatrix, or CompressedMatrix) is in default state if its size is equal to zero. A non-resizable matrix (StaticMatrix and all submatrices) is in default state if all its elements are in default state. For instance, in case the matrix is instantiated for a built-in integral or floating point data type, the function returns true in case all matrix elements are 0 and false in case any vector element is not 0.

isSquare()

If a dense or sparse matrix is a square matrix (i.e. if the number of rows is equal to the number of columns) can be checked via the isSquare() function:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
if( isSquare( A ) ) { ... }

isSymmetric()

Via the isSymmetric() function it is possible to check whether a dense or sparse matrix is symmetric:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isSymmetric( A ) ) { ... }

Note that non-square matrices are never considered to be symmetric!

isUniform()

In order to check if all matrix elements are identical, the isUniform function can be used:

blaze::DynamicMatrix<int> A;
// ... Resizing and initialization
if( isUniform( A ) ) { ... }

Note that in case of a sparse matrix also the zero elements are also taken into account!

isLower()

Via the isLower() function it is possible to check whether a dense or sparse matrix is lower triangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isLower( A ) ) { ... }

Note that non-square matrices are never considered to be lower triangular!

isUniLower()

Via the isUniLower() function it is possible to check whether a dense or sparse matrix is lower unitriangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isUniLower( A ) ) { ... }

Note that non-square matrices are never considered to be lower unitriangular!

isStrictlyLower()

Via the isStrictlyLower() function it is possible to check whether a dense or sparse matrix is strictly lower triangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isStrictlyLower( A ) ) { ... }

Note that non-square matrices are never considered to be strictly lower triangular!

isUpper()

Via the isUpper() function it is possible to check whether a dense or sparse matrix is upper triangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isUpper( A ) ) { ... }

Note that non-square matrices are never considered to be upper triangular!

isUniUpper()

Via the isUniUpper() function it is possible to check whether a dense or sparse matrix is upper unitriangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isUniUpper( A ) ) { ... }

Note that non-square matrices are never considered to be upper unitriangular!

isStrictlyUpper()

Via the isStrictlyUpper() function it is possible to check whether a dense or sparse matrix is strictly upper triangular:

blaze::DynamicMatrix<float> A;
// ... Resizing and initialization
if( isStrictlyUpper( A ) ) { ... }

Note that non-square matrices are never considered to be strictly upper triangular!

isDiagonal()

The isDiagonal() function checks if the given dense or sparse matrix is a diagonal matrix, i.e. if it has only elements on its diagonal and if the non-diagonal elements are default elements:

blaze::CompressedMatrix<float> A;
// ... Resizing and initialization
if( isDiagonal( A ) ) { ... }

Note that non-square matrices are never considered to be diagonal!

isIdentity()

The isIdentity() function checks if the given dense or sparse matrix is an identity matrix, i.e. if all diagonal elements are 1 and all non-diagonal elements are 0:

blaze::CompressedMatrix<float> A;
// ... Resizing and initialization
if( isIdentity( A ) ) { ... }

Note that non-square matrices are never considered to be identity matrices!

det()

The determinant of a square dense matrix can be computed by means of the det() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization
double d = det( A );  // Compute the determinant of A

In case the given dense matrix is not a square matrix, a std::invalid_argument exception is thrown.

Note that the det() function can only be used for dense matrices with float, double, complex<float> or complex<double> element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!

Also note that the function is depending on LAPACK kernels. Thus the function can only be used if a fitting LAPACK library is available and linked to the executable. Otherwise a linker error will be created.

trans()

Matrices can be transposed via the trans() function. Row-major matrices are transposed into a column-major matrix and vice versa:

blaze::DynamicMatrix<int,rowMajor> M1( 5UL, 2UL );
blaze::CompressedMatrix<int,columnMajor> M2( 3UL, 7UL );

M1 = M2;            // Assigning a column-major matrix to a row-major matrix
M1 = trans( M2 );   // Assigning the transpose of M2 (i.e. a row-major matrix) to M1
M1 += trans( M2 );  // Addition assignment of two row-major matrices

ctrans()

The conjugate transpose of a dense or sparse matrix (also called adjoint matrix, Hermetian conjugate, or transjugate) can be computed via the ctrans() function:

blaze::DynamicMatrix< complex<float>, rowMajor > M1( 5UL, 2UL );
blaze::CompressedMatrix< complex<float>, columnMajor > M2( 2UL, 5UL );

M1 = ctrans( M2 );  // Compute the conjugate transpose matrix

Note that the ctrans() function has the same effect as manually applying the conj() and trans() function in any order:

M1 = trans( conj( M2 ) );  // Computing the conjugate transpose matrix
M1 = conj( trans( M2 ) );  // Computing the conjugate transpose matrix

eval() / evaluate()

The evaluate() function forces an evaluation of the given matrix expression and enables an automatic deduction of the correct result type of an operation. The following code example demonstrates its intended use for the multiplication of a lower and a strictly lower dense matrix:

using blaze::DynamicMatrix;
using blaze::LowerMatrix;
using blaze::StrictlyLowerMatrix;

LowerMatrix< DynamicMatrix<double> > A;
StrictlyLowerMatrix< DynamicMatrix<double> > B;
// ... Resizing and initialization

auto C = evaluate( A * B );

In this scenario, the evaluate() function assists in deducing the exact result type of the operation via the auto keyword. Please note that if evaluate() is used in this way, no temporary matrix is created and no copy operation is performed. Instead, the result is directly written to the target matrix due to the return value optimization (RVO). However, if evaluate() is used in combination with an explicit target type, a temporary will be created and a copy operation will be performed if the used type differs from the type returned from the function:

StrictlyLowerMatrix< DynamicMatrix<double> > D( A * B );  // No temporary & no copy operation
LowerMatrix< DynamicMatrix<double> > E( A * B );          // Temporary & copy operation
DynamicMatrix<double> F( A * B );                         // Temporary & copy operation
D = evaluate( A * B );                                    // Temporary & copy operation

Sometimes it might be desirable to explicitly evaluate a sub-expression within a larger expression. However, please note that evaluate() is not intended to be used for this purpose. This task is more elegantly and efficiently handled by the eval() function:

blaze::DynamicMatrix<double> A, B, C, D;

D = A + evaluate( B * C );  // Unnecessary creation of a temporary matrix
D = A + eval( B * C );      // No creation of a temporary matrix

In contrast to the evaluate() function, eval() can take the complete expression into account and therefore can guarantee the most efficient way to evaluate it (see also Intra-Statement Optimization).


Modifying Operations

.resize() / .reserve()

The dimensions of a StaticMatrix are fixed at compile time by the second and third template parameter. In contrast, the number or rows and/or columns of DynamicMatrix and CompressedMatrix can be changed at runtime:

using blaze::DynamicMatrix;
using blaze::CompressedMatrix;

DynamicMatrix<int,rowMajor> M1;
CompressedMatrix<int,columnMajor> M2( 3UL, 2UL );

// Adapting the number of rows and columns via the resize() function. The (optional)
// third parameter specifies whether the existing elements should be preserved. Per
// default, the existing elements are preserved.
M1.resize( 2UL, 2UL );         // Resizing matrix M1 to 2x2 elements. Elements of built-in type remain
                               // uninitialized, elements of class type are default constructed.
M1.resize( 3UL, 1UL, false );  // Resizing M1 to 3x1 elements. The old elements are lost, the
                               // new elements are NOT initialized!
M2.resize( 5UL, 7UL, true );   // Resizing M2 to 5x7 elements. The old elements are preserved.
M2.resize( 3UL, 2UL, false );  // Resizing M2 to 3x2 elements. The old elements are lost.

Note that resizing a matrix invalidates all existing views (see e.g. Submatrices) on the matrix:

using MatrixType = blaze::DynamicMatrix<int,rowMajor>;
using RowType    = blaze::Row<MatrixType>;

MatrixType M1( 10UL, 20UL );    // Creating a 10x20 matrix
RowType row8 = row( M1, 8UL );  // Creating a view on the 8th row of the matrix
M1.resize( 6UL, 20UL );         // Resizing the matrix invalidates the view

When the internal capacity of a matrix is no longer sufficient, the allocation of a larger junk of memory is triggered. In order to avoid frequent reallocations, the reserve() function can be used up front to set the internal capacity:

blaze::DynamicMatrix<int> M1;
M1.reserve( 100 );
M1.rows();      // Returns 0
M1.capacity();  // Returns at least 100

Additionally it is possible to reserve memory in a specific row (for a row-major matrix) or column (for a column-major matrix):

blaze::CompressedMatrix<int> M1( 4UL, 6UL );
M1.reserve( 1, 4 );  // Reserving enough space for four non-zero elements in row 1

.shrinkToFit()

The internal capacity of matrices with dynamic memory is preserved in order to minimize the number of reallocations. For that reason, the resize() and reserve() functions can lead to memory overhead. The shrinkToFit() member function can be used to minimize the internal capacity:

blaze::DynamicMatrix<int> M1( 100UL, 100UL );  // Create a 100x100 integer matrix
M1.resize( 10UL, 10UL );                       // Resize to 10x10, but the capacity is preserved
M1.shrinkToFit();                              // Remove the unused capacity

Please note that due to padding the capacity might not be reduced exactly to rows() times columns(). Please also note that in case a reallocation occurs, all iterators (including end() iterators), all pointers and references to elements of this matrix are invalidated.

reset() / clear()

In order to reset all elements of a dense or sparse matrix, the reset() function can be used. The number of rows and columns of the matrix are preserved:

// Setting up a single precision row-major matrix, whose elements are initialized with 2.0F.
blaze::DynamicMatrix<float> M1( 4UL, 5UL, 2.0F );

// Resetting all elements to 0.0F.
reset( M1 );  // Resetting all elements
M1.rows();    // Returns 4: size and capacity remain unchanged

Alternatively, only a single row or column of the matrix can be resetted:

blaze::DynamicMatrix<int,blaze::rowMajor>    M1( 7UL, 6UL, 5 );  // Setup of a row-major matrix
blaze::DynamicMatrix<int,blaze::columnMajor> M2( 4UL, 5UL, 4 );  // Setup of a column-major matrix

reset( M1, 2UL );  // Resetting the 2nd row of the row-major matrix
reset( M2, 3UL );  // Resetting the 3rd column of the column-major matrix

In order to reset a row of a column-major matrix or a column of a row-major matrix, use a row or column view (see \ref views_rows and views_colums).

In order to return a matrix to its default state (i.e. the state of a default constructed matrix), the clear() function can be used:

// Setting up a single precision row-major matrix, whose elements are initialized with 2.0F.
blaze::DynamicMatrix<float> M1( 4UL, 5UL, 2.0F );

// Resetting all elements to 0.0F.
clear( M1 );  // Resetting the entire matrix
M1.rows();    // Returns 0: size is reset, but capacity remains unchanged

transpose()

In addition to the non-modifying trans() function, matrices can be transposed in-place via the transpose() function:

blaze::DynamicMatrix<int,rowMajor> M( 5UL, 2UL );

transpose( M );  // In-place transpose operation.
M = trans( M );  // Same as above

Note however that the transpose operation fails if ...

  • ... the given matrix has a fixed size and is non-square;
  • ... the given matrix is a triangular matrix;
  • ... the given submatrix affects the restricted parts of a triangular matrix;
  • ... the given submatrix would cause non-deterministic results in a symmetric/Hermitian matrix.

ctranspose()

The ctranspose() function can be used to perform an in-place conjugate transpose operation:

blaze::DynamicMatrix<int,rowMajor> M( 5UL, 2UL );

ctranspose( M );  // In-place conjugate transpose operation.
M = ctrans( M );  // Same as above

Note however that the conjugate transpose operation fails if ...

  • ... the given matrix has a fixed size and is non-square;
  • ... the given matrix is a triangular matrix;
  • ... the given submatrix affects the restricted parts of a triangular matrix;
  • ... the given submatrix would cause non-deterministic results in a symmetric/Hermitian matrix.

swap()

Via the swap() function it is possible to completely swap the contents of two matrices of the same type:

blaze::DynamicMatrix<int,blaze::rowMajor> M1( 10UL, 15UL );
blaze::DynamicMatrix<int,blaze::rowMajor> M2( 20UL, 10UL );

swap( M1, M2 );  // Swapping the contents of M1 and M2

Arithmetic Operations

min() / max()

The min() and max() functions can be used for a single matrix or multiple matrices. If passed a single matrix, the functions return the smallest and largest element of the given dense or sparse matrix, respectively:

using blaze::rowMajor;

blaze::StaticMatrix<int,2UL,3UL,rowMajor> A{ { -5, 2, 7 },
                                             { -4, 0, 1 } };

min( A );  // Returns -5
max( A );  // Returns 7

In case the matrix currently has 0 rows or 0 columns, both functions return 0. Additionally, in case a given sparse matrix is not completely filled, the zero elements are taken into account. For example: the following compressed matrix has only 2 non-zero elements. However, the minimum of this matrix is 0:

blaze::CompressedMatrix<int> B{ { 1, 0, 3 },
                                { 0, 0, 0 } };

min( B );  // Returns 0

If passed two or more dense matrices, the min() and max() functions compute the componentwise minimum or maximum of the given matrices, respectively:

blaze::StaticMatrix<int,2UL,3UL,rowMajor> C{ { -5, 1, -7 }, { 4, 1, 0 } };
blaze::StaticMatrix<int,2UL,3UL,rowMajor> D{ { -5, 3,  0 }, { 2, 2, -2 } };

min( A, C );     // Results in the matrix ( -5, 1, -7 ) ( -4, 0, 0 )
max( A, C, D );  // Results in the matrix ( -5, 3, 7 ) (  4, 2, 1 )

Please note that sparse matrices can only be used in the unary min() and max() functions. Also note that all forms of the min() and max() functions can be used to compute the smallest and largest element of a matrix expression:

min( A + B + C );  // Returns -9, i.e. the smallest value of the resulting matrix
max( A - B - C );  // Returns 11, i.e. the largest value of the resulting matrix

trace()

The trace() function sums the diagonal elements of a square dense or sparse matrix:

blaze::StaticMatrix<int,3UL,3UL> A{ { -1,  2, -3 }
                                  , { -4, -5,  6 }
                                  , {  7, -8, -9 } };

trace( A );  // Returns the sum of the diagonal elements, i.e. -15

In case the given matrix is not a square matrix, a std::invalid_argument exception is thrown.

abs()

The abs() function can be used to compute the absolute values of each element of a matrix. For instance, the following computation

blaze::StaticMatrix<int,2UL,3UL,rowMajor> A{ { -1,  2, -3 },
                                             {  4, -5,  6 } };
blaze::StaticMatrix<int,2UL,3UL,rowMajor> B( abs( A ) );

results in the matrix

B = ( 1 2 3 )
    ( 4 5 6 )

floor() / ceil() / trunc() / round()

The floor(), ceil(), trunc(), and round() functions can be used to round down/up each element of a matrix, respectively:

blaze::StaticMatrix<double,3UL,3UL> A, B;

B = floor( A );  // Rounding down each element of the matrix
B = ceil ( A );  // Rounding up each element of the matrix
B = trunc( A );  // Truncating each element of the matrix
B = round( A );  // Rounding each element of the matrix

conj()

The conj() function can be applied on a dense or sparse matrix to compute the complex conjugate of each element of the matrix:

using blaze::StaticMatrix;

using cplx = std::complex<double>;

// Creating the matrix
//    ( (1,0)  (-2,-1) )
//    ( (1,1)  ( 0, 1) )
StaticMatrix<cplx,2UL,2UL> A{ { cplx( 1.0, 0.0 ), cplx( -2.0, -1.0 ) },
                              { cplx( 1.0, 1.0 ), cplx(  0.0,  1.0 ) } };

// Computing the matrix of conjugate values
//    ( (1, 0)  (-2, 1) )
//    ( (1,-1)  ( 0,-1) )
StaticMatrix<cplx,2UL,2UL> B;
B = conj( A );

Additionally, matrices can be conjugated in-place via the conjugate() function:

blaze::DynamicMatrix<cplx> C( 5UL, 2UL );

conjugate( C );  // In-place conjugate operation.
C = conj( C );   // Same as above

real()

The real() function can be used on a dense or sparse matrix to extract the real part of each element of the matrix:

using blaze::StaticMatrix;

using cplx = std::complex<double>;

// Creating the matrix
//    ( (1,0)  (-2,-1) )
//    ( (1,1)  ( 0, 1) )
StaticMatrix<cplx,2UL,2UL> A{ { cplx( 1.0, 0.0 ), cplx( -2.0, -1.0 ) },
                              { cplx( 1.0, 1.0 ), cplx(  0.0,  1.0 ) } };

// Extracting the real part of each matrix element
//    ( 1 -2 )
//    ( 1  0 )
StaticMatrix<double,2UL,2UL> B;
B = real( A );

imag()

The imag() function can be used on a dense or sparse matrix to extract the imaginary part of each element of the matrix:

using blaze::StaticMatrix;

using cplx = std::complex<double>;

// Creating the matrix
//    ( (1,0)  (-2,-1) )
//    ( (1,1)  ( 0, 1) )
StaticMatrix<cplx,2UL,2UL> A{ { cplx( 1.0, 0.0 ), cplx( -2.0, -1.0 ) },
                              { cplx( 1.0, 1.0 ), cplx(  0.0,  1.0 ) } };

// Extracting the imaginary part of each matrix element
//    ( 0 -1 )
//    ( 1  1 )
StaticMatrix<double,2UL,2UL> B;
B = imag( A );

sqrt() / invsqrt()

Via the sqrt() and invsqrt() functions the (inverse) square root of each element of a matrix can be computed:

blaze::StaticMatrix<double,3UL,3UL> A, B, C;

B = sqrt( A );     // Computes the square root of each element
C = invsqrt( A );  // Computes the inverse square root of each element

Note that in case of sparse matrices only the non-zero elements are taken into account!

cbrt() / invcbrt()

The cbrt() and invcbrt() functions can be used to compute the the (inverse) cubic root of each element of a matrix:

blaze::DynamicMatrix<double> A, B, C;

B = cbrt( A );     // Computes the cubic root of each element
C = invcbrt( A );  // Computes the inverse cubic root of each element

Note that in case of sparse matrices only the non-zero elements are taken into account!

hypot()

The hypot() function can be used to compute the componentwise hypotenous for a pair of dense matrices:

blaze::StaticMatrix<double,3UL,3UL> A, B, C;

C = hypot( A, B );  // Computes the componentwise hypotenuous

clamp()

The clamp() function can be used to restrict all elements of a matrix to a specific range:

blaze::DynamicMatrix<double> A, B;

B = clamp( A, -1.0, 1.0 );  // Restrict all elements to the range [-1..1]

Note that in case of sparse matrices only the non-zero elements are taken into account!

pow()

The pow() function can be used to compute the exponential value of each element of a matrix. If passed a matrix and a numeric exponent, the function computes the exponential value of each element of the matrix using the same exponent. If passed a second matrix, the function computes the componentwise exponential value:

blaze::StaticMatrix<double,3UL,3UL> A, B, C;

C = pow( A, 1.2 );  // Computes the exponential value of each element
C = pow( A, B );    // Computes the componentwise exponential value

exp() / exp2() / exp10()

exp(), exp2() and exp10() compute the base e/2/10 exponential of each element of a matrix, respectively:

blaze::HybridMatrix<double,3UL,3UL> A, B;

B = exp( A );    // Computes the base e exponential of each element
B = exp2( A );   // Computes the base 2 exponential of each element
B = exp10( A );  // Computes the base 10 exponential of each element

Note that in case of sparse matrices only the non-zero elements are taken into account!

log() / log2() / log10()

The log(), log2() and log10() functions can be used to compute the natural, binary and common logarithm of each element of a matrix:

blaze::StaticMatrix<double,3UL,3UL> A, B;

B = log( A );    // Computes the natural logarithm of each element
B = log2( A );   // Computes the binary logarithm of each element
B = log10( A );  // Computes the common logarithm of each element

sin() / cos() / tan() / asin() / acos() / atan()

The following trigonometric functions are available for both dense and sparse matrices:

blaze::DynamicMatrix<double> A, B;

B = sin( A );  // Computes the sine of each element of the matrix
B = cos( A );  // Computes the cosine of each element of the matrix
B = tan( A );  // Computes the tangent of each element of the matrix

B = asin( A );  // Computes the inverse sine of each element of the matrix
B = acos( A );  // Computes the inverse cosine of each element of the matrix
B = atan( A );  // Computes the inverse tangent of each element of the matrix

Note that in case of sparse matrices only the non-zero elements are taken into account!

sinh() / cosh() / tanh() / asinh() / acosh() / atanh()

The following hyperbolic functions are available for both dense and sparse matrices:

blaze::DynamicMatrix<double> A, B;

B = sinh( A );  // Computes the hyperbolic sine of each element of the matrix
B = cosh( A );  // Computes the hyperbolic cosine of each element of the matrix
B = tanh( A );  // Computes the hyperbolic tangent of each element of the matrix

B = asinh( A );  // Computes the inverse hyperbolic sine of each element of the matrix
B = acosh( A );  // Computes the inverse hyperbolic cosine of each element of the matrix
B = atanh( A );  // Computes the inverse hyperbolic tangent of each element of the matrix

atan2()

The multi-valued inverse tangent is available for a pair of dense matrices:

blaze::DynamicMatrix<double> A, B, C;

C = atan2( A, B );  // Computes the componentwise multi-valued inverse tangent

erf() / erfc()

The erf() and erfc() functions compute the (complementary) error function of each element of a matrix:

blaze::StaticMatrix<double,3UL,3UL> A, B;

B = erf( A );   // Computes the error function of each element
B = erfc( A );  // Computes the complementary error function of each element

Note that in case of sparse matrices only the non-zero elements are taken into account!

map() / forEach()

Via the unary and binary map() functions it is possible to execute componentwise custom operations on matrices. The unary map() function can be used to apply a custom operation on each element of a dense or sparse matrix. For instance, the following example demonstrates a custom square root computation via a lambda:

blaze::DynamicMatrix<double> A, B;

B = map( A, []( double d ) { return std::sqrt( d ); } );

The binary map() function can be used to apply an operation pairwise to the elements of two dense matrices. The following example demonstrates the merging of two matrices of double precision values into a matrix of double precision complex numbers:

blaze::DynamicMatrix<double> real{ { 2.1, -4.2 }, { 1.0,  0.6 } };
blaze::DynamicMatrix<double> imag{ { 0.3,  1.4 }, { 2.9, -3.4 } };

blaze::DynamicMatrix< complex<double> > cplx;

// Creating the matrix
//    ( (-2.1,  0.3) (-4.2, -1.4) )
//    ( ( 1.0,  2.9) ( 0.6, -3.4) )
cplx = map( real, imag, []( double r, double i ){ return complex( r, i ); } );

Although the computation can be parallelized it is not vectorized and thus cannot perform at peak performance. However, it is also possible to create vectorized custom operations. See Custom Operations for a detailed overview of the possibilities of custom operations.

Please note that unary custom operations on vectors have been introduced in Blaze 3.0 in form of the forEach() function. With the introduction of binary custom functions, the forEach() function has been renamed to map(). The forEach() function can still be used (even for binary custom operations), but the function might be deprecated in future releases of Blaze.


Norms

norm()

The norm() function computes the L2 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l2 = norm( A );

sqrNorm()

The sqrNorm() function computes the squared L2 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l2 = sqrNorm( A );

l1Norm()

The l1Norm() function computes the squared L1 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l1 = l1Norm( A );

l2Norm()

The l2Norm() function computes the squared L2 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l2 = l2Norm( A );

l3Norm()

The l3Norm() function computes the squared L3 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l3 = l3Norm( A );

l4Norm()

The l4Norm() function computes the squared L4 norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double l4 = l4Norm( A );

lpNorm()

The lpNorm() function computes the general Lp norm of the given dense or sparse matrix, where the norm is specified by either a compile time or a runtime argument:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double lp1 = lpNorm<2>( A );    // Compile time argument
const double lp2 = lpNorm( A, 2.3 );  // Runtime argument

maxNorm()

The maxNorm() function computes the maximum norm of the given dense or sparse matrix:

blaze::DynamicMatrix<double> A;
// ... Resizing and initialization
const double max = maxNorm( A );

Declaration Operations

declsym()

The declsym() operation can be used to explicitly declare any matrix or matrix expression as symmetric:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = declsym( A );

Any matrix or matrix expression that has been declared as symmetric via declsym() will gain all the benefits of a symmetric matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::SymmetricMatrix;

DynamicMatrix<double> A, B, C;
SymmetricMatrix< DynamicMatrix<double> > S;
// ... Resizing and initialization

isSymmetric( declsym( A ) );  // Will always return true without runtime effort

S = declsym( A );  // Omit any runtime check for symmetry

C = declsym( A * B );  // Declare the result of the matrix multiplication as symmetric,
                       // i.e. perform an optimized matrix multiplication

Warning: The declsym() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-symmetric matrix or matrix expression as symmetric via the declsym() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!

declherm()

The declherm() operation can be used to explicitly declare any matrix or matrix expression as Hermitian:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = declherm( A );

Any matrix or matrix expression that has been declared as Hermitian via declherm() will gain all the benefits of an Hermitian matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::HermitianMatrix;

DynamicMatrix<double> A, B, C;
HermitianMatrix< DynamicMatrix<double> > S;
// ... Resizing and initialization

isHermitian( declherm( A ) );  // Will always return true without runtime effort

S = declherm( A );  // Omit any runtime check for Hermitian symmetry

C = declherm( A * B );  // Declare the result of the matrix multiplication as Hermitian,
                        // i.e. perform an optimized matrix multiplication

Warning: The declherm() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-Hermitian matrix or matrix expression as Hermitian via the declherm() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!

decllow()

The decllow() operation can be used to explicitly declare any matrix or matrix expression as lower triangular:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = decllow( A );

Any matrix or matrix expression that has been declared as lower triangular via decllow() will gain all the benefits of a lower triangular matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::LowerMatrix;

DynamicMatrix<double> A, B, C;
LowerMatrix< DynamicMatrix<double> > L;
// ... Resizing and initialization

isLower( decllow( A ) );  // Will always return true without runtime effort

L = decllow( A );  // Omit any runtime check for A being a lower matrix

C = decllow( A * B );  // Declare the result of the matrix multiplication as lower triangular,
                       // i.e. perform an optimized matrix multiplication

Warning: The decllow() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-lower matrix or matrix expression as lower triangular via the decllow() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!

declupp()

The declupp() operation can be used to explicitly declare any matrix or matrix expression as upper triangular:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = declupp( A );

Any matrix or matrix expression that has been declared as upper triangular via declupp() will gain all the benefits of a upper triangular matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::UpperMatrix;

DynamicMatrix<double> A, B, C;
UpperMatrix< DynamicMatrix<double> > U;
// ... Resizing and initialization

isUpper( declupp( A ) );  // Will always return true without runtime effort

U = declupp( A );  // Omit any runtime check for A being a upper matrix

C = declupp( A * B );  // Declare the result of the matrix multiplication as upper triangular,
                       // i.e. perform an optimized matrix multiplication

Warning: The declupp() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-upper matrix or matrix expression as upper triangular via the declupp() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!

decldiag()

The decldiag() operation can be used to explicitly declare any matrix or matrix expression as diagonal:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = decldiag( A );

Any matrix or matrix expression that has been declared as diagonal via decldiag() will gain all the benefits of a diagonal matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::DiagonalMatrix;

DynamicMatrix<double> A, B, C;
DiagonalMatrix< DynamicMatrix<double> > D;
// ... Resizing and initialization

isDiagonal( decldiag( A ) );  // Will always return true without runtime effort

D = decldiag( A );  // Omit any runtime check for A being a diagonal matrix

C = decldiag( A * B );  // Declare the result of the matrix multiplication as diagonal,
                        // i.e. perform an optimized matrix multiplication

Warning: The decldiag() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-diagonal matrix or matrix expression as diagonal via the decldiag() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!

declid()

The declid() operation can be used to explicitly declare any matrix or matrix expression as identity matrix:

blaze::DynamicMatrix<double> A, B;
// ... Resizing and initialization

B = declid( A );

Any matrix or matrix expression that has been declared as identity matrix via declid() will gain all the benefits of an identity matrix, which range from reduced runtime checking to a considerable speed-up in computations:

using blaze::DynamicMatrix;
using blaze::DiagonalMatrix;

DynamicMatrix<double> A, B, C;
DiagonalMatrix< DynamicMatrix<double> > D;
// ... Resizing and initialization

isIdentity( declid( A ) );  // Will always return true without runtime effort

D = declid( A );  // Omit any runtime check for A being a diagonal matrix

C = declid( A ) * B;  // Declare the left operand of the matrix multiplication as an
                      // identity matrix, i.e. perform an optimized matrix multiplication

Warning: The declid() operation has the semantics of a cast: The caller is completely responsible and the system trusts the given information. Declaring a non-identity matrix or matrix expression as identity matrix via the declid() operation leads to undefined behavior (which can be violated invariants or wrong computation results)!


Matrix Inversion

The inverse of a square dense matrix can be computed via the inv() function:

blaze::DynamicMatrix<float,blaze::rowMajor> A, B;
// ... Resizing and initialization
B = inv( A );  // Compute the inverse of A

Alternatively, an in-place inversion of a dense matrix can be performed via the invert() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization
invert( A );  // In-place matrix inversion

Both the inv() and the invert() functions will automatically select the most suited matrix inversion algorithm depending on the size and type of the given matrix. For small matrices of up to 6x6, both functions use manually optimized kernels for maximum performance. For matrices larger than 6x6 the inversion is performed by means of the most suited matrix decomposition method: In case of a general matrix the LU decomposition is used, for symmetric matrices the LDLT decomposition is applied, for Hermitian matrices the LDLH decomposition is performed, and for triangular matrices the inverse is computed via a forward or back substitution.

In case the type of the matrix does not provide additional compile time information about its structure (symmetric, lower, upper, diagonal, ...), the information can be provided manually when calling the invert() function:

using blaze::asGeneral;
using blaze::asSymmetric;
using blaze::asHermitian;
using blaze::asLower;
using blaze::asUniLower;
using blaze::asUpper;
using blaze::asUniUpper;
using blaze::asDiagonal;

invert<asGeneral>  ( A );  // In-place inversion of a general matrix
invert<asSymmetric>( A );  // In-place inversion of a symmetric matrix
invert<asHermitian>( A );  // In-place inversion of a Hermitian matrix
invert<asLower>    ( A );  // In-place inversion of a lower triangular matrix
invert<asUniLower> ( A );  // In-place inversion of a lower unitriangular matrix
invert<asUpper>    ( A );  // In-place inversion of a upper triangular matrix
invert<asUniUpper> ( A );  // In-place inversion of a upper unitriangular matrix
invert<asDiagonal> ( A );  // In-place inversion of a diagonal matrix

Alternatively, via the invert() function it is possible to explicitly specify the inversion algorithm:

using blaze::byLU;
using blaze::byLDLT;
using blaze::byLDLH;
using blaze::byLLH;

// In-place inversion of a general matrix by means of an LU decomposition
invert<byLU>( A );

// In-place inversion of a symmetric indefinite matrix by means of a Bunch-Kaufman decomposition
invert<byLDLT>( A );

// In-place inversion of a Hermitian indefinite matrix by means of a Bunch-Kaufman decomposition
invert<byLDLH>( A );

// In-place inversion of a positive definite matrix by means of a Cholesky decomposition
invert<byLLH>( A );

Whereas the inversion by means of an LU decomposition works for every general square matrix, the inversion by LDLT only works for symmetric indefinite matrices, the inversion by LDLH is restricted to Hermitian indefinite matrices and the Cholesky decomposition (LLH) only works for Hermitian positive definite matrices. Please note that it is in the responsibility of the function caller to guarantee that the selected algorithm is suited for the given matrix. In case this precondition is violated the result can be wrong and might not represent the inverse of the given matrix!

For both the inv() and invert() function the matrix inversion fails if ...

  • ... the given matrix is not a square matrix;
  • ... the given matrix is singular and not invertible.

In all failure cases either a compilation error is created if the failure can be predicted at compile time or a std::invalid_argument exception is thrown.

Note that the matrix inversion can only be used for dense matrices with float, double, complex<float> or complex<double> element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!

Also note that the functions invert the dense matrix by means of LAPACK kernels. Thus the functions can only be used if a fitting LAPACK library is available and linked to the executable. Otherwise a linker error will be created.

Furthermore, it is important to note that it is not possible to use any kind of view on the expression object returned by the inv() function. Also, it is not possible to access individual elements via the function call operator on the expression object:

row( inv( A ), 2UL );  // Compilation error: Views cannot be used on an inv() expression!
inv( A )(1,2);         // Compilation error: It is not possible to access individual elements!

Finally, please note that the inversion functions do not provide any exception safety guarantee, i.e. in case an exception is thrown the matrix may already have been modified.


Matrix Decomposition

Note that all decomposition functions can only be used for dense matrices with float, double, complex<float> or complex<double> element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!

Also note that the functions decompose a dense matrix by means of LAPACK kernels. Thus the functions can only be used if a fitting LAPACK library is available and linked to the executable. Otherwise a linker error will be created.

LU Decomposition

The LU decomposition of a dense matrix can be computed via the lu() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> L, U, P;

lu( A, L, U, P );  // LU decomposition of a row-major matrix

assert( A == L * U * P );
blaze::DynamicMatrix<double,blaze::columnMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::columnMajor> L, U, P;

lu( A, L, U, P );  // LU decomposition of a column-major matrix

assert( A == P * L * U );

The function works for both rowMajor and columnMajor matrices. Note, however, that the three matrices A, L and U are required to have the same storage order. Also, please note that the way the permutation matrix P needs to be applied differs between row-major and column-major matrices, since the algorithm uses column interchanges for row-major matrices and row interchanges for column-major matrices.

Furthermore, lu() can be used with adaptors. For instance, the following example demonstrates the LU decomposition of a symmetric matrix into a lower and upper triangular matrix:

blaze::SymmetricMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::LowerMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > L;
blaze::UpperMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > U;
blaze::DynamicMatrix<double,blaze::columnMajor> P;

lu( A, L, U, P );  // LU decomposition of A

Cholesky Decomposition

The Cholesky (LLH) decomposition of a dense matrix can be computed via the llh() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> L;

llh( A, L );  // LLH decomposition of a row-major matrix

assert( A == L * ctrans( L ) );

The function works for both rowMajor and columnMajor matrices and the two matrices A and L can have any storage order.

Furthermore, llh() can be used with adaptors. For instance, the following example demonstrates the LLH decomposition of a symmetric matrix into a lower triangular matrix:

blaze::SymmetricMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::LowerMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > L;

llh( A, L );  // Cholesky decomposition of A

QR Decomposition

The QR decomposition of a dense matrix can be computed via the qr() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::columnMajor> Q;
blaze::DynamicMatrix<double,blaze::rowMajor> R;

qr( A, Q, R );  // QR decomposition of a row-major matrix

assert( A == Q * R );

The function works for both rowMajor and columnMajor matrices and the three matrices A, Q and R can have any storage order.

Furthermore, qr() can be used with adaptors. For instance, the following example demonstrates the QR decomposition of a symmetric matrix into a general matrix and an upper triangular matrix:

blaze::SymmetricMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> Q;
blaze::UpperMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > R;

qr( A, Q, R );  // QR decomposition of A

RQ Decomposition

Similar to the QR decomposition, the RQ decomposition of a dense matrix can be computed via the rq() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> R;
blaze::DynamicMatrix<double,blaze::columnMajor> Q;

rq( A, R, Q );  // RQ decomposition of a row-major matrix

assert( A == R * Q );

The function works for both rowMajor and columnMajor matrices and the three matrices A, R and Q can have any storage order.

Also the \rq() function can be used in combination with matrix adaptors. For instance, the following example demonstrates the RQ decomposition of an Hermitian matrix into a general matrix and an upper triangular matrix:

blaze::HermitianMatrix< blaze::DynamicMatrix<complex<double>,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::UpperMatrix< blaze::DynamicMatrix<complex<double>,blaze::columnMajor> > R;
blaze::DynamicMatrix<complex<double>,blaze::rowMajor> Q;

rq( A, R, Q );  // RQ decomposition of A

QL Decomposition

The QL decomposition of a dense matrix can be computed via the ql() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> Q;
blaze::DynamicMatrix<double,blaze::columnMajor> L;

ql( A, Q, L );  // QL decomposition of a row-major matrix

assert( A == Q * L );

The function works for both rowMajor and columnMajor matrices and the three matrices A, Q and L can have any storage order.

Also the ql() function can be used in combination with matrix adaptors. For instance, the following example demonstrates the QL decomposition of a symmetric matrix into a general matrix and a lower triangular matrix:

blaze::SymmetricMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> Q;
blaze::LowerMatrix< blaze::DynamicMatrix<double,blaze::columnMajor> > L;

ql( A, Q, L );  // QL decomposition of A

LQ Decomposition

The LQ decomposition of a dense matrix can be computed via the lq() function:

blaze::DynamicMatrix<double,blaze::rowMajor> A;
// ... Resizing and initialization

blaze::DynamicMatrix<double,blaze::rowMajor> L;
blaze::DynamicMatrix<double,blaze::columnMajor> Q;

lq( A, L, Q );  // LQ decomposition of a row-major matrix

assert( A == L * Q );

The function works for both rowMajor and columnMajor matrices and the three matrices A, L and Q can have any storage order.

Furthermore, lq() can be used with adaptors. For instance, the following example demonstrates the LQ decomposition of an Hermitian matrix into a lower triangular matrix and a general matrix:

blaze::HermitianMatrix< blaze::DynamicMatrix<complex<double>,blaze::columnMajor> > A;
// ... Resizing and initialization

blaze::LowerMatrix< blaze::DynamicMatrix<complex<double>,blaze::columnMajor> > L;
blaze::DynamicMatrix<complex<double>,blaze::rowMajor> Q;

lq( A, L, Q );  // LQ decomposition of A

Eigenvalues/Eigenvectors

The eigenvalues and eigenvectors of a dense matrix can be computed via the eigen() functions:

namespace blaze {

template< typename MT, bool SO, typename VT, bool TF >
void eigen( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& w );

template< typename MT1, bool SO1, typename VT, bool TF, typename MT2, bool SO2 >
void eigen( const DenseMatrix<MT1,SO1>& A, DenseVector<VT,TF>& w, DenseMatrix<MT2,SO2>& V );

} // namespace blaze

The first function computes only the eigenvalues of the given n-by-n matrix, the second function additionally computes the eigenvectors. The eigenvalues are returned in the given vector w and the eigenvectors are returned in the given matrix V, which are both resized to the correct dimensions (if possible and necessary).

Depending on the given matrix type, the resulting eigenvalues are either of floating point or complex type: In case the given matrix is either a compile time symmetric matrix with floating point elements or an Hermitian matrix with complex elements, the resulting eigenvalues will be of floating point type and therefore the elements of the given eigenvalue vector are expected to be of floating point type. In all other cases they are expected to be of complex type. Please note that for complex eigenvalues no order of eigenvalues can be assumed, except that complex conjugate pairs of eigenvalues appear consecutively with the eigenvalue having the positive imaginary part first.

In case A is a row-major matrix, the left eigenvectors are returned in the rows of V, in case A is a column-major matrix, the right eigenvectors are returned in the columns of V. In case the given matrix is a compile time symmetric matrix with floating point elements, the resulting eigenvectors will be of floating point type and therefore the elements of the given eigenvector matrix are expected to be of floating point type. In all other cases they are expected to be of complex type.

The following examples give an impression of the computation of eigenvalues and eigenvectors for a general, a symmetric, and an Hermitian matrix:

using blaze::DynamicMatrix;
using blaze::DynamicVector;
using blaze::rowMajor;
using blaze::columnVector;

DynamicMatrix<double,rowMajor> A( 5UL, 5UL );  // The general matrix A
// ... Initialization

DynamicVector<complex<double>,columnVector> w( 5UL );   // The vector for the complex eigenvalues
DynamicMatrix<complex<double>,rowMajor> V( 5UL, 5UL );  // The matrix for the left eigenvectors

eigen( A, w, V );
using blaze::SymmetricMatrix;
using blaze::DynamicMatrix;
using blaze::DynamicVector;
using blaze::rowMajor;
using blaze::columnVector;

SymmetricMatrix< DynamicMatrix<double,rowMajor> > A( 5UL, 5UL );  // The symmetric matrix A
// ... Initialization

DynamicVector<double,columnVector> w( 5UL );       // The vector for the real eigenvalues
DynamicMatrix<double,rowMajor>     V( 5UL, 5UL );  // The matrix for the left eigenvectors

eigen( A, w, V );
using blaze::HermitianMatrix;
using blaze::DynamicMatrix;
using blaze::DynamicVector;
using blaze::rowMajor;
using blaze::columnVector;

HermitianMatrix< DynamicMatrix<complex<double>,rowMajor> > A( 5UL, 5UL );  // The Hermitian matrix A
// ... Initialization

DynamicVector<double,columnVector>      w( 5UL );       // The vector for the real eigenvalues
DynamicMatrix<complex<double>,rowMajor> V( 5UL, 5UL );  // The matrix for the left eigenvectors

eigen( A, w, V );

The functions fail if ...

  • ... the given matrix A is not a square matrix;
  • ... the given vector w is a fixed size vector and the size doesn't match;
  • ... the given matrix V is a fixed size matrix and the dimensions don't match;
  • ... the eigenvalue computation fails.

In all failure cases an exception is thrown.

Note that all eigen() functions can only be used for dense matrices with float, double, complex<float> or complex<double> element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!

Also note that the functions compute the eigenvalues and/or eigenvectors of a dense matrix by means of LAPACK kernels. Thus the functions can only be used if a fitting LAPACK library is available and linked to the executable. Otherwise a linker error will be created.


Singular Values/Singular Vectors

The singular value decomposition (SVD) of a dense matrix can be computed via the svd() functions:

namespace blaze {

template< typename MT, bool SO, typename VT, bool TF >
void svd( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& s );

template< typename MT1, bool SO, typename VT, bool TF, typename MT2, typename MT3 >
void svd( const DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U, DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V );

template< typename MT, bool SO, typename VT, bool TF, typename ST >
size_t svd( const DenseMatrix<MT,SO>& A, DenseVector<VT,TF>& s, ST low, ST upp );

template< typename MT1, bool SO, typename VT, bool TF, typename MT2, typename MT3, typename ST >
size_t svd( const DenseMatrix<MT1,SO>& A, DenseMatrix<MT2,SO>& U, DenseVector<VT,TF>& s, DenseMatrix<MT3,SO>& V, ST low, ST upp );

} // namespace blaze

The first and third function compute only singular values of the given general m-by-n matrix, the second and fourth function additionally compute singular vectors. The resulting singular values are returned in the given vector s, the left singular vectors are returned in the given matrix U, and the right singular vectors are returned in the matrix V. s, U, and V are resized to the correct dimensions (if possible and necessary).

The third and fourth function allow for the specification of a subset of singular values and/or vectors. The number of singular values and vectors to be computed is specified by the lower bound low and the upper bound upp, which either form an integral or a floating point range.

In case low and upp form are of integral type, the function computes all singular values in the index range [low..upp]. The num resulting real and non-negative singular values are stored in descending order in the given vector s, which is either resized (if possible) or expected to be a num-dimensional vector. The resulting left singular vectors are stored in the given matrix U, which is either resized (if possible) or expected to be a m-by-num matrix. The resulting right singular vectors are stored in the given matrix V, which is either resized (if possible) or expected to be a num-by-n matrix.

In case low and upp are of floating point type, the function computes all singular values in the half-open interval (low..upp]. The resulting real and non-negative singular values are stored in descending order in the given vector s, which is either resized (if possible) or expected to be a min(m,n)-dimensional vector. The resulting left singular vectors are stored in the given matrix U, which is either resized (if possible) or expected to be a m-by-min(m,n) matrix. The resulting right singular vectors are stored in the given matrix V, which is either resized (if possible) or expected to be a min(m,n)-by-n matrix.

The functions fail if ...

  • ... the given matrix U is a fixed size matrix and the dimensions don't match;
  • ... the given vector s is a fixed size vector and the size doesn't match;
  • ... the given matrix V is a fixed size matrix and the dimensions don't match;
  • ... the given scalar values don't form a proper range;
  • ... the singular value decomposition fails.

In all failure cases an exception is thrown.

Examples:

using blaze::DynamicMatrix;
using blaze::DynamicVector;
using blaze::rowMajor;
using blaze::columnVector;

DynamicMatrix<double,rowMajor>  A( 5UL, 8UL );  // The general matrix A
// ... Initialization

DynamicMatrix<double,rowMajor>     U;  // The matrix for the left singular vectors
DynamicVector<double,columnVector> s;  // The vector for the singular values
DynamicMatrix<double,rowMajor>     V;  // The matrix for the right singular vectors

svd( A, U, s, V );
using blaze::DynamicMatrix;
using blaze::DynamicVector;
using blaze::rowMajor;
using blaze::columnVector;

DynamicMatrix<complex<double>,rowMajor>  A( 5UL, 8UL );  // The general matrix A
// ... Initialization

DynamicMatrix<complex<double>,rowMajor> U;  // The matrix for the left singular vectors
DynamicVector<double,columnVector>      s;  // The vector for the singular values
DynamicMatrix<complex<double>,rowMajor> V;  // The matrix for the right singular vectors

svd( A, U, s, V, 0, 2 );

Note that all svd() functions can only be used for dense matrices with float, double, complex<float> or complex<double> element type. The attempt to call the function with matrices of any other element type or with a sparse matrix results in a compile time error!

Also note that the functions compute the singular values and/or singular vectors of a dense matrix by means of LAPACK kernels. Thus the functions can only be used if a fitting LAPACK library is available and linked to the executable. Otherwise a linker error will be created.


Previous: Matrix Types ---- Next: Adaptors

Updated