|
|
|
Mat<type>, mat, cx_mat | dense matrix class | |
Col<type>, colvec, vec | dense column vector class | |
Row<type>, rowvec | dense row vector class | |
Cube<type>, cube, cx_cube | dense cube class ("3D matrix") | |
field<object type> | class for storing arbitrary objects in matrix-like or cube-like layouts | |
SpMat<type>, sp_mat, sp_cx_mat | sparse matrix class | |
operators | + − * % / == != <= >= < > && || |
attributes | .n_rows, .n_cols, .n_elem, .n_slices, ... | |
element access | element/object access via (), [] and .at() | |
element initialisation | set elements via initialiser lists | |
.zeros | set all elements to zero | |
.ones | set all elements to one | |
.eye | set elements along main diagonal to one and off-diagonal elements to zero | |
.randu / .randn | set all elements to random values | |
.fill | set all elements to specified value | |
.imbue | imbue (fill) with values provided by functor or lambda function | |
.clean | replace elements below a threshold with zeros | |
.replace | replace specific elements with a new value | |
.clamp | clamp values to lower and upper limits | |
.transform | transform each element via functor or lambda function | |
.for_each | apply a functor or lambda function to each element | |
.set_size | change size without keeping elements (fast) | |
.reshape | change size while keeping elements | |
.resize | change size while keeping elements and preserving layout | |
.copy_size | change size to be same as given object | |
.reset | change size to empty | |
submatrix views | read/write access to contiguous and non-contiguous submatrices | |
subcube views | read/write access to contiguous and non-contiguous subcubes | |
subfield views | read/write access to contiguous subfields | |
.diag | read/write access to matrix diagonals | |
.each_col / .each_row | vector operations applied to each column/row of matrix (aka "broadcasting") | |
.each_slice | matrix operations applied to each slice of cube (aka "broadcasting") | |
.set_imag / .set_real | set imaginary/real part | |
.insert_rows / cols / slices | insert vector/matrix/cube at specified row/column/slice | |
.shed_rows / cols / slices | remove specified rows/columns/slices | |
.swap_rows / cols | swap specified rows or columns | |
.swap | swap contents with given object | |
.memptr | raw pointer to memory | |
.colptr | raw pointer to memory used by specified column | |
iterators (matrices) | iterators and associated member functions for dense matrices and vectors | |
iterators (cubes) | iterators and associated member functions for cubes | |
iterators (sparse matrices) | iterators and associated member functions for sparse matrices | |
iterators (submatrices) | iterators and associated member functions for submatrices & subcubes | |
compat. container functions | compatibility container functions | |
.as_col / .as_row | return flattened matrix as column or row vector | |
.col_as_mat / .row_as_mat | return matrix representation of cube column or cube row | |
.as_dense | return dense vector/matrix representation of sparse matrix expression | |
.t / .st | return matrix transpose | |
.i | return inverse of square matrix | |
.min / .max | return extremum value | |
.index_min / .index_max | return index of extremum value | |
.eval | force evaluation of delayed expression | |
.in_range | check whether given location or span is valid | |
.is_empty | check whether object is empty | |
.is_vec | check whether matrix is a vector | |
.is_sorted | check whether vector or matrix is sorted | |
.is_trimatu / .is_trimatl | check whether matrix is upper/lower triangular | |
.is_diagmat | check whether matrix is diagonal | |
.is_square | check whether matrix is square sized | |
.is_symmetric | check whether matrix is symmetric | |
.is_hermitian | check whether matrix is hermitian | |
.is_sympd | check whether matrix is symmetric/hermitian positive definite | |
.is_zero | check whether all elements are zero | |
.is_finite | check whether all elements are finite | |
.has_inf | check whether any element is ±infinity | |
.has_nan | check whether any element is NaN | |
print object to std::cout or user specified stream | ||
.raw_print | print object without formatting | |
.brief_print | print object in shortened/abridged form | |
.save/.load (matrices & cubes) | save/load matrices and cubes in files or streams | |
.save/.load (fields) | save/load fields in files or streams |
linspace | generate vector with linearly spaced elements | |
logspace | generate vector with logarithmically spaced elements | |
regspace | generate vector with regularly spaced elements | |
randperm | generate vector with random permutation of a sequence of integers | |
eye | generate identity matrix | |
ones | generate object filled with ones | |
zeros | generate object filled with zeros | |
randu | generate object with random values (uniform distribution) | |
randn | generate object with random values (normal distribution) | |
randg | generate object with random values (gamma distribution) | |
randi | generate object with random integer values in specified interval | |
speye | generate sparse identity matrix | |
spones | generate sparse matrix with non-zero elements set to one | |
sprandu / sprandn | generate sparse matrix with non-zero elements set to random values | |
toeplitz | generate Toeplitz matrix |
abs | obtain magnitude of each element | |
accu | accumulate (sum) all elements | |
affmul | affine matrix multiplication | |
all | check whether all elements are non-zero, or satisfy a relational condition | |
any | check whether any element is non-zero, or satisfies a relational condition | |
approx_equal | approximate equality | |
arg | phase angle of each element | |
as_scalar | convert 1x1 matrix to pure scalar | |
clamp | obtain clamped elements according to given limits | |
cond | condition number of matrix | |
conj | obtain complex conjugate of each element | |
conv_to | convert/cast between matrix types | |
cross | cross product | |
cumsum | cumulative sum | |
cumprod | cumulative product | |
det | determinant | |
diagmat | generate diagonal matrix from given matrix or vector | |
diagvec | extract specified diagonal | |
diags / spdiags | generate band matrix from given set of vectors | |
diff | differences between adjacent elements | |
dot / cdot / norm_dot | dot product | |
eps | obtain distance of each element to next largest floating point representation | |
expmat | matrix exponential | |
expmat_sym | matrix exponential of symmetric matrix | |
find | find indices of non-zero elements, or elements satisfying a relational condition | |
find_finite | find indices of finite elements | |
find_nonfinite | find indices of non-finite elements | |
find_nan | find indices of NaN elements | |
find_unique | find indices of unique elements | |
fliplr / flipud | flip matrix left to right or upside down | |
imag / real | extract imaginary/real part | |
ind2sub | convert linear index to subscripts | |
index_min / index_max | indices of extremum values | |
inplace_trans | in-place transpose | |
intersect | find common elements in two vectors/matrices | |
join_rows / join_cols | concatenation of matrices | |
join_slices | concatenation of cubes | |
kron | Kronecker tensor product | |
log_det | log determinant | |
log_det_sympd | log determinant of symmetric positive definite matrix | |
logmat | matrix logarithm | |
logmat_sympd | matrix logarithm of symmetric matrix | |
min / max | return extremum values | |
nonzeros | return non-zero values | |
norm | various norms of vectors and matrices | |
norm2est | fast estimate of the matrix 2-norm | |
normalise | normalise vectors to unit p-norm | |
pow | element-wise power | |
powmat | matrix power | |
prod | product of elements | |
rank | rank of matrix | |
rcond | reciprocal condition number | |
repelem | replicate elements | |
repmat | replicate matrix in block-like fashion | |
reshape | change size while keeping elements | |
resize | change size while keeping elements and preserving layout | |
reverse | reverse order of elements | |
roots | roots of polynomial | |
shift | shift elements | |
shuffle | randomly shuffle elements | |
size | obtain dimensions of given object | |
sort | sort elements | |
sort_index | vector describing sorted order of elements | |
sqrtmat | square root of matrix | |
sqrtmat_sympd | square root of symmetric matrix | |
sum | sum of elements | |
sub2ind | convert subscripts to linear index | |
symmatu / symmatl | generate symmetric matrix from given matrix | |
trace | sum of diagonal elements | |
trans | transpose of matrix | |
trapz | trapezoidal numerical integration | |
trimatu / trimatl | copy upper/lower triangular part | |
trimatu_ind / trimatl_ind | obtain indices of upper/lower triangular part | |
unique | return unique elements | |
vecnorm | obtain vector norm of each row or column of a matrix | |
vectorise | flatten matrix into vector | |
misc functions | miscellaneous element-wise functions: exp, log, sqrt, round, sign, ... | |
trig functions | trigonometric element-wise functions: cos, sin, tan, ... |
chol | Cholesky decomposition | |
eig_sym | eigen decomposition of dense symmetric/hermitian matrix | |
eig_gen | eigen decomposition of dense general square matrix | |
eig_pair | eigen decomposition for pair of general dense square matrices | |
hess | upper Hessenberg decomposition | |
inv | inverse of general square matrix | |
inv_sympd | inverse of symmetric positive definite matrix | |
lu | lower-upper decomposition | |
null | orthonormal basis of null space | |
orth | orthonormal basis of range space | |
pinv | pseudo-inverse / generalised inverse | |
qr | QR decomposition | |
qr_econ | economical QR decomposition | |
qz | generalised Schur decomposition | |
schur | Schur decomposition | |
solve | solve systems of linear equations | |
svd | singular value decomposition | |
svd_econ | economical singular value decomposition | |
syl | Sylvester equation solver |
eigs_sym | limited number of eigenvalues & eigenvectors of sparse symmetric real matrix | |
eigs_gen | limited number of eigenvalues & eigenvectors of sparse general square matrix | |
svds | truncated svd: limited number of singular values & singular vectors of sparse matrix | |
spsolve | solve sparse systems of linear equations | |
spsolve_factoriser | factoriser for solving sparse systems of linear equations |
conv | 1D convolution | |
conv2 | 2D convolution | |
fft / ifft | 1D fast Fourier transform and its inverse | |
fft2 / ifft2 | 2D fast Fourier transform and its inverse | |
interp1 | 1D interpolation | |
interp2 | 2D interpolation | |
polyfit | find polynomial coefficients for data fitting | |
polyval | evaluate polynomial |
stats functions | mean, median, standard deviation, variance | |
cov | covariance | |
cor | correlation | |
hist | histogram of counts | |
histc | histogram of counts with user specified edges | |
quantile | quantiles of a dataset | |
princomp | principal component analysis (PCA) | |
normpdf | probability density function of normal distribution | |
log_normpdf | logarithm version of probability density function of normal distribution | |
normcdf | cumulative distribution function of normal distribution | |
mvnrnd | random vectors from multivariate normal distribution | |
chi2rnd | random numbers from chi-squared distribution | |
wishrnd | random matrix from Wishart distribution | |
iwishrnd | random matrix from inverse Wishart distribution | |
running_stat | running statistics of scalars (one dimensional process/signal) | |
running_stat_vec | running statistics of vectors (multi-dimensional process/signal) | |
kmeans | cluster data into disjoint sets | |
gmm_diag/gmm_full | probabilistic clustering and likelihood calculation via mixture of Gaussians |
constants | pi, inf, NaN, eps, speed of light, ... | |
wall_clock | timer for measuring number of elapsed seconds | |
RNG seed setting | functions for changing RNG seeds | |
output streams | streams for printing warnings and errors | |
uword / sword | shorthand for unsigned and signed integers | |
cx_double / cx_float | shorthand for std::complex<double> and std::complex<float> | |
Matlab/Armadillo syntax differences | examples of Matlab syntax and conceptually corresponding Armadillo syntax | |
example program | short example program | |
config.hpp | configuration options | |
API additions | API stability and list of API additions |
mat() | ||
mat(n_rows, n_cols) | ||
mat(n_rows, n_cols, fill_form) | (elements are initialised according to fill_form) | |
mat(size(X)) | ||
mat(size(X), fill_form) | (elements are initialised according to fill_form) | |
mat(mat) | ||
mat(vec) | ||
mat(rowvec) | ||
mat(initializer_list) | ||
mat(string) | ||
mat(std::vector) | (treated as a column vector) | |
mat(sp_mat) | (for converting a sparse matrix to a dense matrix) | |
cx_mat(mat,mat) | (for constructing a complex matrix out of two real matrices) |
fill::zeros | ↦ | set all elements to 0 (default in Armadillo >= 10.5) |
fill::ones | ↦ | set all elements to 1 |
fill::eye | ↦ | set the elements on the main diagonal to 1 and off-diagonal elements to 0 |
fill::randu | ↦ | set all elements to random values from a uniform distribution in the [0,1] interval |
fill::randn | ↦ | set all elements to random values from a normal/Gaussian distribution with zero mean and unit variance |
fill::value(scalar) | ↦ | set all elements to specified scalar |
fill::none | ↦ | do not initialise the elements (matrix may have garbage values) |
fill::zeros
fill::none
"1 0; 0 1"
;
note that string based initialisation is slower than directly setting the elements or using element initialisation
mat(ptr_aux_mem, n_rows, n_cols, copy_aux_mem = true, strict = false)
mat(const ptr_aux_mem, n_rows, n_cols)
mat::fixed<n_rows, n_cols>
mat::fixed<n_rows, n_cols>(fill_form)
mat::fixed<n_rows, n_cols>(const ptr_aux_mem)
mat A(5, 5, fill::randu); 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.ones(5,6); B.print("B:"); mat::fixed<5,6> F; double aux_mem[24]; mat H(&aux_mem[0], 4, 6, false); // use auxiliary memory
vec() | ||
vec(n_elem) | ||
vec(n_elem, fill_form) | (elements are initialised according to fill_form) | |
vec(size(X)) | ||
vec(size(X), fill_form) | (elements are initialised according to fill_form) | |
vec(vec) | ||
vec(mat) | (std::logic_error exception is thrown if the given matrix has more than one column) | |
vec(initializer_list) | ||
vec(string) | (elements separated by spaces) | |
vec(std::vector) | ||
cx_vec(vec,vec) | (for constructing a complex vector out of two real vectors) |
fill::zeros
fill::none
vec(ptr_aux_mem, number_of_elements, copy_aux_mem = true, strict = false)
vec(const ptr_aux_mem, number_of_elements)
vec::fixed<number_of_elements>
vec::fixed<number_of_elements>(fill_form)
vec::fixed<number_of_elements>(const ptr_aux_mem)
vec x(10); vec y(10, fill::ones); mat A(10, 10, fill::randu); vec z = A.col(5); // extract a column vector
rowvec() | ||
rowvec(n_elem) | ||
rowvec(n_elem, fill_form) | (elements are initialised according to fill_form) | |
rowvec(size(X)) | ||
rowvec(size(X), fill_form) | (elements are initialised according to fill_form) | |
rowvec(rowvec) | ||
rowvec(mat) | (std::logic_error exception is thrown if the given matrix has more than one row) | |
rowvec(initializer_list) | ||
rowvec(string) | (elements separated by spaces) | |
rowvec(std::vector) | ||
cx_rowvec(rowvec,rowvec) | (for constructing a complex row vector out of two real row vectors) |
fill::zeros
fill::none
rowvec(ptr_aux_mem, number_of_elements, copy_aux_mem = true, strict = false)
rowvec(const ptr_aux_mem, number_of_elements)
rowvec::fixed<number_of_elements>
rowvec::fixed<number_of_elements>(fill_form)
rowvec::fixed<number_of_elements>(const ptr_aux_mem)
rowvec x(10); rowvec y(10, fill::ones); mat A(10, 10, fill::randu); rowvec z = A.row(5); // extract a row vector
cube(n_rows, n_cols, n_slices) | ||
cube(n_rows, n_cols, n_slices, fill_form) | (elements are initialised according to fill_form) | |
cube(size(X)) | ||
cube(size(X), fill_form) | (elements are initialised according to fill_form) | |
cube(cube) | ||
cx_cube(cube, cube) | (for constructing a complex cube out of two real cubes) |
fill::zeros | ↦ | set all elements to 0 (default in Armadillo >= 10.5) |
fill::ones | ↦ | set all elements to 1 |
fill::randu | ↦ | set all elements to random values from a uniform distribution in the [0,1] interval |
fill::randn | ↦ | set all elements to random values from a normal/Gaussian distribution with zero mean and unit variance |
fill::value(scalar) | ↦ | set all elements to specified scalar |
fill::none | ↦ | do not initialise the elements (cube may have garbage values) |
fill::zeros
fill::none
cube::fixed<n_rows, n_cols, n_slices>
cube(ptr_aux_mem, n_rows, n_cols, n_slices, copy_aux_mem = true, strict = false)
cube(const ptr_aux_mem, n_rows, n_cols, n_slices)
cube x(1, 2, 3); cube y(4, 5, 6, fill::randu); mat A = y.slice(1); // extract a slice from the cube // (each slice is a matrix) mat B(4, 5, fill::randu); 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();
cube c(5,6,7); c.slice(0) = randu<mat>(10,20); // wrong size
field<object_type>() |
field<object_type>(n_elem) |
field<object_type>(n_rows, n_cols) |
field<object_type>(n_rows, n_cols, n_slices) |
field<object_type>(size(X)) |
field<object_type>(field<object_type>) |
mat A = randn(2,3); mat B = randn(4,5); field<mat> F(2,1); F(0,0) = A; F(1,0) = B; F.print("F:"); F.save("mat_field");
sp_mat() | ||
sp_mat(n_rows, n_cols) | ||
sp_mat(size(X)) | ||
sp_mat(sp_mat) | ||
sp_mat(mat) | (for converting a dense matrix to a sparse matrix) | |
sp_cx_mat(sp_mat,sp_mat) | (for constructing a complex matrix out of two real matrices) |
sp_mat(locations, values, sort_locations = true)
sp_mat(locations, values, n_rows, n_cols, sort_locations = true, check_for_zeros = true)
sp_mat(add_values, locations, values, n_rows, n_cols, sort_locations = true, check_for_zeros = true)
sp_mat(rowind, colptr, values, n_rows, n_cols, check_for_zeros = true)
sp_mat A = sprandu(1000, 2000, 0.01); sp_mat B = sprandu(2000, 1000, 0.01); sp_mat C = 2*B; sp_mat D = A*C; sp_mat E(1000,1000); E(1,2) = 123; // batch insertion of 3 values at // locations (1, 2), (7, 8), (9, 9) umat locations = { { 1, 7, 9 }, { 2, 8, 9 } }; vec values = { 1.0, 2.0, 3.0 }; sp_mat X(locations, values);
+ − * % / == != <= >= < > && ||
==
, !=
, >=
, <=
, >
, <
, &&
, ||
)
each element in the generated object is either 0 or 1, depending on the result of the operation
==
, !=
, >=
, <=
)
are not recommended for matrices of type mat or fmat,
due to the necessarily limited precision of floating-point element types;
consider using approx_equal() instead
+
, −
and %
operators are chained, Armadillo aims to avoid the generation of temporaries;
no temporaries are generated if all given objects are of the same type and size
*
operator is chained, Armadillo aims to find an efficient ordering of the matrix multiplications
mat A(5, 10, fill::randu); mat B(5, 10, fill::randu); mat C(10, 5, fill::randu); 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);
.n_rows
|
number of rows; present in Mat, Col, Row, Cube, field and SpMat | |
.n_cols
|
number of columns; present in Mat, Col, Row, Cube, field and SpMat | |
.n_elem
|
total number of elements; present in Mat, Col, Row, Cube, field and SpMat | |
.n_slices
|
number of slices; present in Cube and field | |
.n_nonzero
|
number of non-zero elements; present in SpMat |
mat X(4,5); cout << "X has " << X.n_cols << " columns" << endl;
(i)
|
|
For vec and rowvec, access the element stored at index i. For mat, cube and field, access the element/object stored at index i under the assumption of a flat layout, with column-major ordering of data (ie. column by column). An exception is thrown if the requested element is out of bounds. | |
.at(i) or [i]
|
As for (i) , but without a bounds check; not recommended; see the caveats below
|
||
(r,c)
|
For mat and 2D field classes, access the element/object stored at row r and column c. An exception is thrown if the requested element is out of bounds. | ||
.at(r,c)
|
As for (r,c) , but without a bounds check; not recommended; see the caveats below
|
||
(r,c,s)
|
For cube and 3D field classes, access the element/object stored at row r, column c, and slice s. An exception is thrown if the requested element is out of bounds. | ||
.at(r,c,s)
|
As for (r,c,s) , but without a bounds check; not recommended; see the caveats below
|
for(uword i=0; i<X.n_elem; ++i) { X(i) = ... }
[r,c]
and [r,c,s]
does not work correctly in C++;
instead use (r,c)
and (r,c,s)
mat M(10, 10, fill::randu); M(9,9) = 123.0; double x = M(1,2); vec v(10, fill::randu); v(9) = 123.0; double y = v(0);
vec v = { 1, 2, 3 }; mat A = { {1, 3, 5}, {2, 4, 6} };
.zeros() | (member function of Mat, Col, Row, SpMat, Cube) | ||
.zeros( n_elem ) | (member function of Col and Row) | ||
.zeros( n_rows, n_cols ) | (member function of Mat and SpMat) | ||
.zeros( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.zeros( size(X) ) | (member function of Mat, Col, Row, Cube, SpMat) |
mat A; A.zeros(5, 10); // or: mat A(5, 10, fill::zeros); mat B; B.zeros( size(A) ); mat C(5, 10, fill::randu); C.zeros();
.ones() | (member function of Mat, Col, Row, Cube) | ||
.ones( n_elem ) | (member function of Col and Row) | ||
.ones( n_rows, n_cols ) | (member function of Mat) | ||
.ones( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.ones( size(X) ) | (member function of Mat, Col, Row, Cube) |
mat A; A.ones(5, 10); // or: mat A(5, 10, fill::ones); mat B; B.ones( size(A) ); mat C(5, 10, fill::randu); C.ones();
mat A; A.eye(5, 5); // or: mat A(5, 5, fill::eye); mat B; B.eye( size(A) ); mat C(5, 5, fill::randu); C.eye();
.randu() | (member function of Mat, Col, Row, Cube) | ||
.randu( n_elem ) | (member function of Col and Row) | ||
.randu( n_rows, n_cols ) | (member function of Mat) | ||
.randu( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.randu( size(X) ) | (member function of Mat, Col, Row, Cube) |
.randn() | (member function of Mat, Col, Row, Cube) | ||
.randn( n_elem ) | (member function of Col and Row) | ||
.randn( n_rows, n_cols ) | (member function of Mat) | ||
.randn( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.randn( size(X) ) | (member function of Mat, Col, Row, Cube) |
mat A; A.randu(5, 10); // or: mat A(5, 10, fill::randu); mat B; B.randu( size(A) ); mat C(5, 10, fill::zeros); C.randu();
mat A(5, 6); A.fill(123.0); // or: mat A(5, 6, fill::value(123.0));
mat A(5, 6, fill::zeros);
std::mt19937 engine; // Mersenne twister random number engine std::uniform_real_distribution<double> distr(0.0, 1.0); mat A(4, 5, fill::none); A.imbue( [&]() { return distr(engine); } );
sp_mat A; A.sprandu(1000, 1000, 0.01); A(12,34) = datum::eps; A(56,78) = -datum::eps; A.clean(datum::eps);
mat A(5, 6, fill::randu); A.diag().fill(datum::nan); A.replace(datum::nan, 0); // replace each NaN with 0
mat A(5, 6, fill::randu); A.clamp(0.2, 0.8);
mat A(4, 5, fill::ones); // add 123 to every element A.transform( [](double val) { return (val + 123.0); } );
// add 123 to each element in a dense matrix mat A(4, 5, fill::ones); A.for_each( [](mat::elem_type& val) { val += 123.0; } ); // NOTE: the '&' is crucial! // add 123 to each non-zero element in a sparse matrix sp_mat S; S.sprandu(1000, 2000, 0.1); S.for_each( [](sp_mat::elem_type& val) { val += 123.0; } ); // NOTE: the '&' is crucial! // set the size of all matrices in field F field<mat> F(2,3); F.for_each( [](mat& X) { X.zeros(4,5); } ); // NOTE: the '&' is crucial!
.set_size( n_elem ) | (member function of Col, Row, field) | ||
.set_size( n_rows, n_cols ) | (member function of Mat, SpMat, field) | ||
.set_size( n_rows, n_cols, n_slices ) | (member function of Cube and field) | ||
.set_size( size(X) ) | (member function of Mat, Col, Row, Cube, SpMat, field) |
mat A; A.set_size(5, 10); // or: mat A(5, 10, fill::none); mat B; B.set_size( size(A) ); // or: mat B(size(A), fill::none); vec v; v.set_size(100); // or: vec v(100, fill::none);
.reshape( n_rows, n_cols ) | (member function of Mat and SpMat) | ||
.reshape( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.reshape( size(X) ) | (member function of Mat, Cube, SpMat) |
mat A(4, 5, fill::randu); A.reshape(5,4);
.resize( n_elem ) | (member function of Col, Row) | ||
.resize( n_rows, n_cols ) | (member function of Mat and SpMat) | ||
.resize( n_rows, n_cols, n_slices ) | (member function of Cube) | ||
.resize( size(X) ) | (member function of Mat, Col, Row, Cube, SpMat) |
mat A(4, 5, fill::randu); A.resize(7,6);
mat A(5, 6, fill::randu); mat B; B.copy_size(A); cout << B.n_rows << endl; cout << B.n_cols << endl;
mat A(5, 5, fill::randu); A.reset();
|
|
mat A(5, 10, fill::zeros); A.submat( 0,1, 2,3 ) = randu<mat>(3,3); A( span(0,2), span(1,3) ) = randu<mat>(3,3); A( 0,1, size(3,3) ) = randu<mat>(3,3); mat B = A.submat( 0,1, 2,3 ); mat C = A( span(0,2), span(1,3) ); mat D = A( 0,1, size(3,3) ); A.col(1) = randu<mat>(5,1); A(span::all, 1) = randu<mat>(5,1); mat X(5, 5, fill::randu); // get all elements of X that are greater than 0.5 vec q = X.elem( find(X > 0.5) ); // add 123 to all elements of X greater than 0.5 X.elem( find(X > 0.5) ) += 123.0; // set four specific elements of X to 1 uvec indices = { 2, 3, 6, 8 }; X.elem(indices) = ones<vec>(4); // add 123 to the last 5 elements of vector a vec a(10, fill::randu); a.tail(5) += 123.0; // add 123 to the first 3 elements of column 2 of X X.col(2).head(3) += 123;
Q.slices( first_slice, last_slice ) Q.row( row_number ) Q.rows( first_row, last_row ) Q.col( col_number ) Q.cols( first_col, last_col ) Q.subcube( first_row, first_col, first_slice, last_row, last_col, last_slice ) Q( span(first_row, last_row), span(first_col, last_col), span(first_slice, last_slice) ) Q( first_row, first_col, first_slice, size(n_rows, n_cols, n_slices) ) Q( first_row, first_col, first_slice, size(R) ) [ R is a cube ] Q.head_slices( number_of_slices ) Q.tail_slices( number_of_slices ) Q.tube( row, col ) Q.tube( first_row, first_col, last_row, last_col ) Q.tube( span(first_row, last_row), span(first_col, last_col) ) Q.tube( first_row, first_col, size(n_rows, n_cols) ) |
Q.elem( vector_of_indices ) Q( vector_of_indices ) Q.slices( vector_of_slice_indices ) |
cube A(2, 3, 4, fill::randu); mat B = A.slice(1); // each slice is a matrix 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); A( 0,0,1, size(2,2,2) ) = randu<cube>(2,2,2); // add 123 to all elements of A greater than 0.5 A.elem( find(A > 0.5) ) += 123.0; cube C = A.head_slices(2); // get first two slices A.head_slices(2) += 123.0;
mat X(5, 5, fill::randu); vec a = X.diag(); vec b = X.diag(1); vec c = X.diag(-2); X.diag() = randu<vec>(5); X.diag() += 6; X.diag().ones(); sp_mat S = sprandu<sp_mat>(10,10,0.1); vec v(S.diag()); // copy sparse diagonal to dense vector
.each_col() | .each_row() | (form 1) | ||
.each_col( vector_of_indices ) | .each_row( vector_of_indices ) | (form 2) | ||
.each_col( lambda_function ) | .each_row( lambda_function ) | (form 3) |
+ | addition | += | in-place addition | |||
− | subtraction | −= | in-place subtraction | |||
% | element-wise multiplication | %= | in-place element-wise multiplication | |||
/ | element-wise division | /= | in-place element-wise division | |||
= | assignment (copy) |
mat X(6, 5, fill::ones); vec v = linspace<vec>(10,15,6); X.each_col() += v; // in-place addition of v to each column vector of X mat Y = X.each_col() + v; // generate Y by adding v to each column vector of X // subtract v from columns 0 through to 3 in X X.cols(0,3).each_col() -= v; uvec indices(2); indices(0) = 2; indices(1) = 4; X.each_col(indices) = v; // copy v to columns 2 and 4 in X X.each_col( [](vec& a){ a.print(); } ); // lambda function with non-const vector const mat& XX = X; XX.each_col( [](const vec& b){ b.print(); } ); // lambda function with const vector
.each_slice() | (form 1) | |
.each_slice( vector_of_indices ) | (form 2) | |
.each_slice( lambda_function ) | (form 3) | |
.each_slice( lambda_function, use_mp ) | (form 4) |
+ | addition | += | in-place addition | |||
− | subtraction | −= | in-place subtraction | |||
% | element-wise multiplication | %= | in-place element-wise multiplication | |||
/ | element-wise division | /= | in-place element-wise division | |||
* | matrix multiplication | *= | in-place matrix multiplication | |||
= | assignment (copy) |
*
and *=
(ie. matrix multiplication)cube C(4, 5, 6, fill::randu); mat M = repmat(linspace<vec>(1,4,4), 1, 5); C.each_slice() += M; // in-place addition of M to each slice of C cube D = C.each_slice() + M; // generate D by adding M to each slice of C uvec indices(2); indices(0) = 2; indices(1) = 4; C.each_slice(indices) = M; // copy M to slices 2 and 4 in C C.each_slice( [](mat& X){ X.print(); } ); // lambda function with non-const matrix const cube& CC = C; CC.each_slice( [](const mat& X){ X.print(); } ); // lambda function with const matrix
mat A(4, 5, fill::randu); mat B(4, 5, fill::randu); cx_mat C(4, 5, fill::zeros); C.set_real(A); C.set_imag(B);
mat A(4, 5, fill::randu); mat B(4, 5, fill::randu); cx_mat C = cx_mat(A,B);
.insert_rows( row_number, X )
.insert_rows( row_number, number_of_rows ) |
(member functions of Mat, Col and Cube) | |
.insert_cols( col_number, X )
.insert_cols( col_number, number_of_cols ) |
(member functions of Mat, Row and Cube) | |
.insert_slices( slice_number, X )
.insert_slices( slice_number, number_of_slices ) |
(member functions of Cube) |
mat A(5, 10, fill::randu); mat B(5, 2, fill::ones ); // 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);
.shed_row( row_number )
.shed_rows( first_row, last_row ) .shed_rows( vector_of_indices ) |
(member function of Mat, Col, SpMat, Cube)
(member function of Mat, Col, SpMat, Cube) (member function of Mat, Col) |
|
.shed_col( column_number )
.shed_cols( first_column, last_column ) .shed_cols( vector_of_indices ) |
(member function of Mat, Row, SpMat, Cube)
(member function of Mat, Row, SpMat, Cube) (member function of Mat, Row) |
|
.shed_slice( slice_number )
.shed_slices( first_slice, last_slice ) .shed_slices( vector_of_indices ) |
(member functions of Cube) |
mat A(5, 10, fill::randu); mat B(5, 10, fill::randu); A.shed_row(2); A.shed_cols(2,4); uvec indices = {4, 6, 8}; B.shed_cols(indices);
mat X(5, 5, fill::randu); X.swap_rows(0,4);
mat A(4, 5, fill::zeros); mat B(6, 7, fill::ones ); A.swap(B);
mat A(5, 5, fill::randu); const mat B(5, 5, fill::randu); double* A_mem = A.memptr(); const double* B_mem = B.memptr();
mat A(5, 5, fill::randu); double* mem = A.colptr(2);
.begin() | |
iterator referring to the first element |
.end() | |
iterator referring to the past-the-end element |
.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 |
.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 |
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 specified columns |
|
||
mat::const_col_iterator
vec::const_col_iterator rowvec::const_col_iterator |
|
random access iterators, for read-only access to the elements of specified columns |
|
||
mat::row_iterator | |
bidirectional iterator, for read/write access to the elements of specified rows |
|
||
mat::const_row_iterator | |
bidirectional iterator, for read-only access to the elements of specified rows |
|
||
vec::row_iterator
rowvec::row_iterator |
|
random access iterators, for read/write access to the elements of specified rows |
|
||
vec::const_row_iterator
rowvec::const_row_iterator |
|
random access iterators, for read-only access to the elements of specified rows |
mat X(5, 6, fill::randu); mat::iterator it = X.begin(); mat::iterator it_end = X.end(); for(; it != it_end; ++it) { cout << (*it) << endl; } mat::col_iterator col_it = X.begin_col(1); // start of column 1 mat::col_iterator col_it_end = X.end_col(3); // end of column 3 for(; col_it != col_it_end; ++col_it) { cout << (*col_it) << endl; (*col_it) = 123.0; }
.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 |
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 iterator, 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 iterator, for read-only access to the elements of a particular slice |
cube X(2, 3, 4, fill::randu); cube::iterator it = X.begin(); cube::iterator it_end = X.end(); for(; it != it_end; ++it) { cout << (*it) << endl; } cube::slice_iterator s_it = X.begin_slice(1); // start of slice 1 cube::slice_iterator s_it_end = X.end_slice(2); // end of slice 2 for(; s_it != s_it_end; ++s_it) { cout << (*s_it) << endl; (*s_it) = 123.0; }
.begin() | |
iterator referring to the first element |
.end() | |
iterator referring to the past-the-end element |
.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 |
.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 |
sp_mat::iterator | |
bidirectional iterator, for read/write access to elements (which are stored column by column) |
sp_mat::const_iterator | |
bidirectional iterator, for read-only access to elements (which are stored column by column) |
|
||
sp_mat::col_iterator | |
bidirectional iterator, for read/write access to the elements of a specific column |
sp_mat::const_col_iterator | |
bidirectional iterator, for read-only access to the elements of a specific column |
|
||
sp_mat::row_iterator | |
bidirectional iterator, for read/write access to the elements of a specific row |
sp_mat::const_row_iterator | |
bidirectional iterator, for read-only access to the elements of a specific row |
|
sp_mat X = sprandu<sp_mat>(1000, 2000, 0.1); sp_mat::const_iterator it = X.begin(); sp_mat::const_iterator it_end = X.end(); for(; it != it_end; ++it) { cout << "val: " << (*it) << endl; cout << "row: " << it.row() << endl; cout << "col: " << it.col() << endl; }
mat X(100, 200, fill::randu); for( double& val : X(span(40,60), span(50,100)) ) { cout << val << endl; val = 123.0; }
.front() | |
access the first element in a vector |
.back() | |
access the last element in a vector |
.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 |
mat A(5, 5, fill::randu); cout << A.size() << endl; A.clear(); cout << A.empty() << endl;
mat X(4, 5, fill::randu); vec v = X.as_col();
cube Q(5, 4, 3, fill::randu); mat A = Q.col_as_mat(2); // size of A: 5x3 mat B = Q.row_as_mat(2); // size of B: 3x4
sp_mat A; A.sprandu(1000, 1000, 0.1); // store the sum of each column of A directly in dense row vector rowvec r = sum(A).as_dense(); // extract column 123 of A directly into dense column vector colvec c = A.col(123).as_dense();
mat A(4, 5, fill::randu); mat B = A.t();
mat A(4, 4, fill::randu); mat X = A.i(); mat Y = (A+A).i();
mat A(5, 5, fill::randu); double max_val = A.max();
mat A(5, 5, fill::randu); uword i = A.index_max(); double max_val = A(i);
cx_mat A( randu<mat>(4,4), randu<mat>(4,4) ); real(A).eval().save("A_real.dat", raw_ascii); imag(A).eval().save("A_imag.dat", raw_ascii);
.in_range( i ) | (member of Mat, Col, Row, Cube, SpMat, field) | ||
.in_range( span(start, end) ) | (member of Mat, Col, Row, Cube, SpMat, field) | ||
.in_range( row, col ) | (member of Mat, Col, Row, SpMat, field) | ||
.in_range( span(start_row, end_row), span(start_col, end_col) ) | (member of Mat, Col, Row, SpMat, field) | ||
.in_range( row, col, slice ) | (member of Cube and field) | ||
.in_range( span(start_row, end_row), span(start_col, end_col), span(start_slice, end_slice) ) | (member of Cube and field) | ||
.in_range( first_row, first_col, size(X) ) (X is a matrix or field) | (member of Mat, Col, Row, SpMat, field) | ||
.in_range( first_row, first_col, size(n_rows, n_cols) ) | (member of Mat, Col, Row, SpMat, field) | ||
.in_range( first_row, first_col, first_slice, size(Q) ) (Q is a cube or field) | (member of Cube and field) | ||
.in_range( first_row, first_col, first_slice, size(n_rows, n_cols n_slices) ) | (member of Cube and field) |
mat A(4, 5, fill::randu); cout << A.in_range(0,0) << endl; // true cout << A.in_range(3,4) << endl; // true cout << A.in_range(4,5) << endl; // false
mat A(5, 5, fill::randu); cout << A.is_empty() << endl; A.reset(); cout << A.is_empty() << endl;
mat A(1, 5, fill::randu); mat B(5, 1, fill::randu); mat C(5, 5, fill::randu); cout << A.is_vec() << endl; cout << B.is_vec() << endl; cout << C.is_vec() << endl;
"ascend" | ↦ | elements are ascending; consecutive elements can be equal; this is the default operation |
"descend" | ↦ | elements are descending; consecutive elements can be equal |
"strictascend" | ↦ | elements are strictly ascending; consecutive elements cannot be equal |
"strictdescend" | ↦ | elements are strictly descending; consecutive elements cannot be equal |
vec a(10, fill::randu); vec b = sort(a); bool check1 = a.is_sorted(); bool check2 = b.is_sorted(); mat A(10, 10, fill::randu); // check whether each column is sorted in descending manner cout << A.is_sorted("descend") << endl; // check whether each row is sorted in ascending manner cout << A.is_sorted("ascend", 1) << endl;
mat A(5, 5, fill::randu); mat B = trimatl(A); cout << A.is_trimatu() << endl; cout << B.is_trimatl() << endl;
mat A(5, 5, fill::randu); mat B = diagmat(A); cout << A.is_diagmat() << endl; cout << B.is_diagmat() << endl;
mat A(5, 5, fill::randu); mat B(6, 7, fill::randu); cout << A.is_square() << endl; cout << B.is_square() << endl;
mat A(5, 5, fill::randu); mat B = A.t() * A; cout << A.is_symmetric() << endl; cout << B.is_symmetric() << endl;
cx_mat A(5, 5, fill::randu); cx_mat B = A.t() * A; cout << A.is_hermitian() << endl; cout << B.is_hermitian() << endl;
mat A(5, 5, fill::randu); mat B = A.t() * A; cout << A.is_sympd() << endl; cout << B.is_sympd() << endl;
mat A(5, 5, fill::zeros); A(0,0) = datum::eps; cout << A.is_zero() << endl; cout << A.is_zero(datum::eps) << endl;
mat A(5, 5, fill::randu); mat B(5, 5, fill::randu); B(1,1) = datum::inf; cout << A.is_finite() << endl; cout << B.is_finite() << endl;
mat A(5, 5, fill::randu); mat B(5, 5, fill::randu); B(1,1) = datum::inf; cout << A.has_inf() << endl; cout << B.has_inf() << endl;
NaN
is not equal to anything, even itself
mat A(5, 5, fill::randu); mat B(5, 5, fill::randu); B(1,1) = datum::nan; cout << A.has_nan() << endl; cout << B.has_nan() << endl;
mat A(5, 5, fill::randu); mat B(6, 6, fill::randu); A.print(); // print a transposed version of A A.t().print(); // "B:" is the optional header line B.print("B:"); cout << A << endl; cout << "B:" << endl; cout << B << endl;
mat A(5, 5, fill::randu); cout.precision(11); cout.setf(ios::fixed); A.raw_print(cout, "A:");
mat A(123, 456, fill::randu); A.brief_print("A:"); // possible output: // // A: // [matrix size: 123x456] // 0.8402 0.7605 0.6218 ... 0.9744 // 0.3944 0.9848 0.0409 ... 0.7799 // 0.7831 0.9350 0.4140 ... 0.8835 // : : : : : // 0.4954 0.1826 0.9848 ... 0.1918
.save( filename )
.save( filename, file_type ) .save( stream ) .save( stream, file_type ) .save( hdf5_name(filename, dataset) ) .save( hdf5_name(filename, dataset, settings) ) .save( csv_name(filename, header) ) .save( csv_name(filename, header, settings) ) |
.load( filename )
.load( filename, file_type ) .load( stream ) .load( stream, file_type ) .load( hdf5_name(filename, dataset) ) .load( hdf5_name(filename, dataset, settings) ) .load( csv_name(filename, header) ) .load( csv_name(filename, header, settings) ) |
auto_detect |
Used only by .load() only: attempt to automatically detect the file type as one of the formats described below;
[ default operation for .load() ] |
|
arma_binary |
Numerical data stored in machine dependent binary format, with a simple header to speed up loading.
The header indicates the type and size of matrix/cube.
[ default operation for .save() ] |
|
arma_ascii |
Numerical data stored in human readable text format, with a simple header to speed up loading.
The header indicates the type and size of matrix/cube.
|
|
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.
|
|
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.
Cubes are loaded as one slice.
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, which is a tuple surrounded by brackets: eg. (1.23,4.56) indicates 1.24 + 4.56i.
|
|
csv_ascii |
Numerical data stored in comma separated value (CSV) text format, without a header.
To save/load with a header, use the csv_name(filename,header) specification instead (more details below).
Handles complex numbers stored in the compound form of 1.24+4.56i.
Applicable to Mat and SpMat.
|
|
coord_ascii |
Numerical data stored as a text file in coordinate list format, without a header.
Only non-zero values are stored.
For real matrices, each line contains information in the following format: row column value
For complex matrices, each line contains information in the following format: row column real_value imag_value
The rows and columns start at zero. Armadillo ≥ 10.3: applicable to Mat and SpMat; Armadillo ≤ 10.2: applicable to SpMat only. Caveat: not supported by auto_detect. |
|
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.
|
|
hdf5_binary |
Numerical data stored in portable HDF5 binary format.
|
-lhdf5
).
HDF5 support can be enabled by defining ARMA_USE_HDF5 before including the armadillo header:
#define ARMA_USE_HDF5
#include <armadillo>
hdf5_opts::trans | save/load the data with columns transposed to rows (and vice versa) | |
hdf5_opts::append | instead of overwriting the file, append the specified dataset to the file; the specified dataset must not already exist in the file | |
hdf5_opts::replace | instead of overwriting the file, replace the specified dataset in the file caveat: HDF5 may not automatically reclaim deleted space; use h5repack to clean HDF5 files |
+
operator; for example: hdf5_opts::trans + hdf5_opts::append
field<std::string>
csv_opts::trans | save/load the data with columns transposed to rows (and vice versa) | |
csv_opts::no_header | assume there is no header line; the header argument is not referenced | |
csv_opts::semicolon | use semicolon (;) instead of comma (,) as the separator character | |
csv_opts::strict | interpret missing values as NaN (not applicable to sparse matrices) |
+
operator; for example: csv_opts::trans + csv_opts::no_header
mat A(5, 6, fill::randu); // default save format is arma_binary A.save("A.bin"); // save in raw_ascii format A.save("A.txt", raw_ascii); // save in CSV format without a header A.save("A.csv", csv_ascii); // save in CSV format with a header field<std::string> header(A.n_cols); header(0) = "foo"; header(1) = "bar"; // etc A.save( csv_name("A.csv", header) ); // save in HDF5 format with internal dataset named as "my_data" A.save(hdf5_name("A.h5", "my_data")); // automatically detect format type while loading mat B; B.load("A.bin"); // force loading in arma_ascii format mat C; C.load("A.txt", arma_ascii); // example of testing for success mat D; bool ok = D.load("A.bin"); if(ok == false) { cout << "problem with loading" << endl; }
.save( name )
.save( name, file_type ) .save( stream ) .save( stream, file_type ) |
.load( name )
.load( name, file_type ) .load( stream ) .load( stream, file_type ) |
auto_detect |
|
|
arma_binary |
|
|
ppm_binary |
|
vec a = linspace(0, 5, 6); rowvec b = linspace<rowvec>(5, 0, 6);
vec a = logspace(0, 5, 6); rowvec b = logspace<rowvec>(5, 0, 6);
vec a = regspace(0, 9); // 0, 1, ..., 9 uvec b = regspace<uvec>(2, 2, 10); // 2, 4, ..., 10 ivec c = regspace<ivec>(0, -1, -10); // 0, -1, ..., -10
uvec X = randperm(10); uvec Y = randperm(10,2);
mat A = eye(5,5); // or: mat A(5,5,fill::eye); fmat B = 123.0 * eye<fmat>(5,5); cx_mat C = eye<cx_mat>( size(B) );
mat A(5, 6, fill::ones)
vec v = ones(10); // or: vec v(10, fill::ones); uvec u = ones<uvec>(10); rowvec r = ones<rowvec>(10); mat A = ones(5,6); // or: mat A(5, 6, fill::ones); fmat B = ones<fmat>(5,6); umat C = ones<umat>(5,6); cube Q = ones(5,6,7); // or: cube Q(5, 6, 7, fill::ones); fcube R = ones<fcube>(5,6,7);
mat A(5, 6, fill::zeros)
vec v = zeros(10); // or: vec v(10, fill::zeros); uvec u = zeros<uvec>(10); rowvec r = zeros<rowvec>(10); mat A = zeros(5,6); // or: mat A(5, 6, fill::zeros); fmat B = zeros<fmat>(5,6); umat C = zeros<umat>(5,6); cube Q = zeros(5,6,7); // or: cube Q(5, 6, 7, fill::zeros); fcube R = zeros<fcube>(5,6,7);
double a = randu(); double b = randu(distr_param(10,20)); vec v1 = randu(5); // or: vec v1(5, fill::randu); vec v2 = randu(5, distr_param(10,20)); rowvec r1 = randu<rowvec>(5); rowvec r2 = randu<rowvec>(5, distr_param(10,20)); mat A1 = randu(5, 6); // or: mat A1(5, 6, fill::randu); mat A2 = randu(5, 6, distr_param(10,20)); fmat B1 = randu<fmat>(5, 6); fmat B2 = randu<fmat>(5, 6, distr_param(10,20));
double a = randn(); double b = randn(distr_param(10,5)); vec v1 = randn(5); // or: vec v1(5, fill::randn); vec v2 = randn(5, distr_param(10,5)); rowvec r1 = randn<rowvec>(5); rowvec r2 = randn<rowvec>(5, distr_param(10,5)); mat A1 = randn(5, 6); // or: mat A1(5, 6, fill::randn); mat A2 = randn(5, 6, distr_param(10,5)); fmat B1 = randn<fmat>(5, 6); fmat B2 = randn<fmat>(5, 6, distr_param(10,5));
x a-1 exp( -x / b ) | ||
p(x | a,b) | = | |
b a Γ(a) |
vec v1 = randg(100); vec v2 = randg(100, distr_param(2,1)); rowvec r1 = randg<rowvec>(100); rowvec r2 = randg<rowvec>(100, distr_param(2,1)); mat A1 = randg(10, 10); mat A2 = randg(10, 10, distr_param(2,1)); fmat B1 = randg<fmat>(10, 10); fmat B2 = randg<fmat>(10, 10, distr_param(2,1));
int a = randi(); int b = randi(distr_param(-10, +20)); imat A1 = randi(5, 6); imat A2 = randi(5, 6, distr_param(-10, +20)); mat B1 = randi<mat>(5, 6); mat B2 = randi<mat>(5, 6, distr_param(-10, +20));
sp_mat A = speye<sp_mat>(5,5);
sp_mat A = sprandu<sp_mat>(100, 200, 0.1); sp_mat B = spones(A);
sp_mat A = sprandu<sp_mat>(100, 200, 0.1);
vec A(5, fill::randu); mat X = toeplitz(A); mat Y = circ_toeplitz(A);
mat A(5, 5, fill::randu); mat B = abs(A); cx_mat X(5, 5, fill::randu); mat Y = abs(X);
mat A(5, 6, fill::randu); mat B(5, 6, fill::randu); double x = accu(A); double y = accu(A % B); // accu(A % B) is a "multiply-and-accumulate" operation // as operator % performs element-wise multiplication
⎡ C0 ⎤ ⎡ A00 A01 A02 ⎤ ⎡ B0 ⎤ ⎢ C1 ⎥ = ⎢ A10 A11 A12 ⎥ x ⎢ B1 ⎥ ⎣ C2 ⎦ ⎣ A20 A21 A22 ⎦ ⎣ 1 ⎦
⎡ C0 ⎤ ⎡ A00 A01 A02 ⎤ ⎡ B0 ⎤ ⎣ C1 ⎦ = ⎣ A10 A11 A12 ⎦ x ⎢ B1 ⎥ ⎣ 1 ⎦
mat A(2, 3, fill::randu); vec B(2, fill::randu); vec C = affmul(A,B);
vec V(10, fill::randu); mat X(5, 5, fill::randu); // status1 will be set to true if vector V has all non-zero elements bool status1 = all(V); // status2 will be set to true if vector V has all elements greater than 0.5 bool status2 = all(V > 0.5); // status3 will be set to true if matrix X has all elements greater than 0.6; // note the use of vectorise() bool status3 = all(vectorise(X) > 0.6); // generate a row vector indicating which columns of X have all elements greater than 0.7 umat A = all(X > 0.7);
vec V(10, fill::randu); mat X(5, 5, fill::randu); // status1 will be set to true if vector V has any non-zero elements bool status1 = any(V); // status2 will be set to true if vector V has any elements greater than 0.5 bool status2 = any(V > 0.5); // status3 will be set to true if matrix X has any elements greater than 0.6; // note the use of vectorise() bool status3 = any(vectorise(X) > 0.6); // generate a row vector indicating which columns of X have elements greater than 0.7 umat A = any(X > 0.7);
"absdiff" | ↦ | scalars x and y are considered equal if |x − y| ≤ tol |
"reldiff" | ↦ | scalars x and y are considered equal if |x − y| / max( |x|, |y| ) ≤ tol |
"both" | ↦ | scalars x and y are considered equal if |x − y| ≤ abs_tol or |x − y| / max( |x|, |y| ) ≤ rel_tol |
mat A(5, 5, fill::randu); mat B = A + 0.001; bool same1 = approx_equal(A, B, "absdiff", 0.002); mat C = 1000 * randu<mat>(5,5); mat D = C + 1; bool same2 = approx_equal(C, D, "reldiff", 0.1); bool same3 = approx_equal(C, D, "both", 2, 0.1);
cx_mat A(5, 5, fill::randu); mat B = arg(A);
rowvec r(5, fill::randu); colvec q(5, fill::randu); mat X(5, 5, fill::randu); // examples of expressions which have optimised implementations 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);
mat A(5, 5, fill::randu); mat B = clamp(A, 0.2, 0.8); mat C = clamp(A, A.min(), 0.8); mat D = clamp(A, 0.2, A.max());
mat A(5, 5, fill::randu); double c = cond(A);
cx_mat X(5, 5, fill::randu); cx_mat Y = conj(X);
mat A(5, 5, fill::randu); 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);
vec a(3, fill::randu); vec b(3, fill::randu); vec c = cross(a,b);
mat A(5, 5, fill::randu); mat B = cumsum(A); mat C = cumsum(A, 1); vec x(10, fill::randu); vec y = cumsum(x);
mat A(5, 5, fill::randu); mat B = cumprod(A); mat C = cumprod(A, 1); vec x(10, fill::randu); vec y = cumprod(x);
val = det( A ) | (form 1) | |
det( val, A ) | (form 2) |
mat A(5, 5, fill::randu); double val1 = det(A); // form 1 double val2; bool success = det(val2, A); // form 2
mat A(5, 5, fill::randu); mat B = diagmat(A); mat C = diagmat(A,1); vec v(5, fill::randu); mat D = diagmat(v); mat E = diagmat(v,1);
mat A(5, 5, fill::randu); vec d = diagvec(A);
diags( V, D, n_rows, n_cols ) |
spdiags( V, D, n_rows, n_cols ) |
mat V(10, 3, fill::randu); ivec D = {-1, 0, +1}; mat X = diags(V, D, 10, 10); sp_mat Y = spdiags(V, D, 10, 10);
vec a = linspace<vec>(1,10,10); vec b = diff(a);
vec a(10, fill::randu); vec b(10, fill::randu); double x = dot(a,b);
mat A(4, 5, fill::randu); mat B = eps(A);
mat A(5, 5, fill::randu); mat B = expmat(A);
mat A(5, 5, fill::randu); mat B = A*A.t(); // make symmetric matrix mat C = expmat_sym(B);
mat A(5, 5, fill::randu); mat B(5, 5, fill::randu); uvec q1 = find(A > B); uvec q2 = find(A > 0.5); uvec q3 = find(A > 0.5, 3, "last"); // change elements of A greater than 0.5 to 1 A.elem( find(A > 0.5) ).ones();
mat A(5, 5, fill::randu); A(1,1) = datum::inf; // accumulate only finite elements double val = accu( A.elem( find_finite(A) ) );
mat A(5, 5, fill::randu); A(1,1) = datum::inf; A(2,2) = datum::nan; // change non-finite elements to zero A.elem( find_nonfinite(A) ).zeros();
mat A(5, 5, fill::randu); A(2,3) = datum::nan; uvec indices = find_nan(A);
true | = | the returned indices are sorted to be ascending (default setting) |
false | = | the returned indices are in arbitrary order (faster operation) |
mat A = { { 2, 2, 4 }, { 4, 6, 6 } }; uvec indices = find_unique(A);
mat A(5, 5, fill::randu); mat B = fliplr(A); mat C = flipud(A);
cx_mat C(5, 5, fill::randu); mat A = imag(C); mat B = real(C);
uvec | sub = ind2sub( size(X), index ) | (form 1) | |
umat | sub = ind2sub( size(X), vector_of_indices ) | (form 2) |
mat M(4, 5, fill::randu); uvec s = ind2sub( size(M), 6 ); cout << "row: " << s(0) << endl; cout << "col: " << s(1) << endl; uvec indices = find(M > 0.5); umat t = ind2sub( size(M), indices ); cube Q(2,3,4); uvec u = ind2sub( size(Q), 8 ); cout << "row: " << u(0) << endl; cout << "col: " << u(1) << endl; cout << "slice: " << u(2) << endl;
index_min( V )
index_min( M ) index_min( M, dim ) index_min( Q ) index_min( Q, dim ) |
index_max( V )
index_max( M ) index_max( M, dim ) index_max( Q ) index_max( Q, dim ) |
vec v(10, fill::randu); uword i = index_max(v); double max_val_in_v = v(i); mat M(5, 6, fill::randu); urowvec ii = index_max(M); ucolvec jj = index_max(M,1); double max_val_in_col_2 = M( ii(2), 2 ); double max_val_in_row_4 = M( 4, jj(4) );
"lowmem"
mat X(4, 5, fill::randu); mat Y(20000, 30000, fill::randu); inplace_trans(X); // use greedy algorithm by default inplace_trans(Y, "lowmem"); // use low-memory (and slow) algorithm
C = intersect( A, B ) | (form 1) | |
intersect( C, iA, iB, A, B ) | (form 2) |
ivec A = regspace<ivec>(4, 1); // 4, 3, 2, 1 ivec B = regspace<ivec>(3, 6); // 3, 4, 5, 6 ivec C = intersect(A,B); // 3, 4 ivec CC; uvec iA; uvec iB; intersect(CC, iA, iB, A, B);
join_rows( A, B )
join_rows( A, B, C ) join_rows( A, B, C, D ) join_cols( A, B ) join_cols( A, B, C ) join_cols( A, B, C, D ) |
join_horiz( A, B )
join_horiz( A, B, C ) join_horiz( A, B, C, D ) join_vert( A, B ) join_vert( A, B, C ) join_vert( A, B, C, D ) |
mat A(4, 5, fill::randu); mat B(4, 6, fill::randu); mat C(6, 5, fill::randu); mat AB = join_rows(A,B); mat AC = join_cols(A,C);
cube C(5, 10, 3, fill::randu); cube D(5, 10, 4, fill::randu); cube E = join_slices(C,D); mat M(10, 20, fill::randu); mat N(10, 20, fill::randu); cube Q = join_slices(M,N); cube R = join_slices(Q,M); cube S = join_slices(M,Q);
mat A(4, 5, fill::randu); mat B(5, 4, fill::randu); mat K = kron(A,B);
complex result = log_det( A ) | (form 1) | |
log_det( val, sign, A ) | (form 2) |
mat A(5, 5, fill::randu); cx_double result = log_det(A); // form 1 double val; double sign; bool ok = log_det(val, sign, A); // form 2
result = log_det_sympd( A ) | (form 1) | |
log_det_sympd( result, A ) | (form 2) |
mat A(5, 5, fill::randu); mat B = A.t() * A; // make symmetric matrix double result1 = log_det_sympd(B); // form 1 double result2; bool success = log_det_sympd(result2, B); // form 2
mat A(5, 5, fill::randu); cx_mat B = logmat(A);
mat A(5, 5, fill::randu); mat B = A*A.t(); // make symmetric matrix mat C = logmat_sympd(B);
min( V )
min( M ) min( M, dim ) min( Q ) min( Q, dim ) min( A, B ) |
max( V )
max( M ) max( M, dim ) max( Q ) max( Q, dim ) max( A, B ) |
colvec v(10, fill::randu); double x = max(v); mat M(10, 10, fill::randu); rowvec a = max(M); rowvec b = max(M,0); colvec c = max(M,1); // element-wise maximum mat X(5, 6, fill::randu); mat Y(5, 6, fill::randu); mat Z = arma::max(X,Y); // use arma:: prefix to distinguish from std::max()
accu(X != 0)
is more efficient
.n_nonzero
attribute is more efficient, eg. X.n_nonzero
sp_mat A = sprandu<sp_mat>(100, 100, 0.1); vec a = nonzeros(A); mat B(100, 100, fill::eye); vec b = nonzeros(B);
"-inf"
, "inf"
, "fro"
"inf"
, "fro"
"-inf"
is the minimum quasi-norm, "inf"
is the maximum norm, "fro"
is the Frobenius norm
accu(X != 0)
vec q(5, fill::randu); double x = norm(q, 2); double y = norm(q, "inf");
mat X(2000, 3000, fill::randu); double x = norm2est(X); double y = norm2est(X, 1e-5); double z = norm2est(X, 1e-4, 10);
vec A(10, fill::randu); vec B = normalise(A); vec C = normalise(A, 1); mat X(5, 6, fill::randu); mat Y = normalise(X); mat Z = normalise(X, 2, 1);
pow( A, scalar ) | (form 1) | |
pow( A, B ) | (form 2) |
mat A(5, 6, fill::randu); mat B(5, 6, fill::randu); mat X = pow(A, 3.45); mat Y = pow(A, B);
mat A(5, 5, fill::randu); mat B = powmat(A, 4); // integer exponent cx_mat C = powmat(A, 4.56); // non-integer exponent
colvec v(10, fill::randu); double x = prod(v); mat M(10, 10, fill::randu); rowvec a = prod(M); rowvec b = prod(M,0); colvec c = prod(M,1);
r = rank( X ) | (form 1) | |
r = rank( X, tolerance ) | ||
rank( r, X ) | (form 2) | |
rank( r, X, tolerance ) |
std::rank
, use the arma::
prefix, ie. arma::rank(X)
mat A(4, 5, fill::randu); uword r1 = rank(A); // form 1 uword r2; bool success = rank(r2, A); // form 2
mat A(5, 5, fill::randu); double r = rcond(A);
n_rows | = | num_copies_per_row | * | A.n_rows |
n_cols | = | num_copies_per_col | * | A.n_cols |
mat A(2, 3, fill::randu); mat B = repelem(A, 4, 5);
n_rows | = | num_copies_per_row | × | A.n_rows |
n_cols | = | num_copies_per_col | × | A.n_cols |
mat A(2, 3, fill::randu); mat B = repmat(A, 4, 5);
mat A(10, 5, fill::randu); mat B = reshape(A, 5, 10);
mat A(4, 5, fill::randu); mat B = resize(A, 7, 6);
vec v(123, fill::randu); vec y = reverse(v); mat A(4, 5, fill::randu); mat B = reverse(A); mat C = reverse(A,1);
vec P(5, fill::randu); cx_vec R = roots(P);
mat A(4, 5, fill::randu); mat B = shift(A, -1); mat C = shift(A, +1);
mat A(4, 5, fill::randu); mat B = shuffle(A);
mat A(5,6); mat B(size(A), fill::zeros); mat C; C.randu(size(A)); mat D = ones<mat>(size(A)); mat E(10,20, fill::ones); E(3,4,size(C)) = C; // access submatrix of E mat F( size(A) + size(E) ); mat G( size(A) * 2 ); cout << "size of A: " << size(A) << endl; bool is_same_size = (size(A) == size(E));
"ascend"
or "descend"
; by default "ascend"
is usedmat A(10, 10, fill::randu); mat B = sort(A);
"ascend"
or "descend"
; by default "ascend"
is usedvec q(10, fill::randu); uvec indices = sort_index(q);
mat A(5, 5, fill::randu); cx_mat B = sqrtmat(A);
mat A(5, 5, fill::randu); mat B = A*A.t(); // make symmetric matrix mat C = sqrtmat_sympd(B);
colvec v(10, fill::randu); double x = sum(v); mat M(10, 10, fill::randu); rowvec a = sum(M); rowvec b = sum(M,0); colvec c = sum(M,1); double y = accu(M); // find the overall sum regardless of object type
uword | index | = sub2ind( size(M), row, col ) | (M is a matrix) |
uvec | indices | = sub2ind( size(M), matrix_of_subscripts ) | |
uword | index | = sub2ind( size(Q), row, col, slice ) | (Q is a cube) |
uvec | indices | = sub2ind( size(Q), matrix_of_subscripts ) |
mat M(4,5); cube Q(4,5,6); uword i = sub2ind( size(M), 2, 3 ); uword j = sub2ind( size(Q), 2, 3, 4 );
mat A(5, 5, fill::randu); mat B = symmatu(A); mat C = symmatl(A);
mat A(5, 5, fill::randu); double x = trace(A);
mat A(5, 10, fill::randu); mat B = trans(A); mat C = A.t(); // equivalent to trans(A), but more compact
vec X = linspace<vec>(0, datum::pi, 1000); vec Y = sin(X); mat Z = trapz(X,Y);
mat A(5, 5, fill::randu); mat U = trimatu(A); mat L = trimatl(A); mat UU = trimatu(A, 1); // omit the main diagonal mat LL = trimatl(A, -1); // omit the main diagonal
mat A(5, 5, fill::randu); uvec upper_indices = trimatu_ind( size(A) ); uvec lower_indices = trimatl_ind( size(A) ); // extract upper/lower triangle into vector vec upper_part = A(upper_indices); vec lower_part = A(lower_indices); // obtain indices without the main diagonal uvec alt_upper_indices = trimatu_ind( size(A), 1); uvec alt_lower_indices = trimatl_ind( size(A), -1);
mat X = { { 1, 2 } { 2, 3 } }; mat Y = unique(X);
"-inf"
(minimum quasi-norm), or "inf"
(maximum norm)
mat X(4, 5, fill::randu); rowvec r = vecnorm(X, 2); colvec c = vecnorm(X, "inf", 1);
mat X(4, 5, fill::randu); vec v = vectorise(X);
exp | log | square | floor | erf | tgamma | sign | ||||||
exp2 | log2 | sqrt | ceil | erfc | lgamma | |||||||
exp10 | log10 | cbrt | round | |||||||||
expm1 | log1p | trunc | ||||||||||
trunc_exp | trunc_log |
exp(A)
|
base-e exponential: e x | |||||||||||||
exp2(A)
|
base-2 exponential: 2 x | |||||||||||||
exp10(A)
|
base-10 exponential: 10 x | |||||||||||||
expm1(A)
|
compute exp(A)-1 accurately for values of A close to zero (only for float and double elements)
|
|||||||||||||
trunc_exp(A)
|
base-e exponential, truncated to avoid infinity (only for float and double elements) | |||||||||||||
log(A)
|
natural log: loge x | |||||||||||||
log2(A)
|
base-2 log: log2 x | |||||||||||||
log10(A)
|
base-10 log: log10 x | |||||||||||||
log1p(A)
|
compute log(1+A) accurately for values of A close to zero (only for float and double elements)
|
|||||||||||||
trunc_log(A)
|
natural log, truncated to avoid ±infinity (only for float and double elements) | |||||||||||||
square(A)
|
square: x 2 | |||||||||||||
sqrt(A)
|
square root: √x | |||||||||||||
cbrt(A)
|
cube root: ∛x | |||||||||||||
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 | |||||||||||||
round(A)
|
round to nearest integer, with halfway cases rounded away from zero | |||||||||||||
trunc(A)
|
round to nearest integer, towards zero | |||||||||||||
erf(A)
|
error function (only for float and double elements) | |||||||||||||
erfc(A)
|
complementary error function (only for float and double elements) | |||||||||||||
tgamma(A)
|
gamma function (only for float and double elements) | |||||||||||||
lgamma(A)
|
natural log of the absolute value of gamma function (only for float and double elements) | |||||||||||||
sign(A)
|
signum function;
for each element a in A, the corresponding element b in B is:
|
mat A(5, 5, fill::randu); mat B = exp(A);
mat X(5, 5, fill::randu); mat Y = cos(X);
R = chol( X ) | (form 1) | |
R = chol( X, layout ) | (form 2) | |
chol( R, X ) | (form 3) | |
chol( R, X, layout ) | (form 4) | |
chol( R, P, X, layout, "vector" ) | (form 5) | |
chol( R, P, X, layout, "matrix" ) | (form 6) |
"upper"
or "lower"
, which specifies whether R is upper or lower triangular
mat A(5, 5, fill::randu); mat X = A.t()*A; mat R1 = chol(X); mat R2 = chol(X, "lower"); mat R3; bool ok = chol(R3, X); mat R; uvec P_vec; umat P_mat; chol(R, P_vec, X, "upper", "vector"); chol(R, P_mat, X, "lower", "matrix");
"dc"
or "std"
"dc"
indicates divide-and-conquer method (default setting)
"std"
indicates standard method
// for matrices with real elements mat A(50, 50, fill::randu); mat B = A.t()*A; // generate a symmetric matrix vec eigval; mat eigvec; eig_sym(eigval, eigvec, B); // for matrices with complex elements cx_mat C(50, 50, fill::randu); cx_mat D = C.t()*C; vec eigval2; cx_mat eigvec2; eig_sym(eigval2, eigvec2, D);
"balance" | ↦ | diagonally scale and permute X to improve conditioning of the eigenvalues |
"nobalance" | ↦ | do not balance X; this is the default operation |
mat A(10, 10, fill::randu); cx_vec eigval; cx_mat eigvec; eig_gen(eigval, eigvec, A); eig_gen(eigval, eigvec, A, "balance");
mat A(10, 10, fill::randu); mat B(10, 10, fill::randu); cx_vec eigval; cx_mat eigvec; eig_pair(eigval, eigvec, A, B);
mat X(20,20, fill::randu); mat U; mat H; hess(U, H, X);
inv_opts::no_ugly | do not provide inverses for poorly conditioned matrices (where rcond < A.n_rows · datum::eps) | |
inv_opts::allow_approx | allow approximate inverses for rank deficient or poorly conditioned matrices |
mat A(5, 5, fill::randu); mat B = inv(A); mat C; bool success = inv(C, A); mat D; double rcond_val; inv(D, rcond_val, A); A.col(1).zeros(); mat E; inv(E, A, inv_opts::allow_approx);
inv_opts::no_ugly | do not provide inverses for poorly conditioned matrices (where rcond < A.n_rows · datum::eps) | |
inv_opts::allow_approx | allow approximate inverses for rank deficient or poorly conditioned symmetric matrices |
mat A(5, 5, fill::randu); mat B = A.t() * A; mat C = inv_sympd(B); mat D; bool success = inv_sympd(D, B); mat E; double rcond_val; inv_sympd(E, rcond_val, B); B.col(1).zeros(); B.row(1).zeros(); mat F; inv_sympd(F, B, inv_opts::allow_approx);
mat A(5, 5, fill::randu); mat L, U, P; lu(L, U, P, A); mat B = P.t()*L*U;
mat A(5, 6, fill::randu); A.row(0).zeros(); A.col(0).zeros(); mat B = null(A);
mat A(5, 6, fill::randu); mat B = orth(A);
"dc"
or "std"
"dc"
indicates divide-and-conquer method (default setting)
"std"
indicates standard method
inv_opts::allow_approx
option is faster
mat A(4, 5, fill::randu); mat B = pinv(A); // use default tolerance mat C = pinv(A, 0.01); // set tolerance to 0.01
qr( Q, R, X ) | (form 1) | |
qr( Q, R, P, X, "vector" ) | (form 2) | |
qr( Q, R, P, X, "matrix" ) | (form 3) |
mat X(5, 5, fill::randu); mat Q; mat R; qr(Q, R, X); uvec P_vec; umat P_mat; qr(Q, R, P_vec, X, "vector"); qr(Q, R, P_mat, X, "matrix");
mat X(6, 5, fill::randu); mat Q; mat R; qr_econ(Q, R, X);
"none" | no ordering (default operation) | |
"lhp" | left-half-plane: eigenvalues with real part < 0 | |
"rhp" | right-half-plane: eigenvalues with real part > 0 | |
"iuc" | inside-unit-circle: eigenvalues with absolute value < 1 | |
"ouc" | outside-unit-circle: eigenvalues with absolute value > 1 |
mat A(10, 10, fill::randu); mat B(10, 10, fill::randu); mat AA; mat BB; mat Q; mat Z; qz(AA, BB, Q, Z, A, B);
mat X(20,20, fill::randu); mat U; mat S; schur(U, S, X);
solve_opts::fast | fast mode: disable determining solution quality via rcond, disable iterative refinement, disable equilibration | |
solve_opts::refine | apply iterative refinement to improve solution quality (matrix A must be square) | |
solve_opts::equilibrate | equilibrate the system before solving (matrix A must be square) | |
solve_opts::likely_sympd | indicate that matrix A is likely symmetric/hermitian positive definite (sympd) | |
solve_opts::allow_ugly | keep solutions of systems that are singular to working precision | |
solve_opts::no_approx | do not find approximate solutions for rank deficient systems | |
solve_opts::force_sym | force use of the symmetric/hermitian solver (not limited to sympd matrices) | |
solve_opts::force_approx | force use of the approximate solver |
+
operator; for example: solve_opts::fast + solve_opts::no_approx
solve_opts::no_approx
option is not enabled, a warning is emitted and an approximate solution is attempted;
#define ARMA_WARN_LEVEL 1
#include <armadillo>
solve_opts::fast
will speed up finding the solution, but for poorly conditioned systems the solution may have lower qualitysolve_opts::likely_sympd
solve_opts::force_approx
is only advised if the system is known to be rank deficient; the approximate solver is considerably slowermat A(5, 5, fill::randu); vec b(5, fill::randu); mat B(5, 5, fill::randu); vec x1 = solve(A, b); vec x2; bool status = solve(x2, A, b); mat X1 = solve(A, B); mat X2 = solve(A, B, solve_opts::fast); // enable fast mode mat X3 = solve(trimatu(A), B); // indicate that A is triangular
"dc"
or "std"
"dc"
indicates divide-and-conquer method (default setting)
"std"
indicates standard method
mat X(5, 5, fill::randu); mat U; vec s; mat V; svd(U,s,V,X);
"both" | = | compute both left and right singular vectors (default operation) | |
"left" | = | compute only left singular vectors | |
"right" | = | compute only right singular vectors |
"dc"
or "std"
"dc"
indicates divide-and-conquer method (default setting)
"std"
indicates standard method
mat X(4, 5, fill::randu); mat U; vec s; mat V; svd_econ(U, s, V, X);
mat A(5, 5, fill::randu); mat B(5, 5, fill::randu); mat C(5, 5, fill::randu); mat X1 = syl(A, B, C); mat X2; syl(X2, A, B, C);
"lm" | = | obtain eigenvalues with largest magnitude (default operation) |
"sm" | = | obtain eigenvalues with smallest magnitude (see the caveats below) |
"la" | = | obtain eigenvalues with largest algebraic value |
"sa" | = | obtain eigenvalues with smallest algebraic value |
struct eigs_opts { double tol; // default: 0 unsigned int maxiter; // default: 1000 unsigned int subdim; // default: max(2*k+1, 20) };
k < subdim ≤ X.n_rows
; recommended value is subdim ≥ 2*k
"sm"
form, use the shift-invert mode with sigma set to 0.0// generate sparse matrix sp_mat A = sprandu<sp_mat>(1000, 1000, 0.1); sp_mat B = A.t()*A; vec eigval; mat eigvec; eigs_sym(eigval, eigvec, B, 5); // find 5 eigenvalues/eigenvectors eigs_opts opts; opts.maxiter = 10000; // increase max iterations to 10000 eigs_sym(eigval, eigvec, B, 5, "lm", opts);
"lm" | = | obtain eigenvalues with largest magnitude (default operation) |
"sm" | = | obtain eigenvalues with smallest magnitude (see the caveats below) |
"lr" | = | obtain eigenvalues with largest real part |
"sr" | = | obtain eigenvalues with smallest real part |
"li" | = | obtain eigenvalues with largest imaginary part |
"si" | = | obtain eigenvalues with smallest imaginary part |
struct eigs_opts { double tol; // default: 0 unsigned int maxiter; // default: 1000 unsigned int subdim; // default: max(2*k+1, 20) };
k + 2 < subdim ≤ X.n_rows
; recommended value is subdim ≥ 2*k + 1
"sm"
form, use the shift-invert mode with sigma set to 0.0// generate sparse matrix sp_mat A = sprandu<sp_mat>(1000, 1000, 0.1); cx_vec eigval; cx_mat eigvec; eigs_gen(eigval, eigvec, A, 5); // find 5 eigenvalues/eigenvectors eigs_opts opts; opts.maxiter = 10000; // increase max iterations to 10000 eigs_gen(eigval, eigvec, A, 5, "lm", opts);
⎡ zeros(X.n_rows, X.n_rows) X ⎤
⎣ X.t() zeros(X.n_cols, X.n_cols) ⎦
sp_mat X = sprandu<sp_mat>(100, 200, 0.1); mat U; vec s; mat V; svds(U, s, V, X, 10);
"superlu"
or "lapack"
; by default "superlu"
is used
"superlu"
, ARMA_USE_SUPERLU must be enabled in config.hpp
"lapack"
, sparse matrix A is converted to a dense matrix before using the LAPACK solver; this considerably increases memory usage
struct superlu_opts { bool allow_ugly; // default: false bool equilibrate; // default: false bool symmetric; // default: false double pivot_thresh; // default: 1.0 permutation_type permutation; // default: superlu_opts::COLAMD refine_type refine; // default: superlu_opts::REF_NONE };
superlu_opts::NATURAL | natural ordering | |
superlu_opts::MMD_ATA | minimum degree ordering on structure of A.t() * A | |
superlu_opts::MMD_AT_PLUS_A | minimum degree ordering on structure of A.t() + A | |
superlu_opts::COLAMD | approximate minimum degree column ordering |
superlu_opts::REF_NONE | no refinement | |
superlu_opts::REF_SINGLE | iterative refinement in single precision | |
superlu_opts::REF_DOUBLE | iterative refinement in double precision | |
superlu_opts::REF_EXTRA | iterative refinement in extra precision |
sp_mat A = sprandu<sp_mat>(1000, 1000, 0.1); vec b(1000, fill::randu); mat B(1000, 5, fill::randu); vec x = spsolve(A, b); // solve one system mat X = spsolve(A, B); // solve several systems bool status = spsolve(x, A, b); // use default solver if(status == false) { cout << "no solution" << endl; } spsolve(x, A, b, "lapack" ); // use LAPACK solver spsolve(x, A, b, "superlu"); // use SuperLU solver superlu_opts opts; opts.allow_ugly = true; opts.equilibrate = true; spsolve(x, A, b, "superlu", opts);
sp_mat A = sprandu<sp_mat>(1000, 1000, 0.1); spsolve_factoriser SF; bool status = SF.factorise(A); if(status == false) { cout << "factorisation failed" << endl; } double rcond_value = SF.rcond(); vec B1(1000, fill::randu); vec B2(1000, fill::randu); vec X1; vec X2; bool solution1_ok = SF.solve(X1,B1); bool solution2_ok = SF.solve(X2,B2); if(solution1_ok == false) { cout << "couldn't find X1" << endl; } if(solution2_ok == false) { cout << "couldn't find X2" << endl; }
"full" | = | return the full convolution (default setting), with the size equal to A.n_elem + B.n_elem - 1 |
"same" | = | return the central part of the convolution, with the same size as vector A |
vec A(256, fill::randu); vec B(16, fill::randu); vec C = conv(A, B); vec D = conv(A, B, "same");
"full" | = | return the full convolution (default setting), with the size equal to size(A) + size(B) - 1 |
"same" | = | return the central part of the convolution, with the same size as matrix A |
mat A(256, 256, fill::randu); mat B(16, 16, fill::randu); mat C = conv2(A, B); mat D = conv2(A, B, "same");
#define ARMA_USE_FFTW3
#include <armadillo>
-lfftw3
)vec X(100, fill::randu); cx_vec Y = fft(X, 128);
mat A(100, 100, fill::randu); cx_mat B = fft2(A); cx_mat C = fft2(A, 128, 128);
"nearest" | = | interpolate using single nearest neighbour |
"linear" | = | linear interpolation between two nearest neighbours (default setting) |
"*nearest" | = | as per "nearest" , but faster by assuming that X and XI are monotonically increasing |
"*linear" | = | as per "linear" , but faster by assuming that X and XI are monotonically increasing |
vec x = linspace<vec>(0, 3, 20); vec y = square(x); vec xx = linspace<vec>(0, 3, 100); vec yy; interp1(x, y, xx, yy); // use linear interpolation by default interp1(x, y, xx, yy, "*linear"); // faster than "linear" interp1(x, y, xx, yy, "nearest");
"nearest" | = | interpolate using nearest neighbours |
"linear" | = | linear interpolation between nearest neighbours (default setting) |
mat Z; Z.load("input_image.pgm", pgm_binary); // load an image in pgm format vec X = regspace(1, Z.n_cols); // X = horizontal spacing vec Y = regspace(1, Z.n_rows); // Y = vertical spacing vec XI = regspace(X.min(), 1.0/2.0, X.max()); // magnify by approx 2 vec YI = regspace(Y.min(), 1.0/3.0, Y.max()); // magnify by approx 3 mat ZI; interp2(X, Y, Z, XI, YI, ZI); // use linear interpolation by default ZI.save("output_image.pgm", pgm_binary);
vec x = linspace<vec>(0,4*datum::pi,100); vec y = cos(x); vec p = polyfit(x,y,10);
vec x1 = linspace<vec>(0,4*datum::pi,100); vec y1 = cos(x1); vec p1 = polyfit(x1,y1,10); vec y2 = polyval(p1,x1);
mean( V )
mean( M ) mean( M, dim ) mean( Q ) mean( Q, dim ) |
⎫ ⎪ ⎬ mean (average value) ⎪ ⎭ |
|
median( V )
median( M ) median( M, dim ) |
⎫ ⎬ median ⎭ |
|
stddev( V )
stddev( V, norm_type ) stddev( M ) stddev( M, norm_type ) stddev( M, norm_type, dim ) |
⎫ ⎪ ⎬ standard deviation ⎪ ⎭ |
|
var( V )
var( V, norm_type ) var( M ) var( M, norm_type ) var( M, norm_type, dim ) |
⎫ ⎪ ⎬ variance ⎪ ⎭ |
|
range( V )
range( M ) range( M, dim ) |
⎫ ⎬ range (difference between max and min) ⎭ |
mat A(5, 5, fill::randu); mat B = mean(A); mat C = var(A); double m = mean(mean(A)); vec v(5, fill::randu); double x = var(v); ivec w = {1, 2, 3, 4}; // integer vector double y = mean( conv_to<vec>::from(w) );
mat X(4, 5, fill::randu); mat Y(4, 5, fill::randu); mat C = cov(X,Y); mat D = cov(X,Y, 1);
mat X(4, 5, fill::randu); mat Y(4, 5, fill::randu); mat R = cor(X,Y); mat S = cor(X,Y, 1);
vec v(1000, fill::randn); // Gaussian distribution uvec h1 = hist(v, 11); uvec h2 = hist(v, linspace<vec>(-2,2,11));
vec v(1000, fill::randn); // Gaussian distribution uvec h = histc(v, linspace<vec>(-2,2,11));
vec V(1000, fill::randn); // Gaussian distribution vec P = { 0.25, 0.50, 0.75 }; vec Q = quantile(V, P);
mat A(5, 4, fill::randu); mat coeff; mat score; vec latent; vec tsquared; princomp(coeff, score, latent, tsquared, A);
1 ⎧ (x-m)2 ⎫
y = ‒‒‒‒ exp ⎪ −0.5 ‒‒‒‒‒‒ ⎪
s·√τ ⎩ s2 ⎭
vec X(10, fill::randu); vec M(10, fill::randu); vec S(10, fill::randu); vec P1 = normpdf(X); vec P2 = normpdf(X, M, S ); vec P3 = normpdf(1.23, M, S ); vec P4 = normpdf(X, 4.56, 7.89); double P5 = normpdf(1.23, 4.56, 7.89);
⎧ 1 ⎧ (x-m)2 ⎫ ⎫
y = log ⎪ ‒‒‒‒ exp ⎪ −0.5 ‒‒‒‒‒‒ ⎪ ⎪
⎩ s·√τ ⎩ s2 ⎭ ⎭
⎧ (x-m)2 ⎫
= −log(s·√τ) + ⎪ −0.5 ‒‒‒‒‒‒ ⎪
⎩ s2 ⎭
vec X(10, fill::randu); vec M(10, fill::randu); vec S(10, fill::randu); vec P1 = log_normpdf(X); vec P2 = log_normpdf(X, M, S ); vec P3 = log_normpdf(1.23, M, S ); vec P4 = log_normpdf(X, 4.56, 7.89); double P5 = log_normpdf(1.23, 4.56, 7.89);
vec X(10, fill::randu); vec M(10, fill::randu); vec S(10, fill::randu); vec P1 = normcdf(X); vec P2 = normcdf(X, M, S ); vec P3 = normcdf(1.23, M, S ); vec P4 = normcdf(X, 4.56, 7.89); double P5 = normcdf(1.23, 4.56, 7.89);
vec M(5, fill::randu); mat B(5, 5, fill::randu); mat C = B.t() * B; mat X = mvnrnd(M, C, 100);
mat X = chi2rnd(2, 5, 6); mat A = randi<mat>(5, 6, distr_param(1, 10)); mat B = chi2rnd(A);
mat X(5, 5, fill::randu); mat S = X.t() * X; mat W = wishrnd(S, 6.7);
mat X(5, 5, fill::randu); mat T = X.t() * X; mat W = iwishrnd(T, 6.7);
X(scalar) | |
update the statistics using the given scalar |
X.min() | |
current minimum value |
X.max() | |
current maximum value |
X.range() | |
current range |
X.mean() | |
current mean or average value |
X.var() and X.var(norm_type) | |
current variance |
X.stddev() and X.stddev(norm_type) | |
current standard deviation |
X.reset() | |
reset all statistics and set the number of samples to zero |
X.count() | |
current number of samples |
running_stat<double> stats; for(uword i=0; i<10000; ++i) { double sample = randn(); stats(sample); } cout << "mean = " << stats.mean() << endl; cout << "var = " << stats.var() << endl; cout << "min = " << stats.min() << endl; cout << "max = " << stats.max() << endl;
X(vector) | |
update the statistics using the given vector |
X.min() | |
vector of current minimum values |
X.max() | |
vector of current maximum values |
X.range() | |
vector of current ranges |
X.mean() | |
vector of current means |
X.var() and X.var(norm_type) | |
vector of current variances |
X.stddev() and X.stddev(norm_type) | |
vector of current standard deviations |
X.cov() and X.cov(norm_type) | |
matrix of current covariances; valid if calc_cov=true during construction of running_stat_vec |
X.reset() | |
reset all statistics and set the number of samples to zero |
X.count() | |
current number of samples |
running_stat_vec<vec> stats; vec sample; for(uword i=0; i<10000; ++i) { sample = randu<vec>(5); stats(sample); } cout << "mean = " << endl << stats.mean() << endl; cout << "var = " << endl << stats.var() << endl; cout << "max = " << endl << stats.max() << endl; // // running_stat_vec<rowvec> more_stats(true); for(uword 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() / (sd.t() * sd);
keep_existing | use the centroids specified in the means matrix as the starting point | |
static_subset | use a subset of the data vectors (repeatable) | |
random_subset | use a subset of the data vectors (random) | |
static_spread | use a maximally spread subset of data vectors (repeatable) | |
random_spread | use a maximally spread subset of data vectors (random start) |
static_spread
and random_spread
can be much more time consuming than with static_subset
and random_subset
uword d = 5; // dimensionality uword N = 10000; // number of vectors mat data(d, N, fill::randu); mat means; bool status = kmeans(means, data, 2, random_subset, 10, true); if(status == false) { cout << "clustering failed" << endl; } means.print("means:");
n_gaus-1 | |||
p(x) | = | ∑ | hg N( x | mg , Cg ) |
g=0 |
M.log_p(V) | |
return a scalar representing the log-likelihood of vector V (of type vec) | |||||||||||||||
M.log_p(V, g) | |
return a scalar representing the log-likelihood of vector V (of type vec), according to Gaussian with index g | |||||||||||||||
|
|
|
|||||||||||||||
M.log_p(X) | |
return a row vector (of type rowvec ) containing log-likelihoods of each column vector in matrix X (of type mat)
|
|||||||||||||||
M.log_p(X, g) | |
return a row vector (of type rowvec ) containing log-likelihoods of each column vector in matrix X (of type mat), according to Gaussian with index g
|
|||||||||||||||
|
|
|
|||||||||||||||
M.sum_log_p(X) | |
return a scalar representing the sum of log-likelihoods for all column vectors in matrix X (of type mat) | |||||||||||||||
M.sum_log_p(X, g) | |
return a scalar representing the sum of log-likelihoods for all column vectors in matrix X (of type mat), according to Gaussian with index g | |||||||||||||||
|
|
|
|||||||||||||||
M.avg_log_p(X) | |
return a scalar representing the average log-likelihood of all column vectors in matrix X (of type mat) | |||||||||||||||
M.avg_log_p(X, g) | |
return a scalar representing the average log-likelihood of all column vectors in matrix X (of type mat), according to Gaussian with index g | |||||||||||||||
|
|
|
|||||||||||||||
M.assign(V, dist_mode) | |
return the index of the closest mean (or Gaussian) to vector V (of type vec); parameter dist_mode is one of:
|
|||||||||||||||
M.assign(X, dist_mode) | |
return a row vector (of type urowvec) containing the indices of the closest means (or Gaussians) to each column vector in matrix X (of type mat); parameter dist_mode is either eucl_dist or prob_dist (as per the .assign() function above)
|
|||||||||||||||
|
|
|
|||||||||||||||
M.raw_hist(X, dist_mode) | |
return a row vector (of type urowvec ) representing the raw histogram of counts;
each entry is the number of counts corresponding to a Gaussian;
each count is the number times the corresponding Gaussian was the closest to each column vector in matrix X;parameter dist_mode is either eucl_dist or prob_dist (as per the .assign() function above)
|
|||||||||||||||
M.norm_hist(X, dist_mode) | |
similar to the .raw_hist() function above;
return a row vector (of type rowvec ) containing normalised counts;
the vector sums to one;parameter dist_mode is either eucl_dist or prob_dist (as per the .assign() function above)
|
|||||||||||||||
|
|
|
|||||||||||||||
M.generate() | |
return a column vector (of type vec) representing a random sample generated according to the model's parameters | |||||||||||||||
M.generate(N) | |
return a matrix (of type mat) containing N column vectors, with each vector representing a random sample generated according to the model's parameters | |||||||||||||||
|
|
|
|||||||||||||||
M.save(filename) | |
save the model to a file and return a bool indicating either success (true) or failure (false) | |||||||||||||||
M.load(filename) | |
load the model from a file and return a bool indicating either success (true) or failure (false) | |||||||||||||||
|
|
|
|||||||||||||||
M.n_gaus() | |
return the number of means/Gaussians in the model | |||||||||||||||
M.n_dims() | |
return the dimensionality of the means/Gaussians in the model | |||||||||||||||
|
|
|
|||||||||||||||
M.reset(n_dims, n_gaus) | |
set the model to have dimensionality n_dims, with n_gaus number of Gaussians; all the means are set to zero, all covariance matrix representations are equivalent to the identity matrix, and all the hefts (weights) are set to be uniform |
|||||||||||||||
|
|
|
|||||||||||||||
M.hefts | |
read-only row vector (of type rowvec) containing the hefts (weights) | |||||||||||||||
M.means | |
read-only matrix (of type mat) containing the means (centroids), stored as column vectors | |||||||||||||||
|
|
|
|||||||||||||||
M.dcovs
[only in gmm_diag] |
|
read-only matrix (of type mat) containing the representation of diagonal covariance matrices, with the set of diagonal covariances for each Gaussian stored as a column vector; applicable only to the gmm_diag class | |||||||||||||||
|
|
|
|||||||||||||||
M.fcovs
[only in gmm_full] |
|
read-only cube containing the full covariance matrices, with each covariance matrix stored as a slice within the cube; applicable only to the gmm_full class | |||||||||||||||
|
|
|
|||||||||||||||
M.set_hefts(V) | |
set the hefts (weights) of the model to be as specified in row vector V (of type rowvec); the number of hefts must match the existing model | |||||||||||||||
M.set_means(X) | |
set the means to be as specified in matrix X (of type mat); the number of means and their dimensionality must match the existing model | |||||||||||||||
|
|
|
|||||||||||||||
M.set_dcovs(X)
[only in gmm_diag] |
|
set the diagonal covariances matrices to be as specified in matrix X (of type mat), with the set of diagonal covariances for each Gaussian stored as a column vector; the number of covariance matrices and their dimensionality must match the existing model; applicable only to the gmm_diag class | |||||||||||||||
|
|
|
|||||||||||||||
M.set_fcovs(X)
[only in gmm_full] |
|
set the full covariances matrices to be as specified in cube X, with each covariance matrix stored as a slice within the cube; the number of covariance matrices and their dimensionality must match the existing model; applicable only to the gmm_full class | |||||||||||||||
|
|
|
|||||||||||||||
M.set_params(means, covs, hefts) | |
set all the parameters at the same time; the type and layout of the parameters is as per the .set_hefts(), .set_means(), .set_dcovs() and .set_fcovs() functions above; the number of Gaussians and dimensionality can be different from the existing model | |||||||||||||||
|
|
|
|||||||||||||||
M.learn(data, n_gaus, dist_mode, seed_mode, km_iter, em_iter, var_floor, print_mode) learn the model parameters via multi-threaded k-means and/or EM algorithms; return a bool value, with true indicating success, and false indicating failure;
the parameters have the following meanings:
|
|||||||||||||||||
data | matrix (of type mat) containing training samples; each sample is stored as a column vector | ||||||||||||||||
n_gaus | set the number of Gaussians to n_gaus; to help convergence, it is recommended that the given data matrix (above) contains at least 10 samples for each Gaussian | ||||||||||||||||
dist_mode |
specifies the distance used during the seeding of initial means and k-means clustering:
|
||||||||||||||||
seed_mode |
specifies how the initial means are seeded prior to running k-means and/or EM algorithms:
caveat: seeding the initial means with static_spread and random_spread
can be much more time consuming than with static_subset and random_subset
|
||||||||||||||||
km_iter | the number of iterations of the k-means algorithm; this is data dependent, but typically 10 iterations are sufficient | ||||||||||||||||
em_iter | the number of iterations of the EM algorithm; this is data dependent, but typically 5 to 10 iterations are sufficient | ||||||||||||||||
var_floor | the variance floor (smallest allowed value) for the diagonal covariances; setting this to a small non-zero value can help with convergence and/or better quality parameter estimates | ||||||||||||||||
print_mode |
either true or false ;
enable or disable printing of progress during the k-means and EM algorithms
|
// create synthetic data with 2 Gaussians uword d = 5; // dimensionality uword N = 10000; // number of vectors mat data(d, N, fill::zeros); vec mean0 = linspace<vec>(1,d,d); vec mean1 = mean0 + 2; uword i = 0; while(i < N) { if(i < N) { data.col(i) = mean0 + randn<vec>(d); ++i; } if(i < N) { data.col(i) = mean0 + randn<vec>(d); ++i; } if(i < N) { data.col(i) = mean1 + randn<vec>(d); ++i; } } // model the data as a diagonal GMM with 2 Gaussians gmm_diag model; bool status = model.learn(data, 2, maha_dist, random_subset, 10, 5, 1e-10, true); if(status == false) { cout << "learning failed" << endl; } model.means.print("means:"); double scalar_likelihood = model.log_p( data.col(0) ); rowvec set_likelihood = model.log_p( data.cols(0,9) ); double overall_likelihood = model.avg_log_p(data); uword gaus_id = model.assign( data.col(0), eucl_dist ); urowvec gaus_ids = model.assign( data.cols(0,9), prob_dist ); urowvec hist1 = model.raw_hist (data, prob_dist); rowvec hist2 = model.norm_hist(data, eucl_dist); model.save("my_model.gmm");
datum::tau
|
τ, the ratio of any circle's circumference to its radius; equivalent to 2π | |
datum::pi
|
π, the ratio of any circle's circumference to its diameter | |
datum::inf
|
∞, infinity | |
datum::nan
|
“not a number” (NaN); caveat: NaN is not equal to anything, even itself | |
datum::eps
|
machine epsilon; approximately 2.2204e-16; difference between 1 and the next representable value | |
datum::e
|
base of the natural logarithm | |
datum::sqrt2
|
square root of 2 | |
datum::log_min
|
log of minimum non-zero value (type and machine dependent) | |
datum::log_max
|
log of maximum value (type and machine dependent) | |
datum::euler
|
Euler's constant, aka Euler-Mascheroni constant | |
datum::gratio
|
golden ratio | |
datum::m_u
|
atomic mass constant (in kg) | |
datum::N_A
|
Avogadro constant | |
datum::k
|
Boltzmann constant (in joules per kelvin) | |
datum::k_evk
|
Boltzmann constant (in eV/K) | |
datum::a_0
|
Bohr radius (in meters) | |
datum::mu_B
|
Bohr magneton | |
datum::Z_0
|
characteristic impedance of vacuum (in ohms) | |
datum::G_0
|
conductance quantum (in siemens) | |
datum::k_e
|
Coulomb's constant (in meters per farad) | |
datum::eps_0
|
electric constant (in farads per meter) | |
datum::m_e
|
electron mass (in kg) | |
datum::eV
|
electron volt (in joules) | |
datum::ec
|
elementary charge (in coulombs) | |
datum::F
|
Faraday constant (in coulombs) | |
datum::alpha
|
fine-structure constant | |
datum::alpha_inv
|
inverse fine-structure constant | |
datum::K_J
|
Josephson constant | |
datum::mu_0
|
magnetic constant (in henries per meter) | |
datum::phi_0
|
magnetic flux quantum (in webers) | |
datum::R
|
molar gas constant (in joules per mole kelvin) | |
datum::G
|
Newtonian constant of gravitation (in newton square meters per kilogram squared) | |
datum::h
|
Planck constant (in joule seconds) | |
datum::h_bar
|
Planck constant over 2 pi, aka reduced Planck constant (in joule seconds) | |
datum::m_p
|
proton mass (in kg) | |
datum::R_inf
|
Rydberg constant (in reciprocal meters) | |
datum::c_0
|
speed of light in vacuum (in meters per second) | |
datum::sigma
|
Stefan-Boltzmann constant | |
datum::R_k
|
von Klitzing constant (in ohms) | |
datum::b
|
Wien wavelength displacement law constant |
datum::nan
is not equal to anything, even itself;
to check whether a scalar x is finite,
use std::isfinite(x)
cout << "speed of light = " << datum::c_0 << endl; cout << "log_max for floats = "; cout << fdatum::log_max << endl; cout << "log_max for doubles = "; cout << datum::log_max << endl;
.tic()
|
|
start the timer |
.toc()
|
|
return the number of seconds since the last call to .tic()
|
wall_clock timer; timer.tic(); // ... do something ... double n = timer.toc(); cout << "number of seconds: " << n << endl;
arma_rng::set_seed( value )
|
|
set the RNG seed to the specified value |
arma_rng::set_seed_random()
|
|
set the RNG seed to a value drawn from std::random_device (if the reported entropy is not zero), or /dev/urandom (on Linux and macOS), or based on the current time (on systems without /dev/urandom) |
#pragma omp parallel { arma_rng::set_seed(123); }
std::atomic<std::size_t> counter(0); #pragma omp parallel { arma_rng::set_seed(123 + counter++); }
std::cout
ARMA_COUT_STREAM
define; see config.hpp
std::cerr
ARMA_CERR_STREAM
define; see config.hpp
ARMA_WARN_LEVEL
define; see config.hpp
#define ARMA_WARN_LEVEL 1 #include <armadillo>
cx_double
|
|
equivalent to std::complex<double>
|
cx_float
|
|
equivalent to std::complex<float>
|
cx_mat X(5, 5, fill::randu); X(1,2) = cx_double(3.4, 5.6); cx_double val = X(2,3);
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( 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( span(p,q), span(r,s), span(t,u) )
|
|||
A'
|
A.t() or trans(A)
|
matrix transpose / Hermitian transpose
(for complex matrices, the conjugate of each element is 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; ]
|
A = { { 1, 2 },
|
element initialisation | ||
X = A(:)
|
X = vectorise(A)
|
|||
X = [ A B ]
|
X = join_horiz(A,B)
|
|||
X = [ A; B ]
|
X = join_vert(A,B)
|
|||
A
|
cout << A << endl;
or A.print("A =");
|
|||
save ‑ascii 'A.txt' A
|
A.save("A.txt", raw_ascii);
|
Matlab/Octave matrices saved as ascii are readable by Armadillo (and vice-versa) | ||
load ‑ascii 'A.txt'
|
A.load("A.txt", raw_ascii);
|
|||
A = randn(2,3);
|
mat A = randn(2,3);
|
fields store arbitrary objects, such as matrices |
#include <iostream> #include <armadillo> using namespace std; using namespace arma; int main() { mat A(4, 5, fill::randu); mat B(4, 5, fill::randu); cout << A*B.t() << endl; return 0; }
g++ example.cpp -o example -std=c++11 -O2 -larmadillo
ARMA_DONT_USE_WRAPPER
|
Disable going through the run-time Armadillo wrapper library (libarmadillo.so) when calling LAPACK, BLAS, ARPACK, and SuperLU functions.
You will need to directly link with BLAS, LAPACK, etc (eg. -lblas -llapack )
|
|||||||||||||
ARMA_USE_LAPACK
|
Enable use of LAPACK, or a high-speed replacement for LAPACK (eg. OpenBLAS, Intel MKL, or the Accelerate framework). Armadillo requires LAPACK for functions such as svd(), inv(), eig_sym(), solve(), etc. | |||||||||||||
ARMA_DONT_USE_LAPACK
|
Disable use of LAPACK; overrides ARMA_USE_LAPACK | |||||||||||||
ARMA_USE_BLAS
|
Enable use of BLAS, or a high-speed replacement for BLAS (eg. OpenBLAS, Intel MKL, or the Accelerate framework). BLAS is used for matrix multiplication. Without BLAS, Armadillo will use a built-in matrix multiplication routine, which might be slower for large matrices. | |||||||||||||
ARMA_DONT_USE_BLAS
|
Disable use of BLAS; overrides ARMA_USE_BLAS | |||||||||||||
ARMA_USE_NEWARP
|
Enable use of NEWARP (built-in alternative to ARPACK). This is used for the eigen decomposition of real (non-complex) sparse matrices, ie. eigs_gen(), eigs_sym() and svds(). Requires ARMA_USE_LAPACK to be enabled. If use of both NEWARP and ARPACK is enabled, NEWARP will be preferred. | |||||||||||||
ARMA_DONT_USE_NEWARP
|
Disable use of NEWARP (built-in alternative to ARPACK); overrides ARMA_USE_NEWARP | |||||||||||||
ARMA_USE_ARPACK
|
Enable use of ARPACK, or a high-speed replacement for ARPACK. Armadillo requires ARPACK for the eigen decomposition of complex sparse matrices, ie. eigs_gen(), eigs_sym() and svds(). If use of NEWARP is disabled, ARPACK will also be used for the eigen decomposition of real sparse matrices. | |||||||||||||
ARMA_DONT_USE_ARPACK
|
Disable use of ARPACK; overrides ARMA_USE_ARPACK | |||||||||||||
ARMA_USE_SUPERLU
|
Enable use of SuperLU, which is used by spsolve() for finding the solutions of sparse systems,
as well as eigs_sym() and eigs_gen() in shift-invert mode.
You will need to link with the superlu library, for example -lsuperlu
|
|||||||||||||
ARMA_DONT_USE_SUPERLU
|
Disable use of SuperLU; overrides ARMA_USE_SUPERLU | |||||||||||||
ARMA_USE_HDF5
|
Enable the ability to save and load matrices stored in the HDF5 format;
the hdf5.h header file must be available on your system and you will need to link with the hdf5 library (eg. -lhdf5 )
|
|||||||||||||
ARMA_DONT_USE_HDF5
|
Disable the use of the HDF5 library; overrides ARMA_USE_HDF5 | |||||||||||||
ARMA_USE_FFTW3
|
Enable use of the FFTW3 library by fft() and ifft();
you will need to link with the FFTW3 library (eg. -lfftw3 )
|
|||||||||||||
ARMA_DONT_USE_FFTW3
|
Disable the use of the FFTW3 library; overrides ARMA_USE_FFTW3 | |||||||||||||
ARMA_DONT_USE_STD_MUTEX
|
Disable use of std::mutex; applicable if your compiler and/or environment doesn't support std::mutex | |||||||||||||
ARMA_DONT_OPTIMISE_BAND
|
Disable automatically optimised handling of band matrices by solve() and chol() | |||||||||||||
ARMA_DONT_OPTIMISE_SYMPD
|
Disable automatically optimised handling of symmetric/hermitian positive definite matrices by solve(), inv(), pinv(), expmat(), logmat(), sqrtmat(), powmat(), rcond() | |||||||||||||
ARMA_USE_OPENMP
|
Use OpenMP for parallelisation of computationally expensive element-wise operations
(such as exp(), log(), cos(), etc).
Automatically enabled when using a compiler which has OpenMP 3.1+ active (eg. the -fopenmp option for gcc and clang).
|
|||||||||||||
ARMA_DONT_USE_OPENMP
|
Disable use of OpenMP for parallelisation of element-wise operations; overrides ARMA_USE_OPENMP | |||||||||||||
ARMA_OPENMP_THRESHOLD
|
The minimum number of elements in a matrix to enable OpenMP based parallelisation of computationally expensive element-wise functions; default value is 320 | |||||||||||||
ARMA_OPENMP_THREADS
|
The maximum number of threads for OpenMP based parallelisation of computationally expensive element-wise functions; default value is 8 | |||||||||||||
ARMA_BLAS_CAPITALS
|
Use capitalised (uppercase) BLAS and LAPACK function names (eg. DGEMM vs dgemm) | |||||||||||||
ARMA_BLAS_UNDERSCORE
|
Append an underscore to BLAS and LAPACK function names (eg. dgemm_ vs dgemm). Enabled by default. | |||||||||||||
ARMA_BLAS_LONG_LONG
|
Use "long long" instead of "int" when calling BLAS and LAPACK functions; the "long long" type is a 64 bit integer type on all platforms | |||||||||||||
ARMA_USE_FORTRAN_HIDDEN_ARGS
|
Use so-called "hidden arguments" when calling BLAS and LAPACK functions. Enabled by default. See Fortran argument passing conventions for more details. | |||||||||||||
ARMA_DONT_USE_FORTRAN_HIDDEN_ARGS
|
Disable use of so-called "hidden arguments" when calling BLAS and LAPACK functions.
May be necessary when using Armadillo in conjunction with broken MKL headers (eg. if you have #include "mkl_lapack.h" in your code).
|
|||||||||||||
ARMA_USE_TBB_ALLOC
|
Use Intel TBB scalable_malloc() and scalable_free() instead of standard malloc() and free() for managing matrix memory | |||||||||||||
ARMA_USE_MKL_ALLOC
|
Use Intel MKL mkl_malloc() and mkl_free() instead of standard malloc() and free() for managing matrix memory | |||||||||||||
ARMA_USE_MKL_TYPES
|
Use Intel MKL types for complex numbers.
You will need to include appropriate MKL headers before the Armadillo header.
You may also need to enable one or more of the following options:
ARMA_BLAS_LONG_LONG , ARMA_DONT_USE_FORTRAN_HIDDEN_ARGS
|
|||||||||||||
ARMA_64BIT_WORD
|
Use 64 bit integers for matrix and vector sizes.
Automatically enabled when using a 64-bit platform, except when using Armadillo in the R environment (via RcppArmadillo).
Useful if matrices/vectors capable of holding more than 4 billion elements are required.
This can also be enabled by adding #define ARMA_64BIT_WORD before each instance of #include <armadillo> .
See also the ARMA_BLAS_LONG_LONG option.
|
|||||||||||||
ARMA_MAT_PREALLOC
|
The number of pre-allocated elements used by matrices and vectors. Must be always enabled and set to an integer that is at least 1. By default set to 16. If you mainly use lots of very small vectors (eg. ≤ 4 elements), change the number to the size of your vectors. | |||||||||||||
ARMA_COUT_STREAM
|
The default stream used for printing matrices and cubes by .print(). Must be always enabled. By default defined to std::cout | |||||||||||||
ARMA_CERR_STREAM
|
The default stream used for printing warnings and errors. Must be always enabled. By default defined to std::cerr | |||||||||||||
ARMA_WARN_LEVEL
|
The level of warning messages printed to ARMA_CERR_STREAM.
Must be an integer ≥ 0. By default defined to 2.
Example usage: #define ARMA_WARN_LEVEL 1 #include <armadillo> |
solve_opts::force_sym
option to solve() to force use of the symmetric/hermitian solver (not limited to sympd matrices)csv_opts::strict
option to loading CSV files to interpret missing values as NaNcheck_for_zeros
option to form 4 of sparse matrix batch constructorsinv_opts::no_ugly
or inv_opts::allow_approx
now use a scaled threshold similar to pinv()set_cout_stream()
and set_cerr_stream()
are now no-ops; instead use the options ARMA_WARN_LEVEL, or ARMA_COUT_STREAM, or ARMA_CERR_STREAMinv_opts::no_ugly
option to inv() and inv_sympd() to disallow inverses of poorly conditioned matricesinv_opts::allow_approx
option in inv() and inv_sympd()inv_opts::allow_approx
option to inv() and inv_sympd() to allow approximate inverses of poorly conditioned matricesX.cols(first_col,last_col)
solve_opts::force_approx
option to force use of the approximate solverfill::value(scalar)
, eg. mat X(4,5,fill::value(123))
fill::none
specifier, eg. mat X(4,5,fill::none)
ARMA_WARN_LEVEL
configuration option, to control the degree of emitted warning messagesARMA_DONT_USE_CXX11_MUTEX
configuration option to disable use of std::mutexsolve_opts::refine
to explicitly enable refinementsolve_opts::likely_sympd
to indicate that the given matrix is likely positive definite#define ARMA_DONT_USE_FORTRAN_HIDDEN_ARGS
before #include <armadillo>
"strictascend"
and "strictdescend"
-std=c++11 -fopenmp
-march=native
in conjunction with -fopenmp
may lead to speed regressions on recent processorscube 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;