ViSP  2.8.0
vpMatrix.cpp
1 /****************************************************************************
2 *
3 * $Id: vpMatrix.cpp 4210 2013-04-16 08:57:46Z fspindle $
4 *
5 * This file is part of the ViSP software.
6 * Copyright (C) 2005 - 2013 by INRIA. All rights reserved.
7 *
8 * This software is free software; you can redistribute it and/or
9 * modify it under the terms of the GNU General Public License
10 * ("GPL") version 2 as published by the Free Software Foundation.
11 * See the file LICENSE.txt at the root directory of this source
12 * distribution for additional information about the GNU GPL.
13 *
14 * For using ViSP with software that can not be combined with the GNU
15 * GPL, please contact INRIA about acquiring a ViSP Professional
16 * Edition License.
17 *
18 * See http://www.irisa.fr/lagadic/visp/visp.html for more information.
19 *
20 * This software was developed at:
21 * INRIA Rennes - Bretagne Atlantique
22 * Campus Universitaire de Beaulieu
23 * 35042 Rennes Cedex
24 * France
25 * http://www.irisa.fr/lagadic
26 *
27 * If you have questions regarding the use of this file, please contact
28 * INRIA at visp@inria.fr
29 *
30 * This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
31 * WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
32 *
33 *
34 * Description:
35 * Matrix manipulation.
36 *
37 * Authors:
38 * Eric Marchand
39 *
40 *****************************************************************************/
41 
42 
43 
49 /*
50 \class vpMatrix
51 
52 \brief Provide simple Matrices computation
53 
54 \author Eric Marchand (Eric.Marchand@irisa.fr) Irisa / Inria Rennes
55 */
56 
57 #include <visp/vpMatrix.h>
58 #include <visp/vpMath.h>
59 #include <visp/vpTranslationVector.h>
60 
61 // Exception
62 #include <visp/vpException.h>
63 #include <visp/vpMatrixException.h>
64 
65 // Debug trace
66 #include <visp/vpDebug.h>
67 
68 #include <stdlib.h>
69 #include <stdio.h>
70 #include <string.h>
71 #include <vector>
72 #include <sstream>
73 #include <algorithm>
74 #include <assert.h>
75 #include <fstream>
76 #include <string>
77 #include <cmath> // std::fabs
78 #include <limits> // numeric_limits
79 
80 #define DEBUG_LEVEL1 0
81 #define DEBUG_LEVEL2 0
82 
83 
84 //Prototypes of specific functions
85 vpMatrix subblock(const vpMatrix &, unsigned int, unsigned int);
86 
92 void
94 {
95 
96  rowNum = 0 ;
97  colNum = 0 ;
98 
99  data = NULL ;
100  rowPtrs = NULL ;
101 
102  dsize = 0 ;
103  trsize =0 ;
104 }
105 
112 {
113  init() ;
114 }
115 
116 
125 vpMatrix::vpMatrix(unsigned int r, unsigned int c)
126 {
127  init() ;
128  resize(r, c);
129 }
130 
135  unsigned int r, unsigned int c,
136  unsigned int nrows, unsigned int ncols)
137 {
138  init() ;
139 
140  if (((r + nrows) > m.rowNum) || ((c + ncols) > m.colNum))
141  {
142  vpERROR_TRACE("\n\t\t SubvpMatrix larger than vpMatrix") ;
144  "\n\t\t SubvpMatrix larger than vpMatrix")) ;
145  }
146 
147  init(m,r,c,nrows,ncols);
148 }
149 
152 {
153  init() ;
154 
155  resize(m.rowNum,m.colNum);
156 
157  memcpy(data,m.data,rowNum*colNum*sizeof(double)) ;
158 }
159 
160 
174 void vpMatrix::resize(const unsigned int nrows, const unsigned int ncols,
175  const bool flagNullify)
176 {
177 
178  if ((nrows == rowNum) && (ncols == colNum))
179  {
180  if (flagNullify)
181  { memset(this->data,0,this->dsize*sizeof(double)) ;}
182  }
183  else
184  {
185  const bool recopyNeeded = (ncols != this ->colNum);
186  double * copyTmp = NULL;
187  unsigned int rowTmp = 0, colTmp=0;
188 
189  vpDEBUG_TRACE (25, "Recopy case per case is required iff number of "
190  "cols has changed (structure of double array is not "
191  "the same in this case.");
192  if (recopyNeeded)
193  {
194  copyTmp = new double[this->dsize];
195  memcpy (copyTmp, this ->data, sizeof(double)*this->dsize);
196  rowTmp=this->rowNum; colTmp=this->colNum;
197  }
198 
199  vpDEBUG_TRACE (25, "Reallocation of this->data array.");
200  this->dsize = nrows*ncols;
201  this->data = (double*)realloc(this->data, this->dsize*sizeof(double));
202  if ((NULL == this->data) && (0 != this->dsize))
203  {
204  if (copyTmp != NULL) delete [] copyTmp;
205  vpERROR_TRACE("\n\t\tMemory allocation error when allocating data") ;
207  "\n\t\t Memory allocation error when "
208  "allocating data")) ;
209  }
210 
211  vpDEBUG_TRACE (25, "Reallocation of this->trsize array.");
212  this->trsize = nrows;
213  this->rowPtrs = (double**)realloc (this->rowPtrs, this->trsize*sizeof(double*));
214  if ((NULL == this->rowPtrs) && (0 != this->dsize))
215  {
216  if (copyTmp != NULL) delete [] copyTmp;
217  vpERROR_TRACE("\n\t\tMemory allocation error when allocating rowPtrs") ;
219  "\n\t\t Memory allocation error when "
220  "allocating rowPtrs")) ;
221  }
222 
223  vpDEBUG_TRACE (25, "Recomputation this->trsize array values.");
224  {
225  double **t= rowPtrs;
226  for (unsigned int i=0; i<dsize; i+=ncols) { *t++ = this->data + i; }
227  }
228 
229  this->rowNum = nrows; this->colNum = ncols;
230 
231  vpDEBUG_TRACE (25, "Recopy of this->data array values or nullify.");
232  if (flagNullify)
233  { memset(this->data,0,this->dsize*sizeof(double)) ;}
234  else
235  {
236  if (recopyNeeded)
237  {
238  vpDEBUG_TRACE (25, "Recopy...");
239  const unsigned int minRow = (this->rowNum<rowTmp)?this->rowNum:rowTmp;
240  const unsigned int minCol = (this->colNum<colTmp)?this->colNum:colTmp;
241  for (unsigned int i=0; i<this->rowNum; ++i)
242  for (unsigned int j=0; j<this->colNum; ++j)
243  {
244  if ((minRow > i) && (minCol > j))
245  {
246  (*this)[i][j] = copyTmp [i*colTmp+j];
247  vpCDEBUG (25) << i << "x" << j << "<- " << i*colTmp+j
248  << "=" << copyTmp [i*colTmp+j] << std::endl;
249  }
250  else {(*this)[i][j] = 0;}
251  }
252  }
253  else { vpDEBUG_TRACE (25,"Nothing to do: already done by realloc.");}
254  }
255 
256  if (copyTmp != NULL) delete [] copyTmp;
257  }
258 
259 }
260 
261 
262 void
263 vpMatrix::init(const vpMatrix &m,unsigned int r, unsigned int c, unsigned int nrows, unsigned int ncols)
264 {
265  try {
266  resize(nrows, ncols) ;
267  }
268  catch(vpException me)
269  {
270  vpERROR_TRACE("Error caught") ;
271  std::cout << me << std::endl ;
272  throw ;
273  }
274 
275  unsigned int rnrows = r+nrows ;
276  unsigned int cncols = c+ncols ;
277  for (unsigned int i=r ; i < rnrows; i++)
278  for (unsigned int j=c ; j < cncols; j++)
279  (*this)[i-r][j-c] = m[i][j] ;
280 }
281 
285 void
287 {
288  if (data != NULL )
289  {
290  free(data);
291  data=NULL;
292  }
293 
294  if (rowPtrs!=NULL)
295  {
296  free(rowPtrs);
297  rowPtrs=NULL ;
298  }
299 }
304 {
305  kill() ;
306 }
307 
308 
315 vpMatrix &
317 {
318  try {
319  resize(B.rowNum, B.colNum) ;
320  // suppress by em 5/12/06
321  // *this = 0;
322  }
323  catch(vpException me)
324  {
325  vpERROR_TRACE("Error caught") ;
326  std::cout << me << std::endl ;
327  throw ;
328  }
329 
330  memcpy(data,B.data,dsize*sizeof(double)) ;
331 
332  return *this;
333 }
334 
336 vpMatrix &
338 {
339 
340  // double t0,t1;
341  //
342  // t0=vpTime::measureTimeMicros();
343  // for (int i=0;i<dsize;i++)
344  // {
345  // *(data+i) = x;
346  // }
347  // t1=vpTime::measureTimeMicros();
348  // std::cout<< t1-t0<< std::endl;
349 
350  // t0=vpTime::measureTimeMicros();
351  for (unsigned int i=0;i<rowNum;i++)
352  for(unsigned int j=0;j<colNum;j++)
353  rowPtrs[i][j] = x;
354 
355  // t1=vpTime::measureTimeMicros();
356  // std::cout<< t1-t0<<" ";
357 
358 
359 
360  return *this;
361 }
362 
363 
367 vpMatrix &
369 {
370 
371  for (unsigned int i=0; i<rowNum; i++) {
372  for (unsigned int j=0; j<colNum; j++) {
373  rowPtrs[i][j] = *x++;
374  }
375  }
376  return *this;
377 }
378 
379 
380 //---------------------------------
381 // Matrix operations.
382 //---------------------------------
383 
393 void vpMatrix::mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
394 {
395  try
396  {
397  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) C.resize(A.rowNum,B.colNum);
398  }
399  catch(vpException me)
400  {
401  vpERROR_TRACE("Error caught") ;
402  std::cout << me << std::endl ;
403  throw ;
404  }
405 
406  if (A.colNum != B.rowNum)
407  {
408  vpERROR_TRACE("\n\t\tvpMatrix mismatch in vpMatrix/vpMatrix multiply") ;
410  "\n\t\tvpMatrix mismatch in "
411  "vpMatrix/vpMatrix multiply")) ;
412  }
413 
414  // 5/12/06 some "very" simple optimization to avoid indexation
415  unsigned int BcolNum = B.colNum;
416  unsigned int BrowNum = B.rowNum;
417  unsigned int i,j,k;
418  double **BrowPtrs = B.rowPtrs;
419  for (i=0;i<A.rowNum;i++)
420  {
421  double *rowptri = A.rowPtrs[i];
422  double *ci = C[i];
423  for (j=0;j<BcolNum;j++)
424  {
425  double s = 0;
426  for (k=0;k<BrowNum;k++) s += rowptri[k] * BrowPtrs[k][j];
427  ci[j] = s;
428  }
429  }
430 }
431 
437 {
438  vpMatrix C;
439 
440  vpMatrix::mult2Matrices(*this,B,C);
441 
442  return C;
443 }
454 void vpMatrix::add2WeightedMatrices(const vpMatrix &A, const double &wA, const vpMatrix &B,const double &wB, vpMatrix &C){
455  try
456  {
457  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) C.resize(A.rowNum,B.colNum);
458  }
459  catch(vpException me)
460  {
461  vpERROR_TRACE("Error caught") ;
462  std::cout << me << std::endl ;
463  throw ;
464  }
465 
466  if ((A.colNum != B.getCols())||(A.rowNum != B.getRows()))
467  {
468  vpERROR_TRACE("\n\t\t vpMatrix mismatch in vpMatrix/vpMatrix addition") ;
470  "\n\t\t vpMatrix mismatch in "
471  "vpMatrix/vpMatrix addition")) ;
472  }
473 
474  double ** ArowPtrs=A.rowPtrs;
475  double ** BrowPtrs=B.rowPtrs;
476  double ** CrowPtrs=C.rowPtrs;
477 
478  for (unsigned int i=0;i<A.rowNum;i++)
479  for(unsigned int j=0;j<A.colNum;j++)
480  CrowPtrs[i][j] = wB*BrowPtrs[i][j]+wA*ArowPtrs[i][j];
481 
482 }
483 
493 void vpMatrix::add2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
494 {
495  try
496  {
497  if ((A.rowNum != C.rowNum) || (B.colNum != C.colNum)) C.resize(A.rowNum,B.colNum);
498  }
499  catch(vpException me)
500  {
501  vpERROR_TRACE("Error caught") ;
502  std::cout << me << std::endl ;
503  throw ;
504  }
505 
506  if ((A.colNum != B.getCols())||(A.rowNum != B.getRows()))
507  {
508  vpERROR_TRACE("\n\t\t vpMatrix mismatch in vpMatrix/vpMatrix addition") ;
510  "\n\t\t vpMatrix mismatch in "
511  "vpMatrix/vpMatrix addition")) ;
512  }
513 
514  double ** ArowPtrs=A.rowPtrs;
515  double ** BrowPtrs=B.rowPtrs;
516  double ** CrowPtrs=C.rowPtrs;
517 
518  // double t0,t1;
519  // t0=vpTime::measureTimeMicros();
520  // for (int i=0;i<A.dsize;i++)
521  // {
522  // *(C.data + i) = *(B.data + i) + *(A.data + i) ;
523  // }
524  // t1=vpTime::measureTimeMicros();
525  // std::cout<< t1-t0<< std::endl;
526 
527  // t0=vpTime::measureTimeMicros();
528  for (unsigned int i=0;i<A.rowNum;i++)
529  for(unsigned int j=0;j<A.colNum;j++)
530  {
531 
532  CrowPtrs[i][j] = BrowPtrs[i][j]+ArowPtrs[i][j];
533  }
534  // t1=vpTime::measureTimeMicros();
535  // std::cout<< t1-t0<<" ";
536 
537 
538 }
539 
545 {
546 
547  vpMatrix C;
548  vpMatrix::add2Matrices(*this,B,C);
549  return C;
550 }
551 
552 
562 void vpMatrix::sub2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
563 {
564  try
565  {
566  if ((A.rowNum != C.rowNum) || (A.colNum != C.colNum)) C.resize(A.rowNum,A.colNum);
567  }
568  catch(vpException me)
569  {
570  vpERROR_TRACE("Error caught") ;
571  std::cout << me << std::endl ;
572  throw ;
573  }
574 
575  if ( (A.colNum != B.getCols())||(A.rowNum != B.getRows()))
576  {
577  vpERROR_TRACE("\n\t\t vpMatrix mismatch in vpMatrix/vpMatrix substraction") ;
579  "\n\t\t vpMatrix mismatch in "
580  "vpMatrix/vpMatrix substraction")) ;
581  }
582 
583 
584 
585  double ** ArowPtrs=A.rowPtrs;
586  double ** BrowPtrs=B.rowPtrs;
587  double ** CrowPtrs=C.rowPtrs;
588 
589  // double t0,t1;
590  //
591  // t0=vpTime::measureTimeMicros();
592  // for (int i=0;i<A.dsize;i++)
593  // {
594  // *(C.data + i) = *(A.data + i) - *(B.data + i) ;
595  // }
596  // t1=vpTime::measureTimeMicros();
597  // std::cout<< t1-t0<< std::endl;
598  //
599  // t0=vpTime::measureTimeMicros();
600  for (unsigned int i=0;i<A.rowNum;i++)
601  for(unsigned int j=0;j<A.colNum;j++)
602  {
603  CrowPtrs[i][j] = ArowPtrs[i][j]-BrowPtrs[i][j];
604  }
605  // t1=vpTime::measureTimeMicros();
606  // std::cout<< t1-t0<<" ";
607 
608 
609 }
610 
616 {
617  vpMatrix C;
618  vpMatrix::sub2Matrices(*this,B,C);
619  return C;
620 }
621 
623 
625 {
626  if ( (colNum != B.getCols())||(rowNum != B.getRows()))
627  {
628  vpERROR_TRACE("\n\t\t vpMatrix mismatch in vpMatrix += addition") ;
630  "\n\t\t vpMatrix mismatch in "
631  "vpMatrix += addition")) ;
632 
633  }
634 
635 
636  double ** BrowPtrs=B.rowPtrs;
637 
638  // double t0,t1;
639  //
640  // t0=vpTime::measureTimeMicros();
641  // for (int i=0;i<dsize;i++)
642  // {
643  // *(data + i) += *(B.data + i) ;
644  // }
645  // t1=vpTime::measureTimeMicros();
646  // std::cout<< t1-t0<< std::endl;
647 
648  // t0=vpTime::measureTimeMicros();
649  for (unsigned int i=0;i<rowNum;i++)
650  for(unsigned int j=0;j<colNum;j++)
651  rowPtrs[i][j] += BrowPtrs[i][j];
652  // t1=vpTime::measureTimeMicros();
653  // std::cout<< t1-t0<<" ";
654 
655 
656 
657 
658 
659  return *this;
660 }
661 
663 
665 {
666  if ( (colNum != B.getCols())||(rowNum != B.getRows()))
667  {
668  vpERROR_TRACE("\n\t\t vpMatrix mismatch in vpMatrix -= substraction") ;
670  "\n\t\t vpMatrix mismatch in "
671  "vpMatrix -= substraction")) ;
672 
673  }
674 
675  double ** BrowPtrs=B.rowPtrs;
676 
677  // double t0,t1;
678 
679  // t0=vpTime::measureTimeMicros();
680  // for (int i=0;i<dsize;i++)
681  // {
682  // *(data + i) -= *(B.data + i) ;
683  // }
684  // t1=vpTime::measureTimeMicros();
685  // std::cout<< t1-t0<< " " ;
686 
687  // t0=vpTime::measureTimeMicros();
688  for (unsigned int i=0;i<rowNum;i++)
689  for(unsigned int j=0;j<colNum;j++)
690  rowPtrs[i][j] -= BrowPtrs[i][j];
691 
692  // t1=vpTime::measureTimeMicros();
693  // std::cout<< t1-t0<<std::endl;
694 
695 
696 
697  return *this;
698 }
699 
710 {
711  try
712  {
713  if ((A.rowNum != C.rowNum) || (A.colNum != C.colNum)) C.resize(A.rowNum,A.colNum);
714  }
715  catch(vpException me)
716  {
717  vpERROR_TRACE("Error caught") ;
718  std::cout << me << std::endl ;
719  throw ;
720  }
721 
722  double ** ArowPtrs=A.rowPtrs;
723  double ** CrowPtrs=C.rowPtrs;
724 
725 
726  // std::cout << "negate" << std::endl;
727  // double t0,t1;
728  //
729  // t0=vpTime::measureTimeMicros();
730  // for (int i=0;i<A.dsize;i++)
731  // {
732  // *(C.data + i) = -*(A.data + i) ;
733  // }
734  // t1=vpTime::measureTimeMicros();
735  // std::cout<< t1-t0<< " ";
736 
737  // t0=vpTime::measureTimeMicros();
738  for (unsigned int i=0;i<A.rowNum;i++)
739  for(unsigned int j=0;j<A.colNum;j++)
740  CrowPtrs[i][j]= -ArowPtrs[i][j];
741  // t1=vpTime::measureTimeMicros();
742  // std::cout<< t1-t0<<std::endl;
743 
744 
745 }
746 
752 {
753  vpMatrix C;
754  vpMatrix::negateMatrix(*this,C);
755  return C;
756 }
757 
759 double
761 {
762  double sum=0.0;
763  double x ;
764 
765  // double t0,t1;
766  //
767  // t0=vpTime::measureTimeMicros();
768  // double *d = data ;
769  // double *n = data+dsize ;
770  // while (d < n )
771  // {
772  // x = *d++ ;
773  // sum += x*x ;
774  // }
775  // t1=vpTime::measureTimeMicros();
776  // std::cout<< t1-t0<<" "<< sum << " ";
777  //
778 
779 
780  // sum= 0.0;
781  // t0=vpTime::measureTimeMicros();
782  for (unsigned int i=0;i<rowNum;i++)
783  for(unsigned int j=0;j<colNum;j++)
784  {
785  x=rowPtrs[i][j];
786  sum+=x*x;
787  }
788  // t1=vpTime::measureTimeMicros();
789  // std::cout<< t1-t0<<" "<< sum << std::endl;
790 
791 
792 
793 
794  return sum;
795 }
796 
797 
798 //---------------------------------
799 // Matrix/vector operations.
800 //---------------------------------
801 
812 {
813  if (A.colNum != b.getRows())
814  {
815  vpERROR_TRACE("vpMatrix mismatch in vpMatrix/vector multiply") ;
817  }
818 
819  try
820  {
821  if (A.rowNum != c.rowNum) c.resize(A.rowNum);
822  }
823  catch(vpException me)
824  {
825  vpERROR_TRACE("Error caught") ;
826  std::cout << me << std::endl ;
827  throw ;
828  }
829 
830  c = 0.0;
831  for (unsigned int j=0;j<A.colNum;j++)
832  {
833  double bj = b[j] ; // optimization em 5/12/2006
834  for (unsigned int i=0;i<A.rowNum;i++)
835  {
836  c[i]+=A.rowPtrs[i][j] * bj;
837  }
838  }
839 }
840 
847 {
848  vpColVector c;
849  vpMatrix::multMatrixVector(*this,b,c);
850  return c;
851 }
852 
856 {
858 
859  if (rowNum != 3 || colNum != 3)
860  {
861  vpERROR_TRACE("vpMatrix mismatch in vpMatrix::operator*(const vpTranslationVector)") ;
863  }
864 
865  for (unsigned int j=0;j<3;j++) c[j]=0 ;
866 
867  for (unsigned int j=0;j<3;j++) {
868  {
869  double bj = b[j] ; // optimization em 5/12/2006
870  for (unsigned int i=0;i<3;i++) {
871  c[i]+=rowPtrs[i][j] * bj;
872  }
873  }
874  }
875  return c ;
876 }
877 
878 //---------------------------------
879 // Matrix/real operations.
880 //---------------------------------
881 
886 vpMatrix operator*(const double &x,const vpMatrix &B)
887 {
888  // Modif EM 13/6/03
889  vpMatrix C ;
890 
891  try {
892  C.resize(B.getRows(),B.getCols());
893  }
894  catch(vpException me)
895  {
896  vpERROR_TRACE("Error caught") ;
897  std::cout << me << std::endl ;
898  throw ;
899  }
900  // double t0,t1;
901 
902  unsigned int Brow = B.getRows() ;
903  unsigned int Bcol = B.getCols() ;
904  // t0=vpTime::measureTimeMicros();
905  //
906  // for (int i=0;i<Brow; i++)
907  // {
908  // double *ci = C[i] ;
909  // double *Bi = B[i] ;
910  // for (int j=0 ; j < Bcol;j++)
911  // ci[j] = Bi[j]*x;
912  // }
913  //
914  // t1=vpTime::measureTimeMicros();
915  // std::cout<< t1-t0<<" ";
916 
917  // t0=vpTime::measureTimeMicros();
918  for (unsigned int i=0;i<Brow;i++)
919  for(unsigned int j=0;j<Bcol;j++)
920  C[i][j]= B[i][j]*x;
921  // t1=vpTime::measureTimeMicros();
922  // std::cout<< t1-t0<<std::endl;
923 
924 
925 
926  return C ;
927 }
928 
931 {
932  vpMatrix C;
933 
934  try {
935  C.resize(rowNum,colNum);
936  }
937  catch(vpException me)
938  {
939  vpERROR_TRACE("Error caught") ;
940  std::cout << me << std::endl ;
941  throw ;
942  }
943  // double t0,t1;
944 
945  double ** CrowPtrs=C.rowPtrs;
946 
947  // t0=vpTime::measureTimeMicros();
948  // for (int i=0;i<dsize;i++)
949  // *(C.data+i) = *(data+i)*x;
950  //
951  // t1=vpTime::measureTimeMicros();
952  // std::cout<< t1-t0<<" ";
953  //
954  // t0=vpTime::measureTimeMicros();
955  for (unsigned int i=0;i<rowNum;i++)
956  for(unsigned int j=0;j<colNum;j++)
957  CrowPtrs[i][j]= rowPtrs[i][j]*x;
958  // t1=vpTime::measureTimeMicros();
959  // std::cout<< t1-t0<<std::endl;
960 
961  return C;
962 }
963 
966 {
967  vpMatrix C;
968 
969  try {
970  C.resize(rowNum,colNum);
971  }
972  catch(vpException me)
973  {
974  vpERROR_TRACE("Error caught") ;
975  vpCERROR << me << std::endl ;
976  throw ;
977  }
978 
979  //if (x == 0) {
980  if (std::fabs(x) <= std::numeric_limits<double>::epsilon()) {
981  vpERROR_TRACE("Divide by zero in method /(double x)") ;
982  throw vpMatrixException(vpMatrixException::divideByZeroError, "Divide by zero in method /(double x)");
983  }
984 
985  double xinv = 1/x ;
986 
987 
988  // double t0,t1;
989  //
990  // t0=vpTime::measureTimeMicros();
991  // for (int i=0;i<dsize;i++)
992  // *(C.data+i) = *(data+i)*xinv ;
993  //
994  // t1=vpTime::measureTimeMicros();
995  // std::cout<< t1-t0<<" ";
996  //
997  // t0=vpTime::measureTimeMicros();
998  for (unsigned int i=0;i<rowNum;i++)
999  for(unsigned int j=0;j<colNum;j++)
1000  C[i][j]=rowPtrs[i][j]*xinv;
1001  // t1=vpTime::measureTimeMicros();
1002  // std::cout<< t1-t0<<std::endl;
1003 
1004  return C;
1005 
1006 }
1007 
1008 
1011 {
1012 
1013  // double t0,t1;
1014  //
1015  // t0=vpTime::measureTimeMicros();
1016  // for (int i=0;i<dsize;i++)
1017  // *(data+i) += x;
1018  //
1019  // t1=vpTime::measureTimeMicros();
1020  // std::cout<< t1-t0<<" ";
1021  //
1022  // t0=vpTime::measureTimeMicros();
1023  for (unsigned int i=0;i<rowNum;i++)
1024  for(unsigned int j=0;j<colNum;j++)
1025  rowPtrs[i][j]+=x;
1026  // t1=vpTime::measureTimeMicros();
1027  // std::cout<< t1-t0<<std::endl;
1028 
1029  return *this;
1030 }
1031 
1032 
1035 {
1036 
1037 
1038  // double t0,t1;
1039  //
1040  // t0=vpTime::measureTimeMicros();
1041  // for (int i=0;i<dsize;i++)
1042  // *(data+i) -= x;
1043  //
1044  // t1=vpTime::measureTimeMicros();
1045  // std::cout<< t1-t0<<" ";
1046  //
1047  // t0=vpTime::measureTimeMicros();
1048  for (unsigned int i=0;i<rowNum;i++)
1049  for(unsigned int j=0;j<colNum;j++)
1050  rowPtrs[i][j]-=x;
1051  // t1=vpTime::measureTimeMicros();
1052  // std::cout<< t1-t0<<std::endl;
1053 
1054  return *this;
1055 }
1056 
1059 {
1060 
1061 
1062  // for (int i=0;i<dsize;i++)
1063  // *(data+i) *= x;
1064 
1065  for (unsigned int i=0;i<rowNum;i++)
1066  for(unsigned int j=0;j<colNum;j++)
1067  rowPtrs[i][j]*=x;
1068 
1069  return *this;
1070 }
1071 
1074 {
1075  //if (x == 0)
1076  if (std::fabs(x) <= std::numeric_limits<double>::epsilon())
1077  throw vpMatrixException(vpMatrixException::divideByZeroError, "Divide by zero in method /=(double x)");
1078 
1079  double xinv = 1/x ;
1080 
1081  // double t0,t1;
1082  //
1083  // t0=vpTime::measureTimeMicros();
1084  // for (int i=0;i<dsize;i++)
1085  // *(data+i) *= xinv;
1086  //
1087  // t1=vpTime::measureTimeMicros();
1088  // std::cout<< t1-t0<<" ";
1089 
1090  // t0=vpTime::measureTimeMicros();
1091  for (unsigned int i=0;i<rowNum;i++)
1092  for(unsigned int j=0;j<colNum;j++)
1093  rowPtrs[i][j]*=xinv;
1094  // t1=vpTime::measureTimeMicros();
1095  // std::cout<< t1-t0<<std::endl;
1096 
1097  return *this;
1098 }
1099 
1100 //----------------------------------------------------------------
1101 // Matrix Operation
1102 //----------------------------------------------------------------
1103 
1104 
1109 void
1110 vpMatrix::setIdentity(const double & val)
1111 {
1112  if (rowNum != colNum)
1113  {
1114  vpERROR_TRACE("non square matrix") ;
1116  }
1117 
1118  unsigned int i,j;
1119  for (i=0;i<rowNum;i++)
1120  for (j=0;j<colNum;j++)
1121  if (i==j) (*this)[i][j] = val ; else (*this)[i][j] = 0;
1122 }
1123 
1131 void
1132 vpMatrix::eye(unsigned int n)
1133 {
1134  try {
1135  eye(n,n) ;
1136  }
1137  catch(vpException me)
1138  {
1139  vpERROR_TRACE("Error caught") ;
1140  vpCERROR << me << std::endl ;
1141  throw ;
1142  }
1143 }
1151 void
1152 vpMatrix::eye(unsigned int m, unsigned int n)
1153 {
1154  try {
1155  resize(m,n) ;
1156  }
1157  catch(vpException me)
1158  {
1159  vpERROR_TRACE("Error caught") ;
1160  vpCERROR << me << std::endl ;
1161  throw ;
1162  }
1163 
1164 
1165  for (unsigned int i=0; i<rowNum; i++)
1166  for (unsigned int j=0; j<colNum; j++)
1167  if (i == j) (*this)[i][j] = 1;
1168  else (*this)[i][j] = 0;
1169 
1170 }
1171 
1172 
1177 {
1178  vpMatrix At ;
1179 
1180  try {
1181  At.resize(colNum,rowNum);
1182  }
1183  catch(vpException me)
1184  {
1185  vpERROR_TRACE("Error caught") ;
1186  vpCERROR << me << std::endl ;
1187  throw ;
1188  }
1189 
1190  unsigned int i,j;
1191  for (i=0;i<rowNum;i++)
1192  {
1193  double *coli = (*this)[i] ;
1194  for (j=0;j<colNum;j++)
1195  At[j][i] = coli[j];
1196  }
1197  return At;
1198 }
1199 
1200 
1207 {
1208  vpMatrix At ;
1209  transpose(At);
1210  return At;
1211 }
1212 
1219 
1220  try {
1221  At.resize(colNum,rowNum);
1222  }
1223  catch(vpException me)
1224  {
1225  vpERROR_TRACE("Error caught") ;
1226  vpCERROR << me << std::endl ;
1227  throw ;
1228  }
1229 
1230  size_t A_step = colNum;
1231  double ** AtRowPtrs = At.rowPtrs;
1232 
1233  for( unsigned int i = 0; i < colNum; i++ )
1234  {
1235  double * row = AtRowPtrs[i];
1236  double * col = rowPtrs[0]+i;
1237  for( unsigned int j = 0; j < rowNum; j++, col+=A_step )
1238  *(row++)=*col;
1239  }
1240 }
1241 
1242 
1249 {
1250  vpMatrix B;
1251 
1252  AAt(B);
1253 
1254  return B;
1255 }
1256 
1268 void vpMatrix::AAt(vpMatrix &B)const {
1269 
1270  try {
1271  if ((B.rowNum != rowNum) || (B.colNum != rowNum)) B.resize(rowNum,rowNum);
1272  }
1273  catch(vpException me)
1274  {
1275  vpERROR_TRACE("Error caught") ;
1276  vpCERROR << me << std::endl ;
1277  throw ;
1278  }
1279 
1280  // compute A*A^T
1281  for(unsigned int i=0;i<rowNum;i++){
1282  for(unsigned int j=i;j<rowNum;j++){
1283  double *pi = rowPtrs[i];// row i
1284  double *pj = rowPtrs[j];// row j
1285 
1286  // sum (row i .* row j)
1287  double ssum=0;
1288  for(unsigned int k=0; k < colNum ;k++)
1289  ssum += *(pi++)* *(pj++);
1290 
1291  B[i][j]=ssum; //upper triangle
1292  if(i!=j)
1293  B[j][i]=ssum; //lower triangle
1294  }
1295  }
1296 }
1297 
1298 
1299 
1311 void vpMatrix::AtA(vpMatrix &B) const
1312 {
1313  try {
1314  if ((B.rowNum != colNum) || (B.colNum != colNum)) B.resize(colNum,colNum);
1315  }
1316  catch(vpException me)
1317  {
1318  vpERROR_TRACE("Error caught") ;
1319  vpCERROR << me << std::endl ;
1320  throw ;
1321  }
1322 
1323  unsigned int i,j,k;
1324  double s;
1325  double *ptr;
1326  double *Bi;
1327  for (i=0;i<colNum;i++)
1328  {
1329  Bi = B[i] ;
1330  for (j=0;j<i;j++)
1331  {
1332  ptr=data;
1333  s = 0 ;
1334  for (k=0;k<rowNum;k++)
1335  {
1336  s +=(*(ptr+i)) * (*(ptr+j));
1337  ptr+=colNum;
1338  }
1339  *Bi++ = s ;
1340  B[j][i] = s;
1341  }
1342  ptr=data;
1343  s = 0 ;
1344  for (k=0;k<rowNum;k++)
1345  {
1346  s +=(*(ptr+i)) * (*(ptr+i));
1347  ptr+=colNum;
1348  }
1349  *Bi = s;
1350  }
1351 }
1352 
1353 
1360 {
1361  vpMatrix B;
1362 
1363  AtA(B);
1364 
1365  return B;
1366 }
1367 
1368 
1374 
1375  try {
1376  if ((out.rowNum != colNum*rowNum) || (out.colNum != 1)) out.resize(rowNum);
1377  }
1378  catch(vpException me)
1379  {
1380  vpERROR_TRACE("Error caught") ;
1381  vpCERROR << me << std::endl ;
1382  throw ;
1383  }
1384 
1385  double *optr=out.data;
1386  for(unsigned int j =0;j<colNum ; j++){
1387  for(unsigned int i =0;i<rowNum ; i++){
1388  *(optr++)=rowPtrs[i][j];
1389  }
1390  }
1391 }
1392 
1398 
1399  vpColVector out(colNum*rowNum);
1400  stackColumns(out);
1401  return out;
1402 }
1403 
1409 
1410  try {
1411  if ((out.rowNum != 1) || (out.colNum != colNum*rowNum)) out.resize(rowNum);
1412  }
1413  catch(vpException me)
1414  {
1415  vpERROR_TRACE("Error caught") ;
1416  vpCERROR << me << std::endl ;
1417  throw ;
1418  }
1419 
1420  double *mdata=data;
1421  double *optr=out.data;
1422  for(unsigned int i =0;i<dsize ; i++){
1423  *(optr++)=*(mdata++);
1424  }
1425 }
1431 
1432  vpRowVector out(colNum*rowNum);
1433  stackRows(out );
1434  return out;
1435 }
1436 
1443 void vpMatrix::kron(const vpMatrix &m1, const vpMatrix &m2 , vpMatrix &out){
1444  unsigned int r1= m1.getRows();
1445  unsigned int c1= m1.getCols();
1446  unsigned int r2= m2.getRows();
1447  unsigned int c2= m2.getCols();
1448 
1449 
1450  if (r1*r2 !=out.rowNum || c1*c2!= out.colNum )
1451  {
1452  vpERROR_TRACE("Kronecker prodect bad dimension of output vpMatrix") ;
1453  throw(vpMatrixException(vpMatrixException::incorrectMatrixSizeError,"Kronecker prodect bad dimension of output vpMatrix"));
1454  }
1455 
1456  for(unsigned int r =0;r<r1 ; r++){
1457  for(unsigned int c =0;c<c1 ; c++){
1458  double alpha = m1[r][c];
1459  double *m2ptr = m2[0];
1460  unsigned int roffset= r*r2;
1461  unsigned int coffset= c*c2;
1462  for(unsigned int rr =0;rr<r2 ; rr++){
1463  for(unsigned int cc =0;cc<c2 ;cc++){
1464  out[roffset+rr][coffset+cc]= alpha* *(m2ptr++);
1465  }
1466  }
1467  }
1468  }
1469 
1470 }
1471 
1477 void vpMatrix::kron(const vpMatrix &m , vpMatrix &out){
1478  kron(*this,m,out);
1479 }
1480 
1487 vpMatrix vpMatrix::kron(const vpMatrix &m1, const vpMatrix &m2 ){
1488 
1489  unsigned int r1= m1.getRows();
1490  unsigned int c1= m1.getCols();
1491  unsigned int r2= m2.getRows();
1492  unsigned int c2= m2.getCols();
1493 
1494  vpMatrix out(r1*r2,c1*c2);
1495 
1496  for(unsigned int r =0;r<r1 ; r++){
1497  for(unsigned int c =0;c<c1 ; c++){
1498  double alpha = m1[r][c];
1499  double *m2ptr = m2[0];
1500  unsigned int roffset= r*r2;
1501  unsigned int coffset= c*c2;
1502  for(unsigned int rr =0;rr<r2 ; rr++){
1503  for(unsigned int cc =0;cc<c2 ;cc++){
1504  out[roffset+rr ][coffset+cc]= alpha* *(m2ptr++);
1505  }
1506  }
1507  }
1508  }
1509  return out;
1510 }
1511 
1512 
1519  return kron(*this,m);
1520 }
1521 
1572 void
1574 {
1575  x = pseudoInverse(1e-6)*b ;
1576 }
1577 
1578 
1629 {
1630  vpColVector X(colNum);
1631 
1632  solveBySVD(B, X);
1633  return X;
1634 }
1635 
1636 
1701 void
1703 {
1704 #if (DEBUG_LEVEL1 == 0) /* no verification */
1705  {
1706  w.resize( this->getCols() );
1707  v.resize( this->getCols(), this->getCols() );
1708 
1709 #if defined (VISP_HAVE_LAPACK)
1710  svdLapack(w,v);
1711 #elif (VISP_HAVE_OPENCV_VERSION >= 0x020101) // Require opencv >= 2.1.1
1712  svdOpenCV(w,v);
1713 #elif defined (VISP_HAVE_GSL) /* be careful of the copy below */
1714  svdGsl(w,v) ;
1715 #else
1716  svdNr(w,v) ;
1717 #endif
1718 
1719  //svdNr(w,v) ;
1720  }
1721 #else /* verification of the SVD */
1722  {
1723  int pb = 0;
1724  unsigned int i,j,k,nrows,ncols;
1725  vpMatrix A, Asvd;
1726 
1727  A = (*this); /* copy because svd is destructive */
1728 
1729  w.resize( this->getCols() );
1730  v.resize( this->getCols(), this->getCols() );
1731 #ifdef VISP_HAVE_GSL /* be careful of the copy above */
1732  svdGsl(w,v) ;
1733 #else
1734  svdNr(w,v) ;
1735 #endif
1736  //svdNr(w,v) ;
1737 
1738  nrows = A.getRows();
1739  ncols = A.getCols();
1740  Asvd.resize(nrows,ncols);
1741 
1742  for (i = 0 ; i < nrows ; i++)
1743  {
1744  for (j = 0 ; j < ncols ; j++)
1745  {
1746  Asvd[i][j] = 0.0;
1747  for (k=0 ; k < ncols ; k++) Asvd[i][j] += (*this)[i][k]*w[k]*v[j][k];
1748  }
1749  }
1750  for (i=0;i<nrows;i++)
1751  {
1752  for (j=0;j<ncols;j++) if (fabs(A[i][j]-Asvd[i][j]) > 1e-6) pb = 1;
1753  }
1754  if (pb == 1)
1755  {
1756  printf("pb in SVD\n");
1757  std::cout << " A : " << std::endl << A << std::endl;
1758  std::cout << " Asvd : " << std::endl << Asvd << std::endl;
1759  }
1760  // else printf("SVD ok ;-)\n"); /* It's so good... */
1761  }
1762 #endif
1763 }
1771 unsigned int
1772 vpMatrix::pseudoInverse(vpMatrix &Ap, double svThreshold) const
1773 {
1774  vpColVector sv ;
1775  return pseudoInverse(Ap, sv, svThreshold) ;
1776 }
1777 
1811 vpMatrix
1812 vpMatrix::pseudoInverse(double svThreshold) const
1813 {
1814  vpMatrix Ap ;
1815  vpColVector sv ;
1816  pseudoInverse(Ap, sv, svThreshold) ;
1817  return Ap ;
1818 }
1819 
1827 unsigned int
1828 vpMatrix::pseudoInverse(vpMatrix &Ap, vpColVector &sv, double svThreshold) const
1829 {
1830  vpMatrix imA, imAt ;
1831  return pseudoInverse(Ap, sv, svThreshold, imA, imAt) ;
1832 }
1833 
1866 unsigned int
1868  vpColVector &sv, double svThreshold,
1869  vpMatrix &imA,
1870  vpMatrix &imAt) const
1871 {
1872 
1873  unsigned int i, j, k ;
1874 
1875  unsigned int nrows, ncols;
1876  unsigned int nrows_orig = getRows() ;
1877  unsigned int ncols_orig = getCols() ;
1878  Ap.resize(ncols_orig,nrows_orig) ;
1879 
1880  if (nrows_orig >= ncols_orig)
1881  {
1882  nrows = nrows_orig;
1883  ncols = ncols_orig;
1884  }
1885  else
1886  {
1887  nrows = ncols_orig;
1888  ncols = nrows_orig;
1889  }
1890 
1891  vpMatrix a(nrows,ncols) ;
1892  vpMatrix a1(ncols,nrows);
1893  vpMatrix v(ncols,ncols) ;
1894  sv.resize(ncols) ;
1895 
1896  if (nrows_orig >= ncols_orig) a = *this;
1897  else a = (*this).t();
1898 
1899  a.svd(sv,v);
1900 
1901  // compute the highest singular value and the rank of h
1902  double maxsv = 0 ;
1903  for (i=0 ; i < ncols ; i++)
1904  if (fabs(sv[i]) > maxsv) maxsv = fabs(sv[i]) ;
1905 
1906  unsigned int rank = 0 ;
1907  for (i=0 ; i < ncols ; i++)
1908  if (fabs(sv[i]) > maxsv*svThreshold) rank++ ;
1909 
1910 
1911 
1912  /*------------------------------------------------------- */
1913  for (i = 0 ; i < ncols ; i++)
1914  {
1915  for (j = 0 ; j < nrows ; j++)
1916  {
1917  a1[i][j] = 0.0;
1918 
1919  for (k=0 ; k < ncols ; k++)
1920  if (fabs(sv[k]) > maxsv*svThreshold)
1921  {
1922  a1[i][j] += v[i][k]*a[j][k]/sv[k];
1923  }
1924  }
1925  }
1926  if (nrows_orig >= ncols_orig) Ap = a1;
1927  else Ap = a1.t();
1928 
1929  if (nrows_orig >= ncols_orig)
1930  {
1931  // compute dim At
1932  imAt.resize(ncols_orig,rank) ;
1933  for (i=0 ; i < ncols_orig ; i++)
1934  for (j=0 ; j < rank ; j++)
1935  imAt[i][j] = v[i][j] ;
1936 
1937  // compute dim A
1938  imA.resize(nrows_orig,rank) ;
1939  for (i=0 ; i < nrows_orig ; i++)
1940  for (j=0 ; j < rank ; j++)
1941  imA[i][j] = a[i][j] ;
1942  }
1943  else
1944  {
1945  // compute dim At
1946  imAt.resize(ncols_orig,rank) ;
1947  for (i=0 ; i < ncols_orig ; i++)
1948  for (j=0 ; j < rank ; j++)
1949  imAt[i][j] = a[i][j] ;
1950 
1951  imA.resize(nrows_orig,rank) ;
1952  for (i=0 ; i < nrows_orig ; i++)
1953  for (j=0 ; j < rank ; j++)
1954  imA[i][j] = v[i][j] ;
1955 
1956  }
1957 
1958 #if (DEBUG_LEVEL1)
1959  {
1960  int pb = 0;
1961  vpMatrix A, ApA, AAp, AApA, ApAAp ;
1962 
1963  nrows = nrows_orig;
1964  ncols = ncols_orig;
1965 
1966  A.resize(nrows,ncols) ;
1967  A = *this ;
1968 
1969  ApA = Ap * A;
1970  AApA = A * ApA;
1971  ApAAp = ApA * Ap;
1972  AAp = A * Ap;
1973 
1974  for (i=0;i<nrows;i++)
1975  {
1976  for (j=0;j<ncols;j++) if (fabs(AApA[i][j]-A[i][j]) > 1e-6) pb = 1;
1977  }
1978  for (i=0;i<ncols;i++)
1979  {
1980  for (j=0;j<nrows;j++) if (fabs(ApAAp[i][j]-Ap[i][j]) > 1e-6) pb = 1;
1981  }
1982  for (i=0;i<nrows;i++)
1983  {
1984  for (j=0;j<nrows;j++) if (fabs(AAp[i][j]-AAp[j][i]) > 1e-6) pb = 1;
1985  }
1986  for (i=0;i<ncols;i++)
1987  {
1988  for (j=0;j<ncols;j++) if (fabs(ApA[i][j]-ApA[j][i]) > 1e-6) pb = 1;
1989  }
1990  if (pb == 1)
1991  {
1992  printf("pb in pseudo inverse\n");
1993  std::cout << " A : " << std::endl << A << std::endl;
1994  std::cout << " Ap : " << std::endl << Ap << std::endl;
1995  std::cout << " A - AApA : " << std::endl << A - AApA << std::endl;
1996  std::cout << " Ap - ApAAp : " << std::endl << Ap - ApAAp << std::endl;
1997  std::cout << " AAp - (AAp)^T : " << std::endl << AAp - AAp.t() << std::endl;
1998  std::cout << " ApA - (ApA)^T : " << std::endl << ApA - ApA.t() << std::endl;
1999  }
2000  // else printf("Ap OK ;-) \n");
2001 
2002  }
2003 #endif
2004 
2005  // std::cout << v << std::endl ;
2006  return rank ;
2007 }
2008 
2009 
2010 
2044 unsigned int
2046  vpColVector &sv, double svThreshold,
2047  vpMatrix &imA,
2048  vpMatrix &imAt,
2049  vpMatrix &kerA) const
2050 {
2051 
2052  unsigned int i, j, k ;
2053 
2054  unsigned int nrows, ncols;
2055  unsigned int nrows_orig = getRows() ;
2056  unsigned int ncols_orig = getCols() ;
2057  Ap.resize(ncols_orig,nrows_orig) ;
2058 
2059  if (nrows_orig >= ncols_orig)
2060  {
2061  nrows = nrows_orig;
2062  ncols = ncols_orig;
2063  }
2064  else
2065  {
2066  nrows = ncols_orig;
2067  ncols = nrows_orig;
2068  }
2069 
2070  vpMatrix a(nrows,ncols) ;
2071  vpMatrix a1(ncols,nrows);
2072  vpMatrix v(ncols,ncols) ;
2073  sv.resize(ncols) ;
2074 
2075  if (nrows_orig >= ncols_orig) a = *this;
2076  else a = (*this).t();
2077 
2078  a.svd(sv,v);
2079 
2080  // compute the highest singular value and the rank of h
2081  double maxsv = 0 ;
2082  for (i=0 ; i < ncols ; i++)
2083  if (fabs(sv[i]) > maxsv) maxsv = fabs(sv[i]) ;
2084 
2085  unsigned int rank = 0 ;
2086  for (i=0 ; i < ncols ; i++)
2087  if (fabs(sv[i]) > maxsv*svThreshold) rank++ ;
2088 
2089 
2090 
2091  /*------------------------------------------------------- */
2092  for (i = 0 ; i < ncols ; i++)
2093  {
2094  for (j = 0 ; j < nrows ; j++)
2095  {
2096  a1[i][j] = 0.0;
2097 
2098  for (k=0 ; k < ncols ; k++)
2099  if (fabs(sv[k]) > maxsv*svThreshold)
2100  {
2101  a1[i][j] += v[i][k]*a[j][k]/sv[k];
2102  }
2103  }
2104  }
2105  if (nrows_orig >= ncols_orig) Ap = a1;
2106  else Ap = a1.t();
2107 
2108  if (nrows_orig >= ncols_orig)
2109  {
2110  // compute dim At
2111  imAt.resize(ncols_orig,rank) ;
2112  for (i=0 ; i < ncols_orig ; i++)
2113  for (j=0 ; j < rank ; j++)
2114  imAt[i][j] = v[i][j] ;
2115 
2116  // compute dim A
2117  imA.resize(nrows_orig,rank) ;
2118  for (i=0 ; i < nrows_orig ; i++)
2119  for (j=0 ; j < rank ; j++)
2120  imA[i][j] = a[i][j] ;
2121  }
2122  else
2123  {
2124  // compute dim At
2125  imAt.resize(ncols_orig,rank) ;
2126  for (i=0 ; i < ncols_orig ; i++)
2127  for (j=0 ; j < rank ; j++)
2128  imAt[i][j] = a[i][j] ;
2129 
2130  imA.resize(nrows_orig,rank) ;
2131  for (i=0 ; i < nrows_orig ; i++)
2132  for (j=0 ; j < rank ; j++)
2133  imA[i][j] = v[i][j] ;
2134 
2135  }
2136 
2137  vpMatrix cons(ncols_orig, ncols_orig);
2138  cons = 0;
2139 
2140  for (j = 0; j < ncols_orig; j++)
2141  {
2142  for (i = 0; i < ncols_orig; i++)
2143  {
2144  if (fabs(sv[i]) <= maxsv*svThreshold)
2145  {
2146  cons[i][j] = v[j][i];
2147  }
2148  }
2149  }
2150 
2151  vpMatrix Ker (ncols_orig-rank, ncols_orig);
2152  k = 0;
2153  for (j = 0; j < ncols_orig ; j++)
2154  {
2155  //if ( cons.row(j+1).sumSquare() != 0)
2156  if ( std::fabs(cons.row(j+1).sumSquare()) > std::numeric_limits<double>::epsilon())
2157  {
2158  for (i = 0; i < cons.getCols(); i++)
2159  Ker[k][i] = cons[j][i];
2160 
2161  k++;
2162  }
2163  }
2164  kerA = Ker;
2165 
2166 #if (DEBUG_LEVEL1)
2167  {
2168  int pb = 0;
2169  vpMatrix A, ApA, AAp, AApA, ApAAp ;
2170 
2171  nrows = nrows_orig;
2172  ncols = ncols_orig;
2173 
2174  A.resize(nrows,ncols) ;
2175  A = *this ;
2176 
2177  ApA = Ap * A;
2178  AApA = A * ApA;
2179  ApAAp = ApA * Ap;
2180  AAp = A * Ap;
2181 
2182  for (i=0;i<nrows;i++)
2183  {
2184  for (j=0;j<ncols;j++) if (fabs(AApA[i][j]-A[i][j]) > 1e-6) pb = 1;
2185  }
2186  for (i=0;i<ncols;i++)
2187  {
2188  for (j=0;j<nrows;j++) if (fabs(ApAAp[i][j]-Ap[i][j]) > 1e-6) pb = 1;
2189  }
2190  for (i=0;i<nrows;i++)
2191  {
2192  for (j=0;j<nrows;j++) if (fabs(AAp[i][j]-AAp[j][i]) > 1e-6) pb = 1;
2193  }
2194  for (i=0;i<ncols;i++)
2195  {
2196  for (j=0;j<ncols;j++) if (fabs(ApA[i][j]-ApA[j][i]) > 1e-6) pb = 1;
2197  }
2198  if (pb == 1)
2199  {
2200  printf("pb in pseudo inverse\n");
2201  std::cout << " A : " << std::endl << A << std::endl;
2202  std::cout << " Ap : " << std::endl << Ap << std::endl;
2203  std::cout << " A - AApA : " << std::endl << A - AApA << std::endl;
2204  std::cout << " Ap - ApAAp : " << std::endl << Ap - ApAAp << std::endl;
2205  std::cout << " AAp - (AAp)^T : " << std::endl << AAp - AAp.t() << std::endl;
2206  std::cout << " ApA - (ApA)^T : " << std::endl << ApA - ApA.t() << std::endl;
2207  std::cout << " KerA : " << std::endl << kerA << std::endl;
2208  }
2209  // else printf("Ap OK ;-) \n");
2210 
2211  }
2212 #endif
2213 
2214  // std::cout << v << std::endl ;
2215  return rank ;
2216 }
2217 
2218 
2225 vpMatrix::row(const unsigned int j)
2226 {
2227  vpRowVector c(getCols()) ;
2228 
2229  for (unsigned int i =0 ; i < getCols() ; i++) c[i] = (*this)[j-1][i] ;
2230  return c ;
2231 }
2232 
2233 
2240 vpMatrix::column(const unsigned int j)
2241 {
2242  vpColVector c(getRows()) ;
2243 
2244  for (unsigned int i =0 ; i < getRows() ; i++) c[i] = (*this)[i][j-1] ;
2245  return c ;
2246 }
2247 
2248 
2249 
2250 
2262 vpMatrix
2264 {
2265  vpMatrix C ;
2266 
2267  try{
2268  stackMatrices(A,B, C) ;
2269  }
2270  catch(vpMatrixException me)
2271  {
2272  vpCERROR << me << std::endl;
2273  throw ;
2274  }
2275 
2276  return C ;
2277 }
2278 
2291 void
2293 {
2294  unsigned int nra = A.getRows() ;
2295  unsigned int nrb = B.getRows() ;
2296 
2297  if (nra !=0)
2298  if (A.getCols() != B.getCols())
2299  {
2300  vpERROR_TRACE("\n\t\t incorrect matrices size") ;
2302  "\n\t\t incorrect matrices size")) ;
2303  }
2304 
2305  try {
2306  C.resize(nra+nrb,B.getCols() ) ;
2307  }
2308  catch(vpException me)
2309  {
2310  vpERROR_TRACE("Error caught") ;
2311  vpCERROR << me << std::endl ;
2312  throw ;
2313  }
2314 
2315  unsigned int i,j ;
2316  for (i=0 ; i < nra ; i++)
2317  for (j=0 ; j < A.getCols() ; j++)
2318  C[i][j] = A[i][j] ;
2319 
2320 
2321  for (i=0 ; i < nrb ; i++)
2322  for (j=0 ; j < B.getCols() ; j++)
2323  {
2324  C[i+nra][j] = B[i][j] ;
2325 
2326  }
2327 
2328 
2329 }
2330 
2331 
2332 
2333 
2334 
2347 vpMatrix
2348 vpMatrix::insert(const vpMatrix &A, const vpMatrix &B,
2349  const unsigned int r, const unsigned int c)
2350 {
2351  vpMatrix C ;
2352 
2353  try{
2354  insert(A,B, C, r, c) ;
2355  }
2356  catch(vpMatrixException me)
2357  {
2358  vpCERROR << me << std::endl;
2359  throw me;
2360  }
2361 
2362  return C ;
2363 }
2364 
2378 void
2379 vpMatrix::insert(const vpMatrix &A, const vpMatrix &B, vpMatrix &C,
2380  const unsigned int r, const unsigned int c)
2381 {
2382  if( ( (r + B.getRows()) <= A.getRows() ) &&
2383  ( (c + B.getCols()) <= A.getRows() ) ){
2384  try {
2385  C.resize(A.getRows(),A.getCols() ) ;
2386  }
2387  catch(vpException me)
2388  {
2389  vpERROR_TRACE("Error caught") ;
2390  vpCERROR << me << std::endl ;
2391  throw ;
2392  }
2393  for(unsigned int i=0; i<A.getRows(); i++){
2394  for(unsigned int j=0; j<A.getCols(); j++){
2395  if(i >= r && i < (r + B.getRows()) && j >= c && j < (c+B.getCols())){
2396  C[i][j] = B[i-r][j-c];
2397  }
2398  else{
2399  C[i][j] = A[i][j];
2400  }
2401  }
2402  }
2403  }
2404  else{
2406  "\n\t\tIncorrect size of the matrix to insert data.");
2407  }
2408 }
2409 
2421 vpMatrix
2423 {
2424  vpMatrix C ;
2425 
2426  try{
2427  juxtaposeMatrices(A,B, C) ;
2428  }
2429  catch(vpMatrixException me)
2430  {
2431  vpCERROR << me << std::endl ;
2432  throw ;
2433  }
2434 
2435  return C ;
2436 }
2437 
2450 void
2452 {
2453  unsigned int nca = A.getCols() ;
2454  unsigned int ncb = B.getCols() ;
2455 
2456  if (nca !=0)
2457  if (A.getRows() != B.getRows())
2458  {
2459  vpERROR_TRACE("\n\t\t incorrect matrices size") ;
2461  "\n\t\t incorrect matrices size")) ;
2462  }
2463 
2464  try {
2465  C.resize(B.getRows(),nca+ncb) ;
2466  }
2467  catch(vpException me)
2468  {
2469  vpERROR_TRACE("Error caught") ;
2470  vpCERROR << me << std::endl ;
2471  throw ;
2472  }
2473 
2474  unsigned int i,j ;
2475  for (i=0 ; i < C.getRows(); i++)
2476  for (j=0 ; j < nca ; j++)
2477  C[i][j] = A[i][j] ;
2478 
2479 
2480  for (i=0 ; i < C.getRows() ; i++)
2481  for (j=0 ; j < ncb ; j++)
2482  {
2483  C[i][nca+j] = B[i][j] ;
2484  }
2485 }
2486 
2522 void
2524 {
2525  unsigned int rows = A.getRows() ;
2526  try {
2527  this->resize(rows,rows) ;
2528  }
2529  catch(vpException me)
2530  {
2531  vpERROR_TRACE("Error caught") ;
2532  vpCERROR << me << std::endl ;
2533  throw ;
2534  }
2535  (*this) = 0 ;
2536  for (unsigned int i=0 ; i< rows ; i++ )
2537  (* this)[i][i] = A[i] ;
2538 }
2550 void
2552 {
2553  unsigned int rows = A.getRows() ;
2554  try {
2555  DA.resize(rows,rows) ;
2556  }
2557  catch(vpException me)
2558  {
2559  vpERROR_TRACE("Error caught") ;
2560  vpCERROR << me << std::endl ;
2561  throw ;
2562  }
2563  DA =0 ;
2564  for (unsigned int i=0 ; i< rows ; i++ )
2565  DA[i][i] = A[i] ;
2566 }
2567 
2568 //--------------------------------------------------------------------
2569 // Output
2570 //--------------------------------------------------------------------
2571 
2572 
2576 std::ostream &operator <<(std::ostream &s,const vpMatrix &m)
2577 {
2578  s.precision(10) ;
2579  for (unsigned int i=0;i<m.getRows();i++) {
2580  for (unsigned int j=0;j<m.getCols();j++){
2581  s << m[i][j] << " ";
2582  }
2583  // We don't add a \n char on the end of the last matrix line
2584  if (i < m.getRows()-1)
2585  s << std::endl;
2586  }
2587 
2588  return s;
2589 }
2590 
2614 int
2615 vpMatrix::print(std::ostream& s, unsigned int length, char const* intro)
2616 {
2617  typedef std::string::size_type size_type;
2618 
2619  unsigned int m = getRows();
2620  unsigned int n = getCols();
2621 
2622  std::vector<std::string> values(m*n);
2623  std::ostringstream oss;
2624  std::ostringstream ossFixed;
2625  // ossFixed <<std::fixed;
2626  ossFixed.setf ( std::ios::fixed, std::ios::floatfield );
2627 
2628  size_type maxBefore=0; // the length of the integral part
2629  size_type maxAfter=0; // number of decimals plus
2630  // one place for the decimal point
2631  for (unsigned int i=0;i<m;++i) {
2632  for (unsigned int j=0;j<n;++j){
2633  oss.str("");
2634  oss << (*this)[i][j];
2635  if (oss.str().find("e")!=std::string::npos){
2636  ossFixed.str("");
2637  ossFixed << (*this)[i][j];
2638  oss.str(ossFixed.str());
2639  }
2640 
2641  values[i*n+j]=oss.str();
2642  size_type thislen=values[i*n+j].size();
2643  size_type p=values[i*n+j].find('.');
2644 
2645  if (p==std::string::npos){
2646  maxBefore=vpMath::maximum(maxBefore, thislen);
2647  // maxAfter remains the same
2648  } else{
2649  maxBefore=vpMath::maximum(maxBefore, p);
2650  maxAfter=vpMath::maximum(maxAfter, thislen-p-1);
2651  }
2652  }
2653  }
2654 
2655  size_type totalLength=length;
2656  // increase totalLength according to maxBefore
2657  totalLength=vpMath::maximum(totalLength,maxBefore);
2658  // decrease maxAfter according to totalLength
2659  maxAfter=std::min(maxAfter, totalLength-maxBefore);
2660  if (maxAfter==1) maxAfter=0;
2661 
2662  // the following line is useful for debugging
2663  std::cerr <<totalLength <<" " <<maxBefore <<" " <<maxAfter <<"\n";
2664 
2665  if (intro) s <<intro;
2666  s <<"["<<m<<","<<n<<"]=\n";
2667 
2668  for (unsigned int i=0;i<m;i++) {
2669  s <<" ";
2670  for (unsigned int j=0;j<n;j++){
2671  size_type p=values[i*n+j].find('.');
2672  s.setf(std::ios::right, std::ios::adjustfield);
2673  s.width((std::streamsize)maxBefore);
2674  s <<values[i*n+j].substr(0,p).c_str();
2675 
2676  if (maxAfter>0){
2677  s.setf(std::ios::left, std::ios::adjustfield);
2678  if (p!=std::string::npos){
2679  s.width((std::streamsize)maxAfter);
2680  s <<values[i*n+j].substr(p,maxAfter).c_str();
2681  } else{
2682  assert(maxAfter>1);
2683  s.width((std::streamsize)maxAfter);
2684  s <<".0";
2685  }
2686  }
2687 
2688  s <<' ';
2689  }
2690  s <<std::endl;
2691  }
2692 
2693  return (int)(maxBefore+maxAfter);
2694 }
2695 
2696 
2705 std::ostream & vpMatrix::
2706 matlabPrint(std::ostream & os)
2707 {
2708 
2709  unsigned int i,j;
2710 
2711  os << "[ ";
2712  for (i=0; i < this->getRows(); ++ i)
2713  {
2714  for (j=0; j < this ->getCols(); ++ j)
2715  {
2716  os << (*this)[i][j] << ", ";
2717  }
2718  if (this ->getRows() != i+1) { os << ";" << std::endl; }
2719  else { os << "]" << std::endl; }
2720  }
2721  return os;
2722 };
2723 
2734 std::ostream & vpMatrix::
2735 maplePrint(std::ostream & os)
2736 {
2737  unsigned int i,j;
2738 
2739  os << "([ " << std::endl;
2740  for (i=0; i < this->getRows(); ++ i)
2741  { os << "[";
2742  for (j=0; j < this->getCols(); ++ j)
2743  {
2744  os << (*this)[i][j] << ", ";
2745  }
2746  os << "]," << std::endl;
2747  }
2748  os << "])" << std::endl;
2749  return os;
2750 };
2751 
2752 
2767 std::ostream & vpMatrix::
2768 cppPrint(std::ostream & os, const char * matrixName, bool octet)
2769 {
2770 
2771  unsigned int i,j;
2772  const char defaultName [] = "A";
2773  if (NULL == matrixName)
2774  {
2775  matrixName = defaultName;
2776  }
2777  os << "vpMatrix " << defaultName
2778  << " (" << this ->getRows ()
2779  << ", " << this ->getCols () << "); " <<std::endl;
2780 
2781  for (i=0; i < this->getRows(); ++ i)
2782  {
2783  for (j=0; j < this ->getCols(); ++ j)
2784  {
2785  if (! octet)
2786  {
2787  os << defaultName << "[" << i << "][" << j
2788  << "] = " << (*this)[i][j] << "; " << std::endl;
2789  }
2790  else
2791  {
2792  for (unsigned int k = 0; k < sizeof(double); ++ k)
2793  {
2794  os << "((unsigned char*)&(" << defaultName
2795  << "[" << i << "][" << j << "]) )[" << k
2796  <<"] = 0x" <<std::hex<<
2797  (unsigned int)((unsigned char*)& ((*this)[i][j])) [k]
2798  << "; " << std::endl;
2799  }
2800  }
2801  }
2802  os << std::endl;
2803  }
2804  return os;
2805 };
2806 
2807 
2815 double
2817 {
2818  double norm=0.0;
2819  double x ;
2820  for (unsigned int i=0;i<dsize;i++) {
2821  x = *(data +i); norm += x*x;
2822  }
2823 
2824  return sqrt(norm);
2825 }
2826 
2827 
2828 
2839 double
2841 {
2842  double norm=0.0;
2843  double x ;
2844  for (unsigned int i=0;i<rowNum;i++){
2845  x = 0;
2846  for (unsigned int j=0; j<colNum;j++){
2847  x += fabs (*(*(rowPtrs + i)+j)) ;
2848  }
2849  if (x > norm) {
2850  norm = x;
2851  }
2852  }
2853  return norm;
2854 }
2855 
2864 double vpMatrix::detByLU() const
2865 {
2866  double det(0);
2867 
2868  // Test wether the matrix is squred
2869  if (rowNum == colNum)
2870  {
2871  // create a temporary matrix that will be modified by LUDcmp
2872  vpMatrix tmp(*this);
2873 
2874  // using th LUdcmp based on NR codes
2875  // it modified the tmp matrix in a special structure of type :
2876  // b11 b12 b13 b14
2877  // a21 b22 b23 b24
2878  // a21 a32 b33 b34
2879  // a31 a42 a43 b44
2880 
2881  unsigned int * perm = new unsigned int[rowNum]; // stores the permutations
2882  int d; // +- 1 fi the number of column interchange is even or odd
2883  tmp.LUDcmp(perm, d);
2884  delete[]perm;
2885 
2886  // compute the determinant that is the product of the eigen values
2887  det = (double) d;
2888  for(unsigned int i=0;i<rowNum;i++)
2889  {
2890  det*=tmp[i][i];
2891  }
2892  }
2893 
2894  else {
2895  vpERROR_TRACE("Determinant Computation : ERR Matrix not squared") ;
2897  "\n\t\tDeterminant Computation : ERR Matrix not squared")) ;
2898 
2899 
2900  }
2901  return det ;
2902 }
2903 
2904 
2905 
2921 {
2922  if(rowNum == 0)
2923  *this = A;
2924  else
2925  *this = vpMatrix::stackMatrices(*this, A);
2926 }
2927 
2928 
2939 void vpMatrix::insert(const vpMatrix&A, const unsigned int r,
2940  const unsigned int c)
2941 {
2942  if( (r + A.getRows() ) <= rowNum && (c + A.getCols() ) <= colNum ){
2943  // recopy matrix A in the current one, does not call static function to avoid initialisation and recopy of matrix
2944  for(unsigned int i=r; i<(r+A.getRows()); i++){
2945  for(unsigned int j=c; j<(c+A.getCols()); j++){
2946  (*this)[i][j] = A[i-r][j-c];
2947  }
2948  }
2949  }
2950  else{
2952  "\n\t\tIncorrect size of the matrix to insert data.");
2953  }
2954 }
2955 
2956 
2997 {
2998  if (rowNum != colNum) {
2999  vpERROR_TRACE("Eigen values computation: ERR Matrix not square") ;
3001  "\n\t\tEigen values computation: ERR Matrix not square")) ;
3002  }
3003 
3004 #ifdef VISP_HAVE_GSL /* be careful of the copy below */
3005  {
3006  // Check if the matrix is symetric: At - A = 0
3007  vpMatrix At_A = (*this).t() - (*this);
3008  for (unsigned int i=0; i < rowNum; i++) {
3009  for (unsigned int j=0; j < rowNum; j++) {
3010  //if (At_A[i][j] != 0) {
3011  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
3012  vpERROR_TRACE("Eigen values computation: ERR Matrix not symmetric") ;
3014  "\n\t\tEigen values computation: "
3015  "ERR Matrix not symmetric")) ;
3016  }
3017  }
3018  }
3019 
3020 
3021  vpColVector evalue(rowNum); // Eigen values
3022 
3023  gsl_vector *eval = gsl_vector_alloc (rowNum);
3024  gsl_matrix *evec = gsl_matrix_alloc (rowNum, colNum);
3025 
3026  gsl_eigen_symmv_workspace * w = gsl_eigen_symmv_alloc (rowNum);
3027  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
3028 
3029  unsigned int Atda = m->tda ;
3030  for (unsigned int i=0 ; i < rowNum ; i++){
3031  unsigned int k = i*Atda ;
3032  for (unsigned int j=0 ; j < colNum ; j++)
3033  m->data[k+j] = (*this)[i][j] ;
3034  }
3035  gsl_eigen_symmv (m, eval, evec, w);
3036 
3037  gsl_eigen_symmv_sort (eval, evec, GSL_EIGEN_SORT_ABS_ASC);
3038 
3039  for (unsigned int i=0; i < rowNum; i++) {
3040  evalue[i] = gsl_vector_get (eval, i);
3041  }
3042 
3043  gsl_eigen_symmv_free (w);
3044  gsl_vector_free (eval);
3045  gsl_matrix_free (m);
3046  gsl_matrix_free (evec);
3047 
3048  return evalue;
3049  }
3050 #else
3051  {
3052  vpERROR_TRACE("Not implemented since the GSL library is not detected.") ;
3054  "\n\t\tEigen values Computation: Not implemented "
3055  "since the GSL library is not detected")) ;
3056  }
3057 #endif
3058 }
3059 
3114 #ifdef VISP_HAVE_GSL /* be careful of the copy below */
3116 #else
3117 void vpMatrix::eigenValues(vpColVector & /* evalue */, vpMatrix & /* evector */)
3118 #endif
3119 {
3120  if (rowNum != colNum) {
3121  vpERROR_TRACE("Eigen values computation: ERR Matrix not square") ;
3123  "\n\t\tEigen values computation: ERR Matrix not square")) ;
3124  }
3125 
3126 #ifdef VISP_HAVE_GSL /* be careful of the copy below */
3127  {
3128  // Check if the matrix is symetric: At - A = 0
3129  vpMatrix At_A = (*this).t() - (*this);
3130  for (unsigned int i=0; i < rowNum; i++) {
3131  for (unsigned int j=0; j < rowNum; j++) {
3132  //if (At_A[i][j] != 0) {
3133  if (std::fabs(At_A[i][j]) > std::numeric_limits<double>::epsilon()) {
3134  vpERROR_TRACE("Eigen values computation: ERR Matrix not symmetric") ;
3136  "\n\t\tEigen values computation: "
3137  "ERR Matrix not symmetric")) ;
3138  }
3139  }
3140  }
3141 
3142  // Resize the output matrices
3143  evalue.resize(rowNum);
3144  evector.resize(rowNum, colNum);
3145 
3146  gsl_vector *eval = gsl_vector_alloc (rowNum);
3147  gsl_matrix *evec = gsl_matrix_alloc (rowNum, colNum);
3148 
3149  gsl_eigen_symmv_workspace * w = gsl_eigen_symmv_alloc (rowNum);
3150  gsl_matrix *m = gsl_matrix_alloc(rowNum, colNum);
3151 
3152  unsigned int Atda = m->tda ;
3153  for (unsigned int i=0 ; i < rowNum ; i++){
3154  unsigned int k = i*Atda ;
3155  for (unsigned int j=0 ; j < colNum ; j++)
3156  m->data[k+j] = (*this)[i][j] ;
3157  }
3158  gsl_eigen_symmv (m, eval, evec, w);
3159 
3160  gsl_eigen_symmv_sort (eval, evec, GSL_EIGEN_SORT_ABS_ASC);
3161 
3162  for (unsigned int i=0; i < rowNum; i++) {
3163  evalue[i] = gsl_vector_get (eval, i);
3164  }
3165  Atda = evec->tda ;
3166  for (unsigned int i=0; i < rowNum; i++) {
3167  unsigned int k = i*Atda ;
3168  for (unsigned int j=0; j < rowNum; j++) {
3169  evector[i][j] = evec->data[k+j];
3170  }
3171  }
3172 
3173 
3174  // {
3175  // int i;
3176 
3177  // for (i = 0; i < rowNum; i++)
3178  // {
3179  // double eval_i
3180  // = gsl_vector_get (eval, i);
3181  // gsl_vector_view evec_i
3182  // = gsl_matrix_column (evec, i);
3183 
3184  // printf ("eigenvalue = %g\n", eval_i);
3185  // printf ("eigenvector = \n");
3186  // gsl_vector_fprintf (stdout,
3187  // &evec_i.vector, "%g");
3188  // }
3189  // }
3190 
3191  gsl_eigen_symmv_free (w);
3192  gsl_vector_free (eval);
3193  gsl_matrix_free (m);
3194  gsl_matrix_free (evec);
3195  }
3196 #else
3197  {
3198  vpERROR_TRACE("Not implemented since the GSL library is not detected.") ;
3200  "\n\t\tEigen values Computation: Not implemented "
3201  "since the GSL library is not detected")) ;
3202  }
3203 #endif
3204 }
3205 
3206 
3216 unsigned int
3217 vpMatrix::kernel(vpMatrix &kerA, double svThreshold)
3218 {
3219 #if (DEBUG_LEVEL1)
3220  std::cout << "begin Kernel" << std::endl ;
3221 #endif
3222  unsigned int i, j ;
3223  unsigned int ncaptor = getRows() ;
3224  unsigned int ddl = getCols() ;
3225  vpMatrix C ;
3226 
3227  if (ncaptor == 0)
3228  std::cout << "Erreur Matrice non initialise" << std::endl ;
3229 
3230 #if (DEBUG_LEVEL2)
3231  {
3232  std::cout << "Interaction matrix L" << std::endl ;
3233  std::cout << *this ;
3234  std::cout << "signaux capteurs : " << ncaptor << std::endl ;
3235  }
3236 #endif
3237 
3238  C.resize(ddl,ncaptor) ;
3239  unsigned int min ;
3240 
3241  if (ncaptor > ddl) min = ddl ; else min = ncaptor ;
3242 
3243  // ! the SVDcmp function inthe matrix lib is destructive
3244 
3245  vpMatrix a1 ;
3246  vpMatrix a2 ;
3247 
3248  vpColVector sv(ddl) ; // singular values
3249  vpMatrix v(ddl,ddl) ;
3250 
3251  if (ncaptor < ddl)
3252  {
3253  a1.resize(ddl,ddl) ;
3254  }
3255  else
3256  {
3257  a1.resize(ncaptor,ddl) ;
3258  }
3259 
3260  a2.resize(ncaptor,ddl) ;
3261 
3262  for (i=0 ; i < ncaptor ; i++)
3263  for (j=0 ; j < ddl ; j++)
3264  {
3265  a1[i][j] = (*this)[i][j] ;
3266  a2[i][j] = (*this)[i][j] ;
3267  }
3268 
3269 
3270  if (ncaptor < ddl)
3271  {
3272  for (i=ncaptor ; i < ddl ; i++)
3273  for (j=0 ; j < ddl ; j++)
3274  {
3275  a1[i][j] = 0 ;
3276  }
3277  a1.svd(sv,v);
3278  }
3279  else
3280  {
3281  a1.svd(sv,v);
3282  }
3283 
3284  // compute the highest singular value and the rank of h
3285 
3286  double maxsv = 0 ;
3287  for (i=0 ; i < ddl ; i++)
3288  if (fabs(sv[i]) > maxsv) maxsv = fabs(sv[i]) ;
3289 
3290 #if (DEBUG_LEVEL2)
3291  {
3292  std::cout << "Singular Value : (" ;
3293  for (i=0 ; i < ddl ; i++) std::cout << sv[i] << " , " ;
3294  std::cout << ")" << std::endl ;
3295  }
3296 #endif
3297 
3298  unsigned int rank = 0 ;
3299  for (i=0 ; i < ddl ; i++)
3300  if (fabs(sv[i]) > maxsv*svThreshold) rank++ ;
3301 
3302 #if (DEBUG_LEVEL2)
3303  {
3304 
3305  for (i = 0 ; i < ddl ; i++)
3306  {
3307  for (j = 0 ; j < ncaptor ; j++)
3308  {
3309  unsigned int k=0 ;
3310  C[i][j] = 0.0;
3311 
3312  // modif le 25 janvier 1999 0.001 <-- maxsv*1.e-ndof
3313  // sinon on peut observer une perte de range de la matrice
3314  // ( d'ou venait ce 0.001 ??? )
3315  for (k=0 ; k < ddl ; k++) if (fabs(sv[k]) > maxsv*svThreshold)
3316  {
3317  C[i][j] += v[i][k]*a1[j][k]/sv[k];
3318  }
3319  }
3320  }
3321 
3322  // cout << v << endl ;
3323  std::cout << C << std::endl ;
3324  }
3325 #endif
3326  /*------------------------------------------------------- */
3327  if (rank != ddl)
3328  {
3329  // Compute the kernel if wanted
3330  if (min < ddl)
3331  {
3332  vpMatrix ch(ddl,ddl) ;
3333  ch = C*a2 ;
3334  ch.svd(sv,v) ;
3335 
3336  }
3337  // if (noyau == 1)
3338  {
3339  double maxsv = 0 ;
3340  for (i=0 ; i < ddl ; i++)
3341  if (fabs(sv[i]) > maxsv) maxsv = fabs(sv[i]) ;
3342  unsigned int rank = 0 ;
3343  for (i=0 ; i < ddl ; i++)
3344  if (fabs(sv[i]) > maxsv*svThreshold) rank++ ;
3345  vpMatrix cons(ddl,ddl) ;
3346 
3347  cons =0 ;
3348  for (j = 0 ; j < ddl ; j++)
3349  {
3350  for (i = 0 ; i < ddl ; i++)
3351  // Change Nicolas Mansard 23/4/04
3352  // was if (fabs(sv[i]) < maxsv*seuilvp)
3353  if (fabs(sv[i]) <= maxsv*svThreshold)
3354  {
3355  cons[i][j] = v[j][i];
3356  }
3357  }
3358 
3359  vpMatrix Ker(ddl-rank,ddl) ;
3360  unsigned int k =0 ;
3361  for (j = 0 ; j < ddl ; j++)
3362  {
3363  // cout << cons.Row(j+1) << " = " << cons.Row(j+1).SumSquare() << endl ;
3364 
3365  //if (cons.row(j+1).sumSquare() !=0)
3366  if (std::fabs(cons.row(j+1).sumSquare()) > std::numeric_limits<double>::epsilon())
3367  {
3368  for (i=0 ; i < cons.getCols() ; i++)
3369  Ker[k][i] = cons[j][i] ;
3370  // Ker.Row(k+1) = cons.Row(j+1) ;
3371  k++;
3372  }
3373  }
3374  kerA = Ker ;
3375  }
3376  }
3377 #if (DEBUG_LEVEL1)
3378  std::cout << "end Kernel" << std::endl ;
3379 #endif
3380  return rank ;
3381 }
3382 
3413 double vpMatrix::det(vpDetMethod method) const
3414 {
3415  double det = 0;
3416 
3417  if ( method == LU_DECOMPOSITION )
3418  {
3419  det = this->detByLU();
3420  }
3421 
3422  return (det);
3423 }
3424 
3425 
3439 bool
3440 vpMatrix::saveMatrix(const char *filename, const vpMatrix &M,
3441  const bool binary, const char *Header)
3442 {
3443  std::fstream file;
3444 
3445  if (!binary)
3446  file.open(filename, std::fstream::out);
3447  else
3448  file.open(filename, std::fstream::out|std::fstream::binary);
3449 
3450  if(!file)
3451  {
3452  file.close();
3453  return false;
3454  }
3455 
3456  else
3457  {
3458  if (!binary)
3459  {
3460  unsigned int i = 0;
3461  file << "# ";
3462  while (Header[i] != '\0')
3463  {
3464  file << Header[i];
3465  if (Header[i] == '\n')
3466  file << "# ";
3467  i++;
3468  }
3469  file << '\0';
3470  file << std::endl;
3471  file << M.getRows() << "\t" << M.getCols() << std::endl;
3472  file << M << std::endl;
3473  }
3474  else
3475  {
3476  int headerSize = 0;
3477  while (Header[headerSize] != '\0') headerSize++;
3478  file.write(Header,headerSize+1);
3479  unsigned int matrixSize;
3480  matrixSize = M.getRows();
3481  file.write((char*)&matrixSize,sizeof(int));
3482  matrixSize = M.getCols();
3483  file.write((char*)&matrixSize,sizeof(int));
3484  double value;
3485  for(unsigned int i = 0; i < M.getRows(); i++)
3486  {
3487  for(unsigned int j = 0; j < M.getCols(); j++)
3488  {
3489  value = M[i][j];
3490  file.write((char*)&value,sizeof(double));
3491  }
3492  }
3493  }
3494  }
3495 
3496  file.close();
3497  return true;
3498 }
3499 
3500 
3512 bool
3513 vpMatrix::loadMatrix(const char *filename, vpMatrix &M, const bool binary,
3514  char *Header)
3515 {
3516  std::fstream file;
3517 
3518  if (!binary)
3519  file.open(filename, std::fstream::in);
3520  else
3521  file.open(filename, std::fstream::in|std::fstream::binary);
3522 
3523  if(!file)
3524  {
3525  file.close();
3526  return false;
3527  }
3528 
3529  else
3530  {
3531  if (!binary)
3532  {
3533  char c='0';
3534  std::string h;
3535  while ((c != '\0') && (c != '\n'))
3536  {
3537  file.read(&c,1);
3538  h+=c;
3539  }
3540  if (Header != NULL)
3541  strncpy(Header, h.c_str(), h.size() + 1);
3542 
3543  unsigned int rows, cols;
3544  file >> rows;
3545  file >> cols;
3546  M.resize(rows,cols);
3547 
3548  double value;
3549  for(unsigned int i = 0; i < rows; i++)
3550  {
3551  for(unsigned int j = 0; j < cols; j++)
3552  {
3553  file >> value;
3554  M[i][j] = value;
3555  }
3556  }
3557  }
3558  else
3559  {
3560  char c='0';
3561  std::string h;
3562  while ((c != '\0') && (c != '\n'))
3563  {
3564  file.read(&c,1);
3565  h+=c;
3566  }
3567  if (Header != NULL)
3568  strncpy(Header, h.c_str(), h.size() + 1);
3569 
3570  unsigned int rows, cols;
3571  file.read((char*)&rows,sizeof(unsigned int));
3572  file.read((char*)&cols,sizeof(unsigned int));
3573  M.resize(rows,cols);
3574 
3575  double value;
3576  for(unsigned int i = 0; i < rows; i++)
3577  {
3578  for(unsigned int j = 0; j < cols; j++)
3579  {
3580  file.read((char*)&value,sizeof(double));
3581  M[i][j] = value;
3582  }
3583  }
3584  }
3585  }
3586 
3587  file.close();
3588  return true;
3589 }
3590 
3591 
3599 vpMatrix
3601 {
3602  if(colNum != rowNum)
3603  {
3604  vpTRACE("The matrix must be square");
3606  "The matrix must be square" ));
3607  }
3608  else
3609  {
3610  vpMatrix _expE(rowNum, colNum);
3611  vpMatrix _expD(rowNum, colNum);
3612  vpMatrix _expX(rowNum, colNum);
3613  vpMatrix _expcX(rowNum, colNum);
3614  vpMatrix _eye(rowNum, colNum);
3615 
3616  _eye.setIdentity();
3617  vpMatrix exp(*this);
3618 
3619  // double f;
3620  int e;
3621  double c = 0.5;
3622  int q = 6;
3623  int p = 1;
3624 
3625  double nA = 0;
3626  for (unsigned int i = 0; i < rowNum;i++)
3627  {
3628  double sum = 0;
3629  for (unsigned int j=0; j < colNum; j++)
3630  {
3631  sum += fabs((*this)[i][j]);
3632  }
3633  if (sum>nA||i==0)
3634  {
3635  nA=sum;
3636  }
3637  }
3638 
3639  /* f = */ frexp(nA, &e);
3640  //double s = (0 > e+1)?0:e+1;
3641  double s = e+1;
3642 
3643  double sca = 1.0 / pow(2.0,s);
3644  exp=sca*exp;
3645  _expX=*this;
3646  _expE=c*exp+_eye;
3647  _expD=-c*exp+_eye;
3648  for (int k=2;k<=q;k++)
3649  {
3650  c = c * ((double)(q-k+1)) / ((double)(k*(2*q-k+1)));
3651  _expcX=exp*_expX;
3652  _expX=_expcX;
3653  _expcX=c*_expX;
3654  _expE=_expE+_expcX;
3655  if (p) _expD=_expD+_expcX;
3656  else _expD=_expD- _expcX;
3657  p = !p;
3658  }
3659  _expX=_expD.inverseByLU();
3660  exp=_expX*_expE;
3661  for (int k=1;k<=s;k++)
3662  {
3663  _expE=exp*exp;
3664  exp=_expE;
3665  }
3666  return exp;
3667  }
3668 }
3669 
3670 double
3672 {
3673  double *dataptr = data;
3674  double min = *dataptr;
3675  dataptr++;
3676  for (unsigned int i = 0; i < dsize-1; i++)
3677  {
3678  if (*dataptr < min) min = *dataptr;
3679  dataptr++;
3680  }
3681  return min;
3682 }
3683 
3684 double
3686 {
3687  double *dataptr = data;
3688  double max = *dataptr;
3689  dataptr++;
3690  for (unsigned int i = 0; i < dsize-1; i++)
3691  {
3692  if (*dataptr > max) max = *dataptr;
3693  dataptr++;
3694  }
3695  return max;
3696 }
3697 
3698 
3699 /**************************************************************************************************************/
3700 /**************************************************************************************************************/
3701 
3702 
3703 //Specific functions
3704 
3705 /*
3706 input:: matrix M(nCols,nRows), nCols > 3, nRows > 3 , nCols == nRows.
3707 
3708 output:: the complement matrix of the element (rowNo,colNo).
3709 This is the matrix obtained from M after elimenating the row rowNo and column colNo
3710 
3711 example:
3712 1 2 3
3713 M = 4 5 6
3714 7 8 9
3715 1 3
3716 subblock(M, 1, 1) give the matrix 7 9
3717 */
3718 vpMatrix subblock(const vpMatrix &M, unsigned int col, unsigned int row)
3719 {
3720  vpMatrix M_comp(M.getRows()-1,M.getCols()-1);
3721 
3722  for ( unsigned int i = 0 ; i < col ; i++)
3723  {
3724  for ( unsigned int j = 0 ; j < row ; j++)
3725  M_comp[i][j]=M[i][j];
3726  for ( unsigned int j = row+1 ; j < M.getRows() ; j++)
3727  M_comp[i][j-1]=M[i][j];
3728  }
3729  for ( unsigned int i = col+1 ; i < M.getCols(); i++)
3730  {
3731  for ( unsigned int j = 0 ; j < row ; j++)
3732  M_comp[i-1][j]=M[i][j];
3733  for ( unsigned int j = row+1 ; j < M.getRows() ; j++)
3734  M_comp[i-1][j-1]=M[i][j];
3735  }
3736  return M_comp;
3737 }
3738 
3739 #undef DEBUG_LEVEL1
3740 #undef DEBUG_LEVEL2
3741 
3742 /*
3743 * Local variables:
3744 * c-basic-offset: 2
3745 * End:
3746 */
vpMatrix operator*(const double &x, const vpMatrix &B)
Definition: vpMatrix.cpp:886
friend VISP_EXPORT std::ostream & operator<<(std::ostream &s, const vpMatrix &m)
std::cout a matrix
Definition: vpMatrix.cpp:2576
#define vpDEBUG_TRACE
Definition: vpDebug.h:454
Definition of the vpMatrix class.
Definition: vpMatrix.h:96
unsigned int kernel(vpMatrix &KerA, double svThreshold=1e-6)
Definition: vpMatrix.cpp:3217
static vpMatrix juxtaposeMatrices(const vpMatrix &A, const vpMatrix &B)
Definition: vpMatrix.cpp:2422
void resize(const unsigned int nrows, const unsigned int ncols, const bool nullify=true)
Definition: vpMatrix.cpp:174
vpRowVector row(const unsigned int i)
Row extraction.
Definition: vpMatrix.cpp:2225
double infinityNorm() const
Definition: vpMatrix.cpp:2840
void init()
Initialization of the object matrix.
Definition: vpMatrix.cpp:93
#define vpERROR_TRACE
Definition: vpDebug.h:379
#define vpTRACE
Definition: vpDebug.h:401
Definition of the row vector class.
Definition: vpRowVector.h:73
static bool saveMatrix(const char *filename, const vpMatrix &M, const bool binary=false, const char *Header="")
Definition: vpMatrix.cpp:3440
static void add2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:493
#define vpCERROR
Definition: vpDebug.h:354
double getMinValue() const
Definition: vpMatrix.cpp:3671
vpMatrix & operator*=(const double x)
Multiply all the element of the matrix by x : Aij = Aij * x.
Definition: vpMatrix.cpp:1058
vpMatrix & operator/=(double x)
Divide all the element of the matrix by x : Aij = Aij / x.
Definition: vpMatrix.cpp:1073
error that can be emited by ViSP classes.
Definition: vpException.h:75
static void sub2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:562
vpMatrix & operator+=(const vpMatrix &B)
Operation A = A + B.
Definition: vpMatrix.cpp:624
double sumSquare() const
return sum of the Aij^2 (for all i, for all j)
Definition: vpMatrix.cpp:760
vpMatrix()
Basic constructor.
Definition: vpMatrix.cpp:111
int print(std::ostream &s, unsigned int length, char const *intro=0)
Definition: vpMatrix.cpp:2615
vpColVector column(const unsigned int j)
Column extraction.
Definition: vpMatrix.cpp:2240
vpMatrix & operator-=(const vpMatrix &B)
Operation A = A - B.
Definition: vpMatrix.cpp:664
static Type maximum(const Type &a, const Type &b)
Definition: vpMath.h:137
double * data
address of the first element of the data array
Definition: vpMatrix.h:116
std::ostream & cppPrint(std::ostream &os, const char *matrixName=NULL, bool octet=false)
Print to be used as part of a C++ code later.
Definition: vpMatrix.cpp:2768
unsigned int trsize
Total row space.
Definition: vpMatrix.h:124
vpColVector eigenValues()
Definition: vpMatrix.cpp:2996
void insert(const vpMatrix &A, const unsigned int r, const unsigned int c)
Definition: vpMatrix.cpp:2939
void solveBySVD(const vpColVector &B, vpColVector &x) const
Definition: vpMatrix.cpp:1573
void setIdentity(const double &val=1.0)
Definition: vpMatrix.cpp:1110
double ** rowPtrs
address of the first element of each rows
Definition: vpMatrix.h:119
void svd(vpColVector &w, vpMatrix &v)
Definition: vpMatrix.cpp:1702
vpMatrix & operator=(const vpMatrix &B)
Copy operator. Allow operation such as A = B.
Definition: vpMatrix.cpp:316
vpMatrix AtA() const
Definition: vpMatrix.cpp:1359
static void multMatrixVector(const vpMatrix &A, const vpColVector &b, vpColVector &c)
Definition: vpMatrix.cpp:811
vpMatrix AAt() const
Definition: vpMatrix.cpp:1248
static bool loadMatrix(const char *filename, vpMatrix &M, const bool binary=false, char *Header=NULL)
Definition: vpMatrix.cpp:3513
static void add2WeightedMatrices(const vpMatrix &A, const double &wA, const vpMatrix &B, const double &wB, vpMatrix &C)
Definition: vpMatrix.cpp:454
static void mult2Matrices(const vpMatrix &A, const vpMatrix &B, vpMatrix &C)
Definition: vpMatrix.cpp:393
vpMatrix operator+(const vpMatrix &B) const
Definition: vpMatrix.cpp:544
vpMatrix transpose() const
Definition: vpMatrix.cpp:1206
void kill()
Destruction of the matrix (Memory de-allocation)
Definition: vpMatrix.cpp:286
#define vpCDEBUG(niv)
Definition: vpDebug.h:478
double det(vpDetMethod method=LU_DECOMPOSITION) const
Definition: vpMatrix.cpp:3413
VISP_EXPORT std::ostream & operator<<(std::ostream &os, const vpImagePoint &ip)
Definition: vpImagePoint.h:529
static vpMatrix stackMatrices(const vpMatrix &A, const vpMatrix &B)
Stack two Matrices C = [ A B ]^T.
Definition: vpMatrix.cpp:2263
void diag(const vpColVector &A)
Definition: vpMatrix.cpp:2523
unsigned int rowNum
number of rows
Definition: vpMatrix.h:110
void eye(unsigned int n)
Definition: vpMatrix.cpp:1132
vpMatrix t() const
Definition: vpMatrix.cpp:1176
double getMaxValue() const
Definition: vpMatrix.cpp:3685
void resize(unsigned int i)
Set the size of the Row vector.
Definition: vpRowVector.h:91
Class that provides a data structure for the column vectors as well as a set of operations on these v...
Definition: vpColVector.h:72
double euclideanNorm() const
Definition: vpMatrix.cpp:2816
vpMatrix inverseByLU() const
unsigned int getCols() const
Return the number of columns of the matrix.
Definition: vpMatrix.h:159
vpRowVector stackRows()
Definition: vpMatrix.cpp:1430
vpMatrix operator/(const double x) const
Cij = Aij / x (A is unchanged)
Definition: vpMatrix.cpp:965
static void createDiagonalMatrix(const vpColVector &A, vpMatrix &DA)
Definition: vpMatrix.cpp:2551
vpMatrix operator*(const vpMatrix &B) const
Definition: vpMatrix.cpp:436
error that can be emited by the vpMatrix class and its derivates
vpColVector stackColumns()
Definition: vpMatrix.cpp:1397
vpMatrix operator-() const
Definition: vpMatrix.cpp:751
std::ostream & maplePrint(std::ostream &os)
Print using MAPLE matrix input format.
Definition: vpMatrix.cpp:2735
vpMatrix pseudoInverse(double svThreshold=1e-6) const
Compute the pseudo inverse of the matrix using the SVD.
Definition: vpMatrix.cpp:1812
unsigned int dsize
Current size (rowNum * colNum)
Definition: vpMatrix.h:122
void kron(const vpMatrix &m1, vpMatrix &out)
Definition: vpMatrix.cpp:1477
unsigned int colNum
number of columns
Definition: vpMatrix.h:112
unsigned int getRows() const
Return the number of rows of the matrix.
Definition: vpMatrix.h:157
static void negateMatrix(const vpMatrix &A, vpMatrix &C)
Definition: vpMatrix.cpp:709
vpMatrix expm()
Definition: vpMatrix.cpp:3600
Class that consider the case of a translation vector.
virtual ~vpMatrix()
Destructor (Memory de-allocation)
Definition: vpMatrix.cpp:303
std::ostream & matlabPrint(std::ostream &os)
Print using matlab syntax, to be put in matlab later.
Definition: vpMatrix.cpp:2706
void resize(const unsigned int i, const bool flagNullify=true)
Definition: vpColVector.h:94