ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TDecompSVD.cxx
Go to the documentation of this file.
1 // @(#)root/matrix:$Id$
2 // Authors: Fons Rademakers, Eddy Offermann Dec 2003
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 ///////////////////////////////////////////////////////////////////////////
13 // //
14 // Single Value Decomposition class //
15 // //
16 // For an (m x n) matrix A with m >= n, the singular value decomposition //
17 // is //
18 // an (m x m) orthogonal matrix fU, an (m x n) diagonal matrix fS, and //
19 // an (n x n) orthogonal matrix fV so that A = U*S*V'. //
20 // //
21 // If the row/column index of A starts at (rowLwb,colLwb) then the //
22 // decomposed matrices/vectors start at : //
23 // fU : (rowLwb,colLwb) //
24 // fV : (colLwb,colLwb) //
25 // fSig : (colLwb) //
26 // //
27 // The diagonal matrix fS is stored in the singular values vector fSig . //
28 // The singular values, fSig[k] = S[k][k], are ordered so that //
29 // fSig[0] >= fSig[1] >= ... >= fSig[n-1]. //
30 // //
31 // The singular value decompostion always exists, so the decomposition //
32 // will (as long as m >=n) never fail. If m < n, the user should add //
33 // sufficient zero rows to A , so that m == n //
34 // //
35 // Here fTol is used to set the threshold on the minimum allowed value //
36 // of the singular values: //
37 // min_singular = fTol*max(fSig[i]) //
38 // //
39 ///////////////////////////////////////////////////////////////////////////
40 
41 #include "TDecompSVD.h"
42 #include "TMath.h"
43 #include "TArrayD.h"
44 
46 
47 ////////////////////////////////////////////////////////////////////////////////
48 /// Constructor for (nrows x ncols) matrix
49 
50 TDecompSVD::TDecompSVD(Int_t nrows,Int_t ncols)
51 {
52  if (nrows < ncols) {
53  Error("TDecompSVD(Int_t,Int_t","matrix rows should be >= columns");
54  return;
55  }
56  fU.ResizeTo(nrows,nrows);
57  fSig.ResizeTo(ncols);
58  fV.ResizeTo(nrows,ncols); // In the end we only need the nColxnCol part
59 }
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 /// Constructor for ([row_lwb..row_upb] x [col_lwb..col_upb]) matrix
63 
64 TDecompSVD::TDecompSVD(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb)
65 {
66  const Int_t nrows = row_upb-row_lwb+1;
67  const Int_t ncols = col_upb-col_lwb+1;
68 
69  if (nrows < ncols) {
70  Error("TDecompSVD(Int_t,Int_t,Int_t,Int_t","matrix rows should be >= columns");
71  return;
72  }
73  fRowLwb = row_lwb;
74  fColLwb = col_lwb;
75  fU.ResizeTo(nrows,nrows);
76  fSig.ResizeTo(ncols);
77  fV.ResizeTo(nrows,ncols); // In the end we only need the nColxnCol part
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Constructor for general matrix A .
82 
84 {
85  R__ASSERT(a.IsValid());
86  if (a.GetNrows() < a.GetNcols()) {
87  Error("TDecompSVD(const TMatrixD &","matrix rows should be >= columns");
88  return;
89  }
90 
92  fTol = a.GetTol();
93  if (tol > 0)
94  fTol = tol;
95 
96  fRowLwb = a.GetRowLwb();
97  fColLwb = a.GetColLwb();
98  const Int_t nRow = a.GetNrows();
99  const Int_t nCol = a.GetNcols();
100 
101  fU.ResizeTo(nRow,nRow);
102  fSig.ResizeTo(nCol);
103  fV.ResizeTo(nRow,nCol); // In the end we only need the nColxnCol part
104 
105  fU.UnitMatrix();
106  memcpy(fV.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
107 }
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Copy constructor
111 
113 {
114  *this = another;
115 }
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// SVD decomposition of matrix
119 /// If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
120 
122 {
123  if (TestBit(kDecomposed)) return kTRUE;
124 
125  if ( !TestBit(kMatrixSet) ) {
126  Error("Decompose()","Matrix has not been set");
127  return kFALSE;
128  }
129 
130  const Int_t nCol = this->GetNcols();
131  const Int_t rowLwb = this->GetRowLwb();
132  const Int_t colLwb = this->GetColLwb();
133 
134  TVectorD offDiag;
135  Double_t work[kWorkMax];
136  if (nCol > kWorkMax) offDiag.ResizeTo(nCol);
137  else offDiag.Use(nCol,work);
138 
139  // step 1: bidiagonalization of A
140  if (!Bidiagonalize(fV,fU,fSig,offDiag))
141  return kFALSE;
142 
143  // step 2: diagonalization of bidiagonal matrix
144  if (!Diagonalize(fV,fU,fSig,offDiag))
145  return kFALSE;
146 
147  // step 3: order singular values and perform permutations
149  fV.ResizeTo(nCol,nCol); fV.Shift(colLwb,colLwb);
150  fSig.Shift(colLwb);
151  fU.Transpose(fU); fU.Shift(rowLwb,colLwb);
153 
154  return kTRUE;
155 }
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Bidiagonalize the (m x n) - matrix a (stored in v) through a series of Householder
159 /// transformations applied to the left (Q^T) and to the right (H) of a ,
160 /// so that A = Q . C . H^T with matrix C bidiagonal. Q and H are orthogonal matrices .
161 ///
162 /// Output:
163 /// v - (n x n) - matrix H in the (n x n) part of v
164 /// u - (m x m) - matrix Q^T
165 /// sDiag - diagonal of the (m x n) C
166 /// oDiag - off-diagonal elements of matrix C
167 ///
168 /// Test code for the output:
169 /// const Int_t nRow = v.GetNrows();
170 /// const Int_t nCol = v.GetNcols();
171 /// TMatrixD H(v); H.ResizeTo(nCol,nCol);
172 /// TMatrixD E1(nCol,nCol); E1.UnitMatrix();
173 /// TMatrixD Ht(TMatrixDBase::kTransposed,H);
174 /// Bool_t ok = kTRUE;
175 /// ok &= VerifyMatrixIdentity(Ht * H,E1,kTRUE,1.0e-13);
176 /// ok &= VerifyMatrixIdentity(H * Ht,E1,kTRUE,1.0e-13);
177 /// TMatrixD E2(nRow,nRow); E2.UnitMatrix();
178 /// TMatrixD Qt(u);
179 /// TMatrixD Q(TMatrixDBase::kTransposed,Qt);
180 /// ok &= VerifyMatrixIdentity(Q * Qt,E2,kTRUE,1.0e-13);
181 /// TMatrixD C(nRow,nCol);
182 /// TMatrixDDiag(C) = sDiag;
183 /// for (Int_t i = 0; i < nCol-1; i++)
184 /// C(i,i+1) = oDiag(i+1);
185 /// TMatrixD A = Q*C*Ht;
186 /// ok &= VerifyMatrixIdentity(A,a,kTRUE,1.0e-13);
187 
189 {
190  const Int_t nRow_v = v.GetNrows();
191  const Int_t nCol_v = v.GetNcols();
192  const Int_t nCol_u = u.GetNcols();
193 
194  TArrayD ups(nCol_v);
195  TArrayD betas(nCol_v);
196 
197  for (Int_t i = 0; i < nCol_v; i++) {
198  // Set up Householder Transformation q(i)
199 
200  Double_t up,beta;
201  if (i < nCol_v-1 || nRow_v > nCol_v) {
202  const TVectorD vc_i = TMatrixDColumn_const(v,i);
203  //if (!DefHouseHolder(vc_i,i,i+1,up,beta))
204  // return kFALSE;
205  DefHouseHolder(vc_i,i,i+1,up,beta);
206 
207  // Apply q(i) to v
208  for (Int_t j = i; j < nCol_v; j++) {
209  TMatrixDColumn vc_j = TMatrixDColumn(v,j);
210  ApplyHouseHolder(vc_i,up,beta,i,i+1,vc_j);
211  }
212 
213  // Apply q(i) to u
214  for (Int_t j = 0; j < nCol_u; j++)
215  {
216  TMatrixDColumn uc_j = TMatrixDColumn(u,j);
217  ApplyHouseHolder(vc_i,up,beta,i,i+1,uc_j);
218  }
219  }
220  if (i < nCol_v-2) {
221  // set up Householder Transformation h(i)
222  const TVectorD vr_i = TMatrixDRow_const(v,i);
223 
224  //if (!DefHouseHolder(vr_i,i+1,i+2,up,beta))
225  // return kFALSE;
226  DefHouseHolder(vr_i,i+1,i+2,up,beta);
227 
228  // save h(i)
229  ups[i] = up;
230  betas[i] = beta;
231 
232  // apply h(i) to v
233  for (Int_t j = i; j < nRow_v; j++) {
234  TMatrixDRow vr_j = TMatrixDRow(v,j);
235  ApplyHouseHolder(vr_i,up,beta,i+1,i+2,vr_j);
236 
237  // save elements i+2,...in row j of matrix v
238  if (j == i) {
239  for (Int_t k = i+2; k < nCol_v; k++)
240  vr_j(k) = vr_i(k);
241  }
242  }
243  }
244  }
245 
246  // copy diagonal of transformed matrix v to sDiag and upper parallel v to oDiag
247  if (nCol_v > 1) {
248  for (Int_t i = 1; i < nCol_v; i++)
249  oDiag(i) = v(i-1,i);
250  }
251  oDiag(0) = 0.;
252  sDiag = TMatrixDDiag(v);
253 
254  // construct product matrix h = h(1)*h(2)*...*h(nCol_v-1), h(nCol_v-1) = I
255 
256  TVectorD vr_i(nCol_v);
257  for (Int_t i = nCol_v-1; i >= 0; i--) {
258  if (i < nCol_v-1)
259  vr_i = TMatrixDRow_const(v,i);
260  TMatrixDRow(v,i) = 0.0;
261  v(i,i) = 1.;
262 
263  if (i < nCol_v-2) {
264  for (Int_t k = i; k < nCol_v; k++) {
265  // householder transformation on k-th column
266  TMatrixDColumn vc_k = TMatrixDColumn(v,k);
267  ApplyHouseHolder(vr_i,ups[i],betas[i],i+1,i+2,vc_k);
268  }
269  }
270  }
271 
272  return kTRUE;
273 }
274 
275 ////////////////////////////////////////////////////////////////////////////////
276 /// Diagonalizes in an iterative fashion the bidiagonal matrix C as described through
277 /// sDiag and oDiag, so that S' = U'^T . C . V' is diagonal. U' and V' are orthogonal
278 /// matrices .
279 ///
280 /// Output:
281 /// v - (n x n) - matrix H . V' in the (n x n) part of v
282 /// u - (m x m) - matrix U'^T . Q^T
283 /// sDiag - diagonal of the (m x n) S'
284 ///
285 /// return convergence flag: 0 -> no convergence
286 /// 1 -> convergence
287 ///
288 /// Test code for the output:
289 /// const Int_t nRow = v.GetNrows();
290 /// const Int_t nCol = v.GetNcols();
291 /// TMatrixD tmp = v; tmp.ResizeTo(nCol,nCol);
292 /// TMatrixD Vprime = Ht*tmp;
293 /// TMatrixD Vprimet(TMatrixDBase::kTransposed,Vprime);
294 /// TMatrixD Uprimet = u*Q;
295 /// TMatrixD Uprime(TMatrixDBase::kTransposed,Uprimet);
296 /// TMatrixD Sprime(nRow,nCol);
297 /// TMatrixDDiag(Sprime) = sDiag;
298 /// ok &= VerifyMatrixIdentity(Uprimet * C * Vprime,Sprime,kTRUE,1.0e-13);
299 /// ok &= VerifyMatrixIdentity(Q*Uprime * Sprime * Vprimet * Ht,a,kTRUE,1.0e-13);
300 
302 {
303  Bool_t ok = kTRUE;
304  Int_t niter = 0;
305  Double_t bmx = sDiag(0);
306 
307  const Int_t nCol_v = v.GetNcols();
308 
309  if (nCol_v > 1) {
310  for (Int_t i = 1; i < nCol_v; i++)
311  bmx = TMath::Max(TMath::Abs(sDiag(i))+TMath::Abs(oDiag(i)),bmx);
312  }
313 
315 
316  const Int_t niterm = 10*nCol_v;
317  for (Int_t k = nCol_v-1; k >= 0; k--) {
318  loop:
319  if (k != 0) {
320  // since sDiag(k) == 0 perform Givens transform with result oDiag[k] = 0
321  if (TMath::Abs(sDiag(k)) < eps*bmx)
322  Diag_1(v,sDiag,oDiag,k);
323 
324  // find l (1 <= l <=k) so that either oDiag(l) = 0 or sDiag(l-1) = 0.
325  // In the latter case transform oDiag(l) to zero. In both cases the matrix
326  // splits and the bottom right minor begins with row l.
327  // If no such l is found set l = 1.
328 
329  Int_t elzero = 0;
330  Int_t l = 0;
331  for (Int_t ll = k; ll >= 0 ; ll--) {
332  l = ll;
333  if (l == 0) {
334  elzero = 0;
335  break;
336  } else if (TMath::Abs(oDiag(l)) < eps*bmx) {
337  elzero = 1;
338  break;
339  } else if (TMath::Abs(sDiag(l-1)) < eps*bmx)
340  elzero = 0;
341  }
342  if (l > 0 && !elzero)
343  Diag_2(sDiag,oDiag,k,l);
344  if (l != k) {
345  // one more QR pass with order k
346  Diag_3(v,u,sDiag,oDiag,k,l);
347  niter++;
348  if (niter <= niterm) goto loop;
349  ::Error("TDecompSVD::Diagonalize","no convergence after %d steps",niter);
350  ok = kFALSE;
351  }
352  }
353 
354  if (sDiag(k) < 0.) {
355  // for negative singular values perform change of sign
356  sDiag(k) = -sDiag(k);
357  TMatrixDColumn(v,k) *= -1.0;
358  }
359  // order is decreased by one in next pass
360  }
361 
362  return ok;
363 }
364 
365 ////////////////////////////////////////////////////////////////////////////////
366 /// Step 1 in the matrix diagonalization
367 
369 {
370  const Int_t nCol_v = v.GetNcols();
371 
372  TMatrixDColumn vc_k = TMatrixDColumn(v,k);
373  for (Int_t i = k-1; i >= 0; i--) {
374  TMatrixDColumn vc_i = TMatrixDColumn(v,i);
375  Double_t h,cs,sn;
376  if (i == k-1)
377  DefAplGivens(sDiag[i],oDiag[i+1],cs,sn);
378  else
379  DefAplGivens(sDiag[i],h,cs,sn);
380  if (i > 0) {
381  h = 0.;
382  ApplyGivens(oDiag[i],h,cs,sn);
383  }
384  for (Int_t j = 0; j < nCol_v; j++)
385  ApplyGivens(vc_i(j),vc_k(j),cs,sn);
386  }
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Step 2 in the matrix diagonalization
391 
393 {
394  for (Int_t i = l; i <= k; i++) {
395  Double_t h,cs,sn;
396  if (i == l)
397  DefAplGivens(sDiag(i),oDiag(i),cs,sn);
398  else
399  DefAplGivens(sDiag(i),h,cs,sn);
400  if (i < k) {
401  h = 0.;
402  ApplyGivens(oDiag(i+1),h,cs,sn);
403  }
404  }
405 }
406 
407 ////////////////////////////////////////////////////////////////////////////////
408 /// Step 3 in the matrix diagonalization
409 
411 {
412  Double_t *pS = sDiag.GetMatrixArray();
413  Double_t *pO = oDiag.GetMatrixArray();
414 
415  // determine shift parameter
416 
417  const Double_t psk1 = pS[k-1];
418  const Double_t psk = pS[k];
419  const Double_t pok1 = pO[k-1];
420  const Double_t pok = pO[k];
421  const Double_t psl = pS[l];
422 
423  Double_t f;
424  if (psl == 0.0 || pok == 0.0 || psk1 == 0.0) {
425  const Double_t b = ((psk1-psk)*(psk1+psk)+pok1*pok1)/2.0;
426  const Double_t c = (psk*pok1)*(psk*pok1);
427 
428  Double_t shift = 0.0;
429  if ((b != 0.0) | (c != 0.0)) {
430  shift = TMath::Sqrt(b*b+c);
431  if (b < 0.0)
432  shift = -shift;
433  shift = c/(b+shift);
434  }
435 
436  f = (psl+psk)*(psl-psk)+shift;
437  } else {
438  f = ((psk1-psk)*(psk1+psk)+(pok1-pok)*(pok1+pok))/(2.*pok*psk1);
439  const Double_t g = TMath::Hypot(1.,f);
440  const Double_t t = (f >= 0.) ? f+g : f-g;
441 
442  f = ((psl-psk)*(psl+psk)+pok*(psk1/t-pok))/psl;
443  }
444 
445  const Int_t nCol_v = v.GetNcols();
446  const Int_t nCol_u = u.GetNcols();
447 
448  Double_t h,cs,sn;
449  Int_t j;
450  for (Int_t i = l; i <= k-1; i++) {
451  if (i == l)
452  // define r[l]
453  DefGivens(f,pO[i+1],cs,sn);
454  else
455  // define r[i]
456  DefAplGivens(pO[i],h,cs,sn);
457 
458  ApplyGivens(pS[i],pO[i+1],cs,sn);
459  h = 0.;
460  ApplyGivens(h,pS[i+1],cs,sn);
461 
462  TMatrixDColumn vc_i = TMatrixDColumn(v,i);
463  TMatrixDColumn vc_i1 = TMatrixDColumn(v,i+1);
464  for (j = 0; j < nCol_v; j++)
465  ApplyGivens(vc_i(j),vc_i1(j),cs,sn);
466  // define t[i]
467  DefAplGivens(pS[i],h,cs,sn);
468  ApplyGivens(pO[i+1],pS[i+1],cs,sn);
469  if (i < k-1) {
470  h = 0.;
471  ApplyGivens(h,pO[i+2],cs,sn);
472  }
473 
474  TMatrixDRow ur_i = TMatrixDRow(u,i);
475  TMatrixDRow ur_i1 = TMatrixDRow(u,i+1);
476  for (j = 0; j < nCol_u; j++)
477  ApplyGivens(ur_i(j),ur_i1(j),cs,sn);
478  }
479 }
480 
481 ////////////////////////////////////////////////////////////////////////////////
482 /// Perform a permutation transformation on the diagonal matrix S', so that
483 /// matrix S'' = U''^T . S' . V'' has diagonal elements ordered such that they
484 /// do not increase.
485 ///
486 /// Output:
487 /// v - (n x n) - matrix H . V' . V'' in the (n x n) part of v
488 /// u - (m x m) - matrix U''^T . U'^T . Q^T
489 /// sDiag - diagonal of the (m x n) S''
490 
492 {
493  const Int_t nCol_v = v.GetNcols();
494  const Int_t nCol_u = u.GetNcols();
495 
496  Double_t *pS = sDiag.GetMatrixArray();
497  Double_t *pV = v.GetMatrixArray();
498  Double_t *pU = u.GetMatrixArray();
499 
500  // order singular values
501 
502  Int_t i,j;
503  if (nCol_v > 1) {
504  while (1) {
505  Bool_t found = kFALSE;
506  i = 1;
507  while (!found && i < nCol_v) {
508  if (pS[i] > pS[i-1])
509  found = kTRUE;
510  else
511  i++;
512  }
513  if (!found) break;
514  for (i = 1; i < nCol_v; i++) {
515  Double_t t = pS[i-1];
516  Int_t k = i-1;
517  for (j = i; j < nCol_v; j++) {
518  if (t < pS[j]) {
519  t = pS[j];
520  k = j;
521  }
522  }
523  if (k != i-1) {
524  // perform permutation on singular values
525  pS[k] = pS[i-1];
526  pS[i-1] = t;
527  // perform permutation on matrix v
528  for (j = 0; j < nCol_v; j++) {
529  const Int_t off_j = j*nCol_v;
530  t = pV[off_j+k];
531  pV[off_j+k] = pV[off_j+i-1];
532  pV[off_j+i-1] = t;
533  }
534  // perform permutation on vector u
535  for (j = 0; j < nCol_u; j++) {
536  const Int_t off_k = k*nCol_u;
537  const Int_t off_i1 = (i-1)*nCol_u;
538  t = pU[off_k+j];
539  pU[off_k+j] = pU[off_i1+j];
540  pU[off_i1+j] = t;
541  }
542  }
543  }
544  }
545  }
546 }
547 
548 ////////////////////////////////////////////////////////////////////////////////
549 /// Reconstruct the original matrix using the decomposition parts
550 
552 {
553  if (TestBit(kSingular)) {
554  Error("GetMatrix()","Matrix is singular");
555  return TMatrixD();
556  }
557  if ( !TestBit(kDecomposed) ) {
558  if (!Decompose()) {
559  Error("GetMatrix()","Decomposition failed");
560  return TMatrixD();
561  }
562  }
563 
564  const Int_t nRows = fU.GetNrows();
565  const Int_t nCols = fV.GetNcols();
566  const Int_t colLwb = this->GetColLwb();
567  TMatrixD s(nRows,nCols); s.Shift(colLwb,colLwb);
568  TMatrixDDiag diag(s); diag = fSig;
570  return fU * s * vt;
571 }
572 
573 ////////////////////////////////////////////////////////////////////////////////
574 /// Set matrix to be decomposed
575 
577 {
578  R__ASSERT(a.IsValid());
579 
580  ResetStatus();
581  if (a.GetNrows() < a.GetNcols()) {
582  Error("TDecompSVD(const TMatrixD &","matrix rows should be >= columns");
583  return;
584  }
585 
587  fCondition = -1.0;
588 
589  fRowLwb = a.GetRowLwb();
590  fColLwb = a.GetColLwb();
591  const Int_t nRow = a.GetNrows();
592  const Int_t nCol = a.GetNcols();
593 
594  fU.ResizeTo(nRow,nRow);
595  fSig.ResizeTo(nCol);
596  fV.ResizeTo(nRow,nCol); // In the end we only need the nColxnCol part
597 
598  fU.UnitMatrix();
599  memcpy(fV.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
600 }
601 
602 ////////////////////////////////////////////////////////////////////////////////
603 /// Solve Ax=b assuming the SVD form of A is stored . Solution returned in b.
604 /// If A is of size (m x n), input vector b should be of size (m), however,
605 /// the solution, returned in b, will be in the first (n) elements .
606 ///
607 /// For m > n , x is the least-squares solution of min(A . x - b)
608 
610 {
611  R__ASSERT(b.IsValid());
612  if (TestBit(kSingular)) {
613  return kFALSE;
614  }
615  if ( !TestBit(kDecomposed) ) {
616  if (!Decompose()) {
617  return kFALSE;
618  }
619  }
620 
621  if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb())
622  {
623  Error("Solve(TVectorD &","vector and matrix incompatible");
624  return kFALSE;
625  }
626 
627  // We start with fU fSig fV^T x = b, and turn it into fV^T x = fSig^-1 fU^T b
628  // Form tmp = fSig^-1 fU^T b but ignore diagonal elements with
629  // fSig(i) < fTol * max(fSig)
630 
631  const Int_t lwb = fU.GetColLwb();
632  const Int_t upb = lwb+fV.GetNcols()-1;
633  const Double_t threshold = fSig(lwb)*fTol;
634 
635  TVectorD tmp(lwb,upb);
636  for (Int_t irow = lwb; irow <= upb; irow++) {
637  Double_t r = 0.0;
638  if (fSig(irow) > threshold) {
639  const TVectorD uc_i = TMatrixDColumn(fU,irow);
640  r = uc_i * b;
641  r /= fSig(irow);
642  }
643  tmp(irow) = r;
644  }
645 
646  if (b.GetNrows() > fV.GetNrows()) {
647  TVectorD tmp2;
648  tmp2.Use(lwb,upb,b.GetMatrixArray());
649  tmp2 = fV*tmp;
650  } else
651  b = fV*tmp;
652 
653  return kTRUE;
654 }
655 
656 ////////////////////////////////////////////////////////////////////////////////
657 /// Solve Ax=b assuming the SVD form of A is stored . Solution returned in the
658 /// matrix column cb b.
659 /// If A is of size (m x n), input vector b should be of size (m), however,
660 /// the solution, returned in b, will be in the first (n) elements .
661 ///
662 /// For m > n , x is the least-squares solution of min(A . x - b)
663 
665 {
666  TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
667  R__ASSERT(b->IsValid());
668  if (TestBit(kSingular)) {
669  return kFALSE;
670  }
671  if ( !TestBit(kDecomposed) ) {
672  if (!Decompose()) {
673  return kFALSE;
674  }
675  }
676 
677  if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb())
678  {
679  Error("Solve(TMatrixDColumn &","vector and matrix incompatible");
680  return kFALSE;
681  }
682 
683  // We start with fU fSig fV^T x = b, and turn it into fV^T x = fSig^-1 fU^T b
684  // Form tmp = fSig^-1 fU^T b but ignore diagonal elements in
685  // fSig(i) < fTol * max(fSig)
686 
687  const Int_t lwb = fU.GetColLwb();
688  const Int_t upb = lwb+fV.GetNcols()-1;
689  const Double_t threshold = fSig(lwb)*fTol;
690 
691  TVectorD tmp(lwb,upb);
692  const TVectorD vb = cb;
693  for (Int_t irow = lwb; irow <= upb; irow++) {
694  Double_t r = 0.0;
695  if (fSig(irow) > threshold) {
696  const TVectorD uc_i = TMatrixDColumn(fU,irow);
697  r = uc_i * vb;
698  r /= fSig(irow);
699  }
700  tmp(irow) = r;
701  }
702 
703  if (b->GetNrows() > fV.GetNrows()) {
704  const TVectorD tmp2 = fV*tmp;
705  TVectorD tmp3(cb);
706  tmp3.SetSub(tmp2.GetLwb(),tmp2);
707  cb = tmp3;
708  } else
709  cb = fV*tmp;
710 
711  return kTRUE;
712 }
713 
714 ////////////////////////////////////////////////////////////////////////////////
715 /// Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b.
716 
718 {
719  R__ASSERT(b.IsValid());
720  if (TestBit(kSingular)) {
721  return kFALSE;
722  }
723  if ( !TestBit(kDecomposed) ) {
724  if (!Decompose()) {
725  return kFALSE;
726  }
727  }
728 
729  if (fU.GetNrows() != fV.GetNrows() || fU.GetRowLwb() != fV.GetRowLwb()) {
730  Error("TransSolve(TVectorD &","matrix should be square");
731  return kFALSE;
732  }
733 
734  if (fV.GetNrows() != b.GetNrows() || fV.GetRowLwb() != b.GetLwb())
735  {
736  Error("TransSolve(TVectorD &","vector and matrix incompatible");
737  return kFALSE;
738  }
739 
740  // We start with fV fSig fU^T x = b, and turn it into fU^T x = fSig^-1 fV^T b
741  // Form tmp = fSig^-1 fV^T b but ignore diagonal elements in
742  // fSig(i) < fTol * max(fSig)
743 
744  const Int_t lwb = fU.GetColLwb();
745  const Int_t upb = lwb+fV.GetNcols()-1;
746  const Double_t threshold = fSig(lwb)*fTol;
747 
748  TVectorD tmp(lwb,upb);
749  for (Int_t i = lwb; i <= upb; i++) {
750  Double_t r = 0.0;
751  if (fSig(i) > threshold) {
752  const TVectorD vc = TMatrixDColumn(fV,i);
753  r = vc * b;
754  r /= fSig(i);
755  }
756  tmp(i) = r;
757  }
758  b = fU*tmp;
759 
760  return kTRUE;
761 }
762 
763 ////////////////////////////////////////////////////////////////////////////////
764 /// Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b.
765 
767 {
768  TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
769  R__ASSERT(b->IsValid());
770  if (TestBit(kSingular)) {
771  return kFALSE;
772  }
773  if ( !TestBit(kDecomposed) ) {
774  if (!Decompose()) {
775  return kFALSE;
776  }
777  }
778 
779  if (fU.GetNrows() != fV.GetNrows() || fU.GetRowLwb() != fV.GetRowLwb()) {
780  Error("TransSolve(TMatrixDColumn &","matrix should be square");
781  return kFALSE;
782  }
783 
784  if (fV.GetNrows() != b->GetNrows() || fV.GetRowLwb() != b->GetRowLwb())
785  {
786  Error("TransSolve(TMatrixDColumn &","vector and matrix incompatible");
787  return kFALSE;
788  }
789 
790  // We start with fV fSig fU^T x = b, and turn it into fU^T x = fSig^-1 fV^T b
791  // Form tmp = fSig^-1 fV^T b but ignore diagonal elements in
792  // fSig(i) < fTol * max(fSig)
793 
794  const Int_t lwb = fU.GetColLwb();
795  const Int_t upb = lwb+fV.GetNcols()-1;
796  const Double_t threshold = fSig(lwb)*fTol;
797 
798  const TVectorD vb = cb;
799  TVectorD tmp(lwb,upb);
800  for (Int_t i = lwb; i <= upb; i++) {
801  Double_t r = 0.0;
802  if (fSig(i) > threshold) {
803  const TVectorD vc = TMatrixDColumn(fV,i);
804  r = vc * vb;
805  r /= fSig(i);
806  }
807  tmp(i) = r;
808  }
809  cb = fU*tmp;
810 
811  return kTRUE;
812 }
813 
814 ////////////////////////////////////////////////////////////////////////////////
815 /// Matrix condition number
816 
818 {
819  if ( !TestBit(kCondition) ) {
820  fCondition = -1;
821  if (TestBit(kSingular))
822  return fCondition;
823  if ( !TestBit(kDecomposed) ) {
824  if (!Decompose())
825  return fCondition;
826  }
827  const Int_t colLwb = GetColLwb();
828  const Int_t nCols = GetNcols();
829  const Double_t max = fSig(colLwb);
830  const Double_t min = fSig(colLwb+nCols-1);
831  fCondition = (min > 0.0) ? max/min : -1.0;
833  }
834  return fCondition;
835 }
836 
837 ////////////////////////////////////////////////////////////////////////////////
838 /// Matrix determinant det = d1*TMath::Power(2.,d2)
839 
841 {
842  if ( !TestBit(kDetermined) ) {
843  if ( !TestBit(kDecomposed) )
844  Decompose();
845  if (TestBit(kSingular)) {
846  fDet1 = 0.0;
847  fDet2 = 0.0;
848  } else {
850  }
852  }
853  d1 = fDet1;
854  d2 = fDet2;
855 }
856 
857 ////////////////////////////////////////////////////////////////////////////////
858 
860 {
861  return fU.GetNrows();
862 }
863 
865 {
866  return fV.GetNcols();
867 }
868 
869 ////////////////////////////////////////////////////////////////////////////////
870 /// For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
871 /// The user should always supply a matrix of size (m x m) !
872 /// If m > n , only the (n x m) part of the returned (pseudo inverse) matrix
873 /// should be used .
874 
876 {
877  const Int_t rowLwb = GetRowLwb();
878  const Int_t colLwb = GetColLwb();
879  const Int_t nRows = fU.GetNrows();
880 
881  if (inv.GetNrows() != nRows || inv.GetNcols() != nRows ||
882  inv.GetRowLwb() != rowLwb || inv.GetColLwb() != colLwb) {
883  Error("Invert(TMatrixD &","Input matrix has wrong shape");
884  return kFALSE;
885  }
886 
887  inv.UnitMatrix();
888  Bool_t status = MultiSolve(inv);
889 
890  return status;
891 }
892 
893 ////////////////////////////////////////////////////////////////////////////////
894 /// For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
895 /// (n x m) Ainv is returned .
896 
898 {
899  const Int_t rowLwb = GetRowLwb();
900  const Int_t colLwb = GetColLwb();
901  const Int_t rowUpb = rowLwb+fU.GetNrows()-1;
902  TMatrixD inv(rowLwb,rowUpb,colLwb,colLwb+fU.GetNrows()-1);
903  inv.UnitMatrix();
904  status = MultiSolve(inv);
905  inv.ResizeTo(rowLwb,rowLwb+fV.GetNcols()-1,colLwb,colLwb+fU.GetNrows()-1);
906 
907  return inv;
908 }
909 
910 ////////////////////////////////////////////////////////////////////////////////
911 /// Print class members
912 
913 void TDecompSVD::Print(Option_t *opt) const
914 {
915  TDecompBase::Print(opt);
916  fU.Print("fU");
917  fV.Print("fV");
918  fSig.Print("fSig");
919 }
920 
921 ////////////////////////////////////////////////////////////////////////////////
922 /// Assignment operator
923 
925 {
926  if (this != &source) {
927  TDecompBase::operator=(source);
928  fU.ResizeTo(source.fU);
929  fU = source.fU;
930  fV.ResizeTo(source.fV);
931  fV = source.fV;
932  fSig.ResizeTo(source.fSig);
933  fSig = source.fSig;
934  }
935  return *this;
936 }
Int_t GetRowLwb() const
Definition: TMatrixTBase.h:132
virtual TMatrixTBase< Element > & UnitMatrix()
Make a unit matrix (matrix need not be a square one).
virtual void Det(Double_t &d1, Double_t &d2)
Matrix determinant det = d1*TMath::Power(2.,d2)
Definition: TDecompSVD.cxx:840
virtual TMatrixTBase< Element > & Shift(Int_t row_shift, Int_t col_shift)
Shift the row index by adding row_shift and the column index by adding col_shift, respectively...
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:291
static Bool_t Bidiagonalize(TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag)
Bidiagonalize the (m x n) - matrix a (stored in v) through a series of Householder transformations ap...
Definition: TDecompSVD.cxx:188
virtual Bool_t MultiSolve(TMatrixD &B)
Solve set of equations with RHS in columns of B.
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
TMatrixT< Element > & Transpose(const TMatrixT< Element > &source)
Transpose matrix source.
Definition: TMatrixT.cxx:1461
static void Diag_3(TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag, Int_t k, Int_t l)
Step 3 in the matrix diagonalization.
Definition: TDecompSVD.cxx:410
Double_t fCondition
Definition: TDecompBase.h:42
virtual Bool_t TransSolve(TVectorD &b)
Solve A^T x=b assuming the SVD form of A is stored . Solution returned in b.
Definition: TDecompSVD.cxx:717
Bool_t IsValid() const
Definition: TVectorT.h:89
Int_t fRowLwb
Definition: TDecompBase.h:43
return c
const char Option_t
Definition: RtypesCore.h:62
void Print(Option_t *opt="") const
Print class members.
Definition: TDecompSVD.cxx:913
Int_t GetNrows() const
Definition: TVectorT.h:81
TH1 * h
Definition: legend2.C:5
void Print(Option_t *opt="") const
Print class members.
#define R__ASSERT(e)
Definition: TError.h:98
double inv(double x)
For comparisons.
Definition: inv.h:58
TMatrixTColumn_const< Double_t > TMatrixDColumn_const
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SetMatrix(const TMatrixD &a)
Set matrix to be decomposed.
Definition: TDecompSVD.cxx:576
TMatrixD fV
Definition: TDecompSVD.h:31
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Definition: TMatrixT.cxx:1202
virtual Int_t GetNcols() const
Definition: TDecompSVD.cxx:864
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TFile * f
Bool_t DefHouseHolder(const TVectorD &vc, Int_t lp, Int_t l, Double_t &up, Double_t &b, Double_t tol=0.0)
Define a Householder-transformation through the parameters up and b .
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
double beta(double x, double y)
Calculates the beta function.
static Bool_t Diagonalize(TMatrixD &v, TMatrixD &u, TVectorD &sDiag, TVectorD &oDiag)
Diagonalizes in an iterative fashion the bidiagonal matrix C as described through sDiag and oDiag...
Definition: TDecompSVD.cxx:301
Double_t fTol
Definition: TDecompBase.h:39
Bool_t IsValid() const
Definition: TMatrixTBase.h:157
Element GetTol() const
Definition: TMatrixTBase.h:139
TVectorT< Element > & Use(Int_t lwb, Int_t upb, Element *data)
Use the array data to fill the vector lwb..upb].
Definition: TVectorT.cxx:346
TDecompBase & operator=(const TDecompBase &source)
Assignment operator.
TVectorT< Element > & SetSub(Int_t row_lwb, const TVectorT< Element > &source)
Insert vector source starting at [row_lwb], thereby overwriting the part [row_lwb..row_lwb+nrows_source];.
Definition: TVectorT.cxx:419
void ApplyGivens(Double_t &z1, Double_t &z2, Double_t c, Double_t s)
Apply a Givens transformation as defined by c and s to the vector compenents v1 and v2 ...
TMatrixTRow_const< Double_t > TMatrixDRow_const
Int_t GetColLwb() const
Definition: TDecompBase.h:77
Int_t fColLwb
Definition: TDecompBase.h:44
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
void Error(const char *location, const char *msgfmt,...)
static void DiagProd(const TVectorD &diag, Double_t tol, Double_t &d1, Double_t &d2)
virtual Bool_t Decompose()
SVD decomposition of matrix If the decomposition succeeds, bit kDecomposed is set ...
Definition: TDecompSVD.cxx:121
TMatrixTRow< Double_t > TMatrixDRow
TVectorT< Element > & Shift(Int_t row_shift)
Definition: TVectorT.h:93
Int_t GetNrows() const
Definition: TMatrixTBase.h:134
Element * GetMatrixArray()
Definition: TVectorT.h:84
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
void Print(Option_t *name="") const
Print the matrix as a table of elements.
void Print(Option_t *option="") const
Print the vector as a list of elements.
Definition: TVectorT.cxx:1360
static void SortSingular(TMatrixD &v, TMatrixD &u, TVectorD &sDiag)
Perform a permutation transformation on the diagonal matrix S', so that matrix S'' = U''^T ...
Definition: TDecompSVD.cxx:491
TThread * t[5]
Definition: threadsh1.C:13
const double tol
ROOT::R::TRInterface & r
Definition: Object.C:4
SVector< double, 2 > v
Definition: Dict.h:5
void ApplyHouseHolder(const TVectorD &vc, Double_t up, Double_t b, Int_t lp, Int_t l, TMatrixDRow &cr)
Apply Householder-transformation.
void DefGivens(Double_t v1, Double_t v2, Double_t &c, Double_t &s)
Defines a Givens-rotation by calculating 2 rotation parameters c and s.
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TLine * l
Definition: textangle.C:4
virtual const Element * GetMatrixArray() const
Definition: TMatrixT.h:222
void ResetStatus()
Definition: TDecompBase.h:46
TMatrixTDiag< Double_t > TMatrixDDiag
virtual Int_t GetNrows() const
Definition: TDecompSVD.cxx:859
TMatrixD fU
Definition: TDecompSVD.h:30
const TMatrixTBase< Element > * GetMatrix() const
Int_t GetColLwb() const
Definition: TMatrixTBase.h:135
REAL epsilon
Definition: triangle.c:617
Double_t fDet2
Definition: TDecompBase.h:41
Double_t fDet1
Definition: TDecompBase.h:40
Int_t GetRowLwb() const
Definition: TDecompBase.h:76
virtual Bool_t Solve(TVectorD &b)
Solve Ax=b assuming the SVD form of A is stored .
Definition: TDecompSVD.cxx:609
static void Diag_1(TMatrixD &v, TVectorD &sDiag, TVectorD &oDiag, Int_t k)
Step 1 in the matrix diagonalization.
Definition: TDecompSVD.cxx:368
Int_t GetLwb() const
Definition: TVectorT.h:79
double Double_t
Definition: RtypesCore.h:55
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
Double_t Hypot(Double_t x, Double_t y)
Definition: TMath.cxx:60
Array of doubles (64 bits per element).
Definition: TArrayD.h:29
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
TMatrixTColumn< Double_t > TMatrixDColumn
const TMatrixD GetMatrix()
Reconstruct the original matrix using the decomposition parts.
Definition: TDecompSVD.cxx:551
TMatrixD Invert()
Definition: TDecompSVD.h:84
void DefAplGivens(Double_t &v1, Double_t &v2, Double_t &c, Double_t &s)
Define and apply a Givens-rotation by calculating 2 rotation parameters c and s.
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
TVectorD fSig
Definition: TDecompSVD.h:32
virtual Double_t Condition()
Matrix condition number.
Definition: TDecompSVD.cxx:817
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
TDecompSVD & operator=(const TDecompSVD &source)
Assignment operator.
Definition: TDecompSVD.cxx:924
const Bool_t kTRUE
Definition: Rtypes.h:91
static void Diag_2(TVectorD &sDiag, TVectorD &oDiag, Int_t k, Int_t l)
Step 2 in the matrix diagonalization.
Definition: TDecompSVD.cxx:392
ClassImp(TDecompSVD) TDecompSVD
Constructor for (nrows x ncols) matrix.
Definition: TDecompSVD.cxx:45