Several mathematical, linear-algebra, or other functions within Shark are not part of any particular class. They are collected here in the doxygen group "shark_globals". More...
Namespaces | |
| namespace | shark |
| AbstractSingleObjectiveOptimizer. | |
Macros | |
| #define | SHARK_LINALG_EIGENSORT_INL |
| Inline implementation of eigenvalue/-vector sorting. More... | |
Functions | |
| template<class T > | |
| T | shark::maxExpInput () |
| Maximum allowed input value for exp. More... | |
| template<class T > | |
| T | shark::minExpInput () |
| Minimum value for exp(x) allowed so that it is not 0. More... | |
| template<class T > | |
| boost::enable_if < boost::is_arithmetic< T >, T > ::type | shark::sqr (const T &x) |
| Calculates x^2. More... | |
| template<class T > | |
| T | shark::cube (const T &x) |
| Calculates x^3. More... | |
| template<class T > | |
| boost::enable_if < boost::is_arithmetic< T >, T > ::type | shark::sigmoid (T x) |
| Logistic function/logistic function. More... | |
| template<class T > | |
| T | shark::safeExp (T x) |
| Thresholded exp function, over- and underflow safe. More... | |
| template<class T > | |
| T | shark::safeLog (T x) |
| Thresholded log function, over- and underflow safe. More... | |
| template<class T > | |
| boost::enable_if < boost::is_arithmetic< T >, T > ::type | shark::softPlus (T x) |
| Numerically stable version of the function log(1+exp(x)). More... | |
| template<class T > | |
| T | shark::copySign (T x, T y) |
Variables | |
| static const double | shark::SQRT_2_PI = boost::math::constants::root_two_pi<double>() |
| Constant for sqrt( 2 * pi ). More... | |
| enum | shark::LabelPosition { shark::FIRST_COLUMN, shark::LAST_COLUMN } |
| Position of the label in a CSV file. More... | |
| template<typename Type > | |
| void | shark::import_csv (Data< Type > &data, std::string fn, std::string separator=",", std::string comment="", std::size_t batchSize=Data< Type >::DefaultBatchSize) |
| Import unlabeled data from a character-separated value file. More... | |
| template<typename Type > | |
| void | shark::export_csv (Data< Type > const &set, std::string fn, std::string separator=",", bool sci=true, unsigned int width=0) |
| Format unlabeled data into a character-separated value file. More... | |
| template<typename InputType , typename LabelType > | |
| void | shark::import_csv (LabeledData< InputType, LabelType > &dataset, std::string fn, LabelPosition lp, std::string separator=",", std::string comment="", bool allowMissingClasses=false, std::map< LabelType, LabelType > const *labelmap=NULL, std::size_t batchSize=LabeledData< InputType, LabelType >::DefaultBatchSize) |
| Import labeled data from a character-separated value file. More... | |
| void | shark::import_csv (LabeledData< RealVector, RealVector > &dataset, std::string fn, LabelPosition lp, std::string separator=",", std::string comment="", std::size_t numberOfOutputs=1) |
| Import regression data from a character-separated value file. More... | |
| template<typename InputType , typename LabelType > | |
| void | shark::export_csv (LabeledData< InputType, LabelType > const &dataset, std::string fn, LabelPosition lp, std::string separator=",", bool sci=true, unsigned int width=0) |
| Format labeled data into a character-separated value file. More... | |
| template<typename InputType > | |
| void | shark::string2data (Data< InputType > &data, const std::string &dataInString, std::size_t batchSize=Data< InputType >::DefaultBatchSize) |
| Construct Shark data from a string. More... | |
| template<typename InputType , typename LabelType > | |
| void | shark::string2data (LabeledData< InputType, LabelType > &dataset, const std::string &dataInString, LabelPosition labelPosition=LAST_COLUMN, bool allowMissingFeatures=false, bool allowMissingClasses=false, std::map< LabelType, LabelType > const *labelmap=NULL, std::size_t batchSize=LabeledData< InputType, LabelType >::DefaultBatchSize) |
| Construct Shark labeled data from a string. More... | |
| template<class I , class L > | |
| CVFolds< LabeledData< I, L > > | shark::createCVIID (LabeledData< I, L > &set, size_t numberOfPartitions, std::size_t batchSize=Data< I >::DefaultBatchSize) |
| Create a partition for cross validation. More... | |
| template<class I , class L > | |
| CVFolds< LabeledData< I, L > > | shark::createCVSameSize (LabeledData< I, L > &set, std::size_t numberOfPartitions, std::size_t batchSize=LabeledData< I, L >::DefaultBatchSize) |
| Create a partition for cross validation. More... | |
| template<class I > | |
| CVFolds< LabeledData< I, unsigned int > > | shark::createCVSameSizeBalanced (LabeledData< I, unsigned int > &set, size_t numberOfPartitions, std::size_t batchSize=Data< I >::DefaultBatchSize) |
| Create a partition for cross validation. More... | |
| template<class I > | |
| CVFolds< LabeledData< I, RealVector > > | shark::createCVSameSizeBalanced (LabeledData< I, RealVector > &set, size_t numberOfPartitions, std::size_t batchSize=Data< I >::DefaultBatchSize) |
| Create a partition for cross validation. More... | |
| template<class I , class L > | |
| CVFolds< LabeledData< I, L > > | shark::createCVIndexed (LabeledData< I, L > &set, size_t numberOfPartitions, std::vector< size_t > indices, std::size_t batchSize=Data< I >::DefaultBatchSize) |
| Create a partition for cross validation from indices. More... | |
| template<class T > | |
| std::ostream & | shark::operator<< (std::ostream &stream, const Data< T > &d) |
| Outstream of elements. More... | |
| template<class Range > | |
| Data< typename boost::range_value< Range > ::type > | shark::createDataFromRange (Range const &inputs, std::size_t maximumBatchSize=0) |
| creates a data object from a range of elements More... | |
| void | shark::import_libsvm (LabeledData< RealVector, unsigned int > &dataset, std::istream &stream, int highestIndex=0) |
| Import data from a LIBSVM file. More... | |
| void | shark::import_libsvm (LabeledData< CompressedRealVector, unsigned int > &dataset, std::istream &stream, int highestIndex=0) |
| Import data from a LIBSVM file. More... | |
| void | shark::import_libsvm (LabeledData< RealVector, unsigned int > &dataset, std::string fn, int highestIndex=0) |
| Import data from a LIBSVM file. More... | |
| void | shark::import_libsvm (LabeledData< CompressedRealVector, unsigned int > &dataset, std::string fn, int highestIndex=0) |
| Import data from a LIBSVM file. More... | |
| template<typename InputType > | |
| void | shark::export_libsvm (LabeledData< InputType, unsigned int > &dataset, const std::string &fn, bool dense=false, bool oneMinusOne=true, bool sortLabels=false) |
| Export data to LIBSVM format. More... | |
| void | shark::detail::writePGM (const char *fileName, const unsigned char *pData, const unsigned int sx, const unsigned int sy) |
| Writes a PGM file. More... | |
| enum | shark::KernelMatrixNormalizationType { shark::NONE, shark::MULTIPLICATIVE_TRACE_ONE, shark::MULTIPLICATIVE_TRACE_N, shark::MULTIPLICATIVE_VARIANCE_ONE, shark::CENTER_ONLY, shark::CENTER_AND_MULTIPLICATIVE_TRACE_ONE } |
| template<typename InputType , typename LabelType > | |
| void | shark::export_kernel_matrix (LabeledData< InputType, LabelType > const &dataset, AbstractKernelFunction< InputType > &kernel, std::ostream &out, KernelMatrixNormalizationType normalizer=NONE, bool scientific=false, unsigned int fieldwidth=0) |
| Write a kernel Gram matrix to stream. More... | |
| template<typename InputType , typename LabelType > | |
| void | shark::export_kernel_matrix (LabeledData< InputType, LabelType > const &dataset, AbstractKernelFunction< InputType > &kernel, std::string fn, KernelMatrixNormalizationType normalizer=NONE, bool sci=false, unsigned int width=0) |
| Write a kernel Gram matrix to file. More... | |
| template<class VectorT , class Function > | |
| void | shark::lnsrch (const VectorT &xold, double fold, VectorT &g, VectorT &p, VectorT &x, double &f, double stpmax, bool &check, Function func) |
| Does a line search, i.e. given a nonlinear function, a starting point and a direction, a new point is calculated where the function has decreased "sufficiently". More... | |
| template<class VectorT , class Function > | |
| void | shark::linmin (VectorT &p, const VectorT &xi, double &fret, Function func, double ax=0.0, double bx=1.0) |
| Minimizes a function of "N" variables. More... | |
| template<class VectorT , class VectorU , class DifferentiableFunction > | |
| void | shark::dlinmin (VectorT &p, const VectorU &xi, double &fret, DifferentiableFunction &func, double ax=0.0, double bx=1.0) |
| template<class Source > | |
| detail::ADLVector< Source & > | shark::init (shark::blas::vector_container< Source > &source) |
| Starting-point for the initialization sequence. More... | |
| template<class Source > | |
| detail::ADLVector< const Source & > | shark::init (const shark::blas::vector_container< Source > &source) |
| Starting-point for the initialization sequence. More... | |
| template<class Source > | |
| detail::ADLVector < shark::blas::vector_range < Source > > | shark::init (const shark::blas::vector_range< Source > &source) |
| Starting-point for the initialization sequence when used for splitting the vector. More... | |
| template<class Source > | |
| detail::ADLVector < shark::blas::matrix_row < Source > > | shark::init (const shark::blas::matrix_row< Source > &source) |
| Specialization for matrix rows. More... | |
| template<class Matrix > | |
| detail::MatrixExpression < const Matrix > | shark::toVector (const shark::blas::matrix_expression< Matrix > &matrix) |
| Linearizes a matrix as a set of row vectors and treats them as a set of vectors for initialization. More... | |
| template<class Matrix > | |
| detail::MatrixExpression< Matrix > | shark::toVector (shark::blas::matrix_expression< Matrix > &matrix) |
| Linearizes a matrix as a set of row vectors and treats them as a set of vectors for initialization. More... | |
| template<class T > | |
| detail::ParameterizableExpression < const T > | shark::parameters (const T &object) |
| Uses the parameters of a parameterizable object for initialization. More... | |
| template<class T > | |
| detail::ParameterizableExpression < T > | shark::parameters (T &object) |
| Uses the parameters of a parameterizable object for initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::const_iterator, detail::VectorExpression < const typename T::value_type & > > | shark::vectorSet (const T &range) |
| Uses a range of vectors for initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::iterator, detail::VectorExpression < typename T::value_type & > > | shark::vectorSet (T &range) |
| Uses a range of vectors for splitting and initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::const_iterator, detail::MatrixExpression < const typename T::value_type > > | shark::matrixSet (const T &range) |
| Uses a range of vectors for initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::iterator, detail::MatrixExpression < typename T::value_type > > | shark::matrixSet (T &range) |
| Uses a range of vectors for splitting and initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::const_iterator, detail::ParameterizableExpression < const typename T::value_type > > | shark::parameterSet (const T &range) |
| Uses a range of parametrizable objects for initialization. More... | |
| template<class T > | |
| detail::InitializerRange < typename T::iterator, detail::ParameterizableExpression < typename T::value_type > > | shark::parameterSet (T &range) |
| Uses a range of parametrizable objects for splitting and initialization. More... | |
| template<class VectorT , class VectorU , class WeightT > | |
| VectorT::value_type | shark::blas::diagonalMahalanobisDistanceSqr (vector_expression< VectorT > const &op1, vector_expression< VectorU > const &op2, vector_expression< WeightT > const &weights) |
| Normalized Euclidian squared distance (squared diagonal Mahalanobis) between two vectors. More... | |
| template<class VectorT , class VectorU > | |
| VectorT::value_type | shark::blas::distanceSqr (vector_expression< VectorT > const &op1, vector_expression< VectorU > const &op2) |
| Squared distance between two vectors. More... | |
| template<class MatrixT , class VectorU , class VectorR > | |
| void | shark::blas::distanceSqr (matrix_expression< MatrixT > const &operands, vector_expression< VectorU > const &op2, vector_expression< VectorR > &distances) |
| Squared distance between a vector and a set of vectors and stores the result in the vector of distances. More... | |
| template<class MatrixT , class VectorU > | |
| vector< typename MatrixT::value_type > | shark::blas::distanceSqr (matrix_expression< MatrixT > const &operands, vector_expression< VectorU > const &op2) |
| Squared distance between a vector and a set of vectors. More... | |
| template<class MatrixT , class VectorU > | |
| vector< typename MatrixT::value_type > | shark::blas::distanceSqr (vector_expression< VectorU > const &op1, matrix_expression< MatrixT > const &operands) |
| Squared distance between a vector and a set of vectors. More... | |
| template<class MatrixT , class MatrixU > | |
| matrix< typename MatrixT::value_type > | shark::blas::distanceSqr (matrix_expression< MatrixT > const &X, matrix_expression< MatrixU > const &Y) |
| Squared distance between the vectors of two sets of vectors. More... | |
| template<class VectorT , class VectorU > | |
| VectorT::value_type | shark::blas::distance (vector_expression< VectorT > const &op1, vector_expression< VectorU > const &op2) |
| Calculates distance between two vectors. More... | |
| template<class VectorT , class VectorU , class WeightT > | |
| VectorT::value_type | shark::blas::diagonalMahalanobisDistance (vector_expression< VectorT > const &op1, vector_expression< VectorU > const &op2, vector_expression< WeightT > const &weights) |
| Normalized euclidian distance (diagonal Mahalanobis) between two vectors. More... | |
| template<class Matrix > | |
| blas::matrix_vector_range < Matrix const > | shark::diag (blas::matrix_expression< Matrix > const &mat) |
| returns the diagonal of a constant square matrix as vector More... | |
| template<class Matrix > | |
| blas::matrix_vector_range< Matrix > | shark::diag (blas::matrix_expression< Matrix > &mat) |
| returns the diagonal of a square matrix as vector More... | |
| template<class Matrix > | |
| void | shark::zero (blas::matrix_expression< Matrix > &mat) |
| Zeros a matrix. If it is sparse, the structure is preserved. More... | |
| template<class Vector > | |
| void | shark::zero (blas::vector_expression< Vector > &vec) |
| Zeros a matrix. If it is sparse, the structure is preserved. More... | |
| template<class Matrix > | |
| void | shark::zero (blas::matrix_range< Matrix > mat) |
| Zeros a subrange of a matrix. If it is sparse, the structure is preserved. More... | |
| template<class Vector > | |
| void | shark::zero (blas::vector_range< Vector > vec) |
| Zeros a subrange of a vector. If it is sparse, the structure is preserved. More... | |
| template<class Vector > | |
| void | shark::zero (blas::matrix_row< Vector > vec) |
| Zeros a row of a matrix. If it is sparse, the structure is preserved. More... | |
| template<class Vector > | |
| void | shark::zero (blas::matrix_column< Vector > vec) |
| Zeros a column of a matrix. If it is sparse, the structure is preserved. More... | |
| template<class V > | |
| std::size_t | shark::nonzeroElements (blas::vector_expression< V > const &vec) |
| template<class Matrix > | |
| void | shark::identity (blas::matrix_expression< Matrix > &mat) |
| Initializes the square matrix A to be the identity matrix. More... | |
| template<class Matrix > | |
| void | shark::ensureSize (blas::matrix_expression< Matrix > &mat, std::size_t rows, std::size_t columns) |
| Ensures that the matrix has the right size. More... | |
| template<class Vector > | |
| void | shark::ensureSize (blas::vector_expression< Vector > &vec, std::size_t size) |
| Ensures that the vector has the right size. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_row< Matrix const > > | shark::triangularRow (blas::matrix_expression< Matrix > const &mat, std::size_t i) |
| Returns the i-th row of an upper triangular matrix excluding the elements right of the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_row< Matrix > > | shark::triangularRow (blas::matrix_expression< Matrix > &mat, std::size_t i) |
| Returns the i-th row of an upper triangular matrix excluding the elements right of the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_row< Matrix const > > | shark::unitTriangularRow (blas::matrix_expression< Matrix > const &mat, std::size_t i) |
| Returns the elements in the i-th row of a lower triangular matrix left of the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_row< Matrix > > | shark::unitTriangularRow (blas::matrix_expression< Matrix > &mat, std::size_t i) |
| Returns the elements in the i-th row of a lower triangular matrix left of the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_column< Matrix const > > | shark::unitTriangularColumn (blas::matrix_expression< Matrix > const &mat, std::size_t i) |
| Returns the elements in the i-th column of the matrix below the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_column< Matrix > > | shark::unitTriangularColumn (blas::matrix_expression< Matrix > &mat, std::size_t i) |
| Returns the elements in the i-th column of the matrix below the diagonal. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_column< Matrix const > > | shark::triangularColumn (blas::matrix_expression< Matrix > const &mat, std::size_t i) |
| Returns the elements in the i-th column of the matrix excluding the zero elements. More... | |
| template<class Matrix > | |
| blas::vector_range < blas::matrix_column< Matrix > > | shark::triangularColumn (blas::matrix_expression< Matrix > &mat, std::size_t i) |
| Returns the elements in the i-th column of the matrix excluding the zero elements. More... | |
| template<class Vector > | |
| VectorRepeater< Vector > | shark::repeat (blas::vector_expression< Vector > const &vector, std::size_t rows) |
| Creates a matrix from a vector by repeating the vector in every row of the matrix. More... | |
| template<class T > | |
| boost::enable_if < boost::is_arithmetic< T > , blas::scalar_vector< T > >::type | shark::repeat (T scalar, std::size_t elements) |
| template<class T > | |
| boost::enable_if < boost::is_arithmetic< T > , blas::scalar_matrix< T > >::type | shark::repeat (T scalar, std::size_t rows, std::size_t columns) |
| template<class Matrix > | |
| blas::matrix_range< Matrix const > | shark::rows (blas::matrix_expression< Matrix > const &mat, std::size_t start, std::size_t end) |
| brief picks a subrange of rows from a matrix. much easier to use than subrange More... | |
| template<class Matrix > | |
| blas::matrix_range< Matrix > | shark::rows (blas::matrix_expression< Matrix > &mat, std::size_t start, std::size_t end) |
| brief picks a subrange of rows from a matrix. much easier to use than subrange More... | |
| template<class Matrix > | |
| blas::matrix_range< Matrix const > | shark::columns (blas::matrix_expression< Matrix > const &mat, std::size_t start, std::size_t end) |
| brief picks a subrange of columns from a matrix. much easier to use than subrange More... | |
| template<class Matrix > | |
| blas::matrix_range< Matrix > | shark::columns (blas::matrix_expression< Matrix > &mat, std::size_t start, std::size_t end) |
| brief picks a subrange of columns from a matrix. much easier to use than subrange More... | |
| template<class MatrixT > | |
| MatrixT::value_type | shark::trace (blas::matrix_expression< MatrixT > const &m) |
| Evaluates the sum of the values at the diagonal of matrix "v". More... | |
| template<class MatrixT , class MatrixL > | |
| void | shark::choleskyDecomposition (blas::matrix_expression< MatrixT > const &A, blas::matrix_expression< MatrixL > &L) |
| Lower triangular Cholesky decomposition. More... | |
| template<class MatrixT , class VectorT > | |
| void | shark::eigensort (MatrixT &vmatA, VectorT &dvecA) |
| Sorts the eigenvalues in vector "dvecA" and the corresponding eigenvectors in matrix "vmatA". More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| void | shark::eigensymmJacobi (MatrixT &amatA, MatrixU &vmatA, VectorT &dvecA) |
| template<class MatrixT , class MatrixU , class VectorT > | |
| void | shark::eigensymmJacobi2 (MatrixT &amatA, MatrixU &vmatA, VectorT &dvecA) |
| template<class MatrixT , class MatrixU , class MatrixV , class VectorT > | |
| void | shark::eigensymm_intermediate (const MatrixT &amatA, MatrixU &hmatA, MatrixV &vmatA, VectorT &dvecA) |
| Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction. More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| void | shark::eigensymm (const MatrixT &A, MatrixU &G, VectorT &l) |
| Used as frontend for eigensymm for calculating the eigenvalues and the normalized eigenvectors of a symmetric matrix 'amatA' using the Givens and Householder reduction. Each time this frontend is called additional memory is allocated for intermediate results. More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| void | shark::eigensymm (const MatrixT &amatA, MatrixU &vmatA, VectorT &dvecA, VectorT &odvecA) |
| Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction without corrupting "amatA" during application. More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| double | shark::eigenerr (const MatrixT &amatA, const MatrixU &vmatA, const VectorT &dvecA, unsigned c) |
| Calculates the relative error of eigenvalue no. "c". More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| unsigned | shark::rank (const MatrixT &amatA, const MatrixU &vmatA, const VectorT &dvecA) |
| Determines the rank of the symmetric matrix "amatA". More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| double | shark::detsymm (MatrixT &amatA, MatrixU &vmatA, VectorT &dvecA) |
| Calculates the determinant of the symmetric matrix "amatA". More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| double | shark::logdetsymm (MatrixT &amatA, MatrixU &vmatA, VectorT &dvecA) |
| Calculates the log of the determinant of the symmetric matrix "amatA". More... | |
| template<class MatrixT , class MatrixU , class MatrixV , class VectorT > | |
| unsigned | shark::rankDecomp (MatrixT &amatA, MatrixU &vmatA, MatrixV &hmatA, VectorT &dvecA) |
| template<class MatrixT > | |
| void | shark::fft (MatrixT &data, int isign) |
| template<class MatrixT > | |
| void | shark::fft (MatrixT &data) |
| Replaces the "data" by its discrete Fourier transform. More... | |
| template<class MatrixT > | |
| void | shark::ifft (MatrixT &data) |
| Replaces the "data" by its inverse discrete Fourier transform. More... | |
| template<class MatrixT > | |
| RealMatrix | shark::invert (const MatrixT &mat) |
| Inverts a matrix with full rank. More... | |
| template<class MatrixT , class MatrixU > | |
| void | shark::invertSymmPositiveDefinite (MatrixT &I, const MatrixU &ArrSymm) |
| Inverts a symmetric positive definite matrix. More... | |
| template<class MatA , class MatU > | |
| void | shark::decomposedGeneralInverse (blas::matrix_expression< MatA > const &matA, blas::matrix_expression< MatU > &matU) |
| For a given square matrix A computes a matrix U where A'=LU^T. More... | |
| template<class MatrixT > | |
| RealMatrix | shark::g_inverse (blas::matrix_expression< MatrixT > const &matrixA) |
| Calculates the generalized inverse matrix of input matrix "matrixA". More... | |
| template<class MatrixT > | |
| void | shark::orthoNormalize (blas::matrix_container< MatrixT > &matrixC) |
| template<class MatrixT , typename RngType > | |
| void | shark::randomRotationMatrix (blas::matrix_container< MatrixT > &matrixC, RngType &rng) |
| Initializes a matrix such that it forms a random rotation matrix. More... | |
| template<class MatrixT > | |
| void | shark::randomRotationMatrix (blas::matrix_container< MatrixT > &matrixC) |
| Initializes a matrix such that it forms a random rotation. More... | |
| template<typename RngType > | |
| RealMatrix | shark::randomRotationMatrix (size_t size, RngType &rng) |
| Creates a random rotation matrix with a certain size using the random number generator rng. More... | |
| RealMatrix | shark::randomRotationMatrix (size_t size) |
| Creates a random rotation matrix with a certain size using the global random number gneerator. More... | |
| template<class X , class R > | |
| X::value_type | shark::createHouseholderReflection (blas::vector_expression< X > const &x, blas::vector_expression< R > &reflection) |
| Generates a Householder reflection from a vector to use with applyHouseholderLeft/Right. More... | |
| template<class Mat , class R , class T > | |
| void | shark::applyHouseholderOnTheRight (blas::matrix_expression< Mat > &matrix, blas::vector_expression< R > const &reflection, T beta) |
| template<class Mat , class R , class T > | |
| void | shark::applyHouseholderOnTheLeft (blas::matrix_expression< Mat > &matrix, blas::vector_expression< R > const &reflection, T const &beta) |
| template<class MatrixT , class Mat > | |
| std::size_t | shark::pivotingRQ (blas::matrix_expression< MatrixT > const &matrixA, blas::matrix_container< Mat > &matrixR, blas::matrix_container< Mat > &matrixQ, blas::permutation_matrix< std::size_t > &permutation) |
| Calculates the RQ-Decomposition of A using pivoting. More... | |
| template<class MatrixT , class MatrixU > | |
| std::size_t | shark::pivotingRQHouseholder (blas::matrix_expression< MatrixT > const &matrixA, blas::matrix_container< MatrixU > &matrixR, blas::matrix_container< MatrixU > &householderTransform, blas::permutation_matrix< std::size_t > &permutation) |
| Determines the RQ Decomposition of the matrix A using pivoting returning the housholder transformation instead of Q. More... | |
| template<class MatrixT , class MatrixU , class VectorT > | |
| unsigned | shark::svdrank (const MatrixT &amatA, MatrixU &umatA, MatrixU &vmatA, VectorT &wvecA) |
| template<class MatrixT , class MatrixU , class VectorT > | |
| void | shark::svd (const MatrixT &amatA, MatrixU &umatA, MatrixU &vmatA, VectorT &wvecA, unsigned maxIterations=200, bool ignoreThreshold=true) |
| Determines the singular value decomposition of a rectangular matrix "amatA". More... | |
| template<class MatrixU , class VectorT > | |
| void | shark::svdsort (MatrixU &umatA, MatrixU &vmatA, VectorT &wvecA) |
| Sorts the singular values in vector "wvecA" by descending order. More... | |
| template<class Vec1T , class Vec2T , class Vec3T > | |
| void | shark::meanvar (Data< Vec1T > const &data, blas::vector_container< Vec2T > &meanVec, blas::vector_container< Vec3T > &varianceVec) |
| Calculates the mean and variance values of a dataset. More... | |
| template<class MatT , class Vec1T , class Vec2T > | |
| void | shark::meanvar (blas::matrix_container< MatT > &data, blas::vector_container< Vec1T > &meanVec, blas::vector_container< Vec2T > &varianceVec) |
| Calculates the mean, variance and covariance values of the input data. More... | |
| template<class Vec1T , class Vec2T , class MatT > | |
| void | shark::meanvar (const Data< Vec1T > &data, blas::vector_container< Vec2T > &meanVec, blas::matrix_container< MatT > &covariance) |
| Calculates the mean and covariance values of a set of data. More... | |
| template<class VectorType > | |
| VectorType | shark::mean (Data< VectorType > const &data) |
| Calculates the mean vector of array "x". More... | |
| template<class MatrixType > | |
| blas::vector< typename MatrixType::value_type > | shark::mean (const blas::matrix_container< MatrixType > &data) |
| Calculates the mean vector of the input vectors. More... | |
| template<class VectorType > | |
| VectorType | shark::variance (const Data< VectorType > &data) |
| Calculates the variance vector of array "x". More... | |
| template<class VectorType > | |
| VectorMatrixTraits< VectorType > ::DenseMatrixType | shark::covariance (const Data< VectorType > &data) |
| Calculates the covariance matrix of the data vectors stored in data. More... | |
| template<class VectorType > | |
| VectorMatrixTraits< VectorType > ::DenseMatrixType | shark::corrcoef (const Data< VectorType > &data) |
| Calculates the coefficient of correlation matrix of the data vectors stored in data. More... | |
| template<class VectorType > | |
| VectorType | shark::mean (UnlabeledData< VectorType > const &data) |
| template<class T > | |
| double | shark::stl_median (std::vector< T > &v) |
| compute median More... | |
| template<class T > | |
| double | shark::stl_percentile (std::vector< T > &v, double p=.25) |
| compute percentilee (Excel way) More... | |
| template<class T > | |
| double | shark::nth (std::vector< T > &v, unsigned n) |
| return nth element after sorting More... | |
| template<class T > | |
| double | shark::stl_mean (std::vector< T > v) |
| compute mean More... | |
| template<class T > | |
| double | shark::stl_correlation (std::vector< T > &v1, std::vector< T > &v2) |
| compute variance More... | |
| template<class T > | |
| double | shark::stl_variance (std::vector< T > &v, bool unbiased=true) |
| compute sample variance More... | |
| template<class InputType , class OutputType > | |
| void | shark::initRandomNormal (AbstractModel< InputType, OutputType > &model, double s) |
| Initialize model parameters normally distributed. More... | |
| template<class InputType , class OutputType > | |
| void | shark::initRandomUniform (AbstractModel< InputType, OutputType > &model, double l, double h) |
| Initialize model parameters uniformly at random. More... | |
Several mathematical, linear-algebra, or other functions within Shark are not part of any particular class. They are collected here in the doxygen group "shark_globals".
| #define SHARK_LINALG_EIGENSORT_INL |
Inline implementation of eigenvalue/-vector sorting.
This file is part of Shark. This library is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3, or (at your option) any later version.
This library is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this library; if not, see http://www.gnu.org/licenses/.
Definition at line 36 of file eigensort.inl.
| Enumerator | |
|---|---|
| NONE | |
| MULTIPLICATIVE_TRACE_ONE | |
| MULTIPLICATIVE_TRACE_N | |
| MULTIPLICATIVE_VARIANCE_ONE | |
| CENTER_ONLY | |
| CENTER_AND_MULTIPLICATIVE_TRACE_ONE | |
Definition at line 46 of file PrecomputedMatrix.h.
| enum shark::LabelPosition |
| void shark::applyHouseholderOnTheLeft | ( | blas::matrix_expression< Mat > & | matrix, |
| blas::vector_expression< R > const & | reflection, | ||
| T const & | beta | ||
| ) |
Definition at line 176 of file rotations.h.
References shark::fast_prod(), shark::blas::noalias(), shark::blas::outer_prod(), shark::size(), SIZE_CHECK, and shark::blas::trans().
| void shark::applyHouseholderOnTheRight | ( | blas::matrix_expression< Mat > & | matrix, |
| blas::vector_expression< R > const & | reflection, | ||
| T | beta | ||
| ) |
Definition at line 149 of file rotations.h.
References shark::fast_prod(), shark::blas::noalias(), shark::blas::outer_prod(), shark::size(), and SIZE_CHECK.
Referenced by shark::pivotingRQHouseholder().
| void shark::choleskyDecomposition | ( | blas::matrix_expression< MatrixT > const & | A, |
| blas::matrix_expression< MatrixL > & | L | ||
| ) |
Lower triangular Cholesky decomposition.
Given an \( m \times m \) symmetric positive definite matrix \(A\), compute the lower triangular matrix \(L\) such that \(A = LL^T \). An exception is thrown if the matrix is not positive definite. If you suspect the matrix to be positive semi-definite, use pivotingCholeskyDecomposition instead
| A | \( m \times m \) matrix, which must be symmetric and positive definite |
| L | \( m \times m \) matrix, which stores the Cholesky factor |
Definition at line 46 of file Cholesky.inl.
References shark::detail::ensureSize(), shark::detail::bindings::potrf(), SHARKEXCEPTION, SIZE_CHECK, and shark::detail::zero().
Referenced by shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::eval(), shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::evalDerivative(), shark::invertSymmPositiveDefinite(), and shark::solveSymmSystemInPlace().
| blas::matrix_range<Matrix const> shark::columns | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | start, | ||
| std::size_t | end | ||
| ) |
brief picks a subrange of columns from a matrix. much easier to use than subrange
Definition at line 418 of file Tools.h.
References shark::blas::subrange().
Referenced by shark::PCA::decoder(), shark::decomposedGeneralInverse(), shark::PCA::encoder(), shark::detail::SubrangeKernelWrapper< InputType >::eval(), shark::repeat(), shark::Centroids::softMembership(), shark::FisherLDA::train(), shark::detail::SubrangeKernelWrapper< InputType >::weightedInputDerivative(), shark::detail::SubrangeKernelWrapper< InputType >::weightedParameterDerivative(), shark::OnlineRNNet::weightedParameterDerivative(), and shark::RNNet::weightedParameterDerivative().
| blas::matrix_range<Matrix> shark::columns | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | start, | ||
| std::size_t | end | ||
| ) |
brief picks a subrange of columns from a matrix. much easier to use than subrange
Definition at line 424 of file Tools.h.
References shark::blas::subrange().
| T shark::copySign | ( | T | x, |
| T | y | ||
| ) |
brief lets x have the same sign as y.
This is the famous well known copysign function from fortran.
Definition at line 159 of file Math.h.
Referenced by shark::svd().
| VectorMatrixTraits< VectorType >::DenseMatrixType shark::corrcoef | ( | const Data< VectorType > & | data | ) |
Calculates the coefficient of correlation matrix of the data vectors stored in data.
Calculates the coefficient of correlation matrix of the data vectors.
Given a matrix \(X = (x_{ij})\) of \(n\) vectors with length \(N\), the function calculates the coefficient of correlation matrix given as
\( r := (r_{kl}) \mbox{,\ } r_{kl} = \frac{c_{kl}}{\Delta x_k \Delta x_l}\mbox{,\ } k,l = 1, \dots, N \)
where \(c_{kl}\) is the entry of the covariance matrix of \(x\) and \(y\) and \(\Delta x_k\) and \(\Delta x_l\) are the standard deviations of \(x_k\) and \(x_l\) respectively.
| data | The \(n \times N\) input matrix. |
Definition at line 233 of file VectorStatistics.inl.
References shark::covariance().
| VectorMatrixTraits< VectorType >::DenseMatrixType shark::covariance | ( | const Data< VectorType > & | data | ) |
Calculates the covariance matrix of the data vectors stored in data.
Calculates the covariance matrix of the data vectors.
Given a Set \(X = (x_{ij})\) of \(n\) vectors with length \(N\), the function calculates the covariance matrix given as
\( Cov = (c_{kl}) \mbox{,\ } c_{kl} = \frac{1}{n - 1} \sum_{i=1}^n (x_{ik} - \overline{x_k})(x_{il} - \overline{x_l})\mbox{,\ } k,l = 1, \dots, N \)
where \(\overline{x_j} = \frac{1}{n} \sum_{i = 1}^n x_{ij}\) is the mean value of \(x_j \mbox{,\ }j = 1, \dots, N\).
| data | The \(n \times N\) input matrix. |
Definition at line 205 of file VectorStatistics.inl.
References shark::mean(), and shark::meanvar().
Referenced by shark::corrcoef(), createData(), shark::meanvar(), shark::NormalizeComponentsWhitening< VectorType >::train(), and shark::LDA::train().
| CVFolds<LabeledData<I,L> > shark::createCVIID | ( | LabeledData< I, L > & | set, |
| size_t | numberOfPartitions, | ||
| std::size_t | batchSize = Data<I>::DefaultBatchSize |
||
| ) |
Create a partition for cross validation.
The subset each training examples belongs to is drawn independently and uniformly distributed. For every partition, all but one subset form the training set, while the remaining one is used for validation. The partitions can be accessed using getCVPartitionName
| set | the input data where the new partitions are created |
| numberOfPartitions | number of partitions to create |
Definition at line 202 of file CVDatasetTools.h.
References shark::createCVIndexed().
Referenced by run_one_trial().
| CVFolds<LabeledData<I,L> > shark::createCVIndexed | ( | LabeledData< I, L > & | set, |
| size_t | numberOfPartitions, | ||
| std::vector< size_t > | indices, | ||
| std::size_t | batchSize = Data<I>::DefaultBatchSize |
||
| ) |
Create a partition for cross validation from indices.
Create a partition from indices. The indices vector for each sample states of what validation partition that sample should become a member. In other words, the index maps a sample to a validation partition, meaning that it will become a part of the training partition for all other folds.
| set | partitions will be subsets of this set |
| numberOfPartitions | number of partitions to create |
| indices | partition indices of the examples in [0, ..., numberOfPartitions[. |
Definition at line 320 of file CVDatasetTools.h.
References shark::LabeledData< InputT, LabelT >::batch(), shark::detail::batchPartitioning(), shark::size(), SIZE_CHECK, shark::subBatch(), and shark::swap().
Referenced by shark::createCVIID().
| CVFolds<LabeledData<I,L> > shark::createCVSameSize | ( | LabeledData< I, L > & | set, |
| std::size_t | numberOfPartitions, | ||
| std::size_t | batchSize = LabeledData<I,L>::DefaultBatchSize |
||
| ) |
Create a partition for cross validation.
Every subset contains (approximately) the same number of elements. For every partition, all but one subset form the training set, while the remaining one is used for validation. The partitions can be accessed using getCVPartitionName
| numberOfPartitions | number of partitions to create |
| set | the input data from which to draw the partitions |
Definition at line 225 of file CVDatasetTools.h.
References shark::detail::batchPartitioning().
Referenced by main(), and shark::CARTTrainer::train().
| CVFolds<LabeledData<I,unsigned int> > shark::createCVSameSizeBalanced | ( | LabeledData< I, unsigned int > & | set, |
| size_t | numberOfPartitions, | ||
| std::size_t | batchSize = Data<I>::DefaultBatchSize |
||
| ) |
Create a partition for cross validation.
Every subset contains (approximately) the same number of elements. For every partition, all but one subset form the training set, while the remaining one is used for validation.
| numberOfPartitions | number of partitions to create |
| set | the input data from which to draw the partitions |
Definition at line 260 of file CVDatasetTools.h.
References shark::detail::createCVSameSizeBalanced(), shark::numberOfClasses(), and shark::DataView< DatasetType >::size().
Referenced by main(), and shark::CARTTrainer::train().
| CVFolds<LabeledData<I,RealVector> > shark::createCVSameSizeBalanced | ( | LabeledData< I, RealVector > & | set, |
| size_t | numberOfPartitions, | ||
| std::size_t | batchSize = Data<I>::DefaultBatchSize |
||
| ) |
Create a partition for cross validation.
Every subset contains (approximately) the same number of elements. For every partition, all but one subset form the training set, while the remaining one is used for validation.
| set | the input data from which to draw the partitions |
| numberOfPartitions | number of partitions to create |
Definition at line 284 of file CVDatasetTools.h.
References shark::DataView< DatasetType >::begin(), shark::detail::createCVSameSizeBalanced(), shark::blas::distance(), shark::DataView< DatasetType >::end(), shark::numberOfClasses(), and shark::DataView< DatasetType >::size().
| Data<typename boost::range_value<Range>::type> shark::createDataFromRange | ( | Range const & | inputs, |
| std::size_t | maximumBatchSize = 0 |
||
| ) |
creates a data object from a range of elements
Definition at line 713 of file Dataset.h.
References shark::Data< Type >::batch(), and shark::size().
Referenced by shark::BarsAndStripes::BarsAndStripes(), createData(), shark::createLabeledDataFromRange(), shark::SvmLogisticInterpretation< InputType >::eval(), shark::SvmLogisticInterpretation< InputType >::evalDerivative(), shark::import_csv(), shark::importHDF5(), shark::importPGMSet(), shark::Centroids::initFromData(), shark::kMeans(), main(), shark::Shifter::Shifter(), and shark::string2data().
| X::value_type shark::createHouseholderReflection | ( | blas::vector_expression< X > const & | x, |
| blas::vector_expression< R > & | reflection | ||
| ) |
Generates a Householder reflection from a vector to use with applyHouseholderLeft/Right.
Given a Vector x=(x0,x1,...,xn), finds a reflection with the property (c, 0,0,...0) = (I-beta v v^t)x and v = (x0-c,x1,x2,...,xn)
Definition at line 116 of file rotations.h.
References shark::blas::noalias(), shark::blas::norm_2(), shark::size(), and SIZE_CHECK.
Referenced by shark::pivotingRQHouseholder().
|
inline |
| void shark::decomposedGeneralInverse | ( | blas::matrix_expression< MatA > const & | matA, |
| blas::matrix_expression< MatU > & | matU | ||
| ) |
For a given square matrix A computes a matrix U where A'=LU^T.
A' is the generalized inverse of A. If A has full rank, the resulting matrix is equivalent to the inverse cholesky decomposition. If it is not, it is not necessarily triangular. This is mostly a helperfunction for g_inverse but also used in other parts of shark
Definition at line 45 of file g_inverse.inl.
References shark::columns(), shark::detail::ensureSize(), shark::identity(), shark::rank(), SIZE_CHECK, shark::swapFullInverted(), shark::symmRankKUpdate(), shark::blas::trans(), and shark::detail::zero().
Referenced by shark::LinearClassifier::importCovarianceMatrix(), and shark::NormalizeComponentsWhitening< VectorType >::train().
| double shark::detsymm | ( | MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the determinant of the symmetric matrix "amatA".
Calculates the determinate of the symmetric matrix "amatA".
Calculates the determinate of matrix amatA by using its n eigenvalues \( x_j \) that first will be calculated. The determinate is then given as:
\( \prod_{j=1}^n x_j \)
| amatA | \( n \times n \) matrix, which is symmetric, so only the bottom triangular matrix must contain values. At the end of the function amatA always contains the full matrix. |
| vmatA | \( n \times n \) matrix, that will contain the scaled eigenvectors at the end of the function. |
| dvecA | n-dimensional vector that will contain the eigenvalues at the end of the function. |
| SharkException | the type of the eception will be "size mismatch" and indicates that amatA is not a square matrix |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 87 of file detsymm.inl.
References shark::eigensymm_intermediate(), and SIZE_CHECK.
| blas::matrix_vector_range<Matrix const> shark::diag | ( | blas::matrix_expression< Matrix > const & | mat | ) |
returns the diagonal of a constant square matrix as vector
given a matrix (1 2 3) A =(4 5 6) (7 8 9)
diag(A) = (1,5,9)
Definition at line 126 of file Tools.h.
References SIZE_CHECK.
Referenced by shark::LinearRegression::train(), and shark::LassoRegression< InputVectorType >::trainInternal().
| blas::matrix_vector_range<Matrix> shark::diag | ( | blas::matrix_expression< Matrix > & | mat | ) |
returns the diagonal of a square matrix as vector
given a matrix (1 2 3) A =(4 5 6) (7 8 9)
diag(A) = (1,5,9)
Definition at line 141 of file Tools.h.
References SIZE_CHECK.
| VectorT::value_type shark::blas::diagonalMahalanobisDistance | ( | vector_expression< VectorT > const & | op1, |
| vector_expression< VectorU > const & | op2, | ||
| vector_expression< WeightT > const & | weights | ||
| ) |
Normalized euclidian distance (diagonal Mahalanobis) between two vectors.
Contrary to some conventions, dimension-wise weights are considered instead of std. deviations: \( d(v) = \left( \sum_i w_i (x_i-z_i)^2 \right)^{1/2} \) nb: the weights themselves are not squared, but multiplied onto the squared components
Definition at line 480 of file Metrics.h.
References shark::blas::diagonalMahalanobisDistanceSqr(), shark::size(), and SIZE_CHECK.
| VectorT::value_type shark::blas::diagonalMahalanobisDistanceSqr | ( | vector_expression< VectorT > const & | op1, |
| vector_expression< VectorU > const & | op2, | ||
| vector_expression< WeightT > const & | weights | ||
| ) |
Normalized Euclidian squared distance (squared diagonal Mahalanobis) between two vectors.
NOTE: The weights themselves are not squared, but multiplied onto the squared components.
Definition at line 347 of file Metrics.h.
References shark::blas::detail::diagonalMahalanobisDistanceSqr(), shark::size(), and SIZE_CHECK.
Referenced by shark::blas::diagonalMahalanobisDistance(), and shark::blas::distanceSqr().
| VectorT::value_type shark::blas::distance | ( | vector_expression< VectorT > const & | op1, |
| vector_expression< VectorU > const & | op2 | ||
| ) |
Calculates distance between two vectors.
Definition at line 464 of file Metrics.h.
References shark::blas::distanceSqr(), shark::size(), and SIZE_CHECK.
Referenced by shark::LCTree< VectorType, CuttingAccuracy >::buildTree(), shark::KHCTree< Container, CuttingAccuracy >::buildTree(), shark::KDTree< InputT >::buildTree(), shark::detail::complement(), shark::FFNet< HiddenNeuron, OutputNeuron >::configure(), shark::createCVSameSizeBalanced(), shark::AbsoluteLoss< VectorType >::eval(), shark::LassoRegression< InputVectorType >::fillData(), shark::SteadyStateIndicatorBasedSelection< Indicator >::leastContributor(), shark::IndicatorBasedSelection< shark::HypervolumeIndicator >::leastContributor(), shark::ApproximatedHypervolumeSelection::leastContributor(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::numberOfParameters(), shark::WilcoxonRankSumTest::operator()(), shark::TournamentSelection::operator()(), shark::RouletteWheelSelection::operator()(), shark::BinaryTournamentSelection< shark::RankShareComparator >::operator()(), shark::EPTournamentSelection< FitnessType >::operator()(), shark::HomogenousNDimFS::operator()(), shark::UniformRankingSelection< FitnessType >::operator()(), shark::LinearRankingSelection< FitnessType >::operator()(), shark::HypervolumeApproximator< Rng >::operator()(), shark::LeastContributorApproximator< Rng, ExactHypervolume >::operator()(), shark::select_mu_komma_lambda(), shark::select_mu_komma_lambda_p(), shark::detail::VectorExpression< Expression >::size(), shark::detail::MatrixExpression< Matrix >::size(), shark::detail::AGE2::step(), shark::detail::AGE::step(), and shark::FrontStore< ResultType, MetaDataType >::writeFront().
| VectorT::value_type shark::blas::distanceSqr | ( | vector_expression< VectorT > const & | op1, |
| vector_expression< VectorU > const & | op2 | ||
| ) |
Squared distance between two vectors.
Definition at line 366 of file Metrics.h.
References shark::blas::diagonalMahalanobisDistanceSqr(), shark::size(), and SIZE_CHECK.
Referenced by shark::LCTree< VectorType, CuttingAccuracy >::calculateNormal(), shark::blas::distance(), shark::blas::distanceSqr(), shark::SquaredLoss< VectorType >::eval(), shark::GaussianRbfKernel< InputType >::eval(), shark::LinearKernel< InputType >::featureDistanceSqr(), shark::JaakkolaHeuristic::JaakkolaHeuristic(), and shark::Centroids::softMembership().
| void shark::blas::distanceSqr | ( | matrix_expression< MatrixT > const & | operands, |
| vector_expression< VectorU > const & | op2, | ||
| vector_expression< VectorR > & | distances | ||
| ) |
Squared distance between a vector and a set of vectors and stores the result in the vector of distances.
The squared distance between the vector and every row-vector of the matrix is calculated. This can be implemented much more efficiently.
Definition at line 383 of file Metrics.h.
References shark::blas::detail::distanceSqrBlockVector(), shark::ensureSize(), shark::size(), and SIZE_CHECK.
| vector<typename MatrixT::value_type> shark::blas::distanceSqr | ( | matrix_expression< MatrixT > const & | operands, |
| vector_expression< VectorU > const & | op2 | ||
| ) |
Squared distance between a vector and a set of vectors.
The squared distance between the vector and every row-vector of the matrix is calculated. This can be implemented much more efficiently.
Definition at line 404 of file Metrics.h.
References shark::blas::distanceSqr(), shark::size(), and SIZE_CHECK.
| vector<typename MatrixT::value_type> shark::blas::distanceSqr | ( | vector_expression< VectorU > const & | op1, |
| matrix_expression< MatrixT > const & | operands | ||
| ) |
Squared distance between a vector and a set of vectors.
The squared distance between the vector and every row-vector of the matrix is calculated. This can be implemented much more efficiently.
Definition at line 421 of file Metrics.h.
References shark::blas::distanceSqr(), shark::size(), and SIZE_CHECK.
| matrix<typename MatrixT::value_type> shark::blas::distanceSqr | ( | matrix_expression< MatrixT > const & | X, |
| matrix_expression< MatrixU > const & | Y | ||
| ) |
Squared distance between the vectors of two sets of vectors.
The squared distance between every row-vector of the first matrix x and every row-vector of the second matrix y is calculated. This can be implemented much more efficiently. The results are returned as a matrix, where the element in the i-th row and the j-th column is distanceSqr(x_i,y_j).
Definition at line 441 of file Metrics.h.
References shark::blas::detail::distanceSqrBlockBlock(), and SIZE_CHECK.
| void shark::dlinmin | ( | VectorT & | p, |
| const VectorU & | xi, | ||
| double & | fret, | ||
| DifferentiableFunction & | func, | ||
| double | ax = 0.0, |
||
| double | bx = 1.0 |
||
| ) |
Minimizes a function of "N" variables by using derivative information.
Performs a minimization of a function of \( N \) variables, i.e. given as input the vectors \( P \) and \( n \) and the function \( f \), the function finds the scalar \( \lambda \) that minimizes \( f(P + \lambda n) \). \( P \) is replaced by \( P + \lambda n \) and \( n \) by \( \lambda n \).
| p | N-dimensional initial starting point for the search, is set to the point where the function takes on a minimum. |
| xi | N-dimensional search direction, is replaced by the actual vector displacement that p was moved. |
| fret | The function value at the new point p. |
| func | The function that will be minimized. |
| ax | guess for the lower bracket |
| bx | guess for the upper bracket |
| SharkException | the type of the exception will be "size mismatch" and indicates that p is not one-dimensional or that p and xi don't have the same length |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 93 of file dlinmin.inl.
References shark::blas::max(), SIGN, and SIZE_CHECK.
Referenced by shark::LineSearch::operator()().
| double shark::eigenerr | ( | const MatrixT & | amatA, |
| const MatrixU & | vmatA, | ||
| const VectorT & | dvecA, | ||
| unsigned | c | ||
| ) |
Calculates the relative error of eigenvalue no. "c".
Calculates the relative error of one eigenvalue with no. "c".
Given a symmetric \( n \times n \) matrix amatA, the matrix vmatA of its eigenvectors and the vector dvecA of the corresponding eigenvalues, this function calculates the relative error of one eigenvalue, denoted by the no. c. If we have a \( n \times n \) matrix A, a matrix x of eigenvectors and a vector \( \lambda \) of corresponding eigenvalues, the relative error of eigenvalue no. c is calculated as
\( \sqrt{\sum_{i=0}^n \left(\sum_{j=0}^n A(i,j) \ast x(j,c) - x(i,c) \ast \lambda(c) \right)^2} \)
| amatA | \( n \times n \) matrix, which has to be symmetric, so only the lower triangular matrix must contain values. The matrix is not changed by the function. |
| vmatA | \( n \times n \) matrix with normalized eigenvectors, each column contains an eigenvector. |
| dvecA | n-dimensional vector with eigenvalues in descending order. |
| c | No. of the considered eigenvalue. |
| SharkException | the type of the exception will be "size mismatch" and indicates that dvecA is not one-dimensional or that amatA or vmatA don't have the same number of rows or columns as dvecA contains number of values |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 92 of file eigenerr.inl.
References SIZE_CHECK.
Referenced by shark::rank().
| void shark::eigensort | ( | MatrixT & | vmatA, |
| VectorT & | dvecA | ||
| ) |
Sorts the eigenvalues in vector "dvecA" and the corresponding eigenvectors in matrix "vmatA".
Sorts the eigenvalues in vector "dvecA" in descending order and the corresponding eigenvectors in matrix "vmatA".
Given the matrix vmatA of eigenvectors and the vector dvecA of corresponding eigenvalues, the values in dvecA will be sorted by descending order and the eigenvectors in vmatA will change their places in a way, that at the end of the function an eigenvalue at position j of vector dvecA will belong to the eigenvector at column j of matrix vmatA. If we've got for example the following result after calling the function:
\( \begin{array}{*{3}{r}} v_{11} & v_{21} & v_{31}\\ v_{12} & v_{22} & v_{32}\\ v_{13} & v_{23} & v_{33}\\ & & \\ v_1 & v_2 & v_3\\ \end{array} \)
then eigenvalue \( v_1 \) has the corresponding eigenvector \( ( v_{11}\ v_{12}\ v_{13} ) \) and \( v_1 > v_2 > v_3 \).
| vmatA | \( n \times n \) matrix with eigenvectors (each column contains an eigenvector, corresponding to one eigenvalue). |
| dvecA | n-dimensional vector with eigenvalues, will contain the eigenvalues in descending order when returning from the function. |
| SharkException | the type of the exception will be "size mismatch" and indicates that dvecA is not one-dimensional or that the number of rows or the number of columns in vmatA is different from the number of values in dvecA |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 92 of file eigensort.inl.
References SIZE_CHECK.
Referenced by shark::eigensymm(), shark::eigensymm_intermediate(), shark::eigensymmJacobi(), and shark::eigensymmJacobi2().
| void shark::eigensymm | ( | const MatrixT & | A, |
| MatrixU & | G, | ||
| VectorT & | l | ||
| ) |
Used as frontend for eigensymm for calculating the eigenvalues and the normalized eigenvectors of a symmetric matrix 'amatA' using the Givens and Householder reduction. Each time this frontend is called additional memory is allocated for intermediate results.
Used as frontend alculating the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction without corrupting "A" during application. Each time this frontend is called additional memory is allocated for intermediate results.
| A | \( n \times n \) matrix, which must be symmetric, so only the bottom triangular matrix must contain values. |
| G | \( n \times n \) matrix with the calculated normalized eigenvectors, each column will contain one eigenvector. |
| l | n-dimensional vector with the calculated eigenvalues in descending order. |
| SharkException |
Definition at line 365 of file eigensymm.inl.
References SIZE_CHECK.
Referenced by shark::Store::operator()(), shark::PCA::setData(), shark::FisherLDA::train(), and shark::detail::TypedMultiVariateNormalDistribution< Rng, RealMatrix, RealVector >::update().
| void shark::eigensymm | ( | const MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| VectorT & | dvecA, | ||
| VectorT & | odvecA | ||
| ) |
Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction without corrupting "amatA" during application.
Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction. This method works without corrupting "amatA" during application by demanding another Array 'odvecA' as an algorithmic buffer instead of using the last row of 'amatA' to store intermediate algorithmic results.
Given a symmetric \( n \times n \) matrix A, this function calculates the eigenvalues \( \lambda \) and the eigenvectors x, defined as
\( Ax = \lambda x \)
where x is a one-column matrix and the matrix multiplication is used for A and x. Here, the Givens reduction as a modification of the Jacobi method is used. Instead of trying to reduce the matrix all the way to diagonal form, we are content to stop when the matrix is tridiagonal. This allows the function to be carried out in a finite number of steps, unlike the Jacobi method, which requires iteration to convergence. So in comparison to the Jacobi method, this function is faster for matrices with an order greater than 10.
| amatA | \( n \times n \) matrix, which must be symmetric, so only the bottom triangular matrix must contain values. |
| vmatA | \( n \times n \) matrix with the calculated normalized eigenvectors, each column will contain an eigenvector. |
| dvecA | n-dimensional vector with the calculated eigenvalues in descending order. |
| odvecA | n-dimensional vector with the calculated offdiagonal of the Householder transformation. |
| SharkException | Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg. |
Definition at line 440 of file eigensymm.inl.
References shark::eigensort(), SHARKEXCEPTION, SIZE_CHECK, and shark::detail::zero().
| void shark::eigensymm_intermediate | ( | const MatrixT & | amatA, |
| MatrixU & | hmatA, | ||
| MatrixV & | vmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction.
Calculates the eigenvalues and the normalized eigenvectors of a symmetric matrix "amatA" using the Givens and Householder reduction, however, "hmatA" contains intermediate results after application.
Given a symmetric \( n \times n \) matrix A, this function calculates the eigenvalues \( \lambda \) and the eigenvectors x, defined as
\( Ax = \lambda x \)
where x is a one-column matrix and the matrix multiplication is used for A and x. Here, the Givens reduction as a modification of the Jacobi method is used. Instead of trying to reduce the matrix all the way to diagonal form, we are content to stop when the matrix is tridiagonal. This allows the function to be carried out in a finite number of steps, unlike the Jacobi method, which requires iteration to convergence. So in comparison to the Jacobi method, this function is faster for matrices with an order greater than 10.
| amatA | \( n \times n \) matrix, which must be symmetric, so only the bottom triangular matrix must contain values. Values below the diagonal will be destroyed. The method uses this matrix as a buffer for intermediate results. |
| hmatA | \( n \times n \) matrix, that is used to store intermediate results (other methods like detsymm use these for further calculations) |
| vmatA | \( n \times n \) matrix with the calculated normalized eigenvectors, each column will contain an eigenvector. |
| dvecA | n-dimensional vector with the calculated eigenvalues in descending order. |
| SharkException |
Definition at line 95 of file eigensymm.inl.
References shark::eigensort(), shark::blas::row(), SHARKEXCEPTION, and SIZE_CHECK.
Referenced by shark::detsymm(), shark::logdetsymm(), and shark::rankDecomp().
| void shark::eigensymmJacobi | ( | MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the eigenvalues and the normalized eigenvectors of the symmetric matrix "amatA" using the Jacobi method.
Given a symmetric \( n \times n \) matrix A, this function calculates the eigenvalues \( \lambda \) and the eigenvectors x, defined as
\( Ax = \lambda x \)
where x is a one-column matrix and the matrix multiplication is used for A and x. For the calculation of the eigenvectors and eigenvalues the so called Jacobi rotation is used to annihilate one of the off-diagonal elements with the basic Jacobi rotation given as matrix of the form
\( P_{pq} = \left( \begin{array}{*{7}{c}} 1 \\ & \dots \\ & & c & \dots & s \\ & & \vdots & 1 & \vdots \\ & & -s & \dots & c \\ & & & & & \dots \\ & & & & & & 1 \\ \end{array} \right) \)
In this matrix all the diagonal elements are unity except for the two elemnts c in rows (and columns) p and q. All off-diagonal elements are zero except the two elements s and - s. The numbers c and s are the cosine of a rotation angle \( \Phi \), so \( c^2 + s^2 = 1\). Successive rotations lead to the off-diagonal elements getting smaller and smaller, until the matrix is diagonal to machine precision. Accumulating the product of the transformations as you go gives the matrix of eigenvectors, while the elements of the final diagonal matrix are the eigenvalues. Use this function for the calculation of the eigenvalues and eigenvectors for matrices amatA with moderate order not greater than 10.
| amatA | \( n \times n \) matrix, which must be symmetric, so only the upper triangular matrix must contain values. Values above the diagonal will be destroyed. |
| vmatA | \( n \times n \) matrix with the calculated normalized eigenvectors, each column will contain an eigenvector. |
| dvecA | n-dimensional vector with the calculated eigenvalues in descending order. |
| SharkException | Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg. |
Definition at line 123 of file eigensymmJacobi.inl.
References shark::eigensort(), SHARKEXCEPTION, and SIZE_CHECK.
| void shark::eigensymmJacobi2 | ( | MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the eigenvalues and the normalized eigenvectors of the symmetric matrix "amatA" using a modified Jacobi method.
Given a symmetric \( n \times n \) matrix A, this function calculates the eigenvalues \( \lambda \) and the eigenvectors x, defined as
\( Ax = \lambda x \)
where x is a one-column matrix and the matrix multiplication is used for A and x. This function uses the Jacobi method as in eigensymmJacobi, but the method is modificated after J. von Neumann to avoid convergence problems when dealing with low machine precision.
| amatA | \( n \times n \) matrix, which must be symmetric, so only the bottom triangular matrix must contain values. Values below the diagonal will be destroyed. |
| vmatA | \( n \times n \) matrix with the calculated normalized eigenvectors, each column will contain one eigenvector. |
| dvecA | n-dimensional vector with the calculated eigenvalues in descending order. |
| SharkException | the type of the exception will be "size mismatch" and indicates that amatA is not a square matrix |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 97 of file eigensymmJacobi2.inl.
References shark::eigensort(), and SIZE_CHECK.
| void shark::ensureSize | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | rows, | ||
| std::size_t | columns | ||
| ) |
Ensures that the matrix has the right size.
Tries to resize mat. If the matrix expression can't be resized a debug assertion is thrown.
Definition at line 217 of file Tools.h.
References shark::detail::ensureSize().
Referenced by shark::approxSolveSymmSystem(), shark::blas::distanceSqr(), shark::RFClassifier::eval(), shark::MissingFeaturesKernelExpansion< InputType >::eval(), shark::DiscreteKernel::eval(), shark::ARDKernelUnconstrained< InputType >::eval(), shark::WeightedSumKernel< InputType >::eval(), shark::DenoisingAutoencoderError< InputType, RngType >::evalDerivative(), shark::Batch< blas::vector< T > >::resize(), shark::Batch< shark::blas::compressed_vector< T > >::resize(), shark::ARDKernelUnconstrained< InputType >::weightedInputDerivative(), shark::ARDKernelUnconstrained< InputType >::weightedParameterDerivative(), shark::NormalizedKernel< InputType >::weightedParameterDerivative(), and shark::WeightedSumKernel< InputType >::weightedParameterDerivative().
| void shark::ensureSize | ( | blas::vector_expression< Vector > & | vec, |
| std::size_t | size | ||
| ) |
Ensures that the vector has the right size.
Tries to resize vec. If the vector expression can't be resized a debug assertion is thrown.
Definition at line 224 of file Tools.h.
References shark::detail::ensureSize().
| void shark::export_csv | ( | Data< Type > const & | set, |
| std::string | fn, | ||
| std::string | separator = ",", |
||
| bool | sci = true, |
||
| unsigned int | width = 0 |
||
| ) |
Format unlabeled data into a character-separated value file.
| set | Container to be exported |
| fn | The file to be written to |
| separator | Separator between entries, typically a comma or a space |
| sci | should the output be in scientific notation? |
| width | argument to std::setw when writing the output |
Definition at line 832 of file Csv.h.
References shark::detail::export_csv().
| void shark::export_csv | ( | LabeledData< InputType, LabelType > const & | dataset, |
| std::string | fn, | ||
| LabelPosition | lp, | ||
| std::string | separator = ",", |
||
| bool | sci = true, |
||
| unsigned int | width = 0 |
||
| ) |
Format labeled data into a character-separated value file.
| dataset | Container to be exported |
| fn | The file to be written to |
| lp | Position of the label in the record, either first or last column |
| separator | Separator between entries, typically a comma or a space |
| sci | should the output be in scientific notation? |
| width | argument to std::setw when writing the output |
Definition at line 903 of file Csv.h.
References shark::Data< Type >::elements(), shark::detail::export_csv(), shark::LabeledData< InputT, LabelT >::inputs(), and shark::LabeledData< InputT, LabelT >::labels().
| void shark::export_kernel_matrix | ( | LabeledData< InputType, LabelType > const & | dataset, |
| AbstractKernelFunction< InputType > & | kernel, | ||
| std::ostream & | out, | ||
| KernelMatrixNormalizationType | normalizer = NONE, |
||
| bool | scientific = false, |
||
| unsigned int | fieldwidth = 0 |
||
| ) |
Write a kernel Gram matrix to stream.
| dataset | data basis for the Gram matrix |
| kernel | pointer to kernel function to be used |
| out | The stream to be written to |
| normalizer | what kind of normalization to apply. see enum declaration for details. |
| sci | should the output be in scientific notation? |
| width | field width for pretty printing |
Definition at line 64 of file PrecomputedMatrix.h.
References shark::CENTER_AND_MULTIPLICATIVE_TRACE_ONE, shark::CENTER_ONLY, shark::Data< Type >::elements(), shark::AbstractKernelFunction< InputTypeT >::eval(), shark::ScaledKernel< InputType >::factor(), shark::LabeledData< InputT, LabelT >::inputs(), shark::LabeledData< InputT, LabelT >::labels(), shark::mean(), shark::MULTIPLICATIVE_TRACE_N, shark::MULTIPLICATIVE_TRACE_ONE, shark::MULTIPLICATIVE_VARIANCE_ONE, shark::NONE, SHARK_ASSERT, SHARKEXCEPTION, shark::size(), shark::DataView< DatasetType >::size(), SIZE_CHECK, shark::trace(), and shark::NormalizeKernelUnitVariance< InputType >::train().
Referenced by shark::export_kernel_matrix().
| void shark::export_kernel_matrix | ( | LabeledData< InputType, LabelType > const & | dataset, |
| AbstractKernelFunction< InputType > & | kernel, | ||
| std::string | fn, | ||
| KernelMatrixNormalizationType | normalizer = NONE, |
||
| bool | sci = false, |
||
| unsigned int | width = 0 |
||
| ) |
Write a kernel Gram matrix to file.
| dataset | data basis for the Gram matrix |
| kernel | pointer to kernel function to be used |
| fn | The filename of the file to be written to |
| normalizer | what kind of normalization to apply. see enum declaration for details. |
| sci | should the output be in scientific notation? |
| width | field width for pretty printing |
Definition at line 227 of file PrecomputedMatrix.h.
References shark::export_kernel_matrix().
| void shark::export_libsvm | ( | LabeledData< InputType, unsigned int > & | dataset, |
| const std::string & | fn, | ||
| bool | dense = false, |
||
| bool | oneMinusOne = true, |
||
| bool | sortLabels = false |
||
| ) |
Export data to LIBSVM format.
| dataset | Container storing the data |
| fn | Output file |
| dense | Flag for using dense output format |
| oneMinusOne | Flag for applying the transformation y<-2y-1 to binary labels |
| sortLabels | Flag for sorting data points according to labels |
Definition at line 119 of file Libsvm.h.
References shark::detail::cmpLabelSortPair(), shark::LabeledData< InputT, LabelT >::element(), shark::inputDimension(), shark::numberOfClasses(), shark::LabeledData< InputT, LabelT >::numberOfElements(), and SHARKEXCEPTION.
| void shark::fft | ( | MatrixT & | data, |
| int | isign | ||
| ) |
Depending on the value of "isign" the "data" is replaced by its discrete Fourier transform or by its inverse discrete Fourier transform.
When isign is set to "-1", then \(H_n \equiv \sum_{k=0}^{N-1} data_k e^{2 \pi i k n / N}\) is calculated.
When isign is set to "1", then \(h_k = \frac{1}{N} \sum_{n=0}^{N-1} data_k e^{-2 \pi i k n / N}\) is calculated.
Please notice, that the number of sample points in data must be an integer power of 2. If there are not enough sample points, create the data array with a size equal to the next higher power of 2 and fill the last empty positions of the array with zero values.
| data | The original sample data (in the time or frequency domain), that will be replaced by its corresponding data of the other domain. |
| isign | Determines the type of transformation. "-1" = a discrete Fourier transform is performed, "1" = an inverse Fourier transform is performed |
Definition at line 118 of file fft.inl.
References shark::size(), and shark::swap().
Referenced by shark::fft(), and shark::ifft().
| void shark::fft | ( | MatrixT & | data | ) |
Replaces the "data" by its discrete Fourier transform.
\(H_n \equiv \sum_{k=0}^{N-1} data_k e^{2 \pi i k n / N}\) is calculated.
Please notice, that the number of sample points in data must be an integer power of 2. If there are not enough sample points, create the data array with a size equal to the next higher power of 2 and fill the last empty positions of the array with zero values.
| data | The original sample data of the time domain, that will be replaced by its corresponding data of the frequency domain. |
Definition at line 260 of file fft.inl.
References shark::fft().
| RealMatrix shark::g_inverse | ( | blas::matrix_expression< MatrixT > const & | matrixA | ) |
Calculates the generalized inverse matrix of input matrix "matrixA".
Given an \( m \times n \) input matrix \( A \) this function uses a Pivoting cholesky decomposition to compute the generalized Inverse with the property \( AA'A = A \)
if m < n the matrix also satisfies AA' = I
If m > n the least squares solution is used.
| matrixA | \( m \times n \) input matrix. |
Referenced by shark::CSvmDerivative< InputType, CacheType >::prepareCSvmParameterDerivative().
| void shark::identity | ( | blas::matrix_expression< Matrix > & | mat | ) |
Initializes the square matrix A to be the identity matrix.
Definition at line 202 of file Tools.h.
References SIZE_CHECK, and shark::zero().
Referenced by shark::decomposedGeneralInverse().
| void shark::ifft | ( | MatrixT & | data | ) |
Replaces the "data" by its inverse discrete Fourier transform.
\(h_k = \frac{1}{N} \sum_{n=0}^{N-1} data_k e^{-2 \pi i k n / N}\) is calculated.
Please notice, that the number of sample points in data must be an integer power of 2. If there are not enough sample points, create the data array with a size equal to the next higher power of 2 and fill the last empty positions of the array with zero values.
| data | The original sample data of the frequency domain, that will be replaced by its corresponding data of the time domain. |
Definition at line 225 of file fft.inl.
References shark::fft(), and shark::size().
| void shark::import_csv | ( | Data< Type > & | data, |
| std::string | fn, | ||
| std::string | separator = ",", |
||
| std::string | comment = "", |
||
| std::size_t | batchSize = Data<Type>::DefaultBatchSize |
||
| ) |
Import unlabeled data from a character-separated value file.
| data | Container storing the loaded data |
| fn | The file to be read from |
| separator | Separator between entries, typically a comma or a space |
| comment | Trailing character indicating comment line |
Definition at line 810 of file Csv.h.
References shark::createDataFromRange(), and shark::detail::import_csv().
Referenced by loadData(), and main().
| void shark::import_csv | ( | LabeledData< InputType, LabelType > & | dataset, |
| std::string | fn, | ||
| LabelPosition | lp, | ||
| std::string | separator = ",", |
||
| std::string | comment = "", |
||
| bool | allowMissingClasses = false, |
||
| std::map< LabelType, LabelType > const * | labelmap = NULL, |
||
| std::size_t | batchSize = LabeledData<InputType, LabelType>::DefaultBatchSize |
||
| ) |
Import labeled data from a character-separated value file.
| dataset | Container storing the loaded data |
| fn | The file to be read from |
| lp | Position of the label in the record, either first or last column |
| separator | Separator between entries, typically a comma or a space |
| comment | Character for indicating a comment, by default empty |
| labelmap | explicit mapping from input data labels to Shark labels (optional) |
Definition at line 853 of file Csv.h.
References shark::createLabeledDataFromRange(), and shark::detail::import_csv().
| void shark::import_csv | ( | LabeledData< RealVector, RealVector > & | dataset, |
| std::string | fn, | ||
| LabelPosition | lp, | ||
| std::string | separator = ",", |
||
| std::string | comment = "", |
||
| std::size_t | numberOfOutputs = 1 |
||
| ) |
Import regression data from a character-separated value file.
| dataset | Container storing the loaded data |
| fn | The file to be read from |
| lp | Position of the label in the record, either first or last column |
| separator | Separator between entries, typically a comma or a space |
| comment | Character for indicating a comment, by default empty |
| numberOfOutputs | Dimensionality of label/output |
Definition at line 879 of file Csv.h.
References shark::createLabeledDataFromRange(), and shark::detail::import_csv_regression().
| void shark::import_libsvm | ( | LabeledData< RealVector, unsigned int > & | dataset, |
| std::istream & | stream, | ||
| int | highestIndex = 0 |
||
| ) |
Import data from a LIBSVM file.
| dataset | container storing the loaded data |
| fn | the file to be read from |
| highestIndex | highest feature index, or 0 for auto-detection |
Definition at line 137 of file LibSVM.cpp.
| void shark::import_libsvm | ( | LabeledData< CompressedRealVector, unsigned int > & | dataset, |
| std::istream & | stream, | ||
| int | highestIndex = 0 |
||
| ) |
Import data from a LIBSVM file.
| dataset | container storing the loaded data |
| fn | the file to be read from |
| highestIndex | highest feature index, or 0 for auto-detection |
Definition at line 145 of file LibSVM.cpp.
| void shark::import_libsvm | ( | LabeledData< RealVector, unsigned int > & | dataset, |
| std::string | fn, | ||
| int | highestIndex = 0 |
||
| ) |
Import data from a LIBSVM file.
| dataset | container storing the loaded data |
| fn | the file to be read from |
| highestIndex | highest feature index, or 0 for auto-detection |
Definition at line 154 of file LibSVM.cpp.
| void shark::import_libsvm | ( | LabeledData< CompressedRealVector, unsigned int > & | dataset, |
| std::string | fn, | ||
| int | highestIndex = 0 |
||
| ) |
Import data from a LIBSVM file.
| dataset | container storing the loaded data |
| fn | the file to be read from |
| highestIndex | highest feature index, or 0 for auto-detection |
Definition at line 163 of file LibSVM.cpp.
| detail::ADLVector<Source&> shark::init | ( | shark::blas::vector_container< Source > & | source | ) |
Starting-point for the initialization sequence.
Usage: init(vector)<<a,b,c where vector is a ublas vector or sub-vector and a,b,c are either scalars or vectors. In debug mode, it is checked that size(vector) == size(a,b,c)
Definition at line 60 of file Initialize.h.
Referenced by shark::blas::axpy_prod(), BOOST_FUSION_ADAPT_TPL_STRUCT(), shark::LooErrorCSvm< InputType, CacheType >::eval(), shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::eval(), shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::evalDerivative(), shark::TypeErasedMultiObjectiveOptimizer< SearchSpace, Optimizer >::init(), shark::blas::opb_prod(), shark::LinearClassifier::parameterVector(), shark::Centroids::parameterVector(), shark::RBFNet::parameterVector(), shark::ProductKernel< InputType >::parameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::parameterVector(), shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::parameterVector(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::parameterVector(), shark::Normalizer< DataType >::parameterVector(), shark::KernelExpansion< InputType >::parameterVector(), shark::AbstractSvmTrainer< InputType, unsigned int >::parameterVector(), shark::AverageEnergyGradient< RBM >::result(), shark::LinearClassifier::setParameterVector(), shark::Centroids::setParameterVector(), shark::RBFNet::setParameterVector(), shark::ProductKernel< InputType >::setParameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::setParameterVector(), shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::setParameterVector(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::setParameterVector(), shark::Normalizer< DataType >::setParameterVector(), shark::KernelExpansion< InputType >::setParameterVector(), shark::AbstractSvmTrainer< InputType, unsigned int >::setParameterVector(), shark::blas::sparse_prod(), shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::weightedDerivatives(), and shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::weightedParameterDerivative().
| detail::ADLVector<const Source&> shark::init | ( | const shark::blas::vector_container< Source > & | source | ) |
Starting-point for the initialization sequence.
Usage: init(vector)<<a,b,c where vector is a ublas vector or sub-vector and a,b,c are either scalars or vectors. In debug mode, it is checked that size(vector) == size(a,b,c)
Definition at line 68 of file Initialize.h.
| detail::ADLVector<shark::blas::vector_range<Source> > shark::init | ( | const shark::blas::vector_range< Source > & | source | ) |
Starting-point for the initialization sequence when used for splitting the vector.
Usage: init(vector)>>a,b,c where vector is a ublas vector or sub-vector and a,b,c are mutable scalars or vectors. In debug mode, it is checked that size(vector) == size(a,b,c) Specialization for ublas vector_range.
Definition at line 84 of file Initialize.h.
| detail::ADLVector<shark::blas::matrix_row<Source> > shark::init | ( | const shark::blas::matrix_row< Source > & | source | ) |
Specialization for matrix rows.
Definition at line 91 of file Initialize.h.
| void shark::initRandomNormal | ( | AbstractModel< InputType, OutputType > & | model, |
| double | s | ||
| ) |
Initialize model parameters normally distributed.
| model,: | model to be initialized |
| s,: | variance of mean-free normal distribution |
Definition at line 316 of file AbstractModel.h.
References shark::IParameterizable::numberOfParameters(), and shark::IParameterizable::setParameterVector().
| void shark::initRandomUniform | ( | AbstractModel< InputType, OutputType > & | model, |
| double | l, | ||
| double | h | ||
| ) |
Initialize model parameters uniformly at random.
| model,: | model to be initialized |
| l,: | lower bound of initialization interval |
| h,: | upper bound of initialization interval |
Definition at line 331 of file AbstractModel.h.
References shark::IParameterizable::numberOfParameters(), and shark::IParameterizable::setParameterVector().
Referenced by experiment(), and main().
| RealMatrix shark::invert | ( | const MatrixT & | mat | ) |
Inverts a matrix with full rank.
Referenced by shark::SpanBoundCSvm< InputType >::eval(), shark::SpanBoundCSvm< InputType >::evalDerivative(), shark::FisherLDA::meanAndScatter(), shark::NegativeAUC< LabelType, OutputType >::NegativeAUC(), shark::NegativeWilcoxonMannWhitneyStatistic< LabelType, OutputType >::NegativeWilcoxonMannWhitneyStatistic(), and shark::KalmanFilter::updateEstimate().
| void shark::invertSymmPositiveDefinite | ( | MatrixT & | I, |
| const MatrixU & | ArrSymm | ||
| ) |
Inverts a symmetric positive definite matrix.
Definition at line 79 of file invert.inl.
References shark::blas::axpy_prod(), shark::choleskyDecomposition(), shark::detail::bindings::potri(), and shark::blas::trans().
| void shark::linmin | ( | VectorT & | p, |
| const VectorT & | xi, | ||
| double & | fret, | ||
| Function | func, | ||
| double | ax = 0.0, |
||
| double | bx = 1.0 |
||
| ) |
Minimizes a function of "N" variables.
Performs a minimization of a function of \( N \) variables, i.e. given as input the vectors \( P \) and \( n \) and the function \( f \), the function finds the scalar \( \lambda \) that minimizes \( f(P + \lambda n) \). \( P \) is replaced by \( P + \lambda n \) and \( n \) by \( \lambda n \).
| p | N-dimensional initial starting point for the search, is set to the point where the function takes on a minimum. |
| xi | N-dimensional search direction, is replaced by the actual vector displacement that p was moved. |
| fret | The function value at the new point p. |
| func | The function that will be minimized. |
| ax | Initial guess for the lower bracket |
| bx | Initial guess for the upper bracket |
| SharkException | the type of the exception will be "size mismatch" and indicates that p is not one-dimensional or that p has not the same size than xi |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 89 of file linmin.inl.
References shark::blas::max(), SIGN, and SIZE_CHECK.
Referenced by shark::LineSearch::operator()().
| void shark::lnsrch | ( | const VectorT & | xold, |
| double | fold, | ||
| VectorT & | g, | ||
| VectorT & | p, | ||
| VectorT & | x, | ||
| double & | f, | ||
| double | stpmax, | ||
| bool & | check, | ||
| Function | func | ||
| ) |
Does a line search, i.e. given a nonlinear function, a starting point and a direction, a new point is calculated where the function has decreased "sufficiently".
Given a nonlinear function, a starting point and a direction, a new point is calculated where the function has decreased "sufficiently".
The line search algorithms are based on the Newton method of approximating root values of monotone decreasing functions. When the derivative \( f' \) of the function \( f \) at a starting point \( x \) on the X-axis can be calculated, the intersection \( x' \) of the tangent at point \( f(x) \) (with gradient \( f'(x) \)) with the X-axis can be used to get a better approximation of the minima at \( x_{min} \).
This function is based on this idea: Given a nonlinear function \( f \), a n-dimensional starting point \( x_{old} \) and a direction \( p \) (known as Newton direction), a new point \( x_{new} \) is calculated as
\( x_{new} = x_{old} + \lambda p, \hspace{2em} 0 < \lambda \leq 1 \)
in a way that \( f(x_{new})\) has decreased sufficiently. Sufficiently means that
\( f(x_{new}) \leq f(x_{old}) + \alpha \nabla f \cdot (x_{new} - x_{old}) \)
where \( 0 < \alpha < 1 \) is a fraction of the initial rate of decrease \( \nabla f \cdot p \).
This function can be used for minimization or solving a set of nonlinear equations of the form \( F(x) = 0 \). Finding the root value (the x-value at which the related function will intersect the X-axis) will then solve the equations. In contrast to cubic line search (cblnsrch.cpp), here the root value can be calculated by only three gradient of function evaluations.
| xold | n-dimensional starting point. |
| fold | The value of function func at point xold. |
| g | The n-dimensional gradient of function func at point xold. |
| p | n-dimensional point to specify the search direction (the Newton direction). |
| x | New n-dimensional point along the direction p from xold where the function func decreases "sufficiently". |
| f | The new function value for point x. |
| stpmax | No. of the steps lnsrch. will do internally, to limit the number avoids evaluating the function in regions where it is undefined or subject to overflow. |
| check | Set to false on a normal exit and to true when x is too close to xold (signals convergence for minimization algorithms). |
| func | Function to decrease. |
| check_exception | the type of the exception will be "size mismatch" and indicates that xold is not one-dimensional or has not the same size than g or p or x |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 129 of file lnsrch.inl.
References shark::blas::max(), SIZE_CHECK, and shark::blas::sum().
| double shark::logdetsymm | ( | MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the log of the determinant of the symmetric matrix "amatA".
Calculates logarithm of the determinant of the symmetric matrix "amatA".
Calculates the logarithm of the determinate of matrix amatA by using its n eigenvalues \( x_j \) that first will be calculated. The determinate is then given as:
\( \prod_{j=1}^n x_j \)
| amatA | \( n \times n \) matrix, which is symmetric, so only the bottom triangular matrix must contain values. At the end of the function amatA always contains the full matrix. |
| vmatA | \( n \times n \) matrix, that will contain the scaled eigenvectors at the end of the function. |
| dvecA | n-dimensional vector that will contain the eigenvalues at the end of the function. |
| SharkException | the type of the eception will be "size mismatch" and indicates that amatA is not a square matrix |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 172 of file detsymm.inl.
References shark::eigensymm_intermediate(), and SIZE_CHECK.
| detail::InitializerRange<typename T::const_iterator,detail::MatrixExpression<const typename T::value_type> > shark::matrixSet | ( | const T & | range | ) |
Uses a range of vectors for initialization.
Sometimes not a single matrix but a set of matrices is needed as argument like std::deque<RealMatrix> set; in this case, matrixSet is needed init(vec)<<vec1,matrixSet(set),vec2;
Definition at line 170 of file Initialize.h.
Referenced by shark::Centroids::parameterVector(), and shark::Centroids::setParameterVector().
| detail::InitializerRange<typename T::iterator,detail::MatrixExpression<typename T::value_type> > shark::matrixSet | ( | T & | range | ) |
Uses a range of vectors for splitting and initialization.
Sometimes not a single matrix but a set of matrices is needed as argument like std::deque<RealMatrix> set; in this case, vectorSet is needed init(vec)>>vec1,matrixSet(set),vec2;
Definition at line 184 of file Initialize.h.
| T shark::maxExpInput | ( | ) |
| VectorType shark::mean | ( | UnlabeledData< VectorType > const & | data | ) |
Definition at line 88 of file VectorStatistics.h.
References shark::mean().
| VectorType shark::mean | ( | Data< VectorType > const & | data | ) |
Calculates the mean vector of array "x".
Calculates the mean vector of the input vectors.
Given a d -dimensional array x with size N1 x ... x Nd, this function calculates the mean vector given as:
\[ mean_j = \frac{1}{N1} \sum_{i=1}^{N1} x_{i,j} \]
Example:
\[ \left( \begin{array}{*{4}{c}} 1 & 2 & 3 & 4\\ 5 & 6 & 7 & 8\\ 9 & 10 & 11 & 12\\ \end{array} \right) \longrightarrow \frac{1}{3} \left( \begin{array}{*{4}{c}} 1+5+9 & 2+6+10 & 3+7+11 & 4+8+12\\ \end{array} \right) \longrightarrow \left( \begin{array}{*{4}{c}} 5 & 6 & 7 & 8\\ \end{array} \right) \]
| data | input data, from which the mean value will be calculated |
Definition at line 135 of file VectorStatistics.inl.
References shark::Data< Type >::batches(), shark::dataDimension(), shark::Data< Type >::empty(), shark::Data< Type >::numberOfElements(), SIZE_CHECK, and shark::sumRows().
Referenced by shark::covariance(), createData(), shark::BaseRng< RNG >::diffGeom(), shark::export_kernel_matrix(), shark::BaseRng< RNG >::gauss(), shark::LinearClassifier::importCovarianceMatrix(), shark::McPegasos< VectorType >::lossGradientADS(), shark::McPegasos< VectorType >::lossGradientATS(), shark::mean(), shark::FisherLDA::meanAndScatter(), shark::meanvar(), shark::Statistics::operator()(), shark::BaseRng< RNG >::poisson(), shark::LinearClassifier::setClassMean(), shark::PCA::setData(), shark::LinearClassifier::setParameterVector(), shark::QpMcDecomp< Matrix >::solveForBias(), shark::NormalizeComponentsWhitening< VectorType >::train(), shark::NormalTrainer::train(), shark::NormalizeComponentsUnitVariance< DataType >::train(), shark::LDA::train(), shark::FisherLDA::train(), shark::QpMcLinearATS::updateWeightVectors(), and shark::QpMcLinearATM::updateWeightVectors().
| blas::vector< typename MatrixType::value_type > shark::mean | ( | blas::matrix_container< MatrixType > const & | data | ) |
Calculates the mean vector of the input vectors.
Definition at line 150 of file VectorStatistics.inl.
References shark::mean(), SIZE_CHECK, and shark::sumRows().
| void shark::meanvar | ( | Data< Vec1T > const & | data, |
| blas::vector_container< Vec2T > & | meanVec, | ||
| blas::vector_container< Vec3T > & | varianceVec | ||
| ) |
Calculates the mean and variance values of a dataset.
Calculates the mean and variance values of the input data.
Given the vector of data, the mean and variance values are calculated as in the functions mean and variance.
| data | Input data. |
| meanVec | Vector of mean values. |
| varianceVec | Vector of variances. |
Definition at line 15 of file VectorStatistics.inl.
References shark::Data< Type >::batches(), shark::dataDimension(), shark::Data< Type >::empty(), shark::mean(), shark::blas::noalias(), shark::Data< Type >::numberOfElements(), shark::repeat(), SIZE_CHECK, shark::sqr(), and shark::sumRows().
Referenced by shark::covariance(), main(), shark::PCA::setData(), shark::NormalizeComponentsWhitening< VectorType >::train(), shark::NormalizeComponentsUnitVariance< DataType >::train(), and shark::variance().
| void shark::meanvar | ( | blas::matrix_container< MatT > & | data, |
| blas::vector_container< Vec1T > & | meanVec, | ||
| blas::vector_container< Vec2T > & | varianceVec | ||
| ) |
Calculates the mean, variance and covariance values of the input data.
Definition at line 41 of file VectorStatistics.inl.
References shark::mean(), shark::blas::noalias(), shark::repeat(), SIZE_CHECK, shark::sqr(), and shark::sumRows().
| void shark::meanvar | ( | const Data< Vec1T > & | data, |
| blas::vector_container< Vec2T > & | meanVec, | ||
| blas::matrix_container< MatT > & | covariance | ||
| ) |
Calculates the mean and covariance values of a set of data.
Calculates the mean, variance and covariance values of the input data.
Given the vector of data, the mean and variance values are calculated as in the functions mean and variance.
| data | Input data. |
| meanVec | Vector of mean values. |
| covariance | Covariance matrix. |
Definition at line 75 of file VectorStatistics.inl.
References shark::Data< Type >::batch(), shark::covariance(), shark::dataDimension(), shark::Data< Type >::empty(), shark::mean(), shark::Data< Type >::numberOfBatches(), shark::Data< Type >::numberOfElements(), shark::repeat(), SIZE_CHECK, shark::symmRankKUpdate(), and shark::blas::trans().
| T shark::minExpInput | ( | ) |
| std::size_t shark::nonzeroElements | ( | blas::vector_expression< V > const & | vec | ) |
Definition at line 185 of file Tools.h.
Referenced by shark::Batch< shark::blas::compressed_vector< T > >::createBatchFromRange().
| double shark::nth | ( | std::vector< T > & | v, |
| unsigned | n | ||
| ) |
return nth element after sorting
Definition at line 155 of file VectorStatistics.h.
References SHARKEXCEPTION.
| std::ostream& shark::operator<< | ( | std::ostream & | stream, |
| const Data< T > & | d | ||
| ) |
| void shark::orthoNormalize | ( | blas::matrix_container< MatrixT > & | matrixC | ) |
Transforms a quadratic matrix such that it forms an orthonormal Basis using Gram-Schmidt Orthonormalisation.
Definition at line 52 of file rotations.h.
References shark::blas::inner_prod(), shark::blas::norm_2(), shark::blas::row(), shark::size(), and SIZE_CHECK.
Referenced by shark::randomRotationMatrix().
| detail::ParameterizableExpression<const T> shark::parameters | ( | const T & | object | ) |
Uses the parameters of a parameterizable object for initialization.
The object doesn't have to be derived from IParameterizable, but needs to offer the methods
Definition at line 119 of file Initialize.h.
Referenced by BOOST_FUSION_ADAPT_TPL_STRUCT(), shark::NestedGridSearch::configure(), shark::PointSearch::configure(), shark::LooErrorCSvm< InputType, CacheType >::eval(), shark::SvmLogisticInterpretation< InputType >::eval(), shark::KernelTargetAlignment< InputType >::evalDerivative(), shark::SvmLogisticInterpretation< InputType >::evalDerivative(), shark::detail::ParameterizableExpression< T >::init(), shark::detail::ParameterizableExpression< T * >::init(), shark::detail::ParameterizableExpression< T *const >::init(), shark::PointSearch::init(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::numberOfParameters(), shark::LinearClassifier::parameterVector(), shark::NearestNeighborClassifier< InputType >::parameterVector(), shark::SoftNearestNeighborClassifier< InputType >::parameterVector(), shark::NearestNeighborRegression< InputType >::parameterVector(), shark::RBFNet::parameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::parameterVector(), shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::parameterVector(), shark::FFNet< HiddenNeuron, OutputNeuron >::parameterVector(), shark::AbstractSvmTrainer< InputType, unsigned int >::parameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::setParameterVector(), shark::detail::ConcatenatedModelWrapper< InputType, IntermediateType, OutputType >::setParameterVector(), shark::AbstractSvmTrainer< InputType, unsigned int >::setParameterVector(), shark::detail::ParameterizableExpression< T >::split(), shark::detail::ParameterizableExpression< T * >::split(), and shark::detail::ParameterizableExpression< T *const >::split().
| detail::ParameterizableExpression<T> shark::parameters | ( | T & | object | ) |
Uses the parameters of a parameterizable object for initialization.
The object doesn't have to be derived from IParameterizable, but needs to offer the methods
Definition at line 126 of file Initialize.h.
| detail::InitializerRange<typename T::const_iterator, detail::ParameterizableExpression<const typename T::value_type> > shark::parameterSet | ( | const T & | range | ) |
Uses a range of parametrizable objects for initialization.
The objects in the set must offer the methods described by IParameterizable. Also pointer to objects are allowed
Definition at line 196 of file Initialize.h.
Referenced by shark::ProductKernel< InputType >::parameterVector(), and shark::ProductKernel< InputType >::setParameterVector().
| detail::InitializerRange<typename T::iterator, detail::ParameterizableExpression<typename T::value_type> > shark::parameterSet | ( | T & | range | ) |
Uses a range of parametrizable objects for splitting and initialization.
The objects in the set must offer the methods described by IParameterizable. Also pointer to objects are allowed
Definition at line 207 of file Initialize.h.
| std::size_t shark::pivotingRQ | ( | blas::matrix_expression< MatrixT > const & | matrixA, |
| blas::matrix_container< Mat > & | matrixR, | ||
| blas::matrix_container< Mat > & | matrixQ, | ||
| blas::permutation_matrix< std::size_t > & | permutation | ||
| ) |
Calculates the RQ-Decomposition of A using pivoting.
Determines the RQ Decomposition of the matrix A using pivoting.
The pivoting RQ-Decomposition finds an orthonormal matrix Q and a lower Triangular matrix R as well as a permuation matrix P such that PA = R*Q. This function is better known as the QR-Decomposition of a transposed matrix B^T = A and B = QR. We depart from the well known algorithm because it is intended to be used with column major matrices. But since shark uses row-major, a QR decomposition is a lot slower.
This Version of the algorithm is based on householder transformations. since it uses pivoting it can be used to determine the rank of a matrix.
Definition at line 150 of file pivotingRQ.inl.
References shark::fast_prod(), shark::pivotingRQHouseholder(), shark::rank(), shark::rows(), shark::symmRankKUpdate(), and shark::blas::trans().
| std::size_t shark::pivotingRQHouseholder | ( | blas::matrix_expression< MatrixT > const & | matrixA, |
| blas::matrix_container< MatrixU > & | matrixR, | ||
| blas::matrix_container< MatrixU > & | householderTransform, | ||
| blas::permutation_matrix< std::size_t > & | permutation | ||
| ) |
Determines the RQ Decomposition of the matrix A using pivoting returning the housholder transformation instead of Q.
The pivoting RQ-Decomposition finds an orthonormal matrix Q and a lower Triangular matrix R as well as a permuation matrix P such that PA = R*Q. Since Q is the multiplication of all householder transformations, It is quite expensive to compute. Often, Q is only an intermediate step in computations which can be carried out more efficiently using the Householder Transformations themselves.
The Matrix format of the householder transform is that the transformations are stored as upper triangular matrix. The first transformation being in the first row and so on.
Definition at line 57 of file pivotingRQ.inl.
References shark::applyHouseholderOnTheRight(), shark::blas::column(), shark::createHouseholderReflection(), shark::blas::min(), shark::blas::noalias(), shark::blas::norm_sqr(), shark::rank(), shark::blas::vector< T, A >::resize(), shark::blas::row(), shark::sqr(), shark::blas::subrange(), shark::sumColumns(), and shark::swap().
Referenced by shark::pivotingRQ().
| void shark::randomRotationMatrix | ( | blas::matrix_container< MatrixT > & | matrixC, |
| RngType & | rng | ||
| ) |
Initializes a matrix such that it forms a random rotation matrix.
The matrix needs to be quadratic and have the proper size (e.g. call matrix::resize before).
Definition at line 72 of file rotations.h.
References shark::orthoNormalize(), shark::size(), and SIZE_CHECK.
Referenced by shark::ELLI1::init(), shark::ELLI2::init(), shark::IHR3::init(), shark::IHR1::init(), shark::IHR4::init(), shark::IHR2::init(), shark::CIGTAB1::init(), shark::IHR6::init(), shark::CIGTAB2::init(), and shark::randomRotationMatrix().
| void shark::randomRotationMatrix | ( | blas::matrix_container< MatrixT > & | matrixC | ) |
Initializes a matrix such that it forms a random rotation.
matrix. The matrix needs to be quadratic and have the proper size (e.g. call matrix::resize before) uses the global RNG.
Definition at line 92 of file rotations.h.
References shark::randomRotationMatrix().
| RealMatrix shark::randomRotationMatrix | ( | size_t | size, |
| RngType & | rng | ||
| ) |
Creates a random rotation matrix with a certain size using the random number generator rng.
Definition at line 98 of file rotations.h.
References shark::randomRotationMatrix().
|
inline |
Creates a random rotation matrix with a certain size using the global random number gneerator.
Definition at line 105 of file rotations.h.
References shark::randomRotationMatrix().
| unsigned shark::rank | ( | const MatrixT & | amatA, |
| const MatrixU & | vmatA, | ||
| const VectorT & | dvecA | ||
| ) |
Determines the rank of the symmetric matrix "amatA".
Determines the numerical rank of the symmetric matrix "amatA".
Given the \( n \times n \) matrix amatA, this function uses the eigenvectors vmatA and the eigenvalues dvecA to calculate the rank of amatA.
| amatA | \( n \times n \) matrix, which is symmetric and is not changed by the function. |
| vmatA | \( n \times n \) matrix with the normalized eigenvectors, each column contains one eigenvector. The matrix is not changed by the function. |
| dvecA | n-dimensional vector of the eigenvalues, given in descending order. The vector is not changed by the function. |
| SharkException | the type of the exception will be "size mismatch" and indicates that dvecA is not one-dimensional or that amatA or vmatA has not the same number of rows or columns than dvecA contains values |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 81 of file rank.inl.
References shark::eigenerr(), and SIZE_CHECK.
Referenced by shark::decomposedGeneralInverse(), shark::WilcoxonRankSumTest::operator()(), shark::SteadyStateIndicatorBasedSelection< Indicator >::operator()(), shark::IndicatorBasedSelection< shark::HypervolumeIndicator >::operator()(), shark::ApproximatedHypervolumeSelection::operator()(), shark::pivotingRQ(), shark::pivotingRQHouseholder(), shark::rankDecomp(), shark::detail::AGE2::step(), and shark::detail::SteadyStateMOCMA< Indicator >::step().
| unsigned shark::rankDecomp | ( | MatrixT & | amatA, |
| MatrixU & | vmatA, | ||
| MatrixV & | hmatA, | ||
| VectorT & | dvecA | ||
| ) |
Calculates the rank of the symmetric matrix "amatA", its eigenvalues and eigenvectors.
Determines the rank of matrix amatA and additionally calculates the eigenvalues and eigenvectors of the matrix. Empty eigenvalues (i.e. eigenvalues equal to zero) are set to the greatest calculated eigenvalue and each eigenvector \( x_j \) is scaled by multiplying it with the scalar value \( \frac{1}{\sqrt{\lambda_j}} \).
| amatA | \( n \times n \) matrix, which is symmetric, so only the bottom triangular matrix must contain values. At the end of the function amatA always contains the full matrix. |
| vmatA | \( n \times n \) matrix, which will contain the scaled eigenvectors. |
| hmatA | \( n \times n \) matrix, that is used to store intermediate results. |
| dvecA | n-dimensional vector with the calculated eigenvalues in descending order. |
| SharkException | the type of the exception will be "size mismatch" and indicates that amatA is not a square matrix |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 89 of file rankDecomp.inl.
References shark::eigensymm_intermediate(), shark::rank(), and SIZE_CHECK.
Referenced by shark::FisherLDA::train().
| VectorRepeater<Vector> shark::repeat | ( | blas::vector_expression< Vector > const & | vector, |
| std::size_t | rows | ||
| ) |
Creates a matrix from a vector by repeating the vector in every row of the matrix.
example: vector = (1,2,3) repeat(vector,3) results in (1,2,3) (1,2,3) (1,2,3)
| vector | the vector which is to be repeated as the rows of the resulting matrix |
| rows | the number of rows of the matrix |
Definition at line 380 of file Tools.h.
References shark::rows().
Referenced by shark::blas::detail::distanceSqrBlockBlock(), shark::SigmoidModel::eval(), shark::MissingFeaturesKernelExpansion< InputType >::eval(), shark::RBFNet::eval(), shark::PolynomialKernel< InputType >::eval(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::eval(), shark::Normalizer< DataType >::eval(), shark::FFNet< HiddenNeuron, OutputNeuron >::eval(), shark::KernelExpansion< InputType >::eval(), shark::PolynomialMutator::init(), shark::SimulatedBinaryCrossover< RealVector >::init(), shark::detail::RealCodedNSGAII< Indicator >::init(), shark::detail::SMSEMOA::init(), shark::meanvar(), shark::PCA::setData(), shark::LinearRegression::train(), and shark::Softmax::weightedInputDerivative().
| boost::enable_if<boost::is_arithmetic<T>, blas::scalar_vector<T> >::type shark::repeat | ( | T | scalar, |
| std::size_t | elements | ||
| ) |
| boost::enable_if<boost::is_arithmetic<T>, blas::scalar_matrix<T> >::type shark::repeat | ( | T | scalar, |
| std::size_t | rows, | ||
| std::size_t | columns | ||
| ) |
brief repeats a single element to form a matrix of size rows x columns
| scalar | the value which is repeated |
| rows | the number of rows of the resulting vector |
| columns | the number of columns of the resulting vector |
Definition at line 401 of file Tools.h.
References shark::columns(), and shark::rows().
| blas::matrix_range<Matrix const> shark::rows | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | start, | ||
| std::size_t | end | ||
| ) |
brief picks a subrange of rows from a matrix. much easier to use than subrange
Definition at line 407 of file Tools.h.
References shark::blas::subrange().
Referenced by shark::detail::SparseFFNetErrorWrapper< HiddenNeuron, OutputNeuron >::eval(), shark::pivotingRQ(), shark::repeat(), shark::PCA::setData(), and shark::FFNet< HiddenNeuron, OutputNeuron >::weightedParameterDerivative().
| blas::matrix_range<Matrix> shark::rows | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | start, | ||
| std::size_t | end | ||
| ) |
brief picks a subrange of rows from a matrix. much easier to use than subrange
Definition at line 412 of file Tools.h.
References shark::blas::subrange().
| T shark::safeExp | ( | T | x | ) |
Thresholded exp function, over- and underflow safe.
Replaces the value of exp(x) for numerical reasons by the a threshold value if it gets too large. Use it only, if there is no other way to get the function stable!
| x | the exponent |
Definition at line 110 of file Math.h.
References shark::blas::max().
| T shark::safeLog | ( | T | x | ) |
Thresholded log function, over- and underflow safe.
Replaces the value of log(x) for numerical reasons by the a threshold value if it gets too low. Use it only, if there is no other way to get the function stable!
| x | the exponent |
Definition at line 125 of file Math.h.
Referenced by shark::NBClassifier< InputType, OutputType >::eval(), shark::AbstractDistribution::logP(), shark::Normal< RngType >::logP(), and shark::SigmoidFitPlatt::train().
| boost::enable_if<boost::is_arithmetic<T>, T>::type shark::sigmoid | ( | T | x | ) |
Logistic function/logistic function.
Calculates the sigmoid function 1/(1+exp(-x)). The type must be arithmetic. For example float,double,long double, int,... but no custom Type.
Definition at line 93 of file Math.h.
Referenced by shark::CrossEntropyIndependent::evalDerivative(), shark::CrossEntropy::evalDerivative(), shark::LogisticNeuron::function(), and shark::BinaryLayer::sufficientStatistics().
| boost::enable_if<boost::is_arithmetic<T>, T>::type shark::softPlus | ( | T | x | ) |
Numerically stable version of the function log(1+exp(x)).
Numerically stable version of the function log(1+exp(x)). This function is the integral of the famous sigmoid function. The type must be arithmetic. For example float,double,long double, int,... but no custom Type.
Definition at line 145 of file Math.h.
Referenced by shark::BinaryLayer::logMarginalize(), and shark::detail::updateLogPartition().
|
inline |
Calculates x^2.
Definition at line 79 of file Math.h.
Referenced by shark::blas::detail::diagonalMahalanobisDistanceSqr(), shark::blas::diagonalMahalanobisNormSqr(), shark::GSP::eval(), shark::LZ9::eval(), shark::Cigar::eval(), shark::CigarDiscus::eval(), shark::Himmelblau::eval(), shark::ZDT2::eval(), shark::LZ1::eval(), shark::LZ7::eval(), shark::LZ4::eval(), shark::LZ5::eval(), shark::LZ6::eval(), shark::LZ8::eval(), shark::LZ3::eval(), shark::Discus::eval(), shark::Fonseca::eval(), shark::LZ2::eval(), shark::IHR2::eval(), shark::IHR4::eval(), shark::CIGTAB1::eval(), shark::Rosenbrock::eval(), shark::IHR6::eval(), shark::Ellipsoid::eval(), shark::DTLZ3::eval(), shark::DTLZ1::eval(), shark::DTLZ4::eval(), shark::DTLZ5::eval(), shark::Rosenbrock::evalDerivative(), shark::FastSigmoidNeuron::functionDerivative(), shark::CARTTrainer::gini(), shark::RFTrainer::gini(), shark::PointSearch::init(), shark::GaussianLayer::logMarginalize(), shark::Normal< RngType >::logP(), shark::meanvar(), shark::RecurrentStructure::neuronDerivative(), shark::cma::Initializer::operator()(), shark::HMGSelectionCriterion::operator()(), shark::Normal< RngType >::p(), shark::Cauchy< RngType >::p(), shark::pivotingRQHouseholder(), shark::LeastContributorApproximator< Rng, ExactHypervolume >::sample(), shark::ARDKernelUnconstrained< InputType >::setParameterVector(), shark::SimpleSigmoidModel::sigmoidDerivative(), shark::KDTree< InputT >::squaredDistanceLowerBound(), shark::KernelMeanClassifier< InputType >::train(), shark::CMA::updateStrategyParameters(), shark::LinearNorm::weightedInputDerivative(), shark::ARDKernelUnconstrained< InputType >::weightedParameterDerivative(), and shark::WeightedSumKernel< InputType >::weightedParameterDerivative().
| double shark::stl_correlation | ( | std::vector< T > & | v1, |
| std::vector< T > & | v2 | ||
| ) |
| double shark::stl_mean | ( | std::vector< T > | v | ) |
compute mean
Definition at line 166 of file VectorStatistics.h.
References SHARKEXCEPTION, and shark::blas::sum().
| double shark::stl_median | ( | std::vector< T > & | v | ) |
| double shark::stl_percentile | ( | std::vector< T > & | v, |
| double | p = .25 |
||
| ) |
compute percentilee (Excel way)
Definition at line 136 of file VectorStatistics.h.
References SHARKEXCEPTION.
| double shark::stl_variance | ( | std::vector< T > & | v, |
| bool | unbiased = true |
||
| ) |
compute sample variance
Definition at line 218 of file VectorStatistics.h.
References SHARKEXCEPTION, and shark::blas::sum().
| void shark::string2data | ( | Data< InputType > & | data, |
| const std::string & | dataInString, | ||
| std::size_t | batchSize = Data<InputType>::DefaultBatchSize |
||
| ) |
Construct Shark data from a string.
| data | Container storing the constructed data |
| dataInString | the string to construct data from |
Definition at line 920 of file Csv.h.
References shark::createDataFromRange(), and shark::detail::import_csv().
| void shark::string2data | ( | LabeledData< InputType, LabelType > & | dataset, |
| const std::string & | dataInString, | ||
| LabelPosition | labelPosition = LAST_COLUMN, |
||
| bool | allowMissingFeatures = false, |
||
| bool | allowMissingClasses = false, |
||
| std::map< LabelType, LabelType > const * | labelmap = NULL, |
||
| std::size_t | batchSize = LabeledData<InputType, LabelType>::DefaultBatchSize |
||
| ) |
Construct Shark labeled data from a string.
| dataset | Container storing the constructed labeled data |
| dataInString | the string to construct data from |
| labelPosition | the position of label. The default value is LAST_COLUMN |
| labelmap | explicit mapping from input data labels to Shark labels (optional) |
Definition at line 938 of file Csv.h.
References shark::createLabeledDataFromRange(), and shark::detail::import_csv().
| void shark::svd | ( | const MatrixT & | amatA, |
| MatrixU & | umatA, | ||
| MatrixU & | vmatA, | ||
| VectorT & | w, | ||
| unsigned | maxIterations = 200, |
||
| bool | ignoreThreshold = true |
||
| ) |
Determines the singular value decomposition of a rectangular matrix "amatA".
Given a \( m \times n \) matrix amatA, this routine computes its singular value decomposition, defined as
\( A = UWV^T \)
where W is an \( n \times n \) diagonal matrix with positive or zero elements, the so-called singular values. The matrices U and V are each orthogonal in the sense that their columns are orthonormal, i.e.
\( UU^T = VV^T = V^TV = 1 \)
| amatA | The input matrix A, with size \( m \times n \) and \( m \geq n \). |
| umatA | The \( m \times n \) column-orthogonal matrix U determined by the function. |
| vmatA | The \( n \times n \) orthogonal matrix V determined by the function. |
| w | n-dimensional vector with the calculated singular values. |
| maxIterations | Number of iterations after which the algorithm gives up, if the solution has still not converged. Default is 200 Iterations. |
| ignoreThreshold | If set to false, the method throws an exception if the threshold maxIterations is exceeded. Otherwise it uses the approximate intermediate results in the further calculations. The default is true. |
| convergence | exception, if the solution has not converged after maxIterations iterations. |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 99 of file svd.inl.
References shark::blas::column(), shark::copySign(), shark::blas::max(), shark::blas::min(), shark::blas::noalias(), shark::blas::norm_1(), and SHARKEXCEPTION.
| unsigned shark::svdrank | ( | const MatrixT & | amatA, |
| MatrixU & | umatA, | ||
| MatrixU & | vmatA, | ||
| VectorT & | wvecA | ||
| ) |
Determines the numerical rank of a rectangular matrix "amatA", when a singular value decomposition for "amatA" has taken place before.
For a singular value decomposition defined as
\( A = UWV^T \)
the resulting orthogonal matrices U and V and the singular values in W sorted by descending order are used to determine the rank of input matrix A.
| amatA | The \( m \times n \) input matrix A, with \( m \geq n \). |
| umatA | The \( m \times n \) column-orthogonal matrix U. |
| vmatA | The \( n \times n \) orthogonal matrix V. |
| wvecA | n-dimensional vector containing the singular values in descending order. |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 91 of file svdrank.inl.
| void shark::svdsort | ( | MatrixU & | umatA, |
| MatrixU & | vmatA, | ||
| VectorT & | wvecA | ||
| ) |
Sorts the singular values in vector "wvecA" by descending order.
For a singular value decomposition defined as
\( A = UWV^T \)
the resulting orthogonal matrices U and V and the singular values in W can be sorted in a way, that the singular values are given in descending order, when leaving the function.
| umatA | The \( m \times n \) matrix U. |
| vmatA | The \( n \times n \) matrix V. |
| wvecA | n-dimensional vector containing the singular values. |
| SharkException | the type of the exception will be "size mismatch" and indicates that wvecA is not one-dimensional |
Please follow the link to view the source code of the example. The example can be executed in the example directory of package LinAlg.
Definition at line 86 of file svdsort.inl.
| detail::MatrixExpression<const Matrix> shark::toVector | ( | const shark::blas::matrix_expression< Matrix > & | matrix | ) |
Linearizes a matrix as a set of row vectors and treats them as a set of vectors for initialization.
Definition at line 105 of file Initialize.h.
Referenced by shark::LinearClassifier::parameterVector(), shark::RBFNet::parameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::parameterVector(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::parameterVector(), shark::KernelExpansion< InputType >::parameterVector(), shark::AverageEnergyGradient< RBM >::result(), shark::LinearClassifier::setParameterVector(), shark::RBFNet::setParameterVector(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::setParameterVector(), shark::detail::LinearModelWrapper< Matrix, InputType, OutputType >::setParameterVector(), and shark::KernelExpansion< InputType >::setParameterVector().
| detail::MatrixExpression<Matrix> shark::toVector | ( | shark::blas::matrix_expression< Matrix > & | matrix | ) |
Linearizes a matrix as a set of row vectors and treats them as a set of vectors for initialization.
Definition at line 110 of file Initialize.h.
| MatrixT::value_type shark::trace | ( | blas::matrix_expression< MatrixT > const & | m | ) |
Evaluates the sum of the values at the diagonal of matrix "v".
Example:
\[ \left( \begin{array}{*{4}{c}} {\bf 1} & 5 & 9 & 13\\ 2 & {\bf 6} & 10 & 14\\ 3 & 7 & {\bf 11} & 15\\ 4 & 8 & 12 & {\bf 16}\\ \end{array} \right) \longrightarrow 1 + 6 + 11 + 16 = 34 \]
| v | square matrix |
Definition at line 451 of file Tools.h.
References SIZE_CHECK.
Referenced by shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::eval(), shark::NegativeGaussianProcessEvidence< InputType, OutputType, LabelType >::evalDerivative(), shark::export_kernel_matrix(), and shark::NormalizeKernelUnitVariance< InputType >::train().
| blas::vector_range<blas::matrix_column<Matrix const> > shark::triangularColumn | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th column of the matrix excluding the zero elements.
Definition at line 352 of file Tools.h.
References shark::blas::column(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_column<Matrix> > shark::triangularColumn | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th column of the matrix excluding the zero elements.
Definition at line 361 of file Tools.h.
References shark::blas::column(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_row<Matrix const> > shark::triangularRow | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | i | ||
| ) |
Returns the i-th row of an upper triangular matrix excluding the elements right of the diagonal.
Definition at line 294 of file Tools.h.
References shark::blas::row(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_row<Matrix> > shark::triangularRow | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | i | ||
| ) |
Returns the i-th row of an upper triangular matrix excluding the elements right of the diagonal.
Definition at line 303 of file Tools.h.
References shark::blas::row(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_column<Matrix const> > shark::unitTriangularColumn | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th column of the matrix below the diagonal.
Definition at line 333 of file Tools.h.
References shark::blas::column(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_column<Matrix> > shark::unitTriangularColumn | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th column of the matrix below the diagonal.
Definition at line 342 of file Tools.h.
References shark::blas::column(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_row<Matrix const> > shark::unitTriangularRow | ( | blas::matrix_expression< Matrix > const & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th row of a lower triangular matrix left of the diagonal.
Definition at line 312 of file Tools.h.
References shark::blas::row(), and shark::blas::subrange().
| blas::vector_range<blas::matrix_row<Matrix> > shark::unitTriangularRow | ( | blas::matrix_expression< Matrix > & | mat, |
| std::size_t | i | ||
| ) |
Returns the elements in the i-th row of a lower triangular matrix left of the diagonal.
Definition at line 321 of file Tools.h.
References shark::blas::row(), and shark::blas::subrange().
| VectorType shark::variance | ( | const Data< VectorType > & | data | ) |
Calculates the variance vector of array "x".
Calculates the variance vector of the input vectors.
Given a d -dimensional array x with size N1 x ... x Nd and mean value vector m, this function calculates the variance vector given as:
\[ variance = \frac{1}{N1} \sum_{i=1}^{N1} (x_i - m_i)^2 \]
| data | input data from which the variance will be calculated |
Definition at line 176 of file VectorStatistics.inl.
References shark::meanvar().
Referenced by shark::Statistics::operator()(), shark::NormalTrainer::train(), and shark::NormalizeComponentsUnitVariance< DataType >::train().
| detail::InitializerRange<typename T::const_iterator,detail::VectorExpression<const typename T::value_type&> > shark::vectorSet | ( | const T & | range | ) |
Uses a range of vectors for initialization.
Sometimes not a single vector but a set of vectors is needed as argument like std::deque<RealVector> set; in this case, vectorSet is needed init(vec)<<vec1,vectorSet(set),vec2;
Definition at line 141 of file Initialize.h.
| detail::InitializerRange<typename T::iterator,detail::VectorExpression<typename T::value_type&> > shark::vectorSet | ( | T & | range | ) |
Uses a range of vectors for splitting and initialization.
Sometimes not a single vector but a set of vectors is needed as argument like std::deque<RealVector> set; in this case, vectorSet is needed init(vec)>>vec1,vectorSet(set),vec2;
Definition at line 155 of file Initialize.h.
| void shark::detail::writePGM | ( | const char * | fileName, |
| const unsigned char * | pData, | ||
| const unsigned int | sx, | ||
| const unsigned int | sy | ||
| ) |
Writes a PGM file.
| fileName | File to write to |
| pData | unsigned char pointer to the data |
| sx | Width of image |
| sy | Height of image |
Definition at line 113 of file Pgm.h.
Referenced by shark::exportPGM().
| void shark::zero | ( | blas::matrix_expression< Matrix > & | mat | ) |
Zeros a matrix. If it is sparse, the structure is preserved.
Definition at line 151 of file Tools.h.
References shark::detail::zero().
Referenced by shark::approxSolveSymmSystem(), shark::calculateKernelMatrixParameterDerivative(), shark::blas::detail::distanceSqrBlockVector(), shark::OnlineRNNet::eval(), shark::RFClassifier::eval(), shark::MissingFeaturesKernelExpansion< InputType >::eval(), shark::SoftNearestNeighborClassifier< InputType >::eval(), shark::NearestNeighborClassifier< InputType >::eval(), shark::NearestNeighborRegression< InputType >::eval(), shark::OneVersusOneClassifier< InputType >::eval(), shark::WeightedSumKernel< InputType >::eval(), shark::DenoisingAutoencoderError< InputType, RngType >::evalDerivative(), shark::KernelTargetAlignment< InputType >::evalDerivative(), shark::SvmLogisticInterpretation< InputType >::evalDerivative(), shark::detail::fast_prod_dense(), shark::detail::fast_prod_detail(), shark::detail::fast_prod_impl(), shark::CARTTrainer::hist(), shark::identity(), shark::FisherLDA::meanAndScatter(), shark::OnlineRNNet::resetInternalState(), shark::PCA::setData(), shark::RBM< VisibleLayerT, HiddenLayerT, RngT >::setStructure(), shark::RecurrentStructure::setStructure(), shark::symmRankKUpdate(), shark::detail::bindings::syrk(), shark::detail::SubrangeKernelWrapper< InputType >::weightedInputDerivative(), shark::ARDKernelUnconstrained< InputType >::weightedInputDerivative(), shark::OnlineRNNet::weightedParameterDerivative(), shark::ARDKernelUnconstrained< InputType >::weightedParameterDerivative(), and shark::FFNet< HiddenNeuron, OutputNeuron >::weightedParameterDerivative().
| void shark::zero | ( | blas::vector_expression< Vector > & | vec | ) |
Zeros a matrix. If it is sparse, the structure is preserved.
Definition at line 156 of file Tools.h.
References shark::detail::zero().
| void shark::zero | ( | blas::matrix_range< Matrix > | mat | ) |
Zeros a subrange of a matrix. If it is sparse, the structure is preserved.
Definition at line 162 of file Tools.h.
References shark::detail::zero().
| void shark::zero | ( | blas::vector_range< Vector > | vec | ) |
Zeros a subrange of a vector. If it is sparse, the structure is preserved.
Definition at line 167 of file Tools.h.
References shark::detail::zero().
| void shark::zero | ( | blas::matrix_row< Vector > | vec | ) |
Zeros a row of a matrix. If it is sparse, the structure is preserved.
Definition at line 172 of file Tools.h.
References shark::detail::zero().
| void shark::zero | ( | blas::matrix_column< Vector > | vec | ) |
Zeros a column of a matrix. If it is sparse, the structure is preserved.
Definition at line 177 of file Tools.h.
References shark::detail::zero().
|
static |
Constant for sqrt( 2 * pi ).
Definition at line 63 of file Math.h.
Referenced by shark::GaussianLayer::logMarginalize(), shark::Normal< RngType >::logP(), and shark::Normal< RngType >::p().