TooN 2.0.0-beta8
|
00001 // -*- c++ -*- 00002 00003 // Copyright (C) 2005,2009 Tom Drummond (twd20@cam.ac.uk), 00004 // Ed Rosten (er258@cam.ac.uk), Gerhard Reitmayr (gr281@cam.ac.uk) 00005 // 00006 // This file is part of the TooN Library. This library is free 00007 // software; you can redistribute it and/or modify it under the 00008 // terms of the GNU General Public License as published by the 00009 // Free Software Foundation; either version 2, or (at your option) 00010 // any later version. 00011 00012 // This library is distributed in the hope that it will be useful, 00013 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00014 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00015 // GNU General Public License for more details. 00016 00017 // You should have received a copy of the GNU General Public License along 00018 // with this library; see the file COPYING. If not, write to the Free 00019 // Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, 00020 // USA. 00021 00022 // As a special exception, you may use this file as part of a free software 00023 // library without restriction. Specifically, if other files instantiate 00024 // templates or use macros or inline functions from this file, or you compile 00025 // this file and link it with other files to produce an executable, this 00026 // file does not by itself cause the resulting executable to be covered by 00027 // the GNU General Public License. This exception does not however 00028 // invalidate any other reasons why the executable file might be covered by 00029 // the GNU General Public License. 00030 00031 00032 #ifndef TOON_INCLUDE_HELPERS_H 00033 #define TOON_INCLUDE_HELPERS_H 00034 00035 #include <TooN/TooN.h> 00036 #include <TooN/gaussian_elimination.h> 00037 #include <cmath> 00038 #include <functional> 00039 #include <utility> 00040 00041 #ifndef M_PI 00042 #define M_PI 3.14159265358979323846 00043 #endif 00044 00045 #ifndef M_SQRT1_2 00046 #define M_SQRT1_2 0.707106781186547524401 00047 #endif 00048 00049 namespace TooN { 00050 00051 ///\deprecated 00052 ///@ingroup gLinAlg 00053 template<int Size, class Precision, class Base> TOON_DEPRECATED void Fill(Vector<Size, Precision, Base>& v, const Precision& p) 00054 { 00055 for(int i=0; i < v.size(); i++) 00056 v[i]= p; 00057 } 00058 00059 ///\deprecated 00060 ///@ingroup gLinAlg 00061 template<int Rows, int Cols, class Precision, class Base> TOON_DEPRECATED void Fill(Matrix<Rows, Cols, Precision, Base>& m, const Precision& p) 00062 { 00063 for(int i=0; i < m.num_rows(); i++) 00064 for(int j=0; j < m.num_cols(); j++) 00065 m[i][j] = p; 00066 } 00067 00068 ///Compute the \f$L_2\f$ norm of \e v 00069 ///@param v \e v 00070 ///@ingroup gLinAlg 00071 template<int Size, class Precision, class Base> inline Precision norm(const Vector<Size, Precision, Base>& v) 00072 { 00073 using std::sqrt; 00074 return sqrt(v*v); 00075 } 00076 00077 ///Compute the \f$L_2^2\f$ norm of \e v 00078 ///@param v \e v 00079 ///@ingroup gLinAlg 00080 template<int Size, class Precision, class Base> inline Precision norm_sq(const Vector<Size, Precision, Base>& v) 00081 { 00082 return v*v; 00083 } 00084 00085 ///Compute the \f$L_1\f$ norm of \e v 00086 ///@param v \e v 00087 ///@ingroup gLinAlg 00088 template<int Size, class Precision, class Base> inline Precision norm_1(const Vector<Size, Precision, Base>& v) 00089 { 00090 using std::abs; 00091 Precision n = 0; 00092 for(int i=0; i < v.size(); i++) 00093 n += abs(v[i]); 00094 return n; 00095 } 00096 00097 ///Compute the \f$L_\infty\f$ norm of \e v 00098 ///@param v \e v 00099 ///@ingroup gLinAlg 00100 template<int Size, class Precision, class Base> inline Precision norm_inf(const Vector<Size, Precision, Base>& v) 00101 { 00102 using std::abs; 00103 using std::max; 00104 Precision n = 0; 00105 n = abs(v[0]); 00106 00107 for(int i=1; i < v.size(); i++) 00108 n = max(n, abs(v[i])); 00109 return n; 00110 } 00111 00112 ///Compute the \f$L_2\f$ norm of \e v. 00113 ///Synonym for norm() 00114 ///@param v \e v 00115 ///@ingroup gLinAlg 00116 template<int Size, class Precision, class Base> inline Precision norm_2(const Vector<Size, Precision, Base>& v) 00117 { 00118 return norm(v); 00119 } 00120 00121 00122 00123 00124 ///Compute a the unit vector \f$\hat{v}\f$. 00125 ///@param v \e v 00126 ///@ingroup gLinAlg 00127 template<int Size, class Precision, class Base> inline Vector<Size, Precision> unit(const Vector<Size, Precision, Base> & v) 00128 { 00129 using std::sqrt; 00130 return TooN::operator*(v,(1/sqrt(v*v))); 00131 } 00132 00133 //Note because of the overload later, this function will ONLY receive sliced vectors. Therefore 00134 //a copy can be made, which is still a slice, so operating on the copy operates on the original 00135 //data. 00136 ///Normalize a vector in place 00137 ///@param v Vector to normalize 00138 ///@ingroup gLinAlg 00139 template<int Size, class Precision, class Base> inline void normalize(Vector<Size, Precision, Base> v) 00140 { 00141 using std::sqrt; 00142 v /= sqrt(v*v); 00143 } 00144 00145 //This overload is required to operate on non-slice vectors 00146 /** 00147 \overload 00148 */ 00149 template<int Size, class Precision> inline void normalize(Vector<Size, Precision> & v) 00150 { 00151 normalize(v.as_slice()); 00152 } 00153 00154 00155 ///For a vector \e v of length \e i, return \f$[v_1, v_2, \cdots, v_{i-1}] / v_i \f$ 00156 ///@param v \e v 00157 ///@ingroup gLinAlg 00158 template<int Size, typename Precision, typename Base> inline Vector<(Size==Dynamic?Dynamic:Size-1), Precision> project( const Vector<Size, Precision, Base> & v){ 00159 static const int Len = (Size==Dynamic?Dynamic:Size-1); 00160 return TooN::operator/(v.template slice<0, Len>(0, v.size()-1),v[v.size() - 1]); 00161 } 00162 00163 //This should probably be done with an operator to prevent an extra new[] for dynamic vectors. 00164 ///For a vector \e v of length \e i, return \f$[v_1, v_2, \cdots, v_{i}, 1]\f$ 00165 ///@param v \e v 00166 ///@ingroup gLinAlg 00167 template<int Size, typename Precision, typename Base> inline Vector<(Size==Dynamic?Dynamic:Size+1), Precision> unproject( const Vector<Size, Precision, Base> & v){ 00168 Vector<(Size==Dynamic?Dynamic:Size+1), Precision> result(v.size()+1); 00169 static const int Len = (Size==Dynamic?Dynamic:Size); 00170 result.template slice<0, Len>(0, v.size()) = v; 00171 result[v.size()] = 1; 00172 return result; 00173 } 00174 00175 /** 00176 \overload 00177 */ 00178 template<int R, int C, typename Precision, typename Base> inline Matrix<R-1, C, Precision> project( const Matrix<R,C, Precision, Base> & m){ 00179 Matrix<R-1, C, Precision> result = m.slice(0,0,R-1,m.num_cols()); 00180 for( int c = 0; c < m.num_cols(); ++c ) { 00181 result.slice(0,c,R-1,1) /= m[R-1][c]; 00182 } 00183 return result; 00184 } 00185 00186 template<int C, typename Precision, typename Base> inline Matrix<-1, C, Precision> project( const Matrix<-1,C, Precision, Base> & m){ 00187 Matrix<-1, C, Precision> result = m.slice(0,0,m.num_rows()-1,m.num_cols()); 00188 for( int c = 0; c < m.num_cols(); ++c ) { 00189 result.slice(0,c,m.num_rows()-1,1) /= m[m.num_rows()-1][c]; 00190 } 00191 return result; 00192 } 00193 00194 /** 00195 \overload 00196 */ 00197 template<int R, int C, typename Precision, typename Base> inline Matrix<R+1, C, Precision> unproject( const Matrix<R, C, Precision, Base> & m){ 00198 Matrix<R+1, C, Precision> result; 00199 result.template slice<0,0,R,C>() = m; 00200 result[R] = Ones; 00201 return result; 00202 } 00203 00204 template<int C, typename Precision, typename Base> inline Matrix<-1, C, Precision> unproject( const Matrix<-1, C, Precision, Base> & m){ 00205 Matrix<-1, C, Precision> result( m.num_rows()+1, m.num_cols() ); 00206 result.slice(0,0,m.num_rows(),m.num_cols()) = m; 00207 result[m.num_rows()] = Ones; 00208 return result; 00209 } 00210 00211 /// Frobenius (root of sum of squares) norm of input matrix \e m 00212 ///@param m \e m 00213 ///@ingroup gLinAlg 00214 template <int R, int C, typename P, typename B> 00215 P inline norm_fro( const Matrix<R,C,P,B> & m ){ 00216 using std::sqrt; 00217 P n = 0; 00218 for(int r = 0; r < m.num_rows(); ++r) 00219 for(int c = 0; c < m.num_cols(); ++c) 00220 n += m[r][c] * m[r][c]; 00221 00222 return sqrt(n); 00223 } 00224 00225 /// \e L<sub>∞</sub> (row sum) norm of input matrix m 00226 /// computes the maximum of the sums of absolute values over rows 00227 ///@ingroup gLinAlg 00228 template <int R, int C, typename P, typename B> 00229 P inline norm_inf( const Matrix<R,C,P,B> & m ){ 00230 using std::abs; 00231 using std::max; 00232 P n = 0; 00233 for(int r = 0; r < m.num_rows(); ++r){ 00234 P s = 0; 00235 for(int c = 0; c < m.num_cols(); ++c) 00236 s += abs(m(r,c)); 00237 n = max(n,s); 00238 } 00239 return n; 00240 } 00241 00242 /// \e L<sub>1</sub> (col sum) norm of input matrix m 00243 /// computes the maximum of the sums of absolute values over columns 00244 ///@ingroup gLinAlg 00245 template <int R, int C, typename P, typename B> 00246 P inline norm_1( const Matrix<R,C,P,B> & m ){ 00247 using std::abs; 00248 using std::max; 00249 P n = 0; 00250 for(int c = 0; c < m.num_cols(); ++c){ 00251 P s = 0; 00252 for(int r = 0; r < m.num_rows(); ++r) 00253 s += abs(m(r,c)); 00254 n = max(n,s); 00255 } 00256 return n; 00257 } 00258 00259 namespace Internal { 00260 ///@internal 00261 ///@brief Exponentiate a matrix using a the Taylor series 00262 ///This will not work if the norm of the matrix is too large. 00263 template <int R, int C, typename P, typename B> 00264 inline Matrix<R, C, P> exp_taylor( const Matrix<R,C,P,B> & m ){ 00265 TooN::SizeMismatch<R, C>::test(m.num_rows(), m.num_cols()); 00266 Matrix<R,C,P> result = TooN::Zeros(m.num_rows(), m.num_cols()); 00267 Matrix<R,C,P> f = TooN::Identity(m.num_rows()); 00268 P k = 1; 00269 while(norm_inf((result+f)-result) > 0){ 00270 result += f; 00271 f = (m * f) / k; 00272 k += 1; 00273 } 00274 return result; 00275 } 00276 00277 ///@internal 00278 ///@brief Logarithm of a matrix using a the Taylor series 00279 ///This will not work if the norm of the matrix is too large. 00280 template <int R, int C, typename P, typename B> 00281 inline Matrix<R, C, P> log_taylor( const Matrix<R,C,P,B> & m ){ 00282 TooN::SizeMismatch<R, C>::test(m.num_rows(), m.num_cols()); 00283 Matrix<R,C,P> X = m - TooN::Identity * 1.0; 00284 Matrix<R,C,P> F = X; 00285 Matrix<R,C,P> sum = TooN::Zeros(m.num_rows(), m.num_cols()); 00286 P k = 1; 00287 while(norm_inf((sum+F/k)-sum) > 0){ 00288 sum += F/k; 00289 F = -F*X; 00290 k += 1; 00291 } 00292 return sum; 00293 } 00294 00295 }; 00296 00297 /// computes the matrix exponential of a matrix m by 00298 /// scaling m by 1/(powers of 2), using Taylor series and 00299 /// squaring again. 00300 /// @param m input matrix, must be square 00301 /// @return result matrix of the same size/type as input 00302 /// @ingroup gLinAlg 00303 template <int R, int C, typename P, typename B> 00304 inline Matrix<R, C, P> exp( const Matrix<R,C,P,B> & m ){ 00305 using std::max; 00306 SizeMismatch<R, C>::test(m.num_rows(), m.num_cols()); 00307 const P l = log2(norm_inf(m)); 00308 const int s = max(0,(int)ceil(l)); 00309 Matrix<R,C,P> result = Internal::exp_taylor(m/(1<<s)); 00310 for(int i = 0; i < s; ++i) 00311 result = result * result; 00312 return result; 00313 } 00314 00315 /// computes a matrix square root of a matrix m by 00316 /// the product form of the Denman and Beavers iteration 00317 /// as given in Chen et al. 'Approximating the logarithm of a matrix to specified accuracy', 00318 /// J. Matrix Anal Appl, 2001. This is used for the matrix 00319 /// logarithm function, but is useable by on its own. 00320 /// @param m input matrix, must be square 00321 /// @return a square root of m of the same size/type as input 00322 /// @ingroup gLinAlg 00323 template <int R, int C, typename P, typename B> 00324 inline Matrix<R, C, P> sqrt( const Matrix<R,C,P,B> & m){ 00325 SizeMismatch<R, C>::test(m.num_rows(), m.num_cols()); 00326 Matrix<R,C,P> M = m; 00327 Matrix<R,C,P> Y = m; 00328 Matrix<R,C,P> M_inv(m.num_rows(), m.num_cols()); 00329 const Matrix<R,C,P> id = Identity(m.num_rows()); 00330 do { 00331 M_inv = gaussian_elimination(M, id); 00332 Y = Y * (id + M_inv) * 0.5; 00333 M = 0.5 * (id + (M + M_inv) * 0.5); 00334 } while(norm_inf(M - M_inv) > 0); 00335 return Y; 00336 } 00337 00338 /// computes the matrix logarithm of a matrix m using the inverse scaling and 00339 /// squaring method. The overall approach is described in 00340 /// Chen et al. 'Approximating the logarithm of a matrix to specified accuracy', 00341 /// J. Matrix Anal Appl, 2001, but this implementation only uses a simple 00342 /// taylor series after the repeated square root operation. 00343 /// @param m input matrix, must be square 00344 /// @return the log of m of the same size/type as input 00345 /// @ingroup gLinAlg 00346 template <int R, int C, typename P, typename B> 00347 inline Matrix<R, C, P> log( const Matrix<R,C,P,B> & m){ 00348 int counter = 0; 00349 Matrix<R,C,P> A = m; 00350 while(norm_inf(A - Identity*1.0) > 0.5){ 00351 ++counter; 00352 A = sqrt(A); 00353 } 00354 return Internal::log_taylor(A) * pow(2.0, counter); 00355 } 00356 00357 /// Returns true if every element is finite 00358 ///@ingroup gLinAlg 00359 template<int S, class P, class B> bool isfinite(const Vector<S, P, B>& v) 00360 { 00361 using std::isfinite; 00362 for(int i=0; i < v.size(); i++) 00363 if(!isfinite(v[i])) 00364 return 0; 00365 return 1; 00366 } 00367 00368 /// Returns true if any element is NaN 00369 ///@ingroup gLinAlg 00370 template<int S, class P, class B> bool isnan(const Vector<S, P, B>& v) 00371 { 00372 using std::isnan; 00373 for(int i=0; i < v.size(); i++) 00374 if(isnan(v[i])) 00375 return 1; 00376 return 0; 00377 } 00378 00379 /// Symmetrize a matrix 00380 ///@param m \e m 00381 ///@return \f$ \frac{m + m^{\mathsf T}}{2} \f$ 00382 ///@ingroup gLinAlg 00383 template<int Rows, int Cols, typename Precision, typename Base> 00384 void Symmetrize(Matrix<Rows,Cols,Precision,Base>& m){ 00385 SizeMismatch<Rows,Cols>::test(m.num_rows(), m.num_cols()); 00386 for(int r=0; r<m.num_rows()-1; r++){ 00387 for(int c=r+1; c<m.num_cols(); c++){ 00388 const Precision temp=(m(r,c)+m(c,r))/2; 00389 m(r,c)=temp; 00390 m(c,r)=temp; 00391 } 00392 } 00393 } 00394 00395 /// computes the trace of a square matrix 00396 ///@ingroup gLinAlg 00397 template<int Rows, int Cols, typename Precision, typename Base> 00398 Precision trace(const Matrix<Rows, Cols, Precision, Base> & m ){ 00399 SizeMismatch<Rows,Cols>::test(m.num_rows(), m.num_cols()); 00400 Precision tr = 0; 00401 for(int i = 0; i < m.num_rows(); ++i) 00402 tr += m(i,i); 00403 return tr; 00404 } 00405 00406 /// creates an returns a cross product matrix M from a 3 vector v, such that for all vectors w, the following holds: v ^ w = M * w 00407 /// @param vec the 3 vector input 00408 /// @return the 3x3 matrix to set to the cross product matrix 00409 ///@ingroup gLinAlg 00410 template<int Size, class P, class B> inline TooN::Matrix<3, 3, P> cross_product_matrix(const Vector<Size, P, B>& vec) 00411 { 00412 SizeMismatch<Size,3>::test(vec.size(), 3); 00413 00414 TooN::Matrix<3> result; 00415 00416 result(0,0) = 0; 00417 result(0,1) = -vec[2]; 00418 result(0,2) = vec[1]; 00419 result(1,0) = vec[2]; 00420 result(1,1) = 0; 00421 result(1,2) = -vec[0]; 00422 result(2,0) = -vec[1]; 00423 result(2,1) = vec[0]; 00424 result(2,2) = 0; 00425 00426 return result; 00427 } 00428 00429 namespace Internal { 00430 template<int Size, typename Precision, typename Base, typename Func, typename Ret> inline Ret accumulate( const Vector<Size, Precision, Base> & v ) { 00431 Func func; 00432 if( v.size() == 0 ) { 00433 return func.null(); // What should we return, exception? 00434 } 00435 func.initialise( v[0], 0 ); 00436 for( int ii = 1; ii < v.size(); ii++ ) { 00437 func( v[ii], ii ); 00438 } 00439 return func.ret(); 00440 } 00441 00442 template<int R, int C, typename Precision, typename Base, typename Func, typename Ret> inline Ret accumulate( const Matrix<R, C, Precision, Base> & m ) { 00443 Func func; 00444 if( m.num_rows() == 0 || m.num_cols() == 0) { 00445 return func.null(); // What should we return, exception? 00446 } 00447 func.initialise( m[0][0], 0, 0 ); 00448 for(int r=0; r<m.num_rows(); r++){ 00449 for(int c=0; c<m.num_cols(); c++){ 00450 func( m[r][c], r, c ); 00451 } 00452 } 00453 return func.ret(); 00454 } 00455 template<int R, int C, typename Precision, typename Base, typename Func, typename Ret> inline Ret accumulate_horizontal( const Matrix<R, C, Precision, Base> & m ) { 00456 Func func( m.num_rows() ); 00457 if( m.num_cols() == 0 || m.num_rows() == 0 ) { 00458 func.null(); // What should we return, exception? 00459 } 00460 for(int r=0; r<m.num_rows(); r++){ 00461 func.initialise( m[r][0], r, 0 ); 00462 for(int c=1; c<m.num_cols(); c++){ 00463 func( m[r][c], r, c ); 00464 } 00465 } 00466 return func.ret(); 00467 } 00468 template<int R, int C, typename Precision, typename Base, typename Func, typename Ret> inline Ret accumulate_vertical( const Matrix<R, C, Precision, Base> & m ) { 00469 Func func( m.num_cols() ); 00470 if( m.num_cols() == 0 || m.num_rows() == 0 ) { 00471 func.null(); // What should we return, exception? 00472 } 00473 for(int c=0; c<m.num_cols(); c++){ 00474 func.initialise( m[0][c], 0, c ); 00475 for(int r=1; r<m.num_rows(); r++){ 00476 func( m[r][c], r, c ); 00477 } 00478 } 00479 return func.ret(); 00480 } 00481 00482 template<typename Precision, typename ComparisonFunctor> 00483 class accumulate_functor_vector { 00484 Precision bestVal; 00485 public: 00486 Precision null() { 00487 return 0; 00488 } 00489 void initialise( Precision initialVal, int ) { 00490 bestVal = initialVal; 00491 } 00492 void operator()( Precision curVal, int ) { 00493 if( ComparisonFunctor()( curVal, bestVal ) ) { 00494 bestVal = curVal; 00495 } 00496 } 00497 Precision ret() { return bestVal; } 00498 }; 00499 template<typename Precision, typename ComparisonFunctor> 00500 class accumulate_element_functor_vector { 00501 Precision bestVal; 00502 int nBestIndex; 00503 public: 00504 std::pair<Precision,int> null() { 00505 return std::pair<Precision,int>( 0, 0 ); 00506 } 00507 void initialise( Precision initialVal, int nIndex ) { 00508 bestVal = initialVal; 00509 nBestIndex = nIndex; 00510 } 00511 void operator()( Precision curVal, int nIndex ) { 00512 if( ComparisonFunctor()( curVal, bestVal ) ) { 00513 bestVal = curVal; 00514 nBestIndex = nIndex; 00515 } 00516 } 00517 std::pair<Precision,int> ret() { 00518 return std::pair<Precision,int>( bestVal, nBestIndex ); 00519 } 00520 }; 00521 template<typename Precision, typename ComparisonFunctor> 00522 class accumulate_functor_matrix { 00523 Precision bestVal; 00524 public: 00525 Precision null() { 00526 return 0; 00527 } 00528 void initialise( Precision initialVal, int, int ) { 00529 bestVal = initialVal; 00530 } 00531 void operator()( Precision curVal, int, int ) { 00532 if( ComparisonFunctor()( curVal, bestVal ) ) { 00533 bestVal = curVal; 00534 } 00535 } 00536 Precision ret() { return bestVal; } 00537 }; 00538 template<typename Precision, typename ComparisonFunctor> 00539 class accumulate_element_functor_matrix { 00540 Precision bestVal; 00541 int nBestRow; 00542 int nBestCol; 00543 public: 00544 std::pair<Precision,std::pair<int,int> > null() { 00545 return std::pair<Precision,std::pair<int,int> >( 0, std::pair<int,int>( 0, 0 ) ); 00546 } 00547 void initialise( Precision initialVal, int nRow, int nCol ) { 00548 bestVal = initialVal; 00549 nBestRow = nRow; 00550 nBestCol = nCol; 00551 } 00552 void operator()( Precision curVal, int nRow, int nCol ) { 00553 if( ComparisonFunctor()( curVal, bestVal ) ) { 00554 bestVal = curVal; 00555 nBestRow = nRow; 00556 nBestCol = nCol; 00557 } 00558 } 00559 std::pair<Precision,std::pair<int,int> > ret() { 00560 return std::pair<Precision,std::pair<int,int> >( bestVal, 00561 std::pair<int,int>( nBestRow, nBestCol ) ); 00562 } 00563 }; 00564 template<typename Precision, typename ComparisonFunctor> 00565 class accumulate_vertical_functor { 00566 Vector<Dynamic,Precision>* bestVal; 00567 public: 00568 accumulate_vertical_functor() { 00569 bestVal = NULL; 00570 } 00571 accumulate_vertical_functor( int nNumCols ) { 00572 bestVal = new Vector<Dynamic,Precision>( nNumCols ); 00573 } 00574 Vector<Dynamic,Precision> null() { 00575 return Vector<Dynamic,Precision>( 0 ); 00576 } 00577 void initialise( Precision initialVal, int, int nCol ) { 00578 (*bestVal)[nCol] = initialVal; 00579 } 00580 void operator()( Precision curVal, int, int nCol ) { 00581 if( ComparisonFunctor()( curVal, (*bestVal)[nCol] ) ) { 00582 (*bestVal)[nCol] = curVal; 00583 } 00584 } 00585 Vector<Dynamic,Precision> ret() { 00586 if( bestVal == NULL ) { 00587 return null(); 00588 } 00589 Vector<Dynamic,Precision> vRet = *bestVal; 00590 delete bestVal; 00591 return vRet; 00592 } 00593 }; 00594 template<typename Precision, typename ComparisonFunctor> 00595 class accumulate_element_vertical_functor { 00596 Vector<Dynamic,Precision>* bestVal; 00597 Vector<Dynamic,Precision>* bestIndices; 00598 public: 00599 accumulate_element_vertical_functor() { 00600 bestVal = NULL; 00601 bestIndices = NULL; 00602 } 00603 accumulate_element_vertical_functor( int nNumCols ) { 00604 bestVal = new Vector<Dynamic,Precision>( nNumCols ); 00605 bestIndices = new Vector<Dynamic,Precision>( nNumCols ); 00606 } 00607 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > null() { 00608 Vector<Dynamic,Precision> vEmpty( 0 ); 00609 return std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> >( vEmpty, vEmpty ); 00610 } 00611 void initialise( Precision initialVal, int nRow, int nCol ) { 00612 (*bestVal)[nCol] = initialVal; 00613 (*bestIndices)[nCol] = nRow; 00614 } 00615 void operator()( Precision curVal, int nRow, int nCol ) { 00616 if( ComparisonFunctor()( curVal, (*bestVal)[nCol] ) ) { 00617 (*bestVal)[nCol] = curVal; 00618 (*bestIndices)[nCol] = nRow; 00619 } 00620 } 00621 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > ret() { 00622 if( bestVal == NULL ) { 00623 return null(); 00624 } 00625 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > vRet = 00626 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > (*bestVal, *bestIndices ); 00627 delete bestVal; bestVal = NULL; 00628 delete bestIndices; bestIndices = NULL; 00629 return vRet; 00630 } 00631 }; 00632 template<typename Precision, typename ComparisonFunctor> 00633 class accumulate_horizontal_functor { 00634 Vector<Dynamic,Precision>* bestVal; 00635 public: 00636 accumulate_horizontal_functor() { 00637 bestVal = NULL; 00638 } 00639 accumulate_horizontal_functor( int nNumRows ) { 00640 bestVal = new Vector<Dynamic,Precision>( nNumRows ); 00641 } 00642 Vector<Dynamic,Precision> null() { 00643 return Vector<Dynamic,Precision>( 0 ); 00644 } 00645 void initialise( Precision initialVal, int nRow, int ) { 00646 (*bestVal)[nRow] = initialVal; 00647 } 00648 void operator()( Precision curVal, int nRow, int ) { 00649 if( ComparisonFunctor()( curVal, (*bestVal)[nRow] ) ) { 00650 (*bestVal)[nRow] = curVal; 00651 } 00652 } 00653 Vector<Dynamic,Precision> ret() { 00654 if( bestVal == NULL ) { 00655 return null(); 00656 } 00657 Vector<Dynamic,Precision> vRet = *bestVal; 00658 delete bestVal; bestVal = NULL; 00659 return vRet; 00660 } 00661 }; 00662 template<typename Precision, typename ComparisonFunctor> 00663 class accumulate_element_horizontal_functor { 00664 Vector<Dynamic,Precision>* bestVal; 00665 Vector<Dynamic,Precision>* bestIndices; 00666 public: 00667 accumulate_element_horizontal_functor() { 00668 bestVal = NULL; 00669 bestIndices = NULL; 00670 } 00671 accumulate_element_horizontal_functor( int nNumRows ) { 00672 bestVal = new Vector<Dynamic,Precision>( nNumRows ); 00673 bestIndices = new Vector<Dynamic,Precision>( nNumRows ); 00674 } 00675 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > null() { 00676 Vector<Dynamic,Precision> vEmpty( 0 ); 00677 return std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> >( vEmpty, vEmpty ); 00678 } 00679 void initialise( Precision initialVal, int nRow, int nCol ) { 00680 (*bestVal)[nRow] = initialVal; 00681 (*bestIndices)[nRow] = nCol; 00682 } 00683 void operator()( Precision curVal, int nRow, int nCol ) { 00684 if( ComparisonFunctor()( curVal, (*bestVal)[nRow] ) ) { 00685 (*bestVal)[nRow] = curVal; 00686 (*bestIndices)[nRow] = nCol; 00687 } 00688 } 00689 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > ret() { 00690 if( bestVal == NULL ) { 00691 return null(); 00692 } 00693 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > vRet = 00694 std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> >( *bestVal, *bestIndices ); 00695 delete bestVal; bestVal = NULL; 00696 delete bestIndices; bestIndices = NULL; 00697 return vRet; 00698 } 00699 }; 00700 } 00701 00702 00703 /// Finds the minimal value of a vector. 00704 /// @param v a vector 00705 /// @return the smallest value of v 00706 template<int Size, typename Precision, typename Base> inline Precision min_value( const Vector<Size, Precision, Base>& v) { 00707 typedef Internal::accumulate_functor_vector<Precision, std::less<Precision> > vector_accumulate_functor; 00708 return Internal::accumulate<Size,Precision,Base, 00709 vector_accumulate_functor, Precision >( v ); 00710 } 00711 /// Finds the largest value of a vector. 00712 /// @param v a vector 00713 /// @return the largest value of v 00714 template<int Size, typename Precision, typename Base> inline Precision max_value( const Vector<Size, Precision, Base>& v) { 00715 typedef Internal::accumulate_functor_vector<Precision, std::greater<Precision> > vector_accumulate_functor; 00716 return Internal::accumulate<Size,Precision,Base, 00717 vector_accumulate_functor, Precision >( v ); 00718 } 00719 00720 /// Finds the smallest value of a matrix. 00721 /// @param m a matrix 00722 /// @return the smallest value of m 00723 template<int R, int C, typename Precision, typename Base> inline Precision min_value( const Matrix<R, C, Precision, Base> & m) { 00724 typedef Internal::accumulate_functor_matrix<Precision, std::less<Precision> > matrix_accumulate_functor; 00725 return Internal::accumulate<R,C,Precision,Base, 00726 matrix_accumulate_functor, Precision>( m ); 00727 } 00728 /// Finds the largest value of a matrix. 00729 /// @param m a matrix 00730 /// @return the largest value of m 00731 template<int R, int C, typename Precision, typename Base> inline Precision max_value( const Matrix<R, C, Precision, Base> & m) { 00732 typedef Internal::accumulate_functor_matrix<Precision, std::greater<Precision> > matrix_accumulate_functor; 00733 return Internal::accumulate<R,C,Precision,Base, 00734 matrix_accumulate_functor, Precision>( m ); 00735 } 00736 /// Finds the smallest values of each column of a matrix. 00737 /// @param m a matrix 00738 /// @return a vector of size C 00739 template<int R, int C, typename Precision, typename Base> inline Vector<Dynamic,Precision> min_value_vertical( const Matrix<R, C, Precision, Base> & m) { 00740 typedef Internal::accumulate_vertical_functor<Precision,std::less<Precision> > matrix_accumulate_vertical_functor; 00741 return Internal::accumulate_vertical<R,C,Precision,Base, 00742 matrix_accumulate_vertical_functor, Vector<Dynamic,Precision> >( m ); 00743 } 00744 /// Finds the largest values of each column of a matrix. 00745 /// @param m a matrix 00746 /// @return a vector of size C 00747 template<int R, int C, typename Precision, typename Base> inline Vector<Dynamic,Precision> max_value_vertical( const Matrix<R, C, Precision, Base> & m) { 00748 typedef Internal::accumulate_vertical_functor<Precision,std::greater<Precision> > matrix_accumulate_vertical_functor; 00749 return Internal::accumulate_vertical<R,C,Precision,Base, 00750 matrix_accumulate_vertical_functor, Vector<Dynamic,Precision> >( m ); 00751 } 00752 /// Finds the smallest values of each row of a matrix. 00753 /// @param m a matrix 00754 /// @return a vector of size R 00755 template<int R, int C, typename Precision, typename Base> inline Vector<Dynamic,Precision> min_value_horizontal( const Matrix<R, C, Precision, Base> & m) { 00756 typedef Internal::accumulate_horizontal_functor<Precision,std::less<Precision> > matrix_accumulate_horizontal_functor; 00757 return Internal::accumulate_horizontal<R,C,Precision,Base, 00758 matrix_accumulate_horizontal_functor, Vector<Dynamic,Precision> >( m ); 00759 } 00760 /// Finds the largest values of each row of a matrix. 00761 /// @param m a matrix 00762 /// @return a vector of size R 00763 template<int R, int C, typename Precision, typename Base> inline Vector<Dynamic,Precision> max_value_horizontal( const Matrix<R, C, Precision, Base> & m) { 00764 typedef Internal::accumulate_horizontal_functor<Precision,std::greater<Precision> > matrix_accumulate_horizontal_functor; 00765 return Internal::accumulate_horizontal<R,C,Precision,Base, 00766 matrix_accumulate_horizontal_functor, Vector<Dynamic,Precision> >( m ); 00767 } 00768 /// Finds the smallest value of a vector and its index. 00769 /// @param v a vector 00770 /// @return a pair containing the smallest value and its index 00771 template<int Size, typename Precision, typename Base> inline std::pair<Precision,int> min_element( const Vector<Size, Precision, Base>& v) { 00772 typedef Internal::accumulate_element_functor_vector<Precision, std::less<Precision> > vector_accumulate_functor; 00773 return Internal::accumulate<Size,Precision,Base, 00774 vector_accumulate_functor, std::pair<Precision,int> >( v ); 00775 } 00776 /// Finds the largest value of a vector and its index. 00777 /// @param v a vector 00778 /// @return a pair containing the largest value and its index 00779 template<int Size, typename Precision, typename Base> inline std::pair<Precision,int> max_element( const Vector<Size, Precision, Base>& v) { 00780 typedef Internal::accumulate_element_functor_vector<Precision, std::greater<Precision> > vector_accumulate_functor; 00781 return Internal::accumulate<Size,Precision,Base, 00782 vector_accumulate_functor, std::pair<Precision,int> >( v ); 00783 } 00784 /// Finds the smallest value of a matrix and its row and column. 00785 /// @param m a matrix 00786 /// @return a pair containing the smallest value and a pair 00787 /// containing its row and column 00788 template<int R, int C, typename Precision, typename Base> inline std::pair<Precision,std::pair<int,int> > min_element( const Matrix<R, C, Precision, Base> & m) { 00789 typedef Internal::accumulate_element_functor_matrix<Precision, std::less<Precision> > matrix_accumulate_functor; 00790 typedef std::pair<Precision,std::pair<int,int> > Ret; 00791 return Internal::accumulate<R,C,Precision,Base, 00792 matrix_accumulate_functor, Ret>( m ); 00793 } 00794 /// Finds the largest value of a matrix and its row and column. 00795 /// @param m a matrix 00796 /// @return a pair containing the largest value and a pair 00797 /// containing its row and column 00798 template<int R, int C, typename Precision, typename Base> inline std::pair<Precision,std::pair<int,int> > max_element( const Matrix<R, C, Precision, Base> & m) { 00799 typedef Internal::accumulate_element_functor_matrix<Precision, std::greater<Precision> > matrix_accumulate_functor; 00800 typedef std::pair<Precision,std::pair<int,int> > Ret; 00801 return Internal::accumulate<R,C,Precision,Base, 00802 matrix_accumulate_functor, Ret>( m ); 00803 } 00804 /// Finds the smallest values of each column of a matrix and their 00805 /// indices. 00806 /// @param m a matrix 00807 /// @return a pair of vectors of size C containg the values and 00808 /// their indices 00809 template<int R, int C, typename Precision, typename Base> inline std::pair<Vector<Dynamic,Precision>,Vector<Dynamic,Precision> > min_element_vertical( const Matrix<R, C, Precision, Base> & m) { 00810 typedef Internal::accumulate_element_vertical_functor<Precision,std::less<Precision> > matrix_accumulate_vertical_functor; 00811 typedef std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > Ret; 00812 return Internal::accumulate_vertical<R,C,Precision,Base, 00813 matrix_accumulate_vertical_functor, Ret >( m ); 00814 } 00815 /// Finds the largest values of each column of a matrix and their 00816 /// indices. 00817 /// @param m a matrix 00818 /// @return a pair of vectors of size C containg the values and 00819 /// their indices 00820 template<int R, int C, typename Precision, typename Base> inline std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > max_element_vertical( const Matrix<R, C, Precision, Base> & m) { 00821 typedef Internal::accumulate_element_vertical_functor<Precision,std::greater<Precision> > matrix_accumulate_vertical_functor; 00822 typedef std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > Ret; 00823 return Internal::accumulate_vertical<R,C,Precision,Base, 00824 matrix_accumulate_vertical_functor, Ret >( m ); 00825 } 00826 /// Finds the smallest values of each row of a matrix and their 00827 /// indices. 00828 /// @param m a matrix 00829 /// @return a pair of vectors of size R containg the values and 00830 /// their indices 00831 template<int R, int C, typename Precision, typename Base> inline std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > min_element_horizontal( const Matrix<R, C, Precision, Base> & m) { 00832 typedef Internal::accumulate_element_horizontal_functor<Precision,std::less<Precision> > matrix_accumulate_vertical_functor; 00833 typedef std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > Ret; 00834 return Internal::accumulate_horizontal<R,C,Precision,Base, 00835 matrix_accumulate_vertical_functor, Ret >( m ); 00836 } 00837 /// Finds the largest values of each row of a matrix and their 00838 /// indices. 00839 /// @param m a matrix 00840 /// @return a pair of vectors of size R containg the values and 00841 /// their indices 00842 template<int R, int C, typename Precision, typename Base> inline std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > max_element_horizontal( const Matrix<R, C, Precision, Base> & m) { 00843 typedef Internal::accumulate_element_horizontal_functor<Precision,std::greater<Precision> > matrix_accumulate_vertical_functor; 00844 typedef std::pair<Vector< Dynamic, Precision >,Vector< Dynamic, Precision > > Ret; 00845 return Internal::accumulate_horizontal<R,C,Precision,Base, 00846 matrix_accumulate_vertical_functor, Ret >( m ); 00847 } 00848 } 00849 #endif