00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041 #include "TDecompSVD.h"
00042 #include "TMath.h"
00043 #include "TArrayD.h"
00044
00045 ClassImp(TDecompSVD)
00046
00047
00048 TDecompSVD::TDecompSVD(Int_t nrows,Int_t ncols)
00049 {
00050
00051
00052 if (nrows < ncols) {
00053 Error("TDecompSVD(Int_t,Int_t","matrix rows should be >= columns");
00054 return;
00055 }
00056 fU.ResizeTo(nrows,nrows);
00057 fSig.ResizeTo(ncols);
00058 fV.ResizeTo(nrows,ncols);
00059 }
00060
00061
00062 TDecompSVD::TDecompSVD(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb)
00063 {
00064
00065
00066 const Int_t nrows = row_upb-row_lwb+1;
00067 const Int_t ncols = col_upb-col_lwb+1;
00068
00069 if (nrows < ncols) {
00070 Error("TDecompSVD(Int_t,Int_t,Int_t,Int_t","matrix rows should be >= columns");
00071 return;
00072 }
00073 fRowLwb = row_lwb;
00074 fColLwb = col_lwb;
00075 fU.ResizeTo(nrows,nrows);
00076 fSig.ResizeTo(ncols);
00077 fV.ResizeTo(nrows,ncols);
00078 }
00079
00080
00081 TDecompSVD::TDecompSVD(const TMatrixD &a,Double_t tol)
00082 {
00083
00084
00085 R__ASSERT(a.IsValid());
00086 if (a.GetNrows() < a.GetNcols()) {
00087 Error("TDecompSVD(const TMatrixD &","matrix rows should be >= columns");
00088 return;
00089 }
00090
00091 SetBit(kMatrixSet);
00092 fTol = a.GetTol();
00093 if (tol > 0)
00094 fTol = tol;
00095
00096 fRowLwb = a.GetRowLwb();
00097 fColLwb = a.GetColLwb();
00098 const Int_t nRow = a.GetNrows();
00099 const Int_t nCol = a.GetNcols();
00100
00101 fU.ResizeTo(nRow,nRow);
00102 fSig.ResizeTo(nCol);
00103 fV.ResizeTo(nRow,nCol);
00104
00105 fU.UnitMatrix();
00106 memcpy(fV.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
00107 }
00108
00109
00110 TDecompSVD::TDecompSVD(const TDecompSVD &another): TDecompBase(another)
00111 {
00112
00113
00114 *this = another;
00115 }
00116
00117
00118 Bool_t TDecompSVD::Decompose()
00119 {
00120
00121
00122
00123 if (TestBit(kDecomposed)) return kTRUE;
00124
00125 if ( !TestBit(kMatrixSet) ) {
00126 Error("Decompose()","Matrix has not been set");
00127 return kFALSE;
00128 }
00129
00130 const Int_t nCol = this->GetNcols();
00131 const Int_t rowLwb = this->GetRowLwb();
00132 const Int_t colLwb = this->GetColLwb();
00133
00134 TVectorD offDiag;
00135 Double_t work[kWorkMax];
00136 if (nCol > kWorkMax) offDiag.ResizeTo(nCol);
00137 else offDiag.Use(nCol,work);
00138
00139
00140 if (!Bidiagonalize(fV,fU,fSig,offDiag))
00141 return kFALSE;
00142
00143
00144 if (!Diagonalize(fV,fU,fSig,offDiag))
00145 return kFALSE;
00146
00147
00148 SortSingular(fV,fU,fSig);
00149 fV.ResizeTo(nCol,nCol); fV.Shift(colLwb,colLwb);
00150 fSig.Shift(colLwb);
00151 fU.Transpose(fU); fU.Shift(rowLwb,colLwb);
00152 SetBit(kDecomposed);
00153
00154 return kTRUE;
00155 }
00156
00157
00158 Bool_t TDecompSVD::Bidiagonalize(TMatrixD &v,TMatrixD &u,TVectorD &sDiag,TVectorD &oDiag)
00159 {
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190 const Int_t nRow_v = v.GetNrows();
00191 const Int_t nCol_v = v.GetNcols();
00192 const Int_t nCol_u = u.GetNcols();
00193
00194 TArrayD ups(nCol_v);
00195 TArrayD betas(nCol_v);
00196
00197 for (Int_t i = 0; i < nCol_v; i++) {
00198
00199
00200 Double_t up,beta;
00201 if (i < nCol_v-1 || nRow_v > nCol_v) {
00202 const TVectorD vc_i = TMatrixDColumn_const(v,i);
00203
00204
00205 DefHouseHolder(vc_i,i,i+1,up,beta);
00206
00207
00208 for (Int_t j = i; j < nCol_v; j++) {
00209 TMatrixDColumn vc_j = TMatrixDColumn(v,j);
00210 ApplyHouseHolder(vc_i,up,beta,i,i+1,vc_j);
00211 }
00212
00213
00214 for (Int_t j = 0; j < nCol_u; j++)
00215 {
00216 TMatrixDColumn uc_j = TMatrixDColumn(u,j);
00217 ApplyHouseHolder(vc_i,up,beta,i,i+1,uc_j);
00218 }
00219 }
00220 if (i < nCol_v-2) {
00221
00222 const TVectorD vr_i = TMatrixDRow_const(v,i);
00223
00224
00225
00226 DefHouseHolder(vr_i,i+1,i+2,up,beta);
00227
00228
00229 ups[i] = up;
00230 betas[i] = beta;
00231
00232
00233 for (Int_t j = i; j < nRow_v; j++) {
00234 TMatrixDRow vr_j = TMatrixDRow(v,j);
00235 ApplyHouseHolder(vr_i,up,beta,i+1,i+2,vr_j);
00236
00237
00238 if (j == i) {
00239 for (Int_t k = i+2; k < nCol_v; k++)
00240 vr_j(k) = vr_i(k);
00241 }
00242 }
00243 }
00244 }
00245
00246
00247 if (nCol_v > 1) {
00248 for (Int_t i = 1; i < nCol_v; i++)
00249 oDiag(i) = v(i-1,i);
00250 }
00251 oDiag(0) = 0.;
00252 sDiag = TMatrixDDiag(v);
00253
00254
00255
00256 TVectorD vr_i(nCol_v);
00257 for (Int_t i = nCol_v-1; i >= 0; i--) {
00258 if (i < nCol_v-1)
00259 vr_i = TMatrixDRow_const(v,i);
00260 TMatrixDRow(v,i) = 0.0;
00261 v(i,i) = 1.;
00262
00263 if (i < nCol_v-2) {
00264 for (Int_t k = i; k < nCol_v; k++) {
00265
00266 TMatrixDColumn vc_k = TMatrixDColumn(v,k);
00267 ApplyHouseHolder(vr_i,ups[i],betas[i],i+1,i+2,vc_k);
00268 }
00269 }
00270 }
00271
00272 return kTRUE;
00273 }
00274
00275
00276 Bool_t TDecompSVD::Diagonalize(TMatrixD &v,TMatrixD &u,TVectorD &sDiag,TVectorD &oDiag)
00277 {
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303 Bool_t ok = kTRUE;
00304 Int_t niter = 0;
00305 Double_t bmx = sDiag(0);
00306
00307 const Int_t nCol_v = v.GetNcols();
00308
00309 if (nCol_v > 1) {
00310 for (Int_t i = 1; i < nCol_v; i++)
00311 bmx = TMath::Max(TMath::Abs(sDiag(i))+TMath::Abs(oDiag(i)),bmx);
00312 }
00313
00314 const Double_t eps = std::numeric_limits<double>::epsilon();
00315
00316 const Int_t niterm = 10*nCol_v;
00317 for (Int_t k = nCol_v-1; k >= 0; k--) {
00318 loop:
00319 if (k != 0) {
00320
00321 if (TMath::Abs(sDiag(k)) < eps*bmx)
00322 Diag_1(v,sDiag,oDiag,k);
00323
00324
00325
00326
00327
00328
00329 Int_t elzero = 0;
00330 Int_t l = 0;
00331 for (Int_t ll = k; ll >= 0 ; ll--) {
00332 l = ll;
00333 if (l == 0) {
00334 elzero = 0;
00335 break;
00336 } else if (TMath::Abs(oDiag(l)) < eps*bmx) {
00337 elzero = 1;
00338 break;
00339 } else if (TMath::Abs(sDiag(l-1)) < eps*bmx)
00340 elzero = 0;
00341 }
00342 if (l > 0 && !elzero)
00343 Diag_2(sDiag,oDiag,k,l);
00344 if (l != k) {
00345
00346 Diag_3(v,u,sDiag,oDiag,k,l);
00347 niter++;
00348 if (niter <= niterm) goto loop;
00349 ::Error("TDecompSVD::Diagonalize","no convergence after %d steps",niter);
00350 ok = kFALSE;
00351 }
00352 }
00353
00354 if (sDiag(k) < 0.) {
00355
00356 sDiag(k) = -sDiag(k);
00357 TMatrixDColumn(v,k) *= -1.0;
00358 }
00359
00360 }
00361
00362 return ok;
00363 }
00364
00365
00366 void TDecompSVD::Diag_1(TMatrixD &v,TVectorD &sDiag,TVectorD &oDiag,Int_t k)
00367 {
00368
00369
00370 const Int_t nCol_v = v.GetNcols();
00371
00372 TMatrixDColumn vc_k = TMatrixDColumn(v,k);
00373 for (Int_t i = k-1; i >= 0; i--) {
00374 TMatrixDColumn vc_i = TMatrixDColumn(v,i);
00375 Double_t h,cs,sn;
00376 if (i == k-1)
00377 DefAplGivens(sDiag[i],oDiag[i+1],cs,sn);
00378 else
00379 DefAplGivens(sDiag[i],h,cs,sn);
00380 if (i > 0) {
00381 h = 0.;
00382 ApplyGivens(oDiag[i],h,cs,sn);
00383 }
00384 for (Int_t j = 0; j < nCol_v; j++)
00385 ApplyGivens(vc_i(j),vc_k(j),cs,sn);
00386 }
00387 }
00388
00389
00390 void TDecompSVD::Diag_2(TVectorD &sDiag,TVectorD &oDiag,Int_t k,Int_t l)
00391 {
00392
00393
00394 for (Int_t i = l; i <= k; i++) {
00395 Double_t h,cs,sn;
00396 if (i == l)
00397 DefAplGivens(sDiag(i),oDiag(i),cs,sn);
00398 else
00399 DefAplGivens(sDiag(i),h,cs,sn);
00400 if (i < k) {
00401 h = 0.;
00402 ApplyGivens(oDiag(i+1),h,cs,sn);
00403 }
00404 }
00405 }
00406
00407
00408 void TDecompSVD::Diag_3(TMatrixD &v,TMatrixD &u,TVectorD &sDiag,TVectorD &oDiag,Int_t k,Int_t l)
00409 {
00410
00411
00412 Double_t *pS = sDiag.GetMatrixArray();
00413 Double_t *pO = oDiag.GetMatrixArray();
00414
00415
00416
00417 const Double_t psk1 = pS[k-1];
00418 const Double_t psk = pS[k];
00419 const Double_t pok1 = pO[k-1];
00420 const Double_t pok = pO[k];
00421 const Double_t psl = pS[l];
00422
00423 Double_t f;
00424 if (psl == 0.0 || pok == 0.0 || psk1 == 0.0) {
00425 const Double_t b = ((psk1-psk)*(psk1+psk)+pok1*pok1)/2.0;
00426 const Double_t c = (psk*pok1)*(psk*pok1);
00427
00428 Double_t shift = 0.0;
00429 if ((b != 0.0) | (c != 0.0)) {
00430 shift = TMath::Sqrt(b*b+c);
00431 if (b < 0.0)
00432 shift = -shift;
00433 shift = c/(b+shift);
00434 }
00435
00436 f = (psl+psk)*(psl-psk)+shift;
00437 } else {
00438 f = ((psk1-psk)*(psk1+psk)+(pok1-pok)*(pok1+pok))/(2.*pok*psk1);
00439 const Double_t g = TMath::Hypot(1.,f);
00440 const Double_t t = (f >= 0.) ? f+g : f-g;
00441
00442 f = ((psl-psk)*(psl+psk)+pok*(psk1/t-pok))/psl;
00443 }
00444
00445 const Int_t nCol_v = v.GetNcols();
00446 const Int_t nCol_u = u.GetNcols();
00447
00448 Double_t h,cs,sn;
00449 Int_t j;
00450 for (Int_t i = l; i <= k-1; i++) {
00451 if (i == l)
00452
00453 DefGivens(f,pO[i+1],cs,sn);
00454 else
00455
00456 DefAplGivens(pO[i],h,cs,sn);
00457
00458 ApplyGivens(pS[i],pO[i+1],cs,sn);
00459 h = 0.;
00460 ApplyGivens(h,pS[i+1],cs,sn);
00461
00462 TMatrixDColumn vc_i = TMatrixDColumn(v,i);
00463 TMatrixDColumn vc_i1 = TMatrixDColumn(v,i+1);
00464 for (j = 0; j < nCol_v; j++)
00465 ApplyGivens(vc_i(j),vc_i1(j),cs,sn);
00466
00467 DefAplGivens(pS[i],h,cs,sn);
00468 ApplyGivens(pO[i+1],pS[i+1],cs,sn);
00469 if (i < k-1) {
00470 h = 0.;
00471 ApplyGivens(h,pO[i+2],cs,sn);
00472 }
00473
00474 TMatrixDRow ur_i = TMatrixDRow(u,i);
00475 TMatrixDRow ur_i1 = TMatrixDRow(u,i+1);
00476 for (j = 0; j < nCol_u; j++)
00477 ApplyGivens(ur_i(j),ur_i1(j),cs,sn);
00478 }
00479 }
00480
00481
00482 void TDecompSVD::SortSingular(TMatrixD &v,TMatrixD &u,TVectorD &sDiag)
00483 {
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493 const Int_t nCol_v = v.GetNcols();
00494 const Int_t nCol_u = u.GetNcols();
00495
00496 Double_t *pS = sDiag.GetMatrixArray();
00497 Double_t *pV = v.GetMatrixArray();
00498 Double_t *pU = u.GetMatrixArray();
00499
00500
00501
00502 Int_t i,j;
00503 if (nCol_v > 1) {
00504 while (1) {
00505 Bool_t found = kFALSE;
00506 i = 1;
00507 while (!found && i < nCol_v) {
00508 if (pS[i] > pS[i-1])
00509 found = kTRUE;
00510 else
00511 i++;
00512 }
00513 if (!found) break;
00514 for (i = 1; i < nCol_v; i++) {
00515 Double_t t = pS[i-1];
00516 Int_t k = i-1;
00517 for (j = i; j < nCol_v; j++) {
00518 if (t < pS[j]) {
00519 t = pS[j];
00520 k = j;
00521 }
00522 }
00523 if (k != i-1) {
00524
00525 pS[k] = pS[i-1];
00526 pS[i-1] = t;
00527
00528 for (j = 0; j < nCol_v; j++) {
00529 const Int_t off_j = j*nCol_v;
00530 t = pV[off_j+k];
00531 pV[off_j+k] = pV[off_j+i-1];
00532 pV[off_j+i-1] = t;
00533 }
00534
00535 for (j = 0; j < nCol_u; j++) {
00536 const Int_t off_k = k*nCol_u;
00537 const Int_t off_i1 = (i-1)*nCol_u;
00538 t = pU[off_k+j];
00539 pU[off_k+j] = pU[off_i1+j];
00540 pU[off_i1+j] = t;
00541 }
00542 }
00543 }
00544 }
00545 }
00546 }
00547
00548
00549 const TMatrixD TDecompSVD::GetMatrix()
00550 {
00551
00552
00553 if (TestBit(kSingular)) {
00554 Error("GetMatrix()","Matrix is singular");
00555 return TMatrixD();
00556 }
00557 if ( !TestBit(kDecomposed) ) {
00558 if (!Decompose()) {
00559 Error("GetMatrix()","Decomposition failed");
00560 return TMatrixD();
00561 }
00562 }
00563
00564 const Int_t nRows = fU.GetNrows();
00565 const Int_t nCols = fV.GetNcols();
00566 const Int_t colLwb = this->GetColLwb();
00567 TMatrixD s(nRows,nCols); s.Shift(colLwb,colLwb);
00568 TMatrixDDiag diag(s); diag = fSig;
00569 const TMatrixD vt(TMatrixD::kTransposed,fV);
00570 return fU * s * vt;
00571 }
00572
00573
00574 void TDecompSVD::SetMatrix(const TMatrixD &a)
00575 {
00576
00577
00578 R__ASSERT(a.IsValid());
00579
00580 ResetStatus();
00581 if (a.GetNrows() < a.GetNcols()) {
00582 Error("TDecompSVD(const TMatrixD &","matrix rows should be >= columns");
00583 return;
00584 }
00585
00586 SetBit(kMatrixSet);
00587 fCondition = -1.0;
00588
00589 fRowLwb = a.GetRowLwb();
00590 fColLwb = a.GetColLwb();
00591 const Int_t nRow = a.GetNrows();
00592 const Int_t nCol = a.GetNcols();
00593
00594 fU.ResizeTo(nRow,nRow);
00595 fSig.ResizeTo(nCol);
00596 fV.ResizeTo(nRow,nCol);
00597
00598 fU.UnitMatrix();
00599 memcpy(fV.GetMatrixArray(),a.GetMatrixArray(),nRow*nCol*sizeof(Double_t));
00600 }
00601
00602
00603 Bool_t TDecompSVD::Solve(TVectorD &b)
00604 {
00605
00606
00607
00608
00609
00610
00611 R__ASSERT(b.IsValid());
00612 if (TestBit(kSingular)) {
00613 return kFALSE;
00614 }
00615 if ( !TestBit(kDecomposed) ) {
00616 if (!Decompose()) {
00617 return kFALSE;
00618 }
00619 }
00620
00621 if (fU.GetNrows() != b.GetNrows() || fU.GetRowLwb() != b.GetLwb())
00622 {
00623 Error("Solve(TVectorD &","vector and matrix incompatible");
00624 return kFALSE;
00625 }
00626
00627
00628
00629
00630
00631 const Int_t lwb = fU.GetColLwb();
00632 const Int_t upb = lwb+fV.GetNcols()-1;
00633 const Double_t threshold = fSig(lwb)*fTol;
00634
00635 TVectorD tmp(lwb,upb);
00636 for (Int_t irow = lwb; irow <= upb; irow++) {
00637 Double_t r = 0.0;
00638 if (fSig(irow) > threshold) {
00639 const TVectorD uc_i = TMatrixDColumn(fU,irow);
00640 r = uc_i * b;
00641 r /= fSig(irow);
00642 }
00643 tmp(irow) = r;
00644 }
00645
00646 if (b.GetNrows() > fV.GetNrows()) {
00647 TVectorD tmp2;
00648 tmp2.Use(lwb,upb,b.GetMatrixArray());
00649 tmp2 = fV*tmp;
00650 } else
00651 b = fV*tmp;
00652
00653 return kTRUE;
00654 }
00655
00656
00657 Bool_t TDecompSVD::Solve(TMatrixDColumn &cb)
00658 {
00659
00660
00661
00662
00663
00664
00665
00666 TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
00667 R__ASSERT(b->IsValid());
00668 if (TestBit(kSingular)) {
00669 return kFALSE;
00670 }
00671 if ( !TestBit(kDecomposed) ) {
00672 if (!Decompose()) {
00673 return kFALSE;
00674 }
00675 }
00676
00677 if (fU.GetNrows() != b->GetNrows() || fU.GetRowLwb() != b->GetRowLwb())
00678 {
00679 Error("Solve(TMatrixDColumn &","vector and matrix incompatible");
00680 return kFALSE;
00681 }
00682
00683
00684
00685
00686
00687 const Int_t lwb = fU.GetColLwb();
00688 const Int_t upb = lwb+fV.GetNcols()-1;
00689 const Double_t threshold = fSig(lwb)*fTol;
00690
00691 TVectorD tmp(lwb,upb);
00692 const TVectorD vb = cb;
00693 for (Int_t irow = lwb; irow <= upb; irow++) {
00694 Double_t r = 0.0;
00695 if (fSig(irow) > threshold) {
00696 const TVectorD uc_i = TMatrixDColumn(fU,irow);
00697 r = uc_i * vb;
00698 r /= fSig(irow);
00699 }
00700 tmp(irow) = r;
00701 }
00702
00703 if (b->GetNrows() > fV.GetNrows()) {
00704 const TVectorD tmp2 = fV*tmp;
00705 TVectorD tmp3(cb);
00706 tmp3.SetSub(tmp2.GetLwb(),tmp2);
00707 cb = tmp3;
00708 } else
00709 cb = fV*tmp;
00710
00711 return kTRUE;
00712 }
00713
00714
00715 Bool_t TDecompSVD::TransSolve(TVectorD &b)
00716 {
00717
00718
00719 R__ASSERT(b.IsValid());
00720 if (TestBit(kSingular)) {
00721 return kFALSE;
00722 }
00723 if ( !TestBit(kDecomposed) ) {
00724 if (!Decompose()) {
00725 return kFALSE;
00726 }
00727 }
00728
00729 if (fU.GetNrows() != fV.GetNrows() || fU.GetRowLwb() != fV.GetRowLwb()) {
00730 Error("TransSolve(TVectorD &","matrix should be square");
00731 return kFALSE;
00732 }
00733
00734 if (fV.GetNrows() != b.GetNrows() || fV.GetRowLwb() != b.GetLwb())
00735 {
00736 Error("TransSolve(TVectorD &","vector and matrix incompatible");
00737 return kFALSE;
00738 }
00739
00740
00741
00742
00743
00744 const Int_t lwb = fU.GetColLwb();
00745 const Int_t upb = lwb+fV.GetNcols()-1;
00746 const Double_t threshold = fSig(lwb)*fTol;
00747
00748 TVectorD tmp(lwb,upb);
00749 for (Int_t i = lwb; i <= upb; i++) {
00750 Double_t r = 0.0;
00751 if (fSig(i) > threshold) {
00752 const TVectorD vc = TMatrixDColumn(fV,i);
00753 r = vc * b;
00754 r /= fSig(i);
00755 }
00756 tmp(i) = r;
00757 }
00758 b = fU*tmp;
00759
00760 return kTRUE;
00761 }
00762
00763
00764 Bool_t TDecompSVD::TransSolve(TMatrixDColumn &cb)
00765 {
00766
00767
00768 TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
00769 R__ASSERT(b->IsValid());
00770 if (TestBit(kSingular)) {
00771 return kFALSE;
00772 }
00773 if ( !TestBit(kDecomposed) ) {
00774 if (!Decompose()) {
00775 return kFALSE;
00776 }
00777 }
00778
00779 if (fU.GetNrows() != fV.GetNrows() || fU.GetRowLwb() != fV.GetRowLwb()) {
00780 Error("TransSolve(TMatrixDColumn &","matrix should be square");
00781 return kFALSE;
00782 }
00783
00784 if (fV.GetNrows() != b->GetNrows() || fV.GetRowLwb() != b->GetRowLwb())
00785 {
00786 Error("TransSolve(TMatrixDColumn &","vector and matrix incompatible");
00787 return kFALSE;
00788 }
00789
00790
00791
00792
00793
00794 const Int_t lwb = fU.GetColLwb();
00795 const Int_t upb = lwb+fV.GetNcols()-1;
00796 const Double_t threshold = fSig(lwb)*fTol;
00797
00798 const TVectorD vb = cb;
00799 TVectorD tmp(lwb,upb);
00800 for (Int_t i = lwb; i <= upb; i++) {
00801 Double_t r = 0.0;
00802 if (fSig(i) > threshold) {
00803 const TVectorD vc = TMatrixDColumn(fV,i);
00804 r = vc * vb;
00805 r /= fSig(i);
00806 }
00807 tmp(i) = r;
00808 }
00809 cb = fU*tmp;
00810
00811 return kTRUE;
00812 }
00813
00814
00815 Double_t TDecompSVD::Condition()
00816 {
00817
00818
00819 if ( !TestBit(kCondition) ) {
00820 fCondition = -1;
00821 if (TestBit(kSingular))
00822 return fCondition;
00823 if ( !TestBit(kDecomposed) ) {
00824 if (!Decompose())
00825 return fCondition;
00826 }
00827 const Int_t colLwb = GetColLwb();
00828 const Int_t nCols = GetNcols();
00829 const Double_t max = fSig(colLwb);
00830 const Double_t min = fSig(colLwb+nCols-1);
00831 fCondition = (min > 0.0) ? max/min : -1.0;
00832 SetBit(kCondition);
00833 }
00834 return fCondition;
00835 }
00836
00837
00838 void TDecompSVD::Det(Double_t &d1,Double_t &d2)
00839 {
00840
00841
00842 if ( !TestBit(kDetermined) ) {
00843 if ( !TestBit(kDecomposed) )
00844 Decompose();
00845 if (TestBit(kSingular)) {
00846 fDet1 = 0.0;
00847 fDet2 = 0.0;
00848 } else {
00849 DiagProd(fSig,fTol,fDet1,fDet2);
00850 }
00851 SetBit(kDetermined);
00852 }
00853 d1 = fDet1;
00854 d2 = fDet2;
00855 }
00856
00857
00858 Bool_t TDecompSVD::Invert(TMatrixD &inv)
00859 {
00860
00861
00862
00863
00864
00865 const Int_t rowLwb = GetRowLwb();
00866 const Int_t colLwb = GetColLwb();
00867 const Int_t nRows = GetNrows();
00868
00869 if (inv.GetNrows() != nRows || inv.GetNcols() != nRows ||
00870 inv.GetRowLwb() != rowLwb || inv.GetColLwb() != colLwb) {
00871 Error("Invert(TMatrixD &","Input matrix has wrong shape");
00872 return kFALSE;
00873 }
00874
00875 inv.UnitMatrix();
00876 Bool_t status = MultiSolve(inv);
00877
00878 return status;
00879 }
00880
00881
00882 TMatrixD TDecompSVD::Invert(Bool_t &status)
00883 {
00884
00885
00886
00887 const Int_t rowLwb = GetRowLwb();
00888 const Int_t colLwb = GetColLwb();
00889 const Int_t rowUpb = rowLwb+GetNrows()-1;
00890 TMatrixD inv(rowLwb,rowUpb,colLwb,colLwb+GetNrows()-1);
00891 inv.UnitMatrix();
00892 status = MultiSolve(inv);
00893 inv.ResizeTo(rowLwb,rowLwb+GetNcols()-1,colLwb,colLwb+GetNrows()-1);
00894
00895 return inv;
00896 }
00897
00898
00899 void TDecompSVD::Print(Option_t *opt) const
00900 {
00901
00902
00903 TDecompBase::Print(opt);
00904 fU.Print("fU");
00905 fV.Print("fV");
00906 fSig.Print("fSig");
00907 }
00908
00909
00910 TDecompSVD &TDecompSVD::operator=(const TDecompSVD &source)
00911 {
00912
00913
00914 if (this != &source) {
00915 TDecompBase::operator=(source);
00916 fU.ResizeTo(source.fU);
00917 fU = source.fU;
00918 fV.ResizeTo(source.fV);
00919 fV = source.fV;
00920 fSig.ResizeTo(source.fSig);
00921 fSig = source.fSig;
00922 }
00923 return *this;
00924 }