iRoCS Toolbox
1.1.0
|
#include <config.hh>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_eigen.h>
#include <gsl/gsl_permute_vector.h>
#include <gsl/gsl_sort_vector.h>
#include <gsl/gsl_sort.h>
#include <gsl/gsl_permutation.h>
#include "TypeTraits.hh"
#include "ATBGSLWrapper.hh"
#include "Quaternion.hh"
#include "ATBLinAlg.icc"
Go to the source code of this file.
Namespaces | |
atb | |
Enumerations | |
enum | atb::SortingType { atb::None, atb::Ascending, atb::Descending, atb::AbsAscending, atb::AbsDescending } |
Functions | |
template<int Dim> | |
void | atb::computeEigenvectors (blitz::TinyMatrix< double, Dim, Dim > const &A, blitz::TinyMatrix< double, Dim, Dim > &U, blitz::TinyVector< double, Dim > const &lambda) |
Computates the eigenvectors of the real symmetric matrix A given the eigenvalues. More... | |
template<int Dim> | |
double | atb::determinant (blitz::TinyMatrix< double, Dim, Dim > const &A) |
Compute the determinant of matrix A. More... | |
template<typename DataT > | |
double | atb::dot (blitz::Array< DataT, 1 > const &vec1, blitz::Array< DataT, 1 > const &vec2) |
Calculation of the inner product (dot product) of the given vectors. More... | |
template<int Dim> | |
void | atb::eigenvalueDecompositionRealSymmetric (blitz::TinyMatrix< double, Dim, Dim > const &A, blitz::TinyVector< double, Dim > &lambda, SortingType sort=Descending) |
Computation of the eigen value decomposition of real symmetric matrices. More... | |
template<int Dim> | |
void | atb::eigenvalueDecompositionRealSymmetric (blitz::TinyMatrix< double, Dim, Dim > const &A, blitz::TinyMatrix< double, Dim, Dim > &U, blitz::TinyVector< double, Dim > &lambda, SortingType sort=Descending) |
Computation of the eigen value decomposition of real symmetric matrices. More... | |
void | atb::eigenvalueDecompositionRealSymmetric (blitz::Array< double, 2 > const &A, blitz::Array< double, 2 > &U, blitz::Array< double, 1 > &lambda, SortingType sort=Descending) |
Computation of the eigen value decomposition of real symmetric matrices. More... | |
template<typename DataT , int Dim> | |
blitz::TinyVector< DataT, Dim+1 > | atb::euclideanToHomogeneous (blitz::TinyVector< DataT, Dim > const &pos) |
Convert the given vector to homogeneous coordinates. More... | |
blitz::TinyVector< double, 3 > | atb::euclideanToSpherical (blitz::TinyVector< double, 3 > const &eucl, blitz::Array< double, 2 > const &trafo) |
Get the spherical (theta, phi, r) coordinates of the given euclidean point. More... | |
blitz::TinyVector< double, 3 > | atb::euclideanToSpherical (blitz::TinyVector< double, 3 > const &eucl) |
Get the spherical (theta, phi, r) coordinates of the given euclidean point. More... | |
double | atb::factorialFraction (int n, int k) |
Compute the value of expression ![]() | |
template<typename DataT , int Dim> | |
blitz::TinyVector< DataT, Dim-1 > | atb::homogeneousToEuclidean (blitz::TinyVector< DataT, Dim > const &pos) |
Convert the given homogeneous vector to euclidean coordinates. More... | |
void | atb::invert (blitz::Array< double, 2 > const &A, blitz::Array< double, 2 > &Ainv) |
Compute the inverse of matrix A. More... | |
blitz::Array< double, 2 > | atb::invert (blitz::Array< double, 2 > const &A) |
Compute the inverse of matrix A. More... | |
template<int Dim> | |
blitz::TinyMatrix< double, Dim, Dim > | atb::invert (blitz::TinyMatrix< double, Dim, Dim > const &A) |
Compute the inverse of matrix A. More... | |
template<typename DataT > | |
blitz::Array< DataT, 2 > | atb::mmMult (blitz::Array< DataT, 2 > const &A, blitz::Array< DataT, 2 > const &B) |
Naive matrix multiplication A * B = C. More... | |
template<typename DataT , int Dim1, int Dim2, int Dim3> | |
blitz::TinyMatrix< DataT, Dim1, Dim3 > | atb::mmMult (blitz::TinyMatrix< DataT, Dim1, Dim2 > const &A, blitz::TinyMatrix< DataT, Dim2, Dim3 > const &B) |
Naive matrix multiplication A * B = C. More... | |
ptrdiff_t | atb::modulo (ptrdiff_t a, ptrdiff_t b) |
Strictly positive modulo, to compute valid indices for Array access modulo Array size. More... | |
template<int Dim> | |
blitz::TinyVector< ptrdiff_t, Dim > | atb::modulo (blitz::TinyVector< ptrdiff_t, Dim > const &a, blitz::TinyVector< ptrdiff_t, Dim > const &b) |
Strictly positive vectorial modulo, to compute valid indices for Array access modulo Array size. More... | |
template<typename MatrixT , typename BaseT , int Dim> | |
blitz::TinyVector< BaseT, Dim > | atb::mvMult (MatrixT const &m, blitz::TinyVector< BaseT, Dim > const &v) |
Quadratic Matrix * Vector product (generic). More... | |
template<typename BaseT , int Dim> | |
blitz::TinyVector< BaseT, Dim > | atb::mvMult (blitz::Array< BaseT, 2 > const &m, blitz::TinyVector< BaseT, Dim > const &v) |
Quadratic Matrix * Vector product. More... | |
template<typename BaseT > | |
blitz::Array< BaseT, 1 > | atb::mvMult (blitz::Array< BaseT, 2 > const &m, blitz::Array< BaseT, 1 > const &v) |
Matrix * Vector product. More... | |
template<typename DataT , int Dim1, int Dim2> | |
blitz::TinyMatrix< DataT, Dim1, Dim2 > | atb::outerProduct (blitz::TinyVector< DataT, Dim1 > const &vec1, blitz::TinyVector< DataT, Dim2 > const &vec2) |
Calculation of the outer product of two vectors. More... | |
template<typename DataT > | |
blitz::TinyVector< DataT, 3 > | atb::rotate (blitz::TinyVector< DataT, 3 > const &vec, blitz::TinyVector< DataT, 3 > const &axis, DataT angle) |
Rotate the given 3-D vector by a specified angle angle around the given axis. More... | |
template<typename DataT > | |
void | atb::sortElements (blitz::TinyVector< DataT, 3 > &v, SortingType sort) |
Sort the values of the given TinyVector according to the given SortingType. More... | |
blitz::TinyVector< double, 3 > | atb::sphericalToEuclidean (blitz::TinyVector< double, 3 > const &spherical, blitz::Array< double, 2 > const &trafo) |
Get the euclidean (z, y, x) coordinates of the given point on the 2-sphere. More... | |
blitz::TinyVector< double, 3 > | atb::sphericalToEuclidean (blitz::TinyVector< double, 3 > const &spherical) |
Get the euclidean (z, y, x) coordinates of the given point on the 2-sphere. More... | |
template<typename DataT , int Dim1, int Dim2> | |
blitz::TinyMatrix< DataT, Dim2, Dim1 > | atb::transpose (blitz::TinyMatrix< DataT, Dim1, Dim2 > const &in) |
Compute the transpose of the given Matrix. More... | |
template<typename DataT > | |
blitz::Array< DataT, 2 > | atb::transpose (blitz::Array< DataT, 2 > const &in) |
Compute the transpose of the given Matrix. More... | |