Reference for Armadillo 2.2.3
(Blue Skies Debauchery)
to Armadillo home page
to NICTA home page



[top] Preamble

 


Matrix, Vector, Cube and Field Classes
Member Functions & Variables
Other Classes
Generated Vectors/Matrices/Cubes
Functions Individually Applied to Each Element of a Matrix/Cube
Scalar Valued Functions of Vectors/Matrices/Cubes
Scalar/Vector Valued Functions of Vectors/Matrices
Vector/Matrix/Cube Valued Functions of Vectors/Matrices/Cubes
Decompositions, Inverses and Equation Solvers
Miscellaneous






Matrix, Vector, Cube and Field Classes



Mat<type>
mat
cx_mat
  • The root template matrix class is Mat<type>, where type can be one of: char, int, float, double, std::complex<double>, etc.

  • For convenience the following typedefs have been defined:
      imat  =  Mat<s32>
      umat  =  Mat<u32>
      fmat  =  Mat<float>
      mat  =  Mat<double>
      cx_fmat  =  Mat<cx_float>
      cx_mat  =  Mat<cx_double>

  • In this documentation the mat type is used for convenience; it is possible to use other types instead, eg. fmat

  • Functions which are wrappers for LAPACK or ATLAS functions (generally matrix decompositions) are only valid for the following types: fmat, mat, cx_fmat, cx_mat

  • Elements are stored with column-major ordering (ie. column by column)

  • Constructors:
    • mat()
    • mat(n_rows, n_cols)
    • mat(mat)
    • mat(vec)
    • mat(rowvec)
    • mat(string)
    • cx_mat(mat,mat)   (for constructing a complex matrix out of two real matrices)

  • The string format for the constructor is elements separated by spaces, and rows denoted by semicolons. For example, the 2x2 identity matrix can be created using the format string "1 0; 0 1". While string based initialisation is compact, directly setting the elements or using element injection is considerably faster.

  • Advanced constructors:

    • mat(aux_mem*, n_rows, n_cols, copy_aux_mem = true, strict = true)

        Create a matrix using data from writeable auxiliary memory. By default the matrix allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the matrix will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

        The strict variable comes into effect only if copy_aux_mem is set to false (ie. the matrix is directly using auxiliary memory). If strict is set to true, the matrix will be bound to the auxiliary memory for its lifetime; the number of elements in the matrix can't be changed (directly or indirectly). If strict is set to false, the matrix will not be bound to the auxiliary memory for its lifetime, ie., the size of the matrix can be changed. If the requested number of elements is different to the size of the auxiliary memory, new memory will be allocated and the auxiliary memory will no longer be used.

    • mat(const aux_mem*, n_rows, n_cols)

        Create a matrix by copying data from read-only auxiliary memory.

    • mat::fixed<n_rows, n_cols>

        Create a fixed size matrix, with the size specified via template arguments. Memory for the matrix is allocated at compile time. This is generally faster than dynamic memory allocation, but the size of the matrix can't be changed afterwards (directly or indirectly).

        For convenience, there are several pre-defined typedefs for each matrix type (where the types are: umat, imat, fmat, mat, cx_fmat, cx_mat). The typedefs specify a square matrix size, ranging from 2x2 to 9x9. The typedefs were defined by simply appending a two digit form of the size to the matrix type -- for example, mat33 is equivalent to mat::fixed<3,3>, while cx_mat44 is equivalent to cx_mat::fixed<4,4>.

    • mat::fixed<n_rows, n_cols>(aux_mem*, copy_aux_mem = true)

        Create a fixed size matrix, with the size specified via template arguments, and using data from writeable auxiliary memory. By default the matrix allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the matrix will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

    • mat::fixed<n_rows, n_cols>(const aux_mem*)

        Create a fixed size matrix, with the size specified via template arguments, and copying data from read-only auxiliary memory.


  • Examples:
      mat A = randu<mat>(5,5);
      double x = A(1,2);
      
      mat B = A + A;
      mat C = A * B;
      mat D = A % B;
      
      cx_mat X(A,B);
      
      B.zeros();
      B.set_size(10,10);
      B.zeros(5,6);
      
      //
      // fixed size matrices:
      
      mat::fixed<5,6> F;
      F.ones();
      
      mat44 G;
      G.randn();
      
      cout << mat22().randu() << endl;
      
      //
      // constructing matrices from
      // auxiliary (external) memory:
      
      double aux_mem[24];
      mat H(aux_mem, 4, 6, false);
      

  • Caveat: For mathematical correctness, scalars are treated as 1x1 matrices during initialisation. As such, the code below will not generate a 5x5 matrix with every element equal to 123.0:
      mat A(5,5);
      A = 123.0;
      
    Use the following code instead:
      mat A(5,5);
      A.fill(123.0);
      
    Or:
      mat A = 123.0 * ones<mat>(5,5);
      

  • See also:



Col<type>
colvec
vec
  • Classes for column vectors (matrices with one column)

  • The Col<type> class is derived from the Mat<type> class and inherits most of the member functions

  • For convenience the following typedefs have been defined:
      ivec, icolvec  =  Col<s32>
      uvec, ucolvec  =  Col<u32>
      fvec, fcolvec  =  Col<float>
      vec, colvec  =  Col<double>
      cx_fvec, cx_fcolvec  =  Col<cx_float>
      cx_vec, cx_colvec  =  Col<cx_double>

  • In this documentation, the vec and colvec types have the same meaning and are used interchangeably

  • In this documentation, the types vec or colvec are used for convenience; it is possible to use other types instead, eg. fvec, fcolvec

  • Functions which take Mat as input can generally also take Col as input. Main exceptions are functions which require square matrices

  • Constructors
    • vec(n_elem=0)
    • vec(vec)
    • vec(mat)   (a std::logic_error exception is thrown if the given matrix has more than one column)
    • vec(string)   (elements separated by spaces)

  • Advanced constructors:

    • vec(aux_mem*, number_of_elements, copy_aux_mem = true, strict = true)

        Create a column vector using data from writeable auxiliary memory. By default the vector allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the vector will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

        The strict variable comes into effect only if copy_aux_mem is set to false (ie. the vector is directly using auxiliary memory). If strict is set to true, the vector will be bound to the auxiliary memory for its lifetime; the number of elements in the vector can't be changed (directly or indirectly). If strict is set to false, the vector will not be bound to the auxiliary memory for its lifetime, ie., the vector's size can be changed. If the requested number of elements is different to the size of the auxiliary memory, new memory will be allocated and the auxiliary memory will no longer be used.

    • vec(const aux_mem*, number_of_elements)

        Create a column vector by copying data from read-only auxiliary memory.

    • vec::fixed<number_of_elements>

        Create a fixed size column vector, with the size specified via the template argument. Memory for the vector is allocated at compile time. This is generally faster than dynamic memory allocation, but the size of the vector can't be changed afterwards (directly or indirectly).

        For convenience, there are several pre-defined typedefs for each vector type (where the types are: uvec, ivec, fvec, vec, cx_fvec, cx_vec as well as the corresponding colvec versions). The pre-defined typedefs specify vector sizes ranging from 2 to 9. The typedefs were defined by simply appending a single digit form of the size to the vector type -- for example, vec3 is equivalent to vec::fixed<3>, while cx_vec4 is equivalent to cx_vec::fixed<4>.

    • vec::fixed<number_of_elements>(aux_mem*, copy_aux_mem = true)

        Create a fixed size column vector, with the size specified via the template argument, and using data from writeable auxiliary memory. By default the vector allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the vector will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

    • vec::fixed<number_of_elements>(const aux_mem*)

        Create a fixed size column vector, with the size specified via the template argument, and copying data from read-only auxiliary memory.


  • Examples:
      vec x(10);
      vec y = zeros<vec>(10,1);
      
      mat A = randu<mat>(10,10);
      vec z = A.col(5); // extract a column vector
      

  • Caveat: For mathematical correctness, scalars are treated as 1x1 matrices during initialisation. As such, the code below will not generate a column vector with every element equal to 123.0:
      vec q(5);
      q = 123.0;
      
    Use the following code instead:
      vec q(5);
      q.fill(123.0);
      
    Or:
      vec q = 123.0 * ones<vec>(5,1);
      

  • See also:



Row<type>
rowvec
  • Classes for row vectors (matrices with one row)

  • The template Row<type> class is derived from the Mat<type> class and inherits most of the member functions

  • For convenience the following typedefs have been defined:
      irowvec  =  Row<s32>
      urowvec  =  Row<u32>
      frowvec  =  Row<float>
      rowvec  =  Row<double>
      cx_frowvec  =  Row<cx_float>
      cx_rowvec  =  Row<cx_double>

  • In this documentation, the rowvec type is used for convenience; it is possible to use other types instead, eg. frowvec

  • Functions which take Mat as input can generally also take Row as input. Main exceptions are functions which require square matrices

  • Constructors
    • rowvec(n_elem=0)
    • rowvec(rowvec)
    • rowvec(mat)   (a std::logic_error exception is thrown if the given matrix has more than one row)
    • rowvec(string)   (elements separated by spaces)

  • Advanced constructors:

    • rowvec(aux_mem*, number_of_elements, copy_aux_mem = true, strict = true)

        Create a row vector using data from writeable auxiliary memory. By default the vector allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the vector will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

        The strict variable comes into effect only if copy_aux_mem is set to false (ie. the vector is directly using auxiliary memory). If strict is set to true, the vector will be bound to the auxiliary memory for its lifetime; the number of elements in the vector can't be changed (directly or indirectly). If strict is set to false, the vector will not be bound to the auxiliary memory for its lifetime, ie., the vector's size can be changed. If the requested number of elements is different to the size of the auxiliary memory, new memory will be allocated and the auxiliary memory will no longer be used.

    • rowvec(const aux_mem*, number_of_elements)

        Create a row vector by copying data from read-only auxiliary memory.

    • rowvec::fixed<number_of_elements>

        Create a fixed size row vector, with the size specified via the template argument. Memory for the vector is allocated at compile time. This is generally faster than dynamic memory allocation, but the size of the vector can't be changed afterwards (directly or indirectly).

        For convenience, there are several pre-defined typedefs for each vector type (where the types are: urowvec, irowvec, frowvec, rowvec, cx_frowvec, cx_rowvec). The pre-defined typedefs specify vector sizes ranging from 2 to 9. The typedefs were defined by simply appending a single digit form of the size to the vector type -- for example, rowvec3 is equivalent to rowvec::fixed<3>, while cx_rowvec4 is equivalent to cx_rowvec::fixed<4>.

    • rowvec::fixed<number_of_elements>(aux_mem*, copy_aux_mem = true)

        Create a fixed size row vector, with the size specified via the template argument, and using data from writeable auxiliary memory. By default the vector allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the vector will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

    • rowvec::fixed<number_of_elements>(const aux_mem*)

        Create a fixed size row vector, with the size specified via the template argument, and copying data from read-only auxiliary memory.


  • Examples:
      rowvec x(10);
      rowvec y = zeros<mat>(1,10);
      
      mat A = randu<mat>(10,10);
      rowvec z = A.row(5); // extract a row vector
      

  • Caveat: For mathematical correctness, scalars are treated as 1x1 matrices during initialisation. As such, the code below will not generate a row vector with every element equal to 123.0:
      rowvec r(5);
      r = 123.0;
      
    Use the following code instead:
      rowvec r(5);
      r.fill(123.0);
      
    Or:
      rowvec r = 123.0 * ones<rowvec>(1,5);
      

  • See also:



Cube<type>
cube
cx_cube
  • Classes for cubes, also known as "3D matrices"

  • The root template cube class is Cube<type>, where type can be one of: char, int, float, double, std::complex<double>, etc

  • For convenience the following typedefs have been defined:
      icube  =  Cube<s32>
      ucube  =  Cube<u32>
      fcube  =  Cube<float>
      cube  =  Cube<double>
      cx_fcube  =  Cube<cx_float>
      cx_cube  =  Cube<cx_double>

  • In this documentation the cube type is used for convenience; it is possible to use other types instead, eg. fcube

  • Cube data is stored as a set of slices (matrices) stored contiguously within memory. Within each slice, elements are stored with column-major ordering (ie. column by column)

  • Each slice can be interpreted as a matrix, hence functions which take Mat as input can generally also take cube slices as input

  • Constructors:
      cube()
      cube(cube)
      cube(n_rows, n_cols, n_slices)
      cx_cube(cube, cube) (for constructing a complex cube out of two real cubes)

  • Advanced constructors:

    • cube::fixed<n_rows, n_cols, n_slices>

        Create a fixed size cube, with the size specified via template arguments. Memory for the cube is allocated at compile time. This is generally faster than dynamic memory allocation, but the size of the cube can't be changed afterwards (directly or indirectly).

    • cube(aux_mem*, n_rows, n_cols, n_slices, copy_aux_mem = true, strict = true)

        Create a cube using data from writeable auxiliary memory. By default the cube allocates its own memory and copies data from the auxiliary memory (for safety). However, if copy_aux_mem is set to false, the cube will instead directly use the auxiliary memory (ie. no copying). This is faster, but can be dangerous unless you know what you're doing!

        The strict variable comes into effect only if copy_aux_mem is set to false (ie. the cube is directly using auxiliary memory). If strict is set to true, the cube will be bound to the auxiliary memory for its lifetime; the number of elements in the cube can't be changed (directly or indirectly). If strict is set to false, the cube will not be bound to the auxiliary memory for its lifetime, ie., the size of the cube can be changed. If the requested number of elements is different to the size of the auxiliary memory, new memory will be allocated and the auxiliary memory will no longer be used.

    • cube(const aux_mem*, n_rows, n_cols, n_slices)

        Create a cube by copying data from read-only auxiliary memory.


  • Examples:
      cube x(1,2,3);
      cube y = randu<cube>(4,5,6);
      
      mat A = y.slice(1);  // extract a slice from the cube
                           // (each slice is a matrix)
      
      mat B = randu<mat>(4,5);
      y.slice(2) = B;     // set a slice in the cube
      
      cube q = y + y;     // cube addition
      cube r = y % y;     // element-wise cube multiplication
      
      cube::fixed<4,5,6> f;
      f.ones();
      

  • Caveats

    • The size of individual slices can't be changed. For example, the following will not work:
        cube c(5,6,7);
        c.slice(0) = randu<mat>(10,20); // wrong size
        
    • For mathematical correctness, scalars are treated as 1x1x1 cubes during initialisation. As such, the code below will not generate a cube with every element equal to 123.0:
        cube c(5,6,7);
        c = 123.0;
        
      Use the following code instead:
        cube c(5,6,7);
        c.fill(123.0);
        
      Or:
        cube c = 123.0 * ones<cube>(5,6,7);
        

  • See also:



field<object type>
  • Class for one and two dimensional fields of arbitrary objects

  • Constructors (where object type is another class, eg. std::string, mat, vec, rowvec, etc):
      field<object type>(n_elem=0)
      field<object type>(n_rows, n_cols)
      field<object type>(field<object type>)

  • Examples:
      // create a field of strings
      field<std::string> S(3,2);
      
      S(0,0) = "hello";
      S(1,0) = "there";
      
      // string fields can be saved as plain text files
      S.save("string_field");
      
      // create a vec field with 3 rows and 2 columns
      field<vec> F(3,2);
      
      // access components of the field
      F(0,0) = vec(5);
      F(1,1) = randu<vec>(6);
      F(2,0).set_size(7);
      
      // access element 1 of the vec stored at 2,0
      double x = F(2,0)(1);
      
      // copy rows
      F.row(0) = F.row(2);
      
      // extract a row of vecs from F
      field<vec> G = F.row(1);
      
      // print the field to the standard output
      G.print("G =");
      
      // save the field to a binary file
      G.save("vec_field");
      

  • See also:





Member Functions & Variables



attributes
    .n_rows   (number of rows)
    .n_cols   (number of columns)
    .n_elem   (total number of elements)
    .n_slices   (number of slices)
    .n_elem_slice   (number of elements per slice)

  • Member variables which are read-only; to change the size, use .set_size(), .copy_size(), .zeros(), .ones(), or .reset()

  • n_rows, n_cols and n_elem are applicable to Mat, Col, Row, Cube and field classes

  • n_slices and n_elem_slice are applicable only to the Cube class

  • For the Col and Row classes, n_elem also indicates vector length

  • The variables are of type u32 (unsigned integer)

  • Examples:
      mat X(4,5);
      cout << "X has " << X.n_cols << " columns" << endl;
      

  • See also:



.colptr(col_number)

.copy_size(A)
  • Member function of Mat, Col, Row, Cube and field

  • Set the size to be the same as object A

  • Object A must be of the same root type as the object being modified (eg. you can't set the size of a matrix by providing a cube)

  • Examples:
      mat A = randu<mat>(5,6);
      mat B;
      B.copy_size(A);
      
      cout << B.n_rows << endl;
      cout << B.n_cols << endl;
      



.diag(k=0)
  • Member function of Mat.

  • Read/write access to the k-th diagonal in a matrix

  • The argument k is optional -- by default the main diagonal is accessed (k=0)

  • For k > 0, the k-th super-diagonal is accessed (top-right corner)

  • For k < 0, the k-th sub-diagonal is accessed (bottom-left corner)

  • An extracted diagonal is interpreted as a column vector

  • Examples:
      mat    X = randu<mat>(5,5);
      
      vec a = X.diag();
      vec b = X.diag(1);
      vec c = X.diag(-2);
      
      X.diag() = randu<vec>(5);
      

  • See also:



element/object access via (), [] and .at()
  • Provide access to individual elements or objects stored in a container object (ie., Mat, Col, Row, Cube, field)

      (n)
       
      For vec and rowvec, access the n-th element. For mat, cube and field, access the n-th element/object under the assumption of a flat layout, with column-major ordering of data (ie. column by column). A std::logic_error exception is thrown if the requested element is out of bounds. The bounds check can be optionally disabled at compile-time (see below).
           
      .at(n) and [n]

      As for (n), but without a bounds check. Not recommended for use unless your code has been thoroughly debugged.
           
      (i,j)

      For mat and field classes, access the element/object stored at the i-th row and j-th column. A std::logic_error exception is thrown if the requested element is out of bounds. The bounds check can be optionally disabled at compile-time (see below).
           
      .at(i,j)

      As for (i,j), but without a bounds check. Not recommended for use unless your code has been thoroughly debugged.
           
      (i,j,k)

      Cube only: access the element stored at the i-th row, j-th column and k-th slice. A std::logic_error exception is thrown if the requested element is out of bounds. The bounds check can be optionally disabled at compile-time (see below).
           
      .at(i,j,k)

      As for (i,j,k), but without a bounds check. Not recommended for use unless your code has been thoroughly debugged.

  • The bounds checks used by the (n), (i,j) and (i,j,k) access forms can be disabled by defining ARMA_NO_DEBUG or NDEBUG macros before including the armadillo header file (eg. #define ARMA_NO_DEBUG). Disabling the bounds checks is not recommended until your code has been thoroughly debugged -- it's better to write correct code first, and then maximise its speed.

  • Examples:
      mat A = randu<mat>(10,10);
      A(9,9) = 123.0;
      double x = A.at(9,9);
      double y = A[99];
      
      vec p = randu<vec>(10,1);
      p(9) = 123.0;
      double z = p[9];
      

  • See also:



element injection via the << stream operator
  • Initialise an instance of Mat, Col, Row or field classes via repeated use of the << operator

  • Special element endr indicates "end of row" (conceptually similar to std::endl)

  • Setting elements via << is a bit slower than directly accessing the elements, but code using << is generally more readable as well as being easier to write

  • Examples:
      mat A;
      
      A << 1 << 2 << 3 << endr
        << 4 << 5 << 6 << endr;
      

  • See also:



.eye()
.eye(n_rows, n_cols)
  • Set the elements along the main diagonal to one and off-diagonal elements set to zero, optionally first resizing to specified dimensions

  • An identity matrix is generated when n_rows = n_cols

  • Examples:
      mat A(5,5);
      A.eye();
      
      mat B;
      B.eye(5,5);
      

  • See also:



.fill(value)
  • Member function of Mat, Col, Row and Cube classes.

  • Sets the elements to a specified value. The type of value must match the type of elements used by the container object (eg. for mat the type is double)

  • Examples:
      mat A(5,5);
      A.fill(123.0);
      

  • See also:



.in_range( i )
(member of Mat, Col, Row, Cube and field)
.in_range( span(start, end) )
(member of Mat, Col, Row, Cube and field)
 
.in_range( row, col )
(member of Mat, Col, Row and field)
.in_range( span(start_row, end_row), span(start_col, end_col) )
(member of Mat, Col, Row and field)
 
.in_range( row, col, slice )
(member of Cube)
.in_range( span(start_row, end_row), span(start_col, end_col), span(start_slice, end_slice) )
(member of Cube)

  • Returns true if the given location or span is currently valid

  • Returns false if the object is empty, the location is out of bounds, or the span is out of bounds

  • Instances of span(a,b) can be replaced by:
    • span() or span::all, to indicate the entire range
    • span(a), to indicate a particular row, column or slice

  • Examples:
      mat A = randu<mat>(4,5);
      
      cout << A.in_range(0,0) << endl;  // true
      cout << A.in_range(3,4) << endl;  // true
      cout << A.in_range(4,5) << endl;  // false
      

  • See also:



.is_empty()
  • Member function of Mat, Col, Row, Cube and field classes

  • Returns true if the object has no elements

  • Returns false if the object has one or more elements

  • Examples:
      mat A = randu<mat>(5,5);
      cout << A.is_empty() << endl;
      
      A.reset();
      cout << A.is_empty() << endl;
      

  • See also:



.is_finite()
  • Member function of Mat, Col, Row and Cube classes

  • Returns true if all elements of the object are finite

  • Returns false if at least one of the elements of the object is non-finite (±infinity or NaN)

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = randu<mat>(5,5);
      
      B(1,1) = math::nan()
      
      cout << A.is_finite() << endl;
      cout << B.is_finite() << endl;
      

  • See also:



.is_square()
  • Member function of the Mat class

  • Returns true if the matrix is square, ie., number of rows is equal to the number of columns

  • Returns false if the matrix is not square

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = randu<mat>(6,7);
      
      cout << A.is_square() << endl;
      cout << B.is_square() << endl;
      

  • See also:



.is_vec()
.is_colvec()
.is_rowvec()
  • Member functions of the Mat class

  • .is_vec():
    • Returns true if the matrix can be interpreted as a vector (either column or row vector)
    • Returns false if the matrix does not have exactly one column or one row

  • .is_colvec():
    • Returns true if the matrix can be interpreted as a column vector
    • Returns false if the matrix does not have exactly one column

  • .is_rowvec():
    • Returns true if the matrix can be interpreted as a row vector
    • Returns false if the matrix does not have exactly one row

  • Caveat: do not assume that the vector has elements if these functions return true -- it is possible to have an empty vector (eg. 0x1)

  • Examples:
      mat A = randu<mat>(1,5);
      mat B = randu<mat>(5,1);
      mat C = randu<mat>(5,5);
      
      cout << A.is_vec() << endl;
      cout << B.is_vec() << endl;
      cout << C.is_vec() << endl;
      

  • See also:



.insert_rows( row_number, X )
.insert_rows( row_number, number_of_rows, set_to_zero = true )

(member functions of Mat and Col)
 
.insert_cols( col_number, X )
.insert_cols( col_number, number_of_cols, set_to_zero = true )

(member functions of Mat and Row)
 
.insert_slices( slice_number, X )
.insert_slices( slice_number, number_of_slices, set_to_zero = true )

(member functions of Cube)

  • Functions with the X argument: insert a copy of X at the specified row/column/slice
    • if inserting rows, X must have the same number of columns as the recipient object
    • if inserting columns, X must have the same number of rows as the recipient object
    • if inserting slices, X must have the same number of rows and columns as the recipient object (ie. all slices must have the same size)

  • Functions with the number_of_... argument: expand the object by creating new rows/columns/slices. By default, the new rows/columns/slices are set to zero. If set_to_zero is false, the memory used by the new rows/columns/slices will not be initialised.

  • Examples:
      mat A = randu<mat>(5,10);
      mat B = ones<mat>(5,2);
      
      // at column 2, insert a copy of B;
      // A will now have 12 columns
      A.insert_cols(2, B);
      
      // at column 1, insert 5 zeroed columns;
      // B will now have 7 columns
      B.insert_cols(1, 5);
      

  • See also:



iterators (matrices & vectors)
  • STL-style iterators and associated member functions of the Mat, Col and Row classes

  • iterator types:

      mat::iterator
      vec::iterator
      rowvec::iterator
       
      random access iterators, for read/write access to elements (which are stored column by column)
         
       
      mat::const_iterator
      vec::const_iterator
      rowvec::const_iterator
       
      random access iterators, for read-only access to elements (which are stored column by column)
         
       
      mat::col_iterator
      vec::col_iterator
      rowvec::col_iterator
       
      random access iterators, for read/write access to the elements of a specific column
         
       
      mat::const_col_iterator
      vec::const_col_iterator
      rowvec::const_col_iterator
       
      random access iterators, for read-only access to the elements of a specific column
         
       
      mat::row_iterator  
      rudimentary forward iterator, for read/write access to the elements of a specific row
         
       
      mat::const_row_iterator  
      rudimentary forward iterator, for read-only access to the elements of a specific row
         
       
      vec::row_iterator
      rowvec::row_iterator
       
      random access iterators, for read/write access to the elements of a specific row
         
       
      vec::const_row_iterator
      rowvec::const_row_iterator
       
      random access iterators, for read-only access to the elements of a specific row


  • Member functions:

      .begin()  
      iterator referring to the first element
      .end()  
      iterator referring to the past-the-end element
       
      .begin_row(row_number)  
      iterator referring to the first element of the specified row
      .end_row(row_number)  
      iterator referring to the past-the-end element of the specified row
       
      .begin_col(col_number)  
      iterator referring to the first element of the specified column
      .end_col(col_number)  
      iterator referring to the past-the-end element of the specified column


  • Examples:
      mat X = randu<mat>(5,5);
      
      
      mat::iterator a = X.begin();
      mat::iterator b = X.end();
      
      for(mat::iterator i=a; i!=b; ++i)
        {
        cout << *i << endl;
        }
      
      
      mat::col_iterator c = X.begin_col(1);  // start of column 1
      mat::col_iterator d = X.end_col(3);    // end of column 3
      
      for(mat::col_iterator i=c; i!=d; ++i)
        {
        cout << *i << endl;
        (*i) = 123.0;
        }
      

  • See also:



iterators (cubes)
  • STL-style iterators and associated member functions of the Cube class

  • iterator types:

      cube::iterator  
      random access iterator, for read/write access to elements; the elements are ordered slice by slice; the elements within each slice are ordered column by column
         
       
      cube::const_iterator  
      random access iterators, for read-only access to elements
         
       
      cube::slice_iterator  
      random access iterator, for read/write access to the elements of a particular slice; the elements are ordered column by column
         
       
      cube::const_slice_iterator  
      random access iterators, for read-only access to the elements of a particular slice


  • Member functions:

      .begin()  
      iterator referring to the first element
      .end()  
      iterator referring to the past-the-end element
       
      .begin_slice(slice_number)  
      iterator referring to the first element of the specified slice
      .end_slice(slice_number)  
      iterator referring to the past-the-end element of the specified slice


  • Examples:
      cube X = randu<cube>(2,3,4);
      
      
      cube::iterator a = X.begin();
      cube::iterator b = X.end();
      
      for(cube::iterator i=a; i!=b; ++i)
        {
        cout << *i << endl;
        }
      
      
      cube::slice_iterator c = X.begin_slice(1);  // start of slice 1
      cube::slice_iterator d = X.end_slice(2);    // end of slice 2
      
      for(cube::slice_iterator i=c; i!=d; ++i)
        {
        cout << *i << endl;
        (*i) = 123.0;
        }
      

  • See also:



.memptr()
  • Member function of Mat, Col, Row and Cube classes

  • Obtain a raw pointer to the memory used for storing elements. Not recommended for use unless you know what you're doing!

  • The function can be used for interfacing with libraries such as FFTW

  • As soon as the size of the matrix/vector/cube is changed, the pointer is no longer valid

  • Data for matrices is stored in a column-by-column order

  • Data for cubes is stored in a slice-by-slice (matrix-by-matrix) order

  • Examples:
            mat A = randu<mat>(5,5);
      const mat B = randu<mat>(5,5);
      
            double* A_mem = A.memptr();
      const double* B_mem = B.memptr();
      

  • See also:



.min()
(member functions of Mat, Col, Row, Cube)
.max()
 
 
 
.min( index_of_min_val )
(member functions of Mat, Col, Row, Cube)
.max( index_of_max_val )
 
 
 
.min( row_of_min_val, col_of_min_val )
(member functions of Mat)
.max( row_of_max_val, col_of_max_val )
 
 
 
.min( row_of_min_val, col_of_min_val, slice_of_min_val )
(member functions of Cube)
.max( row_of_max_val, col_of_max_val, slice_of_max_val )
 

  • Without arguments: return the extremum value of an object

  • With one or more arguments: return the extremum value of an object and store the location of the extremum value in the provided variable(s)

  • The provided variables must be of type u32.

  • Examples:
      vec v = randu<vec>(10);
      
      cout << "min value is " << v.min() << endl;
      
      
      u32    index;
      double min_val = v.min(index);
      
      cout << "index of min value is " << index << endl;
      
      
      mat A = randu<mat>(5,5);
      
      u32    row;
      u32    col;
      double min_val2 = A.max(row,col);
      
      cout << "max value is at " << row << ',' << col << endl;
      

  • See also:



.ones()
.ones(n_elem)
.ones(n_rows, n_cols)
.ones(n_rows, n_cols, n_slices)
  • Set the elements of an object to one, optionally first resizing to specified dimensions

  • .ones() and .ones(n_elem) are member functions of Col and Row

  • .ones() and .ones(n_rows, n_cols) are member functions of Mat

  • .ones() and .ones(n_rows, n_cols, n_slices) are member functions of Cube

  • Examples:
      mat A = randu<mat>(5,10);
      A.ones();      // sets all elements to one
      A.ones(10,20); // sets the size to 10 rows and 20 columns
                     // followed by setting all elements to one
      

  • See also:



operators:   +   -   *   /   %   ==   !=   <=   >=   <   >
  • Overloaded operators for mat, vec, rowvec and cube classes

  • Meanings:

      +    
      Addition of two objects
      -
      Subtraction of one object from another or negation of an object
      /
      Element-wise division of an object by another object or a scalar
      *
      Matrix multiplication of two objects; not applicable to the cube class unless multiplying a cube by a scalar
      %
      Schur product: element-wise multiplication of two objects
      ==
      Element-wise equality evaluation of two objects; generates a matrix of type umat with entries that indicate whether at a given position the two elements from the two objects are equal (1) or not equal (0)
      !=
      Element-wise non-equality evaluation of two objects
      >=
      As for ==, but the check is for "greater than or equal to"
      <=
      As for ==, but the check is for "less than or equal to"
      >
      As for ==, but the check is for "greater than"
      <
      As for ==, but the check is for "less than"

  • A std::logic_error exception is thrown if incompatible object sizes are used

  • If the +, - and % operators are chained, Armadillo will try to avoid the generation of temporaries; no temporaries are generated if all given objects are of the same type and size

  • If the * operator is chained, Armadillo will try to find an efficient ordering of the matrix multiplications

  • Caveat: operators involving an equality comparison (ie., ==, !=, >=, <=) may not work as expected for floating point element types (ie., float, double) due to the necessarily limited precision of these types; in other words, these operators are (in general) not recommended for matrices of type mat or fmat


  • Examples:
      mat A = randu<mat>(5,10);
      mat B = randu<mat>(5,10);
      mat C = randu<mat>(10,5);
      
      mat P = A + B;
      mat Q = A - B;
      mat R = -B;
      mat S = A / 123.0;
      mat T = A % B;
      mat U = A * C;
      
      // V is constructed without temporaries
      mat V = A + B + A + B;
      
      imat AA = "1 2 3; 4 5 6; 7 8 9;";
      imat BB = "3 2 1; 6 5 4; 9 8 7;";
      
      // compare elements
      umat ZZ = (AA >= BB);
      

  • See also:



.print(header="")
.print(stream, header="")
  • Member function of Mat, Col, Row, Cube and field

  • The first form prints the contents of an object to the std::cout stream, with an optional header line

  • The second form prints to a user specified stream

  • It's also possible to print objects using the << stream operator

  • Elements of a field can only be printed if there is an associated operator<< function defined

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = randu<mat>(6,6);
      
      A.print();
      
      // "B:" is the optional header line
      B.print("B:");
      
      cout << A << endl;
      cout << "B:" << endl << B << endl;
      

  • See also:



.print_trans(header="")
.print_trans(stream, header="")
  • Member function of Mat, Col and Row

  • Similar to .print(), with the difference being that a transposed version of the object is printed



.raw_print(header="")
.raw_print(stream, header="")
  • Member function of Mat, Col, Row and Cube

  • Similar to the .print() member function, with the difference being that no formatting of the output is done -- ie. the user can set the stream's parameters such as precision, cell width, etc.

  • If the cell width is set to zero, a space is printed between the elements

  • Examples:
      mat A = randu<mat>(5,5);
      
      cout.precision(11);
      cout.setf(ios::fixed);
      
      A.raw_print(cout, "A =");
      



.raw_print_trans(header="")
.raw_print_trans(stream, header="")
  • Member function of Mat, Col and Row

  • Similar to .raw_print(), with the difference being that a transposed version of the object is printed



.randu()
.randu(n_elem)
.randu(n_rows, n_cols)
.randu(n_rows, n_cols, n_slices)

.randn()
.randn(n_elem)
.randn(n_rows, n_cols)
.randn(n_rows, n_cols, n_slices)
  • Fill an object with random values, optionally first resizing to specified dimensions

  • .randu() uses a uniform distribution in the [0,1] interval

  • .randn() uses a normal/Gaussian distribution with zero mean and unit variance

  • To change the seed, use the std::srand() function

  • Examples:
      mat A(4,5);
      A.randu();
      
      mat B;
      B.randu(6,7);
      

  • See also:



.reset()
  • Member function of Mat, Col, Row, Cube and field

  • Causes an object to have no elements

  • Examples:
      mat A = randu<mat>(5, 5);
      A.reset();
      

  • See also:



.reshape(n_rows, n_cols, dim=0)
(member function of Mat, Col, Row)
.reshape(n_rows, n_cols, n_slices, dim=0)
(member function of Cube)

  • Reshape the object to the specified size, while preserving elements; the order of preservation of the elements can be either column-wise (dim=0) or row-wise (dim=1)

  • The new total number of elements (according to the specified size) doesn't have to be the same as the current total number of elements in the object

  • If the total number of elements in the object is less than the specified size, the remaining elements in the reshaped object are set to zero

  • If the total number of elements in the given object is greater than the specified size, only a subset of the elements is preserved

  • Caveat: .reshape() is slower than .set_size(), which doesn't preserve data

  • Examples:
      mat A = randu<mat>(4,5);
      A.reshape(5,4);
      

  • See also:



.save(name, file_type = arma_binary)
.save(stream, file_type = arma_binary)

.load(name, file_type = auto_detect)
.load(stream, file_type = auto_detect)

.quiet_save(name, file_type = arma_binary)
.quiet_save(stream, file_type = arma_binary)

.quiet_load(name, file_type = auto_detect)
.quiet_load(stream, file_type = auto_detect)
  • Member functions of Mat, Col, Row and Cube classes

  • Store/retrieve data in files or streams

  • On success, save(), load(), quiet_save(), and quite_load() will return a bool set to true

  • save() and quiet_save() will return a bool set to false if the saving process fails

  • load() and quiet_load() will return a bool set to false if the loading process fails; additionally, the object will be reset so it has no elements

  • load() and save() will print warning messages if any problems are encountered

  • quiet_load() and quiet_save() do not print any error messages

  • The following file formats are supported:

      auto_detect
      for load() and quiet_load(): try to automatically detect the file type as one of the formats described below. This is the default operation.

      raw_ascii
      Numerical data stored in raw ASCII format, without a header. The numbers are separated by whitespace. The number of columns must be the same in each row. Data which was saved in Matlab/Octave using the -ascii option can be read in Armadillo, except for complex numbers. Complex numbers are stored in standard C++ notation (a tuple surrounded by brackets: eg. (1.23,4.56) indicates 1.24 + 4.56i). Cubes are loaded as one slice.

      raw_binary
      Numerical data stored in machine dependent raw binary format, without a header. Matrices are loaded to have one column, while cubes are loaded to have one slice with one column. The .reshape() function can be used to alter the size of the loaded matrix/cube without losing data.

      arma_ascii
      Numerical data stored in human readable text format, with a simple header to speed up loading. The header indicates the type of matrix as well as the number of rows and columns. For cubes, the header additionally specifies the number of slices.

      arma_binary
      Numerical data stored in machine dependent binary format, with a simple header to speed up loading. The header indicates the type of matrix as well as the number of rows and columns. For cubes, the header additionally specifies the number of slices.

      csv_ascii
      Numerical data stored in comma separated value (CSV) text format, without a header. Applicable to Mat only.

      pgm_binary
      Image data stored in Portable Gray Map (PGM) format. Applicable to Mat only. Saving int, float or double matrices is a lossy operation, as each element is copied and converted to an 8 bit representation. As such the matrix should have values in the [0,255] interval, otherwise the resulting image may not display correctly.

      ppm_binary
      Image data stored in Portable Pixel Map (PPM) format. Applicable to Cube only. Saving int, float or double matrices is a lossy operation, as each element is copied and converted to an 8 bit representation. As such the cube/field should have values in the [0,255] interval, otherwise the resulting image may not display correctly.


  • Examples:
      mat A = randu<mat>(5,5);
      
      A.save("A1.mat");  // default save format is arma_binary
      A.save("A2.mat", arma_ascii);
      
      mat B;
      // automatically detect format type
      B.load("A1.mat");
      
      mat C;
      // force loading in the arma_ascii format
      C.load("A2.mat", arma_ascii);
      
      
      // example of saving/loading using a stream
      std::stringstream s;
      A.save(s);
      
      mat D;
      D.load(s);
      
      
      // example of testing for success
      mat E;
      bool status = E.load("A2.mat");
      
      if(status == true)
        {
        cout << "loaded okay" << endl;
        }
      else
        {
        cout << "problem with loading" << endl;
        }
      

  • See also:



.save(name, file_type = arma_binary)
.save(stream, file_type = arma_binary)

.load(name, file_type = auto_detect)
.load(stream, file_type = auto_detect)

.quiet_save(name, file_type = arma_binary)
.quiet_save(stream, file_type = arma_binary)

.quiet_load(name, file_type = auto_detect)
.quiet_load(stream, file_type = auto_detect)
  • Member functions of the field class

  • Store/retrieve fields in files or stream

  • On success, save(), load(), quiet_save(), and quite_load() will return a bool set to true

  • save() and quiet_save() will return a bool set to false if the saving process fails

  • load() and quiet_load() will return a bool set to false if the loading process fails; additionally, the field will be reset so it has no elements

  • load() and save() will print warning messages if any problems are encountered

  • quiet_load() and quiet_save() do not print any error messages

  • Fields with objects of type std::string are saved and loaded as raw text files. The text files do not have a header. Each string is separated by a whitespace. load() and quiet_load() will only accept text files that have the same number of strings on each line. The strings can have variable lengths.

  • Other than storing string fields as text files, the following file formats are supported:

      auto_detect

    • load(): try to automatically detect the field format type as one of the formats described below. This is the default operation.

    • arma_binary

    • Objects are stored in machine dependent binary format.
    • Default type for fields of type Mat, Col or Row.
    • Only applicable to fields of type Mat, Col or Row.

    • ppm_binary

    • Image data stored in Portable Pixmap Map (PPM) format.
    • Only applicable to fields of type Mat, Col or Row.
    • .load(): Loads the specified image and stores the red, green and blue components as three separate matrices. The resulting field is comprised of the three matrices, with the red, green and blue components in the first, second and third matrix, respectively.
    • .save(): Saves a field with exactly three matrices of equal size as an image. It is assumed that the red, green and blue components are stored in the first, second and third matrix, respectively. Saving int, float or double matrices is a lossy operation, as each matrix element is copied and converted to an 8 bit representation.

  • See also:



.set_imag(X)
.set_real(X)
  • Member functions of Mat, Col, Row and Cube

  • Set the imaginary/real part of an object

  • X must have the same size as the recipient object

  • Examples:
         mat A = randu<mat>(4,5);
         mat B = randu<mat>(4,5);
      
      cx_mat C = zeros<mat>(4,5);
      
      C.set_real(A);
      C.set_imag(B);
      

  • Caveat: if you want to directly construct a complex matrix out of two real matrices, the following code is faster:
         mat A = randu<mat>(4,5);
         mat B = randu<mat>(4,5);
      
      cx_mat C = cx_mat(A,B);
      

  • See also:



.set_size(n_elem)
(member function of Col, Row, and field)
.set_size(n_rows, n_cols)
(member function of Mat and field)
.set_size(n_rows, n_cols, n_slices)
(member function of Cube)

  • Changes the size of an object

  • If the requested number of elements is equal to the old number of elements, existing memory is reused

  • If the requested number of elements is not equal to the old number of elements, new memory is used

  • If you need to explicitly preserve data, use .reshape()

  • Examples:
      mat A;
      A.set_size(5,10);
      
      vec q;
      q.set_size(100);
      

  • See also:



.shed_row( row_number )
.shed_rows( first_row, last_row )

(member functions of Mat and Col)
 
.shed_col( column_number )
.shed_cols( first_column, last_column )

(member functions of Mat and Row)
 
.shed_slice( slice_number )
.shed_slices( first_slice, last_slice )

(member functions of Cube)



STL container functions
  • Mat, Col and Row classes provide the following member functions that mimic the containers in the C++ Standard Template Library:

    .clear()  
    causes an object to have no elements
       
     
    .empty()  
    returns true if the object has no elements; returns false if the object has one or more elements
       
     
    .size()  
    returns the total number of elements

  • Examples:
      mat A = randu<mat>(5,5);
      cout << A.size() << endl;
      
      A.clear();
      cout << A.empty() << endl;
      

  • See also:



submatrix views
  • A collection of member functions of Mat, Col and Row classes that provide submatrix views

  • For a matrix or vector X, the subviews are accessed as:

    • X.col( col_number )
      X( span::all, col_number )
      X( span(first_row, last_row), col_number )

      X.unsafe_col( col_number )

      X.row( row_number )
      X( row_number, span::all )
      X( row_number, span(first_col, last_col) )

      X.cols( first_col, last_col )
      X.rows( first_row, last_row )

      X.submat( first_row, first_col, last_row, last_col )
      X.submat( span(first_row, last_row), span(first_col, last_col) )

      X( span(first_row, last_row), span(first_col, last_col) )

      X.elem( vector_of_indices )

  • For a vector V (column or row vector), there is an additional method:

    • V.subvec( first_index, last_index )

  • Instances of span::all, to indicate an entire range, can be replaced by span(), where no number is specified

  • In the function X.elem(vector_of_indices), elements specified in vector_of_indices are accessed. X is interpreted as one long vector, with column-by-column ordering of the elements of X. The vector_of_indices must evaluate to be a vector of type uvec (eg., generated by find()). The aggregate set of the specified elements is treated as a column vector (eg., the output of X.elem() is always a column vector).

  • The function .unsafe_col() is provided for speed reasons and should be used only if you know what you're doing. The function creates a seemingly independent Col vector object (eg. vec), but the vector actually uses memory from the existing matrix object. As such, the created Col vector is currently not alias safe and does not take into account that the parent matrix object could be deleted. If deleted memory is accessed through the created Col vector, it will cause memory corruption and/or a crash.

  • Examples:
      mat A = zeros<mat>(5,10);
      
      A.submat(0,1,2,3) = randu<mat>(3,3);
      
      // the following three statements
      // access the same part of A
      mat B = A.submat(0,1,2,3);
      mat C = A.submat( span(0,2), span(1,3) );
      mat D = A( span(0,2), span(1,3) );
      
      // the following two statements
      // access the same part of A
      A.col(1)        = randu<mat>(5,1);
      A(span::all, 1) = randu<mat>(5,1);
      
      mat X = randu<mat>(5,5);
      
      // get all elements of X that are greater than 0.5
      vec q = X.elem( find(X > 0.5) );
      
      // set four specific elements of X to 1
      uvec indices;
      indices << 2 << 3 << 6 << 8;
      
      X.elem(indices) = ones<vec>(4);
      

  • See also:



subcube views and slices
  • A collection of member functions of the Cube class that provide subcube views

  • For a cube Q, the subviews are accessed as:

    • Q.slice( slice_number )
      Q.slices( first_slice, last_slice )

      Q.subcube( first_row, first_col, first_slice, last_row, last_col, last_slice )
      Q.subcube( span(first_row, last_row), span(first_col, last_col), span(first_slice, last_slice) )

      Q( span(first_row, last_row), span(first_col, last_col), span(first_slice, last_slice) )

  • Instances of span(a,b) can be replaced by:
    • span() or span::all, to indicate the entire range
    • span(a), to indicate a particular row, column or slice

  • An individual slice, accessed via .slice(), is an instance of the Mat class (a reference to a matrix is provided)

  • Examples:
      cube A = randu<cube>(2,3,4);
      mat B = A.slice(1);
      
      A.slice(0) = randu<mat>(2,3);
      A.slice(0)(1,2) = 99.0;
      
      A.subcube(0,0,1,  1,1,2) = randu<cube>(2,2,2);
      A( span(0,1), span(0,1), span(1,2) ) = randu<cube>(2,2,2);
      
      

  • See also:



subfield views
  • A collection of member functions of the field class that provide subfield views

  • For a field F, the subfields are accessed as:

    • F.row( row_number )
      F.col( col_number )

      F.rows( first_row, last_row )
      F.cols( first_col, last_col )

      F.subfield( first_row, first_col, last_row, last_col )
      F.subfield( span(first_row, last_row), span(first_col, last_col) )

      F( span(first_row, last_row), span(first_col, last_col) )

  • Instances of span(a,b) can be replaced by:
    • span() or span::all, to indicate the entire range
    • span(a), to indicate a particular row or column

  • See also:



.swap_rows(row1, row2)
.swap_cols(col1, col2)
  • Member functions of Mat, Col and Row classes

  • Swap the contents of specified rows or columns

  • Examples:
      mat X = randu<mat>(5,5);
      X.swap_rows(0,4);
      

  • See also:



.zeros()
.zeros(n_elem)
.zeros(n_rows, n_cols)
.zeros(n_rows, n_cols, n_slices)
  • Set the elements of an object to zero, optionally first resizing to specified dimensions

  • .zeros() and .zeros(n_elem) are member function of Col and Row

  • .zeros() and .zeros(n_rows, n_cols) are member functions of Mat

  • .zeros() and .zeros(n_rows, n_cols, n_slices) are member functions of Cube

  • Examples:
      mat A = randu<mat>(5,10);
      A.zeros();      // sets all elements to zero
      A.zeros(10,20); // sets the size to 10 rows and 20 columns
                      // followed by setting all elements to zero
      

  • See also:





Other Classes



running_stat<type>
  • Class for keeping statistics of a continuously sampled one dimensional process/signal. Useful if the storage of individual samples (scalars) is not necessary or desired. Also useful if the number of samples is not known beforehand or exceeds available memory.

  • type should be one of: float, double, cx_float, cx_double

  • Member functions:
      .operator()(scalar)  
      update the statistics so far using the given scalar
      .mean()  
      get the mean or average value so far
      .var(norm_type=0)  
      get the variance so far
      .stddev(norm_type=0)  
      get the standard deviation so far
      .min()  
      get the minimum value so far
      .max()  
      get the maximum value so far
      .reset()  
      reset all statistics and set the number of samples to zero
      .count()  
      get the number of samples so far

  • For the .var() and .stddev() functions, the default norm_type=0 performs normalisation using N-1 (where N is the number of samples so far), providing the best unbiased estimator. Using norm_type=1 causes normalisation to be done using N, which provides the second moment around the mean.

  • The return type of .count() depends on the underlying form of type: it is either float or double.

  • Examples:
      running_stat<double> stats;
      
      for(u32 i=0; i<10000; ++i)
        {
        double sample = double(rand())/RAND_MAX;
        stats(sample);
        }
      
      cout << "mean = " << stats.mean() << endl;
      cout << "var  = " << stats.var()  << endl;
      cout << "min  = " << stats.min()  << endl;
      cout << "max  = " << stats.max()  << endl;
      

  • See also:



running_stat_vec<type>(calc_cov = false)
  • Class for keeping statistics of a continuously sampled multi-dimensional process/signal. Useful if the storage of individual samples is not necessary or desired. Also useful if the number of samples is not known beforehand or exceeds available memory.

  • This class is similar to running_stat, with the difference being that vectors are processed instead of single values.

  • type must match the element type used by the sample vectors (ie. it must be one of float, double, cx_float, cx_double)

  • type can be inferred via the use of the elem_type member typedef of a vector class. For example, fvec::elem_type will be interpreted as float, while cx_vec::elem_type will be interpreted as cx_double.

  • Member functions:
      .operator()(vector)  
      update the statistics so far using the given vector
      .mean()  
      get the mean vector so far
      .var(norm_type=0)  
      get the vector of variances so far
      .stddev(norm_type=0)  
      get the vector of standard deviations so far
      .min()  
      get the vector of minimum values so far
      .max()  
      get the vector of maximum values so far
      .cov(norm_type=0)  
      get the covariance matrix so far
      NOTE: this only works if calc_cov is set to true during the construction of the class
      .reset()  
      reset all statistics and set the number of samples to zero
      .count()  
      get the number of samples so far

  • For the .var() and .stddev() functions, the default norm_type=0 performs normalisation using N-1 (where N is the number of samples so far), providing the best unbiased estimator. Using norm_type=1 causes normalisation to be done using N, which provides the second moment around the mean.

  • The return type of .count() depends on the underlying form of type: it is either float or double.

  • Examples:
      running_stat_vec<rowvec::elem_type> stats;
      
      rowvec sample;
      
      for(u32 i=0; i<10000; ++i)
        {
        sample = randu<rowvec>(5);
        stats(sample);
        }
      
      cout << "mean = " << stats.mean() << endl;
      cout << "var  = " << stats.var()  << endl;
      cout << "min  = " << stats.min()  << endl;
      cout << "max  = " << stats.max()  << endl;
      
      //
      //
      
      running_stat_vec<rowvec::elem_type> more_stats(true);
      
      for(u32 i=0; i<20; ++i)
        {
        sample = randu<rowvec>(3);
        
        sample(1) -= sample(0);
        sample(2) += sample(1);
        
        more_stats(sample);
        }
      
      cout << "covariance matrix = " << endl;
      cout << more_stats.cov() << endl;
      
      rowvec sd = more_stats.stddev();
      
      cout << "correlations = " << endl;
      cout << more_stats.cov() / (trans(sd) * sd);
      

  • See also:



wall_clock
  • Simple wall clock timer class, for measuring the number of elapsed seconds between two intervals

  • Examples:
      wall_clock timer;
      
      mat A = randu<mat>(4,4);
      mat B = randu<mat>(4,4);
      mat C;
      
      timer.tic();
      for(u32 i=0; i<100000; ++i)
        C = A + B + A + B;
      
      double n_secs = timer.toc();
      cout << "took " << n_secs << " seconds" << endl;
      





Generated Vectors/Matrices



eye(n_rows, n_cols)
  • Generate a matrix with the elements along the main diagonal set to one and off-diagonal elements set to zero

  • An identity matrix is generated when n_rows = n_cols

  • Usage:
    • matrix_type X = eye<matrix_type>(n_rows, n_cols)

  • Examples:
      mat A = eye<mat>(5,5);
      mat B = 123.0 * eye<mat>(5,5);
      

  • See also:



linspace(start, end, N=100)
  • Generate a vector with N elements; the values of the elements linearly increase from start upto (and including) end

  • Usage:
    • vector_type v = linspace<vector_type>(start, end, N)
    • matrix_type X = linspace<matrix_type>(start, end, N)

  • If a matrix_type is specified, the resultant matrix will have one column

  • Examples:
      vec v = linspace<vec>(10, 20, 5);
      mat X = linspace<mat>(10, 20, 5);
      



ones(n_elem)
ones(n_rows, n_cols)
ones(n_rows, n_cols, n_slices)
  • Generate a vector, matrix or cube with all elements set to one

  • Usage:
    • vector_type v = ones<vector_type>(n_elem)
    • matrix_type X = ones<matrix_type>(n_rows, n_cols)
    • cube_type Q = ones<cube_type>(n_rows, n_cols, n_slices)

  • Examples:
      vec  v = ones<vec>(10);
      mat  A = ones<mat>(5,6);
      cube Q = ones<cube>(5,6,7);
      
      mat  B = 123.0 * ones<mat>(5,6);
      

  • See also:



randu(n_elem)
randu(n_rows, n_cols)
randu(n_rows, n_cols, n_slices)

randn(n_elem)
randn(n_rows, n_cols)
randn(n_rows, n_cols, n_slices)
  • Generate a vector, matrix or cube with the elements set to random values

  • randu() uses a uniform distribution in the [0,1] interval

  • randn() uses a normal/Gaussian distribution with zero mean and unit variance

  • Usage:
    • vector_type v = randu<vector_type>(n_elem)
    • matrix_type X = randu<matrix_type>(n_rows, n_cols)
    • cube_type Q = randu<cube_type>(n_rows, n_cols, n_slices)

  • To change the seed, use the std::srand() function.

  • Examples:
      vec  v = randu<vec>(5);
      mat  A = randu<mat>(5,6);
      cube Q = randu<cube>(5,6,7);
      
  • See also:



repmat(A, num_copies_per_row, num_copies_per_col)
  • Generate a matrix by replicating matrix A in a block-like fashion

  • The generated matrix has the following size:
      rows = num_copies_per_row * A.n_rows
      cols = num_copies_per_col * A.n_cols

  • Examples:
      mat A = randu<mat>(2, 3);
      
      mat B = repmat(A, 4, 5);
      



toeplitz(A)
toeplitz(A,B)
circ_toeplitz(A)

zeros(n_elem)
zeros(n_rows, n_cols)
zeros(n_rows, n_cols, n_slices)
  • Generate a vector, matrix or cube with the elements set to zero

  • Usage:
    • vector_type v = zeros<vector_type>(n_elem)
    • matrix_type X = zeros<matrix_type>(n_rows, n_cols)
    • cube_type X = zeros<cube_type>(n_rows, n_cols, n_slices)

  • Examples:
      vec  v = zeros<vec>(5);
      mat  A = zeros<mat>(5,6);
      cube Q = zeros<cube>(5,6,7);
      

  • See also:
    • .zeros() (member function of Mat, Col, Row and Cube)
    • .ones() (member function of Mat, Col, Row and Cube)
    • ones()





Functions Individually Applied to Each Element of a Matrix/Cube



abs(mat)
abs(cube)
abs(cx_mat)
abs(cx_cube)
  • Obtain the magnitude of each element

  • Usage for non-complex matrices:
    • matrix_type Y = abs(X)
    • X and Y must have the same matrix_type

  • Usage for non-complex cubes:
    • cube_type Y = abs(X)
    • X and Y must have the same cube_type

  • Usage for complex matrices:
    • non_complex_matrix_type Y = abs(X)
    • X must be a have complex matrix type, eg., cx_mat or cx_fmat
    • The type of Y must be related to the type of X, eg., if X has the type cx_mat, then the type of Y must be mat

  • Usage for complex cubes:
    • non_complex_cube_type Y = abs(X)
    • X must be a have complex cube type, eg., cx_cube or cx_fcube
    • The type of Y must be related to the type of X, eg., if X has the type cx_cube, then the type of Y must be cube

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = abs(A); 
      
      cx_mat X = randu<cx_mat>(5,5);
      mat    Y = abs(X);
      



eps(X)


miscellaneous functions:
  exp, exp2, exp10, trunc_exp,
  log, log2, log10, trunc_log,
  pow, sqrt, square
  floor, ceil

  • Apply a function to each element

  • Usage:
    • matrix_type B = misc_fn(A)
    • cube_type B = misc_fn(A)
    • A and B must have the same matrix_type/cube_type
    • misc_fn(A) is one of:
        exp(A)    base-e exponential, ex
        exp2(A)    base-2 exponential, 2x
        exp10(A)    base-10 exponential, 10x
        trunc_exp(A)   base-e exponential, truncated to avoid infinity
        (only for elements with type float or double)
        log(A)    natural log, loge x
        log2(A)    base-2 log, log2 x
        log10(A)    base-10 log, log10 x
        trunc_log(A)   natural log, truncated to avoid ±infinity
        (only for elements with type float or double)
        pow(A, p)    raise to the power of p, xp
        sqrt(A)    square root, x½
        square(A)    square, x2
        floor(A)   largest integral value that is not greater than the input value
        ceil(A)   smallest integral value that is not less than the input value

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = exp(A);
      



trigonometric functions (cos, sin, tan, ...)
  • Apply a trigonometric function to each element

  • Usage:
    • matrix_type Y = trig_fn(X)
    • cube_type Y = trig_fn(X)
    • X and Y must have the same matrix_type/cube_type
    • trig_fn is one of:
      • cos family: cos, acos, cosh, acosh
      • sin family: sin, asin, sinh, asinh
      • tan family: tan, atan, tanh, atanh

  • Examples:
      mat X = randu<mat>(5,5);
      mat Y = cos(X);
      





Scalar Valued Functions of Vectors/Matrices/Cubes



accu(mat)
accu(cube)
  • Accumulate (sum) all elements

  • Examples:
      mat A = randu<mat>(5,5);
      double x = accu(A);
      
      mat B = randu<mat>(5,5);
      double y = accu(A % B);
      
      // operator % performs element-wise multiplication,
      // hence accu(A % B) is a "multiply-and-accumulate"
      // operation
      

  • See also:



as_scalar(expression)
  • Evaluate an expression that results in a 1x1 matrix, followed by converting the 1x1 matrix to a pure scalar

  • If a binary or trinary expression is given (ie. 2 or 3 terms), the function will try to exploit the fact that the result is a 1x1 matrix by using optimised expression evaluations

  • Examples:
      rowvec r = randu<rowvec>(5);
      colvec q = randu<colvec>(5);
      mat    X = randu<mat>(5,5);
      
      // examples of some expressions
      // for which optimised implementations exist
      
      double a = as_scalar(r*q);
      double b = as_scalar(r*X*q);
      double c = as_scalar(r*diagmat(X)*q);
      double d = as_scalar(r*inv(diagmat(X))*q);
      

  • See also:



det(A, slow = false)
  • Determinant of square matrix A

  • If A is not square, a std::logic_error exception is thrown

  • Caveat: for large matrices you may want to use log_det() instead

  • For matrix sizes ≤ 4x4, a fast algorithm is used by default. In rare instances, the fast algorithm might be less precise than the standard algorithm. To force the use of the standard algorithm, set the slow argument to true

  • Examples:
      mat    A = randu<mat>(5,5);
      double x = det(A);
      
      mat44  B = randu<mat>(4,4);
      
      double y = det(B);       // use fast algorithm by default
      double z = det(B, true); // use slow algorithm
      

  • See also:



dot(A, B)
cdot(A, B)
norm_dot(A, B)
  • dot(A,B): dot product of A and B, under the assumption that A and B are vectors with the same number of elements

  • cdot(A,B): as per dot(A,B), but the complex conjugate of A is used

  • norm_dot(A,B): normalised version of dot(A,B)

  • Examples:
      vec a = randu<vec>(10);
      vec b = randu<vec>(10);
      
      double x = dot(a,b);
      

  • See also:



log_det(val, sign, A)
  • Log determinant of square matrix A, such that the determinant is equal to exp(val)*sign

  • If A is not square, a std::logic_error exception is thrown

  • Examples:
      mat A = randu<mat>(5,5);
      
      double val;
      double sign;
      
      log_det(val, sign, A);
      

  • See also:



norm(X, p)
  • Compute the p-norm of X, where X can be a vector or a matrix

  • For vectors, p is an integer ≥1, or one of: "-inf", "inf", "fro"

  • For matrices, p is one of: 1, 2, "inf", "fro"

  • "-inf" is the minimum norm, "inf" is the maximum norm, while "fro" is the Frobenius norm

  • To obtain the zero norm or Hamming norm (ie. the number of non-zero elements), you may want to use this expression: accu(X != 0).

  • Examples:
      vec    q = randu<vec>(5);
      double x = norm(q, 2);
      double y = norm(q, "inf");
      

  • See also:



rank(X, tolerance = default)
  • Returns the rank of matrix X

  • Any singular values less than default tolerance are treated as zero

  • The default tolerance is max(X.n_rows, X.n_cols)*eps(sigma), where sigma is the largest singular value of X

  • The computation is based on singular value decomposition; if the decomposition fails, a std::runtime_error exception is thrown

  • Examples:
      mat A = randu<mat>(4,5);
      u32 r = rank(A);
      

  • See also:



trace(mat)
  • Sum of the diagonal elements of a square matrix

  • A std::logic_error exception is thrown if the given matrix is not square

  • Examples:
      mat A = randu<mat>(5,5);
      double x = trace(A);
      

  • See also:





Scalar/Vector Valued Functions of Vectors/Matrices



diagvec(A, k=0)
  • Extract the k-th diagonal from matrix A

  • The argument k is optional -- by default the main diagonal is extracted (k=0)

  • For k > 0, the k-th super-diagonal is extracted (top-right corner)

  • For k < 0, the k-th sub-diagonal is extracted (bottom-left corner)

  • An extracted a diagonal is interpreted as a column vector

  • Examples:
      mat A = randu<mat>(5,5);
      vec d = diagvec(A);
      

  • See also:



min(mat, dim=0)
min(rowvec)
min(colvec)

max(mat, dim=0)
max(rowvec)
max(colvec)
  • For a matrix argument, return the extremum value for each column (dim=0), or each row (dim=1)

  • For a vector argument, return the extremum value

  • Examples:
      colvec q = randu<colvec>(10,1);
      double x = max(q);
      
      mat    A = randu<mat>(10,10);
      rowvec b = max(A);
      
      // same result as max(A)
      // the 0 explicitly indicates
      // "traverse across rows"
      rowvec c = max(A,0); 
      
      // the 1 explicitly indicates
      // "traverse across columns"
      colvec d = max(A,1);
      
      // find the overall maximum value
      double y = max(max(A));
      

  • See also:



prod(mat, dim=0)
prod(rowvec)
prod(colvec)
  • For a matrix argument, return the product of elements in each column (dim=0), or each row (dim=1)

  • For a vector argument, return the product of all elements

  • Examples:
      colvec q = randu<colvec>(10,1);
      double x = prod(q);
      
      mat    A = randu<mat>(10,10);
      rowvec b = prod(A);
      
      // same result as prod(A)
      // the 0 explicitly indicates
      // "traverse across rows"
      rowvec c = prod(A,0);
      
      // the 1 explicitly indicates
      // "traverse across columns"
      colvec d = prod(A,1);
      
      // find the overall product
      double y = prod(prod(A));
      

  • See also:



sum(mat, dim=0)
sum(rowvec)
sum(colvec)
  • For a matrix argument, return the sum of elements in each column (dim=0), or each row (dim=1)

  • For a vector argument, return the sum of all elements

  • To get a sum of all the elements regardless of the argument type (ie. matrix or vector), you may wish to use accu() instead

  • Examples:
      colvec q = randu<colvec>(10,1);
      double x = sum(q);
      
      mat    A = randu<mat>(10,10);
      rowvec b = sum(A);
      
      // same result as sum(A)
      // the 0 explicitly indicates
      // "traverse across rows"
      rowvec c = sum(A,0);
      
      // the 1 explicitly indicates
      // "traverse across columns"
      colvec d = sum(A,1);
      
      // find the overall sum
      double y = sum(sum(A));
      

  • See also:



statistics: mean, median, stddev, var
    mean(mat, dim=0)
    mean(colvec)
    mean(rowvec)

      mean (average value)
    median(mat, dim=0)
    median(colvec)
    median(rowvec)

      median
    stddev(mat, norm_type=0, dim=0)
    stddev(colvec, norm_type=0)
    stddev(rowvec, norm_type=0)

      standard deviation
    var(mat, norm_type=0, dim=0)
    var(colvec, norm_type=0)
    var(rowvec, norm_type=0)

      variance

  • For a matrix argument, find a particular statistic for each column (dim=0), or each row (dim=1)

  • For a vector argument, return a particular statistic calculated using all the elements of the vector

  • For the var() and stddev() functions, the default norm_type=0 performs normalisation using N-1 (where N is the number of samples), providing the best unbiased estimator. Using norm_type=1 causes normalisation to be done using N, which provides the second moment around the mean

  • Examples:
      mat A    = randu<mat>(5,5);
      mat B    = mean(A);
      mat C    = var(A);
      double m = mean(mean(A));
      
      vec    q = randu<vec>(5);
      double v = var(q);
      

  • See also:





Vector/Matrix/Cube Valued Functions of Vectors/Matrices/Cubes



C = conv(A, B)
  • Convolution of vectors A and B.

  • If A and B are polynomial coefficient vectors, convolving them is equivalent to multiplying the two polynomials

  • The convolution operation is also equivalent to FIR filtering

  • The orientation of the result vector is the same as the orientation of A (ie. column or row vector)

  • Examples:
      vec A = randu<vec>(128) - 0.5;
      vec B = randu<vec>(128) - 0.5;
      
      vec C = conv(A,B);
      

  • See also:



conv_to<type>::from(X)
  • A form of casting

  • Convert between matrix/vector types (eg. mat to fmat), as well as cube types (eg. cube to fcube)

  • Conversion between std::vector and Armadillo matrices/vectors is also possible

  • Conversion of a mat object into colvec, rowvec or std::vector is possible if the object can be interpreted as a vector

  • Examples:
      mat  A = randu<mat>(5,5);
      fmat B = conv_to<fmat>::from(A);
      
      typedef std::vector<double> stdvec;
      
      stdvec x(3);
      x[0] = 0.0; x[1] = 1.0;  x[2] = 2.0;
      
      colvec y = conv_to< colvec >::from(x);
      stdvec z = conv_to< stdvec >::from(y); 
      

  • See also:



conj(cx_mat)
conj(cx_cube)
  • Obtain the complex conjugate of each element in a complex matrix/cube

  • Examples:
      cx_mat X = randu<cx_mat>(5,5);
      cx_mat Y = conj(X);
      

  • See also:



cor(X, Y, norm_type=0)
cor(X, norm_type=0)
  • For two matrix arguments X and Y, if each row of X and Y is an observation and each column is a variable, the (i,j)-th entry of cor(X,Y) is the correlation coefficient between the i-th variable in X and the j-th variable in Y

  • For vector arguments, the type of vector is ignored and each element in the vector is treated as an observation

  • For matrices, X and Y must have the same dimensions

  • For vectors, X and Y must have the same number of elements

  • cor(X) is equivalent to cor(X, X), also called autocorrelation

  • The default norm_type=0 performs normalisation of the correlation matrix using N-1 (where N is the number of observations). Using norm_type=1 causes normalisation to be done using N

  • Examples:
      mat X = randu<mat>(4,5);
      mat Y = randu<mat>(4,5);
      
      mat R = cor(X,Y);
      

  • See also:



cov(X, Y, norm_type=0)
cov(X, norm_type=0)
  • For two matrix arguments X and Y, if each row of X and Y is an observation and each column is a variable, the (i,j)-th entry of cov(X,Y) is the covariance between the i-th variable in X and the j-th variable in Y

  • For vector arguments, the type of vector is ignored and each element in the vector is treated as an observation

  • For matrices, X and Y must have the same dimensions

  • For vectors, X and Y must have the same number of elements

  • cov(X) is equivalent to cov(X, X)

  • The default norm_type=0 performs normalisation using N-1 (where N is the number of observations), providing the best unbiased estimation of the covariance matrix (if the observations are from a normal distribution). Using norm_type=1 causes normalisation to be done using N, which provides the second moment matrix of the observations about their mean

  • Examples:
      mat X = randu<mat>(4,5);
      mat Y = randu<mat>(4,5);
      
      mat C = cov(X,Y);
      

  • See also:



cross(A, B)


cumsum(mat, dim=0)
cumsum(rowvec)
cumsum(colvec)
  • For a matrix argument, return a matrix containing the cumulative sum of elements in each column (dim=0, the default), or each row (dim=1)

  • For a vector argument, return a vector of the same orientation, containing the cumulative sum of elements

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = cumsum(A);
      
      vec x = randu<vec>(10);
      vec y = cumsum(x);
      

  • See also:



diagmat(mat)
diagmat(rowvec)
diagmat(colvec)
  • Interpret a matrix or vector as a diagonal matrix

  • For mat, given matrix must be square; the main diagonal is copied and all other elements in the generated matrix are set to zero

  • For colvec and rowvec, elements of the vector are placed on the main diagonal in the generated matrix and all other elements are set to zero

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = diagmat(A);
      mat C = A*diagmat(A);
      
      rowvec q = randu<rowvec>(5);
      colvec r = randu<colvec>(5);
      mat    X = diagmat(q)*diagmat(r);
      

  • See also:



find(X, k=0, s="first")
  • Return a column vector of the indices of non-zero elements of X

  • The output vector must have the type uvec or umat (ie. the indices are stored as unsigned integers of type u32)

  • The input matrix X is interpreted as a vector, with column-by-column ordering of the elements of X

  • Relational operators can be used instead of X, eg. A > 0.5

  • If k=0 (default), return the indices of all non-zero elements, otherwise return at most k of their indices

  • If s="first" (default), return at most the first k indices of the non-zero elements

  • If s="last", return at most the last k indices of the non-zero elements

  • Examples:
      mat  A  = randu<mat>(5,5);
      mat  B  = randu<mat>(5,5);
      
      uvec q1 = find(A > B);
      uvec q2 = find(A > 0.5);
      uvec q3 = find(A > 0.5, 3, "last");
      

  • See also:



fliplr(mat)
flipud(mat)
  • fliplr(): generate a copy of the input matrix, with the order of the columns reversed

  • flipud(): generate a copy of the input matrix, with the order of the rows reversed

  • Examples:
      mat A = randu<mat>(5,5);
      
      mat B = fliplr(A);
      mat C = flipud(A);
      

  • See also:



imag(cx_mat)
imag(cx_cube)

real(cx_mat)
real(cx_cube)
  • Extract the imaginary/real part of a complex matrix/cube

  • Examples:
      cx_mat C = randu<cx_mat>(5,5);
      
      mat    A = imag(C);
      mat    B = real(C);
      

  • See also:



join_rows(mat A, mat B)
join_cols(mat A, mat B)
join_slices(cube A, cube B)
  • join_rows(): for two matrices A and B, append each row of B to its respective row of A; matrices A and B must have the same number of rows

  • join_cols(): for two matrices A and B, append each column of B to its respective column of A; matrices A and B must have the same number of columns

  • join_slices(): for two cubes A and B, append the slices of B to the slices of A; cubes A and B have the same number of rows and columns (ie. all slices must have the same size)

  • Examples:
      mat A = randu<mat>(4,5);
      mat B = randu<mat>(4,6);
      mat C = randu<mat>(6,5);
      
      mat X = join_rows(A,B);
      mat Y = join_cols(A,C);
      

  • See also:



kron(A,B)
  • Kronecker tensor product.

  • Using matrix A (with n rows and p columns) and matrix B (with m rows and q columns), kron(A,B) returns a matrix (with nm rows and pq columns) which denotes the tensor product of A and B

  • Examples:
      mat A = randu<mat>(4,5);
      mat B = randu<mat>(5,4);
      
      mat K = kron(A,B);
      

  • See also:



reshape(mat, n_rows, n_cols, dim=0)
reshape(cube, n_rows, n_cols, n_slices, dim=0)
  • Generate a matrix/cube sized according to given size specifications, whose elements are taken from the given matrix/cube, either column-wise (dim=0) or row-wise (dim=1)

  • The total number of elements in the generated matrix/cube doesn't have to be the same as the total number of elements in the given matrix/cube

  • If the total number of elements in the given matrix/cube is less than the specified size, the remaining elements in the generated matrix/cube are set to zero

  • If the total number of elements in the given matrix/cube is greater than the specified size, only a subset of elements is taken from the given matrix/cube

  • Examples:
      mat A = randu<mat>(10, 5);
      mat B = reshape(A, 5, 10);
      

  • See also:



shuffle(mat, dim=0)
shuffle(rowvec, dim=0)
shuffle(colvec, dim=0)

  • Shuffle the rows (dim=0) or columns (dim=1) of a matrix or vector

  • Examples:
      mat A = randu<mat>(4,5);
      mat B = shuffle(A);
      



sort(mat, sort_type=0, dim=0)
sort(rowvec, sort_type=0)
sort(colvec, sort_type=0)
  • For a matrix argument, return a matrix with the elements of the input matrix sorted in each column (dim=0), or each row (dim=1)

  • sort_type=0 (default) indicates an ascending sort

  • sort_type=1 indicates a descending sort

  • For a vector argument, return a vector which is a sorted version of the input vector

  • Examples:
      mat A = randu<mat>(10,10);
      mat B = sort(A);
      



sort_index(colvec, sort_type=0)
sort_index(rowvec, sort_type=0)
  • Return a vector which describes the sorted order of the given vector's elements (ie. it contains the indices of the given vector's elements)

  • The output vector must have the type uvec or umat (ie. the indices are stored as unsigned integers of type u32)

  • sort_type=0 (default) indicates an ascending sort

  • sort_type=1 indicates a descending sort

  • Examples:
      vec  q       = randu<vec>(10);
      uvec indices = sort_index(q);
      

  • See also:



symmatu(A)
symmatl(A)
  • symmatu(A): interpret square matrix A as symmetric, reflecting the upper triangle to the lower triangle

  • symmatl(A): interpret square matrix A as symmetric, reflecting the lower triangle to the upper triangle

  • If A is non-square, a std::logic_error exception is thrown

  • Examples:
      mat A = randu<mat>(5,5);
      
      mat B = symmatu(A);
      mat C = symmatl(A);
      

  • See also:



strans(mat)
strans(colvec)
strans(rowvec)
  • Simple matrix transpose, without taking the conjugate of the elements (complex matrices)

  • Use trans() instead, unless you explicitly need to take the transpose of a complex matrix without taking the conjugate of the elements

  • See also:



trans(mat)
trans(colvec)
trans(rowvec)
  • Matrix transpose / Hermitian transpose

  • If a given object has real elements, a normal transpose is done

  • If a given object has complex elements, a Hermitian transpose is done (ie. the conjugate of the elements is taken during the transpose operation)

  • Caveat: for complex matrices, the functionality of trans() has changed in version 2.0:
    • in version 1.2.x and earlier, trans() does not take the conjugate of complex elements
    • in version 1.2.x and earlier, the deprecated htrans() function is used for the Hermitian transpose

  • Examples:
      mat A = randu<mat>(5,10);
      mat B = trans(A);
      

  • See also:



trimatu(A)
trimatl(A)




Decompositions, Inverses and Equation Solvers



R = chol(X)
chol(R, X)
  • Cholesky decomposition of X, such that trans(R)*R = X

  • X must be a symmetric, positive-definite matrix

  • If the decomposition fails, R is reset and:
    • chol(X) throws a std::runtime_error exception
    • chol(R,X) returns a bool set to false

  • Examples:
      mat X = randu<mat>(5,5);
      mat Y = trans(X)*X;
      
      mat R = chol(Y);
      

  • See also:



vec eigval = eig_sym(mat X)
vec eigval = eig_sym(cx_mat X)

eig_sym(vec eigval, mat X)
eig_sym(vec eigval, cx_mat X)

eig_sym(vec eigval, mat eigvec, mat X)
eig_sym(vec eigval, cx_mat eigvec, cx_mat X)
  • Eigen decomposition of symmetric/hermitian matrix X

  • The eigenvalues and corresponding eigenvectors are stored in eigval and eigvec, respectively

  • The eigenvalues are in ascending order

  • If X is not square, a std::logic_error exception is thrown

  • If the decomposition fails, the output objects are reset and:
    • eig_sym(X) throws a std::runtime_error exception
    • eig_sym(eigval, X) and eig_sym(eigval, eigvec, X) return a bool set to false

  • There is currently no check whether X is symmetric

  • Examples:
      mat A = randu<mat>(10,10);
      mat B = trans(A)*A;  // generate a symmetric matrix
      
      vec eigval;
      mat eigvec;
      
      eig_sym(eigval, eigvec, B);
      

  • See also:



eig_gen(cx_vec eigval, cx_mat eigvec, mat X, side='r')
eig_gen(cx_vec eigval, cx_mat eigvec, cx_mat X, side='r')

eig_gen(cx_vec eigval, mat l_eigvec, mat r_eigvec, mat X)
eig_gen(cx_vec eigval, cx_mat l_eigvec, cx_mat r_eigvec, cx_mat X)
  • Eigen decomposition of general (non-symmetric/non-hermitian) square matrix X

  • The eigenvalues and corresponding eigenvectors are stored in eigval and eigvec, respectively

  • For the first two forms, side='r' (default) specifies that right eigenvectors are computed, while side='l' specifies that left eigenvectors are computed

  • For the last two forms, both left and right eigenvectors are computed

  • The eigenvalues are not guaranteed to be ordered

  • If X is not square, a std::logic_error exception is thrown

  • If the decomposition fails, the output objects are reset and eig_gen() returns a bool set to false

  • Examples:
      mat A = randu<mat>(10,10);
      
      cx_vec eigval;
      cx_mat eigvec;
      
      eig_gen(eigval, eigvec, A);
      

  • See also:



B = inv(A, slow = false)
inv(B, A, slow = false)
  • Inverse of square matrix A

  • If A is known to be a triangular matrix, the inverse can be computed faster by explicitly marking the matrix as triangular through trimatu() or trimatl()

  • If A is known to be a positive-definite symmetric matrix, the inverse can be computed faster by explicitly marking the matrix using sympd()

  • If A is not square, a std::logic_error exception is thrown

  • If A appears to be singular, B is reset and:
    • inv(A) throws a std::runtime_error exception
    • inv(B,A) returns a bool set to false

  • If you want to solve a system of linear equations, eg., X = inv(A)*B, the solve() function is generally more efficient

  • For matrix sizes ≤ 4x4, a fast inverse algorithm is used by default. In rare instances, the fast algorithm might be less precise than the standard algorithm. To force the use of the standard algorithm, set the slow argument to true

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = inv(A);
      
      
      // Diagonal elements in C are set to the
      // reciprocal of the corresponding elements in A.
      // Off-diagonal elements in C are set to zero.
      mat C = inv( diagmat(A) );
      
      
      // tell inv() to look only at the upper triangular part of A
      mat D = inv( trimatu(A) );
      
      
      // tell inv() that AA is a symmetric positive definite matrix
      mat AA = A*trans(A);
      mat E  = inv( sympd(AA) );
      
      
      mat44 F = randu<mat>(4,4);
      
      mat G = inv(F);       // use fast algorithm by default
      mat H = inv(F, true); // use slow algorithm
      

  • See also:



lu(mat L, mat U, mat P, mat X)
lu(mat L, mat U, mat X)
  • Lower-upper decomposition (with partial pivoting) of matrix X

  • The first form provides a lower-triangular matrix L, an upper-triangular matrix U, and a permutation matrix P, such that trans(P)*L*U = X

  • The second form provides permuted L and U, such that L*U = X. Note that in this case L is generally not lower-triangular

  • If the decomposition fails, the output objects are reset and lu() returns a bool set to false

  • Examples:
      mat A = randu<mat>(5,5);
      
      mat L, U, P;
      
      lu(L, U, P, A);
      
      mat B = trans(P)*L*U;
      

  • See also:



B = pinv(A, tolerance = default)
pinv(B, A, tolerance = default)
  • Moore-Penrose pseudo-inverse of matrix A

  • The computation is based on singular value decomposition; if the decomposition fails, B is reset and:
    • pinv(A) throws a std::runtime_error exception
    • pinv(B,A) returns a bool set to false

  • Any singular values less than tol are treated as zero

  • For matrix A with m rows and n columns, the default tolerance is max(m,n)*norm(A)*math::eps(), where math::eps() denotes the difference between 1 and the least value greater than 1 that is representable

  • Examples:
      mat A = randu<mat>(4,5);
      mat B = pinv(A);
      

  • See also:



mat coeff = princomp(mat X)
cx_mat coeff = princomp(cx_mat X)

princomp(mat coeff, mat X)
princomp(cx_mat coeff, cx_mat X)

princomp(mat coeff, mat score, mat X)
princomp(cx_mat coeff, cx_mat score, cx_mat X)

princomp(mat coeff, mat score, vec latent, mat X)
princomp(cx_mat coeff, cx_mat score, vec latent, cx_mat X)

princomp(mat coeff, mat score, vec latent, vec tsquared, mat X)
princomp(cx_mat coeff, cx_mat score, vec latent, cx_vec tsquared, cx_mat X)

  • Principal component analysis of matrix X

  • Each row of X is an observation and each column is a variable

  • output objects:
    • coeff: principal component coefficients
    • score: projected data
    • latent: eigenvalues of the covariance matrix of X
    • tsquared: Hotteling's statistic for each sample

  • The computation is based on singular value decomposition; if the decomposition fails, the output objects are reset and:
    • princomp(X) throws a std::runtime_error exception
    • remaining forms of princomp() return a bool set to false

  • Examples:
      mat A = randu<mat>(5,4);
      
      mat coeff;
      mat score;
      vec latent;
      vec tsquared;
      
      princomp(coeff, score, latent, tsquared, A);
      

  • See also:



qr(Q,R,X)
  • Decomposition of matrix X into an orthogonal (Q) and a right triangular matrix (R), such that Q*R = X

  • If the decomposition fails, Q and R are reset and the function returns a bool set to false

  • Examples:
      mat X = randu<mat>(5,5);
      mat Q, R;
      
      qr(Q,R,X);
      

  • See also:



X = solve(A, B, slow = false)
solve(X, A, B, slow = false)
  • Solve a system of linear equations, ie., A*X = B, where X is unknown

  • For a square matrix A, this function is conceptually the same as X = inv(A)*B, but is more efficient

  • Similar functionality to the "\" (left division operator) operator in Matlab/Octave, ie. X = A \ B

  • The number of rows in A and B must be the same

  • If A is known to be a triangular matrix, the solution can be computed faster by explicitly marking the matrix as triangular through trimatu() or trimatl()

  • If A is non-square (and hence also non-triangular), solve() will also try to provide approximate solutions to under-determined as well as over-determined systems

  • If no solution is found, X is reset and:
    • solve(A,B) throws a std::runtime_error exception
    • solve(X,A,B) returns a bool set to false

  • For matrix sizes ≤ 4x4, a fast algorithm is used by default. In rare instances, the fast algorithm might be less precise than the standard algorithm. To force the use of the standard algorithm, set the slow argument to true

  • NOTE: Old versions of the ATLAS library (eg. 3.6) can corrupt memory and crash your program; the standard LAPACK library and later versions of ATLAS (eg. 3.8) work without problems

  • Examples:
      mat A = randu<mat>(5,5);
      vec b = randu<vec>(5);
      mat B = randu<mat>(5,5);
      
      vec x = solve(A, b);
      mat X = solve(A, B);
      
      vec x2;
      solve(x2, A, b);
      
      // tell solve() to look only at the upper triangular part of A
      mat Y = solve( trimatu(A), B );
      
      
      mat44 C = randu<mat>(4,4);
      mat44 D = randu<mat>(4,4);
      
      mat E = solve(C, D);       // use fast algorithm by default
      mat F = solve(C, D, true); // use slow algorithm
      
      

  • See also:



vec s = svd(mat X)
vec s = svd(cx_mat X)

svd(vec s, mat X),
svd(vec s, cx_mat X)

svd(mat U, vec s, mat V, mat X)
svd(cx_mat U, vec s, cx_mat V, cx_mat X)
  • The single and two argument versions compute the singular values of X

  • The four argument version computes the full singular value decomposition of X

  • If X is square, it can be reconstructed using X = U*diagmat(s)*trans(V)

  • The singular values are in descending order

  • If the decomposition fails, the output objects are reset and:
    • svd(X) throws a std::runtime_error exception
    • svd(s,X) and svd(U,s,V,X) return a bool set to false

  • NOTE: Old versions of the ATLAS library (eg. 3.6) can corrupt memory and crash your program; the standard LAPACK library and later versions of ATLAS (eg. 3.8) work without problems

  • Examples:
      mat X = randu<mat>(5,5);
      
      mat U;
      vec s;
      mat V;
      svd(U,s,V,X);
      

  • See also:



svd_econ(mat U, vec s, mat V, mat X, mode = 'b')
svd_econ(cx_mat U, vec s, cx_mat V, cx_mat X, mode = 'b')
  • Economical singular value decomposition of X

  • mode is one of:
    • 'l': compute only left singular vectors
    • 'r': compute only right singular vectors
    • 'b': compute both left and right singular vectors (default)

  • The singular values are in descending order

  • If the decomposition fails, the output objects are reset and bool set to false is returned

  • Examples:
      mat X = randu<mat>(4,5);
      
      mat U;
      vec s;
      mat V;
      svd_econ(U, s, V, X, 'l');
      

  • See also:



X = syl(A, B, C)
syl(X, A, B, C)
  • Solve the Sylvester equation, ie., AX + XB + C = 0, where X is unknown.

  • Matrices A, B and C must be square sized.

  • If no solution is found, X is reset and:
    • syl(A,B,C) throws a std::runtime_error exception
    • svd(X,A,B,C) returns a bool set to false

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = randu<mat>(5,5);
      mat C = randu<mat>(5,5);
      
      mat X1 = syl(A, B, C);
      
      mat X2;
      syl(X2, A, B, C);
      

  • See also:





Miscellaneous



is_finite(X)
  • Returns true if all elements in X are finite

  • Returns false if at least one element in X is non-finite (±infinity or NaN)

  • X can be a scalar (eg. double), vector, matrix or cube

  • Examples:
      mat A = randu<mat>(5,5);
      mat B = randu<mat>(5,5);
      
      B(1,1) = math::nan();
      
      cout << is_finite(A) << endl;
      cout << is_finite(B) << endl;
      
      cout << is_finite( 0.123456789 ) << endl;
      cout << is_finite( math::nan() ) << endl;
      cout << is_finite( math::inf() ) << endl;
      

  • See also:



logging of warnings and errors

set_stream_err1(user_stream)
set_stream_err2(user_stream)

std::ostream& x = get_stream_err1()
std::ostream& x = get_stream_err2()
  • By default, Armadillo prints warnings and messages associated with std::logic_error and std::runtime_error exceptions to the std::cout stream

  • set_stream_err1(): change the stream for messages associated with std::logic_error exceptions (eg. out of bounds accesses)

  • set_stream_err2(): change the stream for warnings and messages associated with std::runtime_error exceptions (eg. failed decompositions)

  • get_stream_err1(): get a reference to the stream for messages associated with std::logic_error exceptions

  • get_stream_err2(): get a reference to the stream for warnings and messages associated with std::runtime_error exceptions

  • Examples:
      // print "hello" to the current err1 stream
      get_stream_err1() << "hello" << endl;
      
      // change the err2 stream to be a file
      ofstream f("my_log.txt");
      set_stream_err2(f);
      
      // trying to invert a singular matrix
      // will print a message to the err2 stream
      // and throw an exception
      mat X = zeros<mat>(5,5);
      mat Y = inv(X);
      
      // disable messages being printed to the err2 stream
      std::ostream nullstream(0);
      set_stream_err2(nullstream);
      

  • Caveat: set_stream_err1() and set_stream_err2() will not change the stream used by .print()

  • See also:


math constants (pi, e, euler, gratio, sqrt2, eps, log_min, log_max, nan, inf)
  • Collection of constants, with their precision and/or value dependant on the numerical type and/or machine used.

  • The constants are stored as static functions in the Math<type> class, where type is either float or double.

  • For convenience, Math<float> has been typedefed as fmath, while Math<double> has been typedefed as math.

  • Meaning of the constants:
      math::pi()   π, the ratio of any circle's circumference to its diameter
      math::e()   base of the natural logarithm
      math::euler()   Euler's constant, aka Euler-Mascheroni constant
      math::gratio()   golden ratio
      math::sqrt2()   square root of 2
      math::eps()   the difference between 1 and the least value greater than 1 that is representable
      math::log_min()   log of minimum non-zero value
      math::log_max()   log of maximum value
      math::nan()   “not a number” (NaN)
      math::inf()   infinity

  • Caveat: nan() is not equal to anything, even itself; if you wish to check whether a given number x is finite, use is_finite(x).

  • Examples:
      cout << "2.0 * pi = " << 2.0 * math::pi() << endl;
      
      cout << "log_max for floats = ";
      cout << fmath::log_max() << endl;
      
      cout << "log_max for doubles = ";
      cout << math::log_max() << endl;
      
  • See also:



physical constants (speed of light, etc)
  • Collection of fundamental physical constants, mainly taken from NIST and some from WolframAlpha on 2009-06-23.

  • Constants from NIST are in turn sourced from the 2006 CODATA values.

  • The constants are stored as static functions in the Phy<type> class, where type is either float or double.

  • For convenience, Phy<float> has been typedefed as fphy, while Phy<double> has been typedefed as phy.

  • Meaning of the constants:
      phy::m_u()   atomic mass constant (in kg)
      phy::N_A()   Avogadro constant
      phy::k()   Boltzmann constant (in joules per kelvin)
      phy::k_evk()   Boltzmann constant (in eV/K)
      phy::a_0()   Bohr radius (in meters)
      phy::mu_B()   Bohr magneton
      phy::Z_0()   characteristic impedance of vacuum (in ohms)
      phy::G_0()   conductance quantum (in siemens)
      phy::k_e()   Coulomb's constant (in meters per farad)
      phy::eps_0()   electric constant (in farads per meter)
      phy::m_e()   electron mass (in kg)
      phy::eV()   electron volt (in joules)
      phy::e()   elementary charge (in coulombs)
      phy::F()   Faraday constant (in coulombs)
      phy::alpha()   fine-structure constant
      phy::alpha_inv()   inverse fine-structure constant
      phy::K_J()   Josephson constant
      phy::mu_0()   magnetic constant (in henries per meter)
      phy::phi_0()   magnetic flux quantum (in webers)
      phy::R()   molar gas constant (in joules per mole kelvin)
      phy::G()   Newtonian constant of gravitation (in newton square meters per kilogram squared)
      phy::h()   Planck constant (in joule seconds)
      phy::h_bar()   Planck constant over 2 pi, aka reduced Planck constant (in joule seconds)
      phy::m_p()   proton mass (in kg)
      phy::R_inf()   Rydberg constant (in reciprocal meters)
      phy::c_0()   speed of light in vacuum (in meters per second)
      phy::sigma()   Stefan-Boltzmann constant
      phy::R_k()   von Klitzing constant (in ohms)
      phy::b()   Wien wavelength displacement law constant

  • Examples:
      cout << "speed of light = " << phy::c_0() << endl;
      

  • See also:



log_add(log_a, log_b)
  • Safe replacement for log(exp(log_a) + exp(log_b))

  • Usage:
    • scalar_type log_c = log_add(log_a, log_b)
    • scalar_type is either float or double
    • log_a, log_b and log_c must have the same type



s32, u32
  • s32 is a typedef for a signed integer with a minimum width of 32 bits; this is typically signed int, but is machine dependent

  • u32 is a typedef for an unsigned integer with a minimum width of 32 bits; this is typically unsigned int, but is machine dependent

  • See also:



cx_float, cx_double

Examples of Matlab/Octave syntax and conceptually corresponding Armadillo syntax

    Matlab/Octave   Armadillo   Notes
             
    A(1, 1)   A(0, 0)   indexing in Armadillo starts at 0
    A(k, k)   A(k-1, k-1)    
             
    size(A,1)   A.n_rows   read only
    size(A,2)   A.n_cols    
    size(Q,3)   Q.n_slices   Q is a cube (3D array)
    numel(A)   A.n_elem    
             
    A(:, k)   A.col(k)   this is a conceptual example only; exact conversion from Matlab/Octave to Armadillo syntax will require taking into account that indexing starts at 0
    A(k, :)   A.row(k)    
    A(:, p:q)   A.cols(p, q)    
    A(p:q, :)   A.rows(p, q)    
             
    A(p:q, r:s)   A.submat(p, r, q, s)   A.submat(first_row, first_col, last_row, last_col)
        or    
        A( span(p,q), span(r,s) )   A( span(first_row, last_row), span(first_col, last_col) )
             
    Q(:, :, k)   Q.slice(k)   Q is a cube (3D array)
    Q(:, :, t:u)   Q.slices(t, u)    
             
    Q(p:q, r:s, t:u)   Q.subcube(p, r, t, q, s, u)   .subcube(first_row, first_col, first_slice, last_row, last_col, last_slice)
        or    
        Q( span(p,q), span(r,s), span(t,u) )    
             
    A'   trans(A)   matrix transpose / Hermitian transpose
    (for complex matrices, the conjugate of each element is taken)
    A.'   strans(A)   simple matrix transpose
    (for complex matrices, the conjugate of each element is not taken)
             
    A = zeros(size(A))   A.zeros()    
    A = ones(size(A))   A.ones()    
    A = zeros(k)   A = zeros<mat>(k,k)    
    A = ones(k)   A = ones<mat>(k,k)    
             
    C = complex(A,B)   cx_mat C = cx_mat(A,B)    
             
    A .* B   A % B   element-wise multiplication
    A ./ B   A / B   element-wise division
    A \ B   solve(A,B)   conceptually similar to inv(A)*B, but more efficient
    A = A + 1;   A++    
    A = A - 1;   A--    
             
    A = [ 1 2; 3 4; ]   << 1 << 2 << endr
       << 3 << 4 << endr;
      initialisation via element injection; special element endr indicates end of row
             
    X = [ A  B ]   X = join_rows(A,B)    
    X = [ A; B ]   X = join_cols(A,B)    
             
    A   cout << A << endl;
    or
    A.print("A =");
       
             
    save -ascii 'A.dat' A   A.save("A.dat", raw_ascii);   Matlab/Octave matrices saved as ascii are readable by Armadillo (and vice-versa)
    load -ascii 'A.dat'   A.load("A.dat", raw_ascii);    
             
    S = { 'abc'; 'def' }   field<std::string> S(2);
    S(0) = "abc";
    S(1) = "def";
      fields can store arbitrary objects, in a 1D or 2D layout



example program

  • If you save the program below as example.cpp, under Linux you can compile it using:
    g++ example.cpp -o example -O1 -larmadillo
    • #include <iostream>
      #include <armadillo>
      
      using namespace std;
      using namespace arma;
      
      int main(int argc, char** argv)
        {
        mat A = randu<mat>(4,5);
        mat B = randu<mat>(4,5);
        
        cout << A*trans(B) << endl;
        
        return 0;
        }
      
  • You may also want to have a look at the example programs that come with the Armadillo archive.

  • As Armadillo is a template library, we strongly recommended to have optimisation enabled when compiling programs (eg. when compiling with GCC, use the -O1 or -O2 options).



API Changes & Additions

Armadillo's version number is X.Y.Z, where X is a major version, Y is a minor version, and Z is the patch level (indicating bug fixes).

Within each major version (eg. 1.x), minor versions with an even number (eg. 1.2) are backwards compatible with earlier even minor versions (eg. 1.0). For example, code written for version 1.0 will work with version 1.2. However, as each minor version may have more features (ie. API extensions) than earlier versions, code written for version 1.2 doesn't necessarily work with 1.0.

An odd minor version number (eg. 0.9, 1.3, 1.9) indicates an experimental version. Experimental versions are generally faster and have more functionality, but their APIs have not been finalised yet.

In general, we don't like changes to existing APIs and prefer not to break any user software. However, to allow evolution, we reserve the right to change the APIs in future major versions of Armadillo, while remaining backwards compatible wherever possible (eg. 2.x may have different APIs than 1.x). Also, in a rare instance the user API may need to be altered if a bug fix absolutely requires it.

Below is a list of additions and changes since prior versions of Armadillo:
  • Added in 2.2:
  • Changed in 2.0:
  • Added in 2.0:
  • Added in 1.2:
    • .min() & .max() member functions of Mat and Cube
    • floor() and ceil()
    • representation of “not a number”: math::nan()
    • representation of infinity: math::inf()
    • standalone is_finite()
    • .in_range() can use span() arguments
    • fixed size matrices and vectors can use auxiliary (external) memory
    • submatrices and subfields can be accessed via X( span(a,b)span(c,d) )
    • subcubes can be accessed via X( span(a,b)span(c,d)span(e,f) )
    • the two argument version of span can be replaced by span::all or span(), to indicate an entire range
    • for cubes, the two argument version of span can be replaced by a single argument version, span(a), to indicate a single column, row or slice
    • arbitrary "flat" subcubes can be interpreted as matrices; for example:
        cube Q = randu<cube>(5,3,4);
        mat  A = Q( span(1), span(1,2), span::all );
        // A has a size of 2x4
        
        vec v = ones<vec>(4);
        Q( span(1), span(1), span::all ) = v;
        
    • interpretation of matrices as triangular through trimatu() / trimatl()
    • explicit handling of triangular matrices by solve() and inv()
    • extended syntax for submatrices, including access to elements whose indices are specified in a vector
    • ability to change the stream used for logging of errors and warnings
    • ability to save/load matrices in raw binary format
    • cumulative sum function: cumsum()

  • Changed in 1.0 (compared to earlier 0.x development versions):
    • the 3 argument version of lu(), eg. lu(L,U,X), provides L and U which should be the same as produced by Octave 3.2 (this was not the case in versions prior to 0.9.90)

    • rand() has been replaced by randu(); this has been done to avoid confusion with std::rand(), which generates random numbers in a different interval

    • In versions earlier than 0.9.0, some multiplication operations directly converted result matrices with a size of 1x1 into scalars. This is no longer the case. If you know the result of an expression will be a 1x1 matrix and wish to treat it as a pure scalar, use the as_scalar() wrapping function

    • Almost all functions have been placed in the delayed operations framework (for speed purposes). This may affect code which assumed that the output of some functions was a pure matrix. The solution is easy, as explained below.

      In general, Armadillo queues operations before executing them. As such, the direct output of an operation or function cannot be assumed to be a directly accessible matrix. The queued operations are executed when the output needs to be stored in a matrix, eg. mat B = trans(A) or mat B(trans(A)). If you need to force the execution of the delayed operations, place the operation or function inside the corresponding Mat constructor. For example, if your code assumed that the output of some functions was a pure matrix, eg. chol(m).diag(), change the code to mat(chol(m)).diag(). Similarly, if you need to pass the result of an operation such as A+B to one of your own functions, use my_function( mat(A+B) ).