TDecompLU.cxx

Go to the documentation of this file.
00001 // @(#)root/matrix:$Id: TDecompLU.cxx 35904 2010-09-30 10:43:35Z brun $
00002 // Authors: Fons Rademakers, Eddy Offermann  Dec 2003
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "TDecompLU.h"
00013 #include "TMath.h"
00014 
00015 ClassImp(TDecompLU)
00016 
00017 ///////////////////////////////////////////////////////////////////////////
00018 //                                                                       //
00019 // LU Decomposition class                                                //
00020 //                                                                       //
00021 // Decompose  a general n x n matrix A into P A = L U                    //
00022 //                                                                       //
00023 // where P is a permutation matrix, L is unit lower triangular and U     //
00024 // is upper triangular.                                                  //
00025 // L is stored in the strict lower triangular part of the matrix fLU.    //
00026 // The diagonal elements of L are unity and are not stored.              //
00027 // U is stored in the diagonal and upper triangular part of the matrix   //
00028 // fU.                                                                   //
00029 // P is stored in the index array fIndex : j = fIndex[i] indicates that  //
00030 // row j and row i should be swapped .                                   //
00031 //                                                                       //
00032 // fSign gives the sign of the permutation, (-1)^n, where n is the       //
00033 // number of interchanges in the permutation.                            //
00034 //                                                                       //
00035 // fLU has the same indexing range as matrix A .                         //
00036 //                                                                       //
00037 // The decomposition fails if a diagonal element of abs(fLU) is == 0,    //
00038 // The matrix fUL is made invalid .                                      //
00039 //                                                                       //
00040 ///////////////////////////////////////////////////////////////////////////
00041 
00042 //______________________________________________________________________________
00043 TDecompLU::TDecompLU()
00044 {
00045 // Default constructor
00046 
00047    fSign = 0.0;
00048    fNIndex = 0;
00049    fIndex = 0;
00050    fImplicitPivot = 0;
00051 }
00052 
00053 //______________________________________________________________________________
00054 TDecompLU::TDecompLU(Int_t nrows)
00055 {
00056 // Constructor for (nrows x nrows) matrix
00057 
00058    fSign = 1.0;
00059    fNIndex = nrows;
00060    fIndex = new Int_t[fNIndex];
00061    memset(fIndex,0,fNIndex*sizeof(Int_t));
00062    fImplicitPivot = 0;
00063    fLU.ResizeTo(nrows,nrows);
00064 }
00065 
00066 //______________________________________________________________________________
00067 TDecompLU::TDecompLU(Int_t row_lwb,Int_t row_upb)
00068 {
00069 // Constructor for ([row_lwb..row_upb] x [row_lwb..row_upb]) matrix
00070 
00071    const Int_t nrows = row_upb-row_lwb+1;
00072    fSign = 1.0;
00073    fNIndex = nrows;
00074    fIndex = new Int_t[fNIndex];
00075    memset(fIndex,0,fNIndex*sizeof(Int_t));
00076    fRowLwb = row_lwb;
00077    fColLwb = row_lwb;
00078    fImplicitPivot = 0;
00079    fLU.ResizeTo(row_lwb,row_lwb+nrows-1,row_lwb,row_lwb+nrows-1);
00080 }
00081 
00082 //______________________________________________________________________________
00083 TDecompLU::TDecompLU(const TMatrixD &a,Double_t tol,Int_t implicit)
00084 {
00085 // Constructor for matrix a
00086 
00087    R__ASSERT(a.IsValid());
00088 
00089    if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
00090       Error("TDecompLU(const TMatrixD &","matrix should be square");
00091       return;
00092    }
00093 
00094    SetBit(kMatrixSet);
00095    fCondition = a.Norm1();
00096    fImplicitPivot = implicit;
00097    fTol = a.GetTol();
00098    if (tol > 0)
00099       fTol = tol;
00100 
00101    fSign = 1.0;
00102    fNIndex = a.GetNcols();
00103    fIndex = new Int_t[fNIndex];
00104    memset(fIndex,0,fNIndex*sizeof(Int_t));
00105 
00106    fRowLwb = a.GetRowLwb();
00107    fColLwb = a.GetColLwb();
00108    fLU.ResizeTo(a);
00109    fLU = a;
00110 }
00111 
00112 //______________________________________________________________________________
00113 TDecompLU::TDecompLU(const TDecompLU &another) : TDecompBase(another)
00114 {
00115 // Copy constructor
00116 
00117    fNIndex = 0;
00118    fIndex  = 0;
00119    *this = another;
00120 }
00121 
00122 //______________________________________________________________________________
00123 Bool_t TDecompLU::Decompose()
00124 {
00125 // Matrix A is decomposed in components U and L so that P * A = U * L
00126 // If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular
00127 
00128    if (TestBit(kDecomposed)) return kTRUE;
00129 
00130    if ( !TestBit(kMatrixSet) ) {
00131       Error("Decompose()","Matrix has not been set");
00132       return kFALSE;
00133    }
00134 
00135    Int_t nrZeros = 0;
00136    Bool_t ok;
00137    if (fImplicitPivot)
00138       ok = DecomposeLUCrout(fLU,fIndex,fSign,fTol,nrZeros);
00139    else
00140       ok = DecomposeLUGauss(fLU,fIndex,fSign,fTol,nrZeros);
00141 
00142    if (!ok) SetBit(kSingular);
00143    else     SetBit(kDecomposed);
00144 
00145    return ok;
00146 }
00147 
00148 //______________________________________________________________________________
00149 const TMatrixD TDecompLU::GetMatrix()
00150 {
00151 // Reconstruct the original matrix using the decomposition parts
00152 
00153    if (TestBit(kSingular)) {
00154       Error("GetMatrix()","Matrix is singular");
00155       return TMatrixD();
00156    }
00157    if ( !TestBit(kDecomposed) ) {
00158       if (!Decompose()) {
00159          Error("GetMatrix()","Decomposition failed");
00160          return TMatrixD();
00161       }
00162    }
00163 
00164    TMatrixD mL = fLU;
00165    TMatrixD mU = fLU;
00166    Double_t * const pU = mU.GetMatrixArray();
00167    Double_t * const pL = mL.GetMatrixArray();
00168    const Int_t n = fLU.GetNcols();
00169    for (Int_t irow = 0; irow < n; irow++) {
00170       const Int_t off_row = irow*n;
00171          for (Int_t icol = 0; icol < n; icol++) {
00172             if (icol > irow)      pL[off_row+icol] = 0.;
00173             else if (icol < irow) pU[off_row+icol] = 0.;
00174             else                  pL[off_row+icol] = 1.;
00175          }
00176       }
00177 
00178    TMatrixD a = mL*mU;
00179 
00180    // swap rows
00181 
00182    Double_t * const pA = a.GetMatrixArray();
00183    for (Int_t i = n-1; i >= 0; i--) {
00184       const Int_t j = fIndex[i];
00185       if (j != i) {
00186          const Int_t off_j = j*n;
00187          const Int_t off_i = i*n;
00188          for (Int_t k = 0; k < n; k++) {
00189             const Double_t tmp = pA[off_j+k];
00190             pA[off_j+k] = pA[off_i+k];
00191             pA[off_i+k] = tmp;
00192          }
00193       }
00194    }
00195 
00196    return a;
00197 }
00198 
00199 //______________________________________________________________________________
00200 void TDecompLU::SetMatrix(const TMatrixD &a)
00201 {
00202 // Set matrix to be decomposed
00203 
00204    R__ASSERT(a.IsValid());
00205 
00206    ResetStatus();
00207    if (a.GetNrows() != a.GetNcols() || a.GetRowLwb() != a.GetColLwb()) {
00208       Error("TDecompLU(const TMatrixD &","matrix should be square");
00209       return;
00210    }
00211 
00212    SetBit(kMatrixSet);
00213    fCondition = a.Norm1();
00214 
00215    fSign = 1.0;
00216    if (fNIndex != a.GetNcols()) {
00217       fNIndex = a.GetNcols();
00218       delete [] fIndex;
00219       fIndex = new Int_t[fNIndex];
00220       memset(fIndex,0,fNIndex*sizeof(Int_t));
00221    }
00222 
00223    fRowLwb = a.GetRowLwb();
00224    fColLwb = a.GetColLwb();
00225    fLU.ResizeTo(a);
00226    fLU = a;
00227 }
00228 
00229 //______________________________________________________________________________
00230 Bool_t TDecompLU::Solve(TVectorD &b)
00231 {
00232 // Solve Ax=b assuming the LU form of A is stored in fLU, but assume b has *not*
00233 // been transformed.  Solution returned in b.
00234 
00235    R__ASSERT(b.IsValid());
00236    if (TestBit(kSingular)) {
00237       Error("Solve()","Matrix is singular");
00238       return kFALSE;
00239    }
00240    if ( !TestBit(kDecomposed) ) {
00241       if (!Decompose()) {
00242          Error("Solve()","Decomposition failed");
00243          return kFALSE;
00244       }
00245    }
00246 
00247    if (fLU.GetNrows() != b.GetNrows() || fLU.GetRowLwb() != b.GetLwb()) {
00248       Error("Solve(TVectorD &","vector and matrix incompatible");
00249       return kFALSE;
00250    }
00251 
00252    const Int_t n = fLU.GetNrows();
00253 
00254    const Double_t *pLU = fLU.GetMatrixArray();
00255    Double_t *pb  = b.GetMatrixArray();
00256 
00257    Int_t i;
00258 
00259    // Check for zero diagonals
00260    for (i = 0; i < n ; i++) {
00261       const Int_t off_i = i*n;
00262       if (TMath::Abs(pLU[off_i+i]) < fTol) {
00263          Error("Solve(TVectorD &b)","LU[%d,%d]=%.4e < %.4e",i,i,pLU[off_i+i],fTol);
00264          return kFALSE;
00265       }
00266    }
00267 
00268    // Transform b allowing for leading zeros
00269    Int_t nonzero = -1;
00270    for (i = 0; i < n; i++) {
00271       const Int_t off_i = i*n;
00272       const Int_t iperm = fIndex[i];
00273       Double_t r = pb[iperm];
00274       pb[iperm] = pb[i];
00275       if (nonzero >= 0)
00276          for (Int_t j = nonzero; j < i; j++)
00277             r -= pLU[off_i+j]*pb[j];
00278       else if (r != 0.0)
00279          nonzero = i;
00280       pb[i] = r;
00281    }
00282 
00283    // Backward substitution
00284    for (i = n-1; i >= 0; i--) {
00285       const Int_t off_i = i*n;
00286       Double_t r = pb[i];
00287       for (Int_t j = i+1; j < n; j++)
00288          r -= pLU[off_i+j]*pb[j];
00289       pb[i] = r/pLU[off_i+i];
00290    }
00291 
00292    return kTRUE;
00293 }
00294 
00295 //______________________________________________________________________________
00296 Bool_t TDecompLU::Solve(TMatrixDColumn &cb)
00297 {
00298 // Solve Ax=b assuming the LU form of A is stored in fLU, but assume b has *not*
00299 // been transformed.  Solution returned in b.
00300 
00301    TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
00302    R__ASSERT(b->IsValid());
00303    if (TestBit(kSingular)) {
00304       Error("Solve()","Matrix is singular");
00305       return kFALSE;
00306    }
00307    if ( !TestBit(kDecomposed) ) {
00308       if (!Decompose()) {
00309          Error("Solve()","Decomposition failed");
00310          return kFALSE;
00311       }
00312    }
00313 
00314    if (fLU.GetNrows() != b->GetNrows() || fLU.GetRowLwb() != b->GetRowLwb()) {
00315       Error("Solve(TMatrixDColumn &","vector and matrix incompatible");
00316       return kFALSE;
00317    }
00318 
00319    const Int_t     n   = fLU.GetNrows();
00320    const Double_t *pLU = fLU.GetMatrixArray();
00321  
00322    Int_t i;
00323 
00324    // Check for zero diagonals
00325    for (i = 0; i < n ; i++) {
00326       const Int_t off_i = i*n;
00327       if (TMath::Abs(pLU[off_i+i]) < fTol) {
00328          Error("Solve(TMatrixDColumn &cb)","LU[%d,%d]=%.4e < %.4e",i,i,pLU[off_i+i],fTol);
00329          return kFALSE;
00330       }
00331    }
00332 
00333    Double_t *pcb = cb.GetPtr();
00334    const Int_t     inc = cb.GetInc();
00335 
00336    // Transform b allowing for leading zeros
00337    Int_t nonzero = -1;
00338    for (i = 0; i < n; i++) {
00339       const Int_t off_i  = i*n;
00340       const Int_t off_i2 = i*inc;
00341       const Int_t iperm = fIndex[i];
00342       const Int_t off_iperm = iperm*inc;
00343       Double_t r = pcb[off_iperm];
00344       pcb[off_iperm] = pcb[off_i2];
00345       if (nonzero >=0)
00346          for (Int_t j = nonzero; j <= i-1; j++)
00347             r -= pLU[off_i+j]*pcb[j*inc];
00348       else if (r != 0.0)
00349          nonzero = i;
00350       pcb[off_i2] = r;
00351    }
00352 
00353    // Backward substitution
00354    for (i = n-1; i >= 0; i--) {
00355       const Int_t off_i  = i*n;
00356       const Int_t off_i2 = i*inc;
00357       Double_t r = pcb[off_i2];
00358       for (Int_t j = i+1; j < n; j++)
00359          r -= pLU[off_i+j]*pcb[j*inc];
00360       pcb[off_i2] = r/pLU[off_i+i];
00361    }
00362 
00363    return kTRUE;
00364 }
00365 
00366 
00367 //______________________________________________________________________________
00368 Bool_t TDecompLU::TransSolve(TVectorD &b)
00369 {
00370 // Solve A^T x=b assuming the LU form of A^T is stored in fLU, but assume b has *not*
00371 // been transformed.  Solution returned in b.
00372 
00373    R__ASSERT(b.IsValid());
00374    if (TestBit(kSingular)) {
00375       Error("TransSolve()","Matrix is singular");
00376       return kFALSE;
00377    }
00378    if ( !TestBit(kDecomposed) ) {
00379       if (!Decompose()) {
00380          Error("TransSolve()","Decomposition failed");
00381          return kFALSE;
00382       }
00383    }
00384 
00385    if (fLU.GetNrows() != b.GetNrows() || fLU.GetRowLwb() != b.GetLwb()) {
00386       Error("TransSolve(TVectorD &","vector and matrix incompatible");
00387       return kFALSE;
00388    }
00389 
00390    const Int_t n = fLU.GetNrows();
00391 
00392    const Double_t *pLU = fLU.GetMatrixArray();
00393    Double_t *pb  = b.GetMatrixArray();
00394 
00395    Int_t i;
00396 
00397    // Check for zero diagonals
00398    for (i = 0; i < n ; i++) {
00399       const Int_t off_i = i*n;
00400       if (TMath::Abs(pLU[off_i+i]) < fTol) {
00401          Error("TransSolve(TVectorD &b)","LU[%d,%d]=%.4e < %.4e",i,i,pLU[off_i+i],fTol);
00402          return kFALSE;
00403       }
00404    }
00405 
00406    // Forward Substitution
00407    for (i = 0; i < n; i++) {
00408       const Int_t off_i = i*n;
00409       Double_t r = pb[i];
00410       for (Int_t j = 0; j < i ; j++) {
00411          const Int_t off_j = j*n;
00412          r -= pLU[off_j+i]*pb[j];
00413       }
00414       pb[i] = r/pLU[off_i+i];
00415    }
00416 
00417    // Transform b allowing for leading zeros
00418    Int_t nonzero = -1;
00419    for (i = n-1 ; i >= 0; i--) {
00420       Double_t r = pb[i];
00421       if (nonzero >= 0) {
00422          for (Int_t j = i+1; j <= nonzero; j++) {
00423             const Int_t off_j = j*n;
00424             r -= pLU[off_j+i]*pb[j];
00425          }
00426       } else if (r != 0.0)
00427          nonzero = i;
00428       const Int_t iperm = fIndex[i];
00429       pb[i]     = pb[iperm];
00430       pb[iperm] = r;
00431    }
00432 
00433    return kTRUE;
00434 }
00435 
00436 //______________________________________________________________________________
00437 Bool_t TDecompLU::TransSolve(TMatrixDColumn &cb)
00438 {
00439 // Solve A^T x=b assuming the LU form of A^T is stored in fLU, but assume b has *not*
00440 // been transformed.  Solution returned in b.
00441 
00442    TMatrixDBase *b = const_cast<TMatrixDBase *>(cb.GetMatrix());
00443    R__ASSERT(b->IsValid());
00444    if (TestBit(kSingular)) {
00445       Error("TransSolve()","Matrix is singular");
00446       return kFALSE;
00447    }
00448    if ( !TestBit(kDecomposed) ) {
00449       if (!Decompose()) {
00450          Error("TransSolve()","Decomposition failed");
00451          return kFALSE;
00452       }
00453    }
00454 
00455    if (fLU.GetNrows() != b->GetNrows() || fLU.GetRowLwb() != b->GetRowLwb()) {
00456       Error("TransSolve(TMatrixDColumn &","vector and matrix incompatible");
00457       return kFALSE;
00458    }
00459 
00460    const Int_t n   = fLU.GetNrows();
00461    const Int_t lwb = fLU.GetRowLwb();
00462  
00463    const Double_t *pLU = fLU.GetMatrixArray();
00464 
00465    Int_t i;
00466 
00467    // Check for zero diagonals
00468    for (i = 0; i < n ; i++) {
00469       const Int_t off_i = i*n;
00470       if (TMath::Abs(pLU[off_i+i]) < fTol) {
00471          Error("TransSolve(TMatrixDColumn &cb)","LU[%d,%d]=%.4e < %.4e",i,i,pLU[off_i+i],fTol);
00472          return kFALSE;
00473       }
00474    }
00475 
00476    // Forward Substitution
00477    for (i = 0; i < n; i++) {
00478       const Int_t off_i = i*n;
00479       Double_t r = cb(i+lwb);
00480       for (Int_t j = 0; j < i ; j++) {
00481          const Int_t off_j = j*n;
00482          r -= pLU[off_j+i]*cb(j+lwb);
00483       }
00484       cb(i+lwb) = r/pLU[off_i+i];
00485    }
00486 
00487    // Transform b allowing for leading zeros
00488    Int_t nonzero = -1;
00489    for (i = n-1 ; i >= 0; i--) {
00490       Double_t r = cb(i+lwb);
00491       if (nonzero >= 0) {
00492          for (Int_t j = i+1; j <= nonzero; j++) {
00493             const Int_t off_j = j*n;
00494             r -= pLU[off_j+i]*cb(j+lwb);
00495          }
00496       } else if (r != 0.0)
00497          nonzero = i;
00498       const Int_t iperm = fIndex[i];
00499       cb(i+lwb)     = cb(iperm+lwb);
00500       cb(iperm+lwb) = r;
00501    }
00502 
00503    return kTRUE;
00504 }
00505 
00506 //______________________________________________________________________________
00507 void TDecompLU::Det(Double_t &d1,Double_t &d2)
00508 {
00509 // Calculate determinant det = d1*TMath::Power(2.,d2) 
00510 
00511    if ( !TestBit(kDetermined) ) {
00512       if ( !TestBit(kDecomposed) )
00513          Decompose();
00514       TDecompBase::Det(d1,d2);
00515       fDet1 *= fSign;
00516       SetBit(kDetermined);
00517    }
00518    d1 = fDet1;
00519    d2 = fDet2;
00520 }
00521 
00522 //______________________________________________________________________________
00523 Bool_t TDecompLU::Invert(TMatrixD &inv)
00524 {
00525 // For a matrix A(m,m), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
00526 // (m x m) Ainv is returned .
00527 
00528    if (inv.GetNrows()  != GetNrows()  || inv.GetNcols()  != GetNcols() ||
00529         inv.GetRowLwb() != GetRowLwb() || inv.GetColLwb() != GetColLwb()) {
00530       Error("Invert(TMatrixD &","Input matrix has wrong shape");
00531       return kFALSE;
00532    }
00533 
00534    inv.UnitMatrix();
00535    const Bool_t status = MultiSolve(inv);
00536 
00537    return status;
00538 }
00539 
00540 //______________________________________________________________________________
00541 TMatrixD TDecompLU::Invert(Bool_t &status)
00542 {
00543 // For a matrix A(m,n), its inverse A_inv is defined as A * A_inv = A_inv * A = unit
00544 // (n x m) Ainv is returned .
00545 
00546    const Int_t rowLwb = GetRowLwb();
00547    const Int_t rowUpb = rowLwb+GetNrows()-1;
00548 
00549    TMatrixD inv(rowLwb,rowUpb,rowLwb,rowUpb);
00550    inv.UnitMatrix();
00551    status = MultiSolve(inv);
00552 
00553    return inv;
00554 }
00555 
00556 //______________________________________________________________________________
00557 void TDecompLU::Print(Option_t *opt) const
00558 {
00559    //Print internals of this object
00560    TDecompBase::Print(opt);
00561    printf("fImplicitPivot = %d\n",fImplicitPivot);
00562    printf("fSign          = %f\n",fSign);
00563    printf("fIndex:\n");
00564    for (Int_t i = 0; i < fNIndex; i++)
00565       printf("[%d] = %d\n",i,fIndex[i]);
00566    fLU.Print("fLU");
00567 }
00568 
00569 //______________________________________________________________________________
00570 TDecompLU &TDecompLU::operator=(const TDecompLU &source)
00571 {
00572    //assignement operator
00573    if (this != &source) {
00574       TDecompBase::operator=(source);
00575       fLU.ResizeTo(source.fLU);
00576       fLU   = source.fLU;
00577       fSign = source.fSign;
00578       fImplicitPivot = source.fImplicitPivot;
00579       if (fNIndex != source.fNIndex) {
00580          if (fIndex)
00581             delete [] fIndex;
00582          fNIndex = source.fNIndex;
00583          fIndex = new Int_t[fNIndex];
00584       }
00585       if (fIndex) memcpy(fIndex,source.fIndex,fNIndex*sizeof(Int_t));
00586    }
00587    return *this;
00588 }
00589 
00590 //______________________________________________________________________________
00591 Bool_t TDecompLU::DecomposeLUCrout(TMatrixD &lu,Int_t *index,Double_t &sign,
00592                                    Double_t tol,Int_t &nrZeros)
00593 {
00594 // Crout/Doolittle algorithm of LU decomposing a square matrix, with implicit partial
00595 // pivoting.  The decomposition is stored in fLU: U is explicit in the upper triag
00596 // and L is in multiplier form in the subdiagionals .
00597 // Row permutations are mapped out in fIndex. fSign, used for calculating the
00598 // determinant, is +/- 1 for even/odd row permutations. .
00599 
00600    const Int_t     n     = lu.GetNcols();
00601    Double_t *pLU   = lu.GetMatrixArray();
00602 
00603    Double_t work[kWorkMax];
00604    Bool_t isAllocated = kFALSE;
00605    Double_t *scale = work;
00606    if (n > kWorkMax) {
00607       isAllocated = kTRUE;
00608       scale = new Double_t[n];
00609    }
00610 
00611    sign    = 1.0;
00612    nrZeros = 0;
00613    // Find implicit scaling factors for each row
00614    for (Int_t i = 0; i < n ; i++) {
00615       const Int_t off_i = i*n;
00616       Double_t max = 0.0;
00617       for (Int_t j = 0; j < n; j++) {
00618          const Double_t tmp = TMath::Abs(pLU[off_i+j]);
00619          if (tmp > max)
00620             max = tmp;
00621       }
00622       scale[i] = (max == 0.0 ? 0.0 : 1.0/max);
00623    }
00624 
00625    for (Int_t j = 0; j < n; j++) {
00626       const Int_t off_j = j*n;
00627       // Run down jth column from top to diag, to form the elements of U.
00628       for (Int_t i = 0; i < j; i++) {
00629          const Int_t off_i = i*n;
00630          Double_t r = pLU[off_i+j];
00631          for (Int_t k = 0; k < i; k++) {
00632             const Int_t off_k = k*n;
00633             r -= pLU[off_i+k]*pLU[off_k+j];
00634          }
00635          pLU[off_i+j] = r;
00636       }
00637 
00638       // Run down jth subdiag to form the residuals after the elimination of
00639       // the first j-1 subdiags.  These residuals divided by the appropriate
00640       // diagonal term will become the multipliers in the elimination of the jth.
00641       // subdiag. Find fIndex of largest scaled term in imax.
00642 
00643       Double_t max = 0.0;
00644       Int_t imax = 0;
00645       for (Int_t i = j; i < n; i++) {
00646          const Int_t off_i = i*n;
00647          Double_t r = pLU[off_i+j];
00648          for (Int_t k = 0; k < j; k++) {
00649             const Int_t off_k = k*n;
00650             r -= pLU[off_i+k]*pLU[off_k+j];
00651          }
00652          pLU[off_i+j] = r;
00653          const Double_t tmp = scale[i]*TMath::Abs(r);
00654          if (tmp >= max) {
00655             max = tmp;
00656             imax = i;
00657          }
00658       }
00659 
00660       // Permute current row with imax
00661       if (j != imax) {
00662          const Int_t off_imax = imax*n;
00663          for (Int_t k = 0; k < n; k++ ) {
00664             const Double_t tmp = pLU[off_imax+k];
00665             pLU[off_imax+k] = pLU[off_j+k];
00666             pLU[off_j+k]    = tmp;
00667          }
00668          sign = -sign;
00669          scale[imax] = scale[j];
00670       }
00671       index[j] = imax;
00672 
00673       // If diag term is not zero divide subdiag to form multipliers.
00674       if (pLU[off_j+j] != 0.0) {
00675          if (TMath::Abs(pLU[off_j+j]) < tol)
00676             nrZeros++;
00677          if (j != n-1) {
00678             const Double_t tmp = 1.0/pLU[off_j+j];
00679             for (Int_t i = j+1; i < n; i++) {
00680                const Int_t off_i = i*n;
00681                pLU[off_i+j] *= tmp;
00682             }
00683          }
00684       } else {
00685          ::Error("TDecompLU::DecomposeLUCrout","matrix is singular");
00686          if (isAllocated)  delete [] scale;
00687          return kFALSE;
00688       }
00689    }
00690 
00691    if (isAllocated)
00692       delete [] scale;
00693 
00694    return kTRUE;
00695 }
00696 
00697 //______________________________________________________________________________
00698 Bool_t TDecompLU::DecomposeLUGauss(TMatrixD &lu,Int_t *index,Double_t &sign,
00699                                    Double_t tol,Int_t &nrZeros)
00700 {
00701 // LU decomposition using Gaussain Elimination with partial pivoting (See Golub &
00702 // Van Loan, Matrix Computations, Algorithm 3.4.1) of a square matrix .
00703 // The decomposition is stored in fLU: U is explicit in the upper triag and L is in
00704 // multiplier form in the subdiagionals . Row permutations are mapped out in fIndex.
00705 // fSign, used for calculating the determinant, is +/- 1 for even/odd row permutations.
00706 // Since this algorithm uses partial pivoting without scaling like in Crout/Doolitle.
00707 // it is somewhat faster but less precise .
00708 
00709    const Int_t     n   = lu.GetNcols();
00710    Double_t *pLU = lu.GetMatrixArray();
00711 
00712    sign    = 1.0;
00713    nrZeros = 0;
00714 
00715    index[n-1] = n-1;
00716    for (Int_t j = 0; j < n-1; j++) {
00717       const Int_t off_j = j*n;
00718 
00719       // Find maximum in the j-th column
00720 
00721       Double_t max = TMath::Abs(pLU[off_j+j]);
00722       Int_t i_pivot = j;
00723 
00724       for (Int_t i = j+1; i < n; i++) {
00725          const Int_t off_i = i*n;
00726          const Double_t mLUij = TMath::Abs(pLU[off_i+j]);
00727 
00728          if (mLUij > max) {
00729             max = mLUij;
00730             i_pivot = i;
00731          }
00732       }
00733   
00734       if (i_pivot != j) {
00735          const Int_t off_ipov = i_pivot*n;
00736          for (Int_t k = 0; k < n; k++ ) {
00737             const Double_t tmp = pLU[off_ipov+k];
00738             pLU[off_ipov+k] = pLU[off_j+k];
00739             pLU[off_j+k]    = tmp;
00740          }
00741          sign = -sign;
00742       }
00743       index[j] = i_pivot;
00744 
00745       const Double_t mLUjj = pLU[off_j+j];
00746 
00747       if (mLUjj != 0.0) {
00748          if (TMath::Abs(mLUjj) < tol)
00749             nrZeros++;
00750          for (Int_t i = j+1; i < n; i++) {
00751             const Int_t off_i = i*n;
00752             const Double_t mLUij = pLU[off_i+j]/mLUjj;
00753             pLU[off_i+j] = mLUij;
00754 
00755             for (Int_t k = j+1; k < n; k++) {
00756                const Double_t mLUik = pLU[off_i+k];
00757                const Double_t mLUjk = pLU[off_j+k];
00758                pLU[off_i+k] = mLUik-mLUij*mLUjk;
00759             }
00760          }
00761       } else {
00762          ::Error("TDecompLU::DecomposeLUGauss","matrix is singular");
00763          return kFALSE;
00764       }
00765    }
00766 
00767    return kTRUE;
00768 }
00769 
00770 //______________________________________________________________________________
00771 Bool_t TDecompLU::InvertLU(TMatrixD &lu,Double_t tol,Double_t *det)
00772 {
00773    // Calculate matrix inversion through in place forward/backward substitution
00774 
00775    if (det)
00776       *det = 0.0;
00777 
00778    if (lu.GetNrows() != lu.GetNcols() || lu.GetRowLwb() != lu.GetColLwb()) {
00779       ::Error("TDecompLU::InvertLU","matrix should be square");
00780       return kFALSE;
00781    }
00782 
00783    const Int_t     n   = lu.GetNcols();
00784    Double_t *pLU = lu.GetMatrixArray();
00785 
00786    Int_t worki[kWorkMax];
00787    Bool_t isAllocatedI = kFALSE;
00788    Int_t *index = worki;
00789    if (n > kWorkMax) {
00790       isAllocatedI = kTRUE;
00791       index = new Int_t[n];
00792    }
00793 
00794    Double_t sign = 1.0;
00795    Int_t nrZeros = 0;
00796    if (!DecomposeLUCrout(lu,index,sign,tol,nrZeros) || nrZeros > 0) {
00797       if (isAllocatedI)
00798          delete [] index;
00799       ::Error("TDecompLU::InvertLU","matrix is singular, %d diag elements < tolerance of %.4e",nrZeros,tol);
00800       return kFALSE;
00801    }
00802 
00803    if (det) {
00804       Double_t d1;
00805       Double_t d2;
00806       const TVectorD diagv = TMatrixDDiag_const(lu);
00807       DiagProd(diagv,tol,d1,d2);
00808       d1 *= sign;
00809       *det = d1*TMath::Power(2.0,d2);
00810    }
00811 
00812    //  Form inv(U).
00813 
00814    Int_t j;
00815 
00816    for (j = 0; j < n; j++) {
00817       const Int_t off_j = j*n;
00818 
00819       pLU[off_j+j] = 1./pLU[off_j+j];
00820       const Double_t mLU_jj = -pLU[off_j+j];
00821 
00822 //    Compute elements 0:j-1 of j-th column.
00823 
00824       Double_t *pX = pLU+j;
00825       Int_t k;
00826       for (k = 0; k <= j-1; k++) {
00827          const Int_t off_k = k*n;
00828          if ( pX[off_k] != 0.0 ) {
00829             const Double_t tmp = pX[off_k];
00830             for (Int_t i = 0; i <= k-1; i++) {
00831                const Int_t off_i = i*n;
00832                pX[off_i] += tmp*pLU[off_i+k];
00833             }
00834             pX[off_k] *= pLU[off_k+k];
00835          }
00836       }
00837       for (k = 0; k <= j-1; k++) {
00838          const Int_t off_k = k*n;
00839          pX[off_k] *= mLU_jj;
00840       }
00841    }
00842 
00843    // Solve the equation inv(A)*L = inv(U) for inv(A).
00844 
00845    Double_t workd[kWorkMax];
00846    Bool_t isAllocatedD = kFALSE;
00847    Double_t *pWorkd = workd;
00848    if (n > kWorkMax) {
00849       isAllocatedD = kTRUE;
00850       pWorkd = new Double_t[n];
00851    }
00852 
00853    for (j = n-1; j >= 0; j--) {
00854 
00855       // Copy current column j of L to WORK and replace with zeros.
00856       for (Int_t i = j+1; i < n; i++) {
00857          const Int_t off_i = i*n;
00858          pWorkd[i] = pLU[off_i+j];
00859          pLU[off_i+j] = 0.0;
00860       }
00861 
00862       // Compute current column of inv(A).
00863 
00864       if (j < n-1) {
00865          const Double_t *mp = pLU+j+1;  // Matrix row ptr
00866          Double_t *tp = pLU+j;          // Target vector ptr
00867 
00868          for (Int_t irow = 0; irow < n; irow++) {
00869             Double_t sum = 0.;
00870             const Double_t *sp = pWorkd+j+1; // Source vector ptr
00871             for (Int_t icol = 0; icol < n-1-j ; icol++)
00872                sum += *mp++ * *sp++;
00873             *tp = -sum + *tp;
00874             mp += j+1;
00875             tp += n;
00876          }
00877       }
00878    }
00879 
00880    if (isAllocatedD)
00881       delete [] pWorkd;
00882 
00883    // Apply column interchanges.
00884    for (j = n-1; j >= 0; j--) {
00885       const Int_t jperm = index[j];
00886       if (jperm != j) {
00887          for (Int_t i = 0; i < n; i++) {
00888             const Int_t off_i = i*n;
00889             const Double_t tmp = pLU[off_i+jperm];
00890             pLU[off_i+jperm] = pLU[off_i+j];
00891             pLU[off_i+j]     = tmp;
00892          }
00893       }
00894    }
00895 
00896    if (isAllocatedI)
00897       delete [] index;
00898 
00899    return kTRUE;
00900 }

Generated on Tue Jul 5 14:36:27 2011 for ROOT_528-00b_version by  doxygen 1.5.1