TMatrixTSparse.cxx

Go to the documentation of this file.
00001 // @(#)root/matrix:$Id: TMatrixTSparse.cxx 35649 2010-09-23 13:18:11Z moneta $
00002 // Authors: Fons Rademakers, Eddy Offermann   Feb 2004
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 //////////////////////////////////////////////////////////////////////////
00013 //                                                                      //
00014 // TMatrixTSparse                                                       //
00015 //                                                                      //
00016 // Template class of a general sparse matrix in the Harwell-Boeing      //
00017 // format                                                               //
00018 //                                                                      //
00019 // Besides the usual shape/size decsriptors of a matrix like fNrows,    //
00020 // fRowLwb,fNcols and fColLwb, we also store a row index, fRowIndex and //
00021 // column index, fColIndex only for those elements unequal zero:        //
00022 //                                                                      //
00023 // fRowIndex[0,..,fNrows]:    Stores for each row the index range of    //
00024 //                            the elements in the data and column array //
00025 // fColIndex[0,..,fNelems-1]: Stores the column number for each data    //
00026 //                            element != 0                              //
00027 //                                                                      //
00028 // As an example how to access all sparse data elements:                //
00029 //                                                                      //
00030 // for (Int_t irow = 0; irow < this->fNrows; irow++) {                  //
00031 //   const Int_t sIndex = fRowIndex[irow];                              //
00032 //   const Int_t eIndex = fRowIndex[irow+1];                            //
00033 //   for (Int_t index = sIndex; index < eIndex; index++) {              //
00034 //     const Int_t icol = fColIndex[index];                             //
00035 //     const Element data = fElements[index];                           //
00036 //     printf("data(%d,%d) = %.4e\n",irow+this->fRowLwb,icol+           //
00037 //                                               this->fColLwb,data);   //
00038 //   }                                                                  //
00039 // }                                                                    //
00040 //                                                                      //
00041 // When checking whether sparse matrices are compatible (like in an     //
00042 // assigment !), not only the shape parameters are compared but also    //
00043 // the sparse structure through fRowIndex and fColIndex .               //
00044 //                                                                      //
00045 // Several methods exist to fill a sparse matrix with data entries.     //
00046 // Most are the same like for dense matrices but some care has to be    //
00047 // taken with regard to performance. In the constructor, always the     //
00048 // shape of the matrix has to be specified in some form . Data can be   //
00049 // entered through the following methods :                              //
00050 // 1. constructor                                                       //
00051 //    TMatrixTSparse(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,         //
00052 //                   Int_t col_upb,Int_t nr_nonzeros,                   //
00053 //                   Int_t *row, Int_t *col,Element *data);            //
00054 //    It uses SetMatrixArray(..), see below                             //
00055 // 2. copy constructors                                                 //
00056 // 3. SetMatrixArray(Int_t nr,Int_t *irow,Int_t *icol,Element *data)   //
00057 //    where it is expected that the irow,icol and data array contain    //
00058 //    nr entries . Only the entries with non-zero data[i] value are     //
00059 //    inserted !                                                        //
00060 // 4. TMatrixTSparse a(n,m); for(....) { a(i,j) = ....                  //
00061 //    This is a very flexible method but expensive :                    //
00062 //    - if no entry for slot (i,j) is found in the sparse index table   //
00063 //      it will be entered, which involves some memory management !     //
00064 //    - before invoking this method in a loop it is smart to first      //
00065 //      set the index table through a call to SetSparseIndex(..)        //
00066 // 5. SetSub(Int_t row_lwb,Int_t col_lwb,const TMatrixTBase &source)    //
00067 //    the matrix to be inserted at position (row_lwb,col_lwb) can be    //
00068 //    both dense or sparse .                                            //
00069 //                                                                      //
00070 //////////////////////////////////////////////////////////////////////////
00071 
00072 #include "TMatrixTSparse.h"
00073 #include "TMatrixT.h"
00074 #include "TMath.h"
00075 
00076 #ifndef R__ALPHA
00077 templateClassImp(TMatrixTSparse)
00078 #endif
00079 
00080 //______________________________________________________________________________
00081 template<class Element>
00082 TMatrixTSparse<Element>::TMatrixTSparse(Int_t no_rows,Int_t no_cols)
00083 {
00084   // Space is allocated for row/column indices and data, but the sparse structure
00085   // information has still to be set !
00086 
00087    Allocate(no_rows,no_cols,0,0,1);
00088 }
00089 
00090 //______________________________________________________________________________
00091 template<class Element>
00092 TMatrixTSparse<Element>::TMatrixTSparse(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb)
00093 {
00094   // Space is allocated for row/column indices and data, but the sparse structure
00095   // information has still to be set !
00096 
00097    Allocate(row_upb-row_lwb+1,col_upb-col_lwb+1,row_lwb,col_lwb,1);
00098 }
00099 
00100 //______________________________________________________________________________
00101 template<class Element>
00102 TMatrixTSparse<Element>::TMatrixTSparse(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
00103                                         Int_t nr,Int_t *row, Int_t *col,Element *data)
00104 {
00105   // Space is allocated for row/column indices and data. Sparse row/column index
00106   // structure together with data is coming from the arrays, row, col and data, resp .
00107 
00108    const Int_t irowmin = TMath::LocMin(nr,row);
00109    const Int_t irowmax = TMath::LocMax(nr,row);
00110    const Int_t icolmin = TMath::LocMin(nr,col);
00111    const Int_t icolmax = TMath::LocMax(nr,col);
00112 
00113    if (row[irowmin] < row_lwb || row[irowmax] > row_upb) {
00114       Error("TMatrixTSparse","Inconsistency between row index and its range");
00115       if (row[irowmin] < row_lwb) {
00116          Info("TMatrixTSparse","row index lower bound adjusted to %d",row[irowmin]);
00117          row_lwb = row[irowmin];
00118       }
00119       if (row[irowmax] > row_upb) {
00120          Info("TMatrixTSparse","row index upper bound adjusted to %d",row[irowmax]);
00121          col_lwb = col[irowmax];
00122       }
00123    }
00124    if (col[icolmin] < col_lwb || col[icolmax] > col_upb) {
00125       Error("TMatrixTSparse","Inconsistency between column index and its range");
00126       if (col[icolmin] < col_lwb) {
00127          Info("TMatrixTSparse","column index lower bound adjusted to %d",col[icolmin]);
00128          col_lwb = col[icolmin];
00129       }
00130       if (col[icolmax] > col_upb) {
00131          Info("TMatrixTSparse","column index upper bound adjusted to %d",col[icolmax]);
00132          col_upb = col[icolmax];
00133       }
00134    }
00135 
00136    Allocate(row_upb-row_lwb+1,col_upb-col_lwb+1,row_lwb,col_lwb,1,nr);
00137 
00138    SetMatrixArray(nr,row,col,data);
00139 }
00140 
00141 //______________________________________________________________________________
00142 template<class Element>
00143 TMatrixTSparse<Element>::TMatrixTSparse(const TMatrixTSparse<Element> &another) : TMatrixTBase<Element>(another)
00144 {
00145    Allocate(another.GetNrows(),another.GetNcols(),another.GetRowLwb(),another.GetColLwb(),1,
00146             another.GetNoElements());
00147    memcpy(fRowIndex,another.GetRowIndexArray(),this->fNrowIndex*sizeof(Int_t));
00148    memcpy(fColIndex,another.GetColIndexArray(),this->fNelems*sizeof(Int_t));
00149 
00150    *this = another;
00151 }
00152 
00153 //______________________________________________________________________________
00154 template<class Element>
00155 TMatrixTSparse<Element>::TMatrixTSparse(const TMatrixT<Element> &another) : TMatrixTBase<Element>(another)
00156 {
00157    const Int_t nr_nonzeros = another.NonZeros();
00158    Allocate(another.GetNrows(),another.GetNcols(),another.GetRowLwb(),another.GetColLwb(),1,nr_nonzeros);
00159    SetSparseIndex(another);
00160    *this = another;
00161 }
00162 
00163 //______________________________________________________________________________
00164 template<class Element>
00165 TMatrixTSparse<Element>::TMatrixTSparse(EMatrixCreatorsOp1 op,const TMatrixTSparse<Element> &prototype)
00166 {
00167   // Create a matrix applying a specific operation to the prototype.
00168   // Supported operations are: kZero, kUnit, kTransposed and kAtA
00169 
00170    R__ASSERT(prototype.IsValid());
00171 
00172    Int_t nr_nonzeros = 0;
00173 
00174    switch(op) {
00175       case kZero:
00176       {
00177          Allocate(prototype.GetNrows(),prototype.GetNcols(),
00178                   prototype.GetRowLwb(),prototype.GetColLwb(),1,nr_nonzeros);
00179          break;
00180       }
00181       case kUnit:
00182       {
00183          const Int_t nrows  = prototype.GetNrows();
00184          const Int_t ncols  = prototype.GetNcols();
00185          const Int_t rowLwb = prototype.GetRowLwb();
00186          const Int_t colLwb = prototype.GetColLwb();
00187          for (Int_t i = rowLwb; i <= rowLwb+nrows-1; i++)
00188             for (Int_t j = colLwb; j <= colLwb+ncols-1; j++)
00189                if (i==j) nr_nonzeros++;
00190             Allocate(nrows,ncols,rowLwb,colLwb,1,nr_nonzeros);
00191             UnitMatrix();
00192             break;
00193       }
00194       case kTransposed:
00195       {
00196           Allocate(prototype.GetNcols(), prototype.GetNrows(),
00197                    prototype.GetColLwb(),prototype.GetRowLwb(),1,prototype.GetNoElements());
00198           Transpose(prototype);
00199           break;
00200       }
00201       case kAtA:
00202       {
00203          const TMatrixTSparse<Element> at(TMatrixTSparse<Element>::kTransposed,prototype);
00204          AMultBt(at,at,1);
00205          break;
00206       }
00207 
00208       default:
00209          Error("TMatrixTSparse(EMatrixCreatorOp1)","operation %d not yet implemented", op);
00210    }
00211 }
00212 
00213 //______________________________________________________________________________
00214 template<class Element>
00215 TMatrixTSparse<Element>::TMatrixTSparse(const TMatrixTSparse<Element> &a,EMatrixCreatorsOp2 op,const TMatrixTSparse<Element> &b)
00216 {
00217   // Create a matrix applying a specific operation to two prototypes.
00218   // Supported operations are: kMult (a*b), kMultTranspose (a*b'), kPlus (a+b), kMinus (a-b)
00219 
00220    R__ASSERT(a.IsValid());
00221    R__ASSERT(b.IsValid());
00222 
00223    switch(op) {
00224       case kMult:
00225          AMultB(a,b,1);
00226          break;
00227 
00228       case kMultTranspose:
00229          AMultBt(a,b,1);
00230          break;
00231 
00232       case kPlus:
00233          APlusB(a,b,1);
00234          break;
00235 
00236       case kMinus:
00237          AMinusB(a,b,1);
00238          break;
00239 
00240       default:
00241          Error("TMatrixTSparse(EMatrixCreatorOp2)", "operation %d not yet implemented",op);
00242    }
00243 }
00244 
00245 //______________________________________________________________________________
00246 template<class Element>
00247 TMatrixTSparse<Element>::TMatrixTSparse(const TMatrixTSparse<Element> &a,EMatrixCreatorsOp2 op,const TMatrixT<Element> &b)
00248 {
00249   // Create a matrix applying a specific operation to two prototypes.
00250   // Supported operations are: kMult (a*b), kMultTranspose (a*b'), kPlus (a+b), kMinus (a-b)
00251 
00252    R__ASSERT(a.IsValid());
00253    R__ASSERT(b.IsValid());
00254 
00255    switch(op) {
00256       case kMult:
00257          AMultB(a,b,1);
00258          break;
00259 
00260       case kMultTranspose:
00261          AMultBt(a,b,1);
00262          break;
00263 
00264       case kPlus:
00265          APlusB(a,b,1);
00266          break;
00267 
00268       case kMinus:
00269          AMinusB(a,b,1);
00270          break;
00271 
00272       default:
00273          Error("TMatrixTSparse(EMatrixCreatorOp2)", "operation %d not yet implemented",op);
00274    }
00275 }
00276 
00277 //______________________________________________________________________________
00278 template<class Element>
00279 TMatrixTSparse<Element>::TMatrixTSparse(const TMatrixT<Element> &a,EMatrixCreatorsOp2 op,const TMatrixTSparse<Element> &b)
00280 {
00281   // Create a matrix applying a specific operation to two prototypes.
00282   // Supported operations are: kMult (a*b), kMultTranspose (a*b'), kPlus (a+b), kMinus (a-b)
00283 
00284    R__ASSERT(a.IsValid());
00285    R__ASSERT(b.IsValid());
00286 
00287    switch(op) {
00288       case kMult:
00289          AMultB(a,b,1);
00290          break;
00291 
00292       case kMultTranspose:
00293          AMultBt(a,b,1);
00294          break;
00295 
00296       case kPlus:
00297          APlusB(a,b,1);
00298          break;
00299 
00300       case kMinus:
00301          AMinusB(a,b,1);
00302          break;
00303 
00304       default:
00305          Error("TMatrixTSparse(EMatrixCreatorOp2)", "operation %d not yet implemented",op);
00306    }
00307 }
00308 
00309 //______________________________________________________________________________
00310 template<class Element>
00311 void TMatrixTSparse<Element>::Allocate(Int_t no_rows,Int_t no_cols,Int_t row_lwb,Int_t col_lwb,
00312                               Int_t init,Int_t nr_nonzeros)
00313 {
00314   // Allocate new matrix. Arguments are number of rows, columns, row lowerbound (0 default)
00315   // and column lowerbound (0 default), 0 initialization flag and number of non-zero
00316   // elements (only relevant for sparse format).
00317 
00318    if ( (nr_nonzeros > 0 && (no_rows == 0 || no_cols == 0)) ||
00319        (no_rows < 0 || no_cols < 0 || nr_nonzeros < 0) )
00320    {
00321       Error("Allocate","no_rows=%d no_cols=%d non_zeros=%d",no_rows,no_cols,nr_nonzeros);
00322       this->Invalidate();
00323       return;
00324    }
00325 
00326    this->MakeValid();
00327    this->fNrows     = no_rows;
00328    this->fNcols     = no_cols;
00329    this->fRowLwb    = row_lwb;
00330    this->fColLwb    = col_lwb;
00331    this->fNrowIndex = this->fNrows+1;
00332    this->fNelems    = nr_nonzeros;
00333    this->fIsOwner   = kTRUE;
00334    this->fTol       = std::numeric_limits<Element>::epsilon();
00335 
00336    fRowIndex = new Int_t[this->fNrowIndex];
00337    if (init)
00338       memset(fRowIndex,0,this->fNrowIndex*sizeof(Int_t));
00339 
00340    if (this->fNelems > 0) {
00341       fElements = new Element[this->fNelems];
00342       fColIndex = new Int_t   [this->fNelems];
00343       if (init) {
00344          memset(fElements,0,this->fNelems*sizeof(Element));
00345          memset(fColIndex,0,this->fNelems*sizeof(Int_t));
00346       }
00347    } else {
00348       fElements = 0;
00349       fColIndex = 0;
00350    }
00351 }
00352 
00353 //______________________________________________________________________________
00354 template<class Element>
00355 TMatrixTBase<Element> &TMatrixTSparse<Element>::InsertRow(Int_t rown,Int_t coln,const Element *v,Int_t n)
00356 {
00357   // Insert in row rown, n elements of array v at column coln
00358 
00359    const Int_t arown = rown-this->fRowLwb;
00360    const Int_t acoln = coln-this->fColLwb;
00361    const Int_t nr = (n > 0) ? n : this->fNcols;
00362 
00363    if (gMatrixCheck) {
00364       if (arown >= this->fNrows || arown < 0) {
00365          Error("InsertRow","row %d out of matrix range",rown);
00366          return *this;
00367       }
00368 
00369       if (acoln >= this->fNcols || acoln < 0) {
00370          Error("InsertRow","column %d out of matrix range",coln);
00371          return *this;
00372       }
00373 
00374       if (acoln+nr > this->fNcols || nr < 0) {
00375          Error("InsertRow","row length %d out of range",nr);
00376          return *this;
00377       }
00378    }
00379 
00380    const Int_t sIndex = fRowIndex[arown];
00381    const Int_t eIndex = fRowIndex[arown+1];
00382 
00383    // check first how many slots are available from [acoln,..,acoln+nr-1]
00384    // also note lIndex and rIndex so that [sIndex..lIndex] and [rIndex..eIndex-1]
00385    // contain the row entries except for the region to be inserted
00386 
00387    Int_t nslots = 0;
00388    Int_t lIndex = sIndex-1;
00389    Int_t rIndex = sIndex-1;
00390    Int_t index;
00391    for (index = sIndex; index < eIndex; index++) {
00392       const Int_t icol = fColIndex[index];
00393       rIndex++;
00394       if (icol >= acoln+nr) break;
00395       if (icol >= acoln) nslots++;
00396       else               lIndex++;
00397    }
00398 
00399    const Int_t nelems_old = this->fNelems;
00400    const Int_t    *colIndex_old = fColIndex;
00401    const Element *elements_old = fElements;
00402 
00403    const Int_t ndiff = nr-nslots;
00404    this->fNelems += ndiff;
00405    fColIndex = new Int_t[this->fNelems];
00406    fElements = new Element[this->fNelems];
00407 
00408    for (Int_t irow = arown+1; irow < this->fNrows+1; irow++)
00409       fRowIndex[irow] += ndiff;
00410 
00411    if (lIndex+1 > 0) {
00412       memmove(fColIndex,colIndex_old,(lIndex+1)*sizeof(Int_t));
00413       memmove(fElements,elements_old,(lIndex+1)*sizeof(Element));
00414    }
00415 
00416    if (nelems_old > 0 && nelems_old-rIndex > 0) {
00417       memmove(fColIndex+rIndex+ndiff,colIndex_old+rIndex,(nelems_old-rIndex)*sizeof(Int_t));
00418       memmove(fElements+rIndex+ndiff,elements_old+rIndex,(nelems_old-rIndex)*sizeof(Element));
00419    }
00420 
00421    index = lIndex+1;
00422    for (Int_t i = 0; i < nr; i++) {
00423       fColIndex[index] = acoln+i;
00424       fElements[index] = v[i];
00425       index++;
00426    }
00427 
00428    if (colIndex_old) delete [] (Int_t*)    colIndex_old;
00429    if (elements_old) delete [] (Element*) elements_old;
00430 
00431    R__ASSERT(this->fNelems == fRowIndex[this->fNrowIndex-1]);
00432 
00433    return *this;
00434 }
00435 
00436 //______________________________________________________________________________
00437 template<class Element>
00438 void TMatrixTSparse<Element>::ExtractRow(Int_t rown, Int_t coln, Element *v,Int_t n) const
00439 {
00440   // Store in array v, n matrix elements of row rown starting at column coln
00441 
00442    const Int_t arown = rown-this->fRowLwb;
00443    const Int_t acoln = coln-this->fColLwb;
00444    const Int_t nr = (n > 0) ? n : this->fNcols;
00445 
00446    if (gMatrixCheck) {
00447       if (arown >= this->fNrows || arown < 0) {
00448          Error("ExtractRow","row %d out of matrix range",rown);
00449          return;
00450       }
00451 
00452       if (acoln >= this->fNcols || acoln < 0) {
00453          Error("ExtractRow","column %d out of matrix range",coln);
00454          return;
00455       }
00456 
00457       if (acoln+nr > this->fNcols || nr < 0) {
00458          Error("ExtractRow","row length %d out of range",nr);
00459          return;
00460       }
00461    }
00462 
00463    const Int_t sIndex = fRowIndex[arown];
00464    const Int_t eIndex = fRowIndex[arown+1];
00465 
00466    memset(v,0,nr*sizeof(Element));
00467    const Int_t   * const pColIndex = GetColIndexArray();
00468    const Element * const pData     = GetMatrixArray();
00469    for (Int_t index = sIndex; index < eIndex; index++) {
00470       const Int_t icol = pColIndex[index];
00471       if (icol < acoln || icol >= acoln+nr) continue;
00472        v[icol-acoln] = pData[index];
00473    }
00474 }
00475 
00476 //______________________________________________________________________________
00477 template<class Element>
00478 void TMatrixTSparse<Element>::AMultBt(const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr)
00479 {
00480   // General matrix multiplication. Create a matrix C such that C = A * B'.
00481   // Note, matrix C is allocated for constr=1.
00482 
00483    if (gMatrixCheck) {
00484       R__ASSERT(a.IsValid());
00485       R__ASSERT(b.IsValid());
00486 
00487       if (a.GetNcols() != b.GetNcols() || a.GetColLwb() != b.GetColLwb()) {
00488          Error("AMultBt","A and B columns incompatible");
00489          return;
00490       }
00491 
00492       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00493          Error("AMultB","this = &a");
00494          return;
00495       }
00496 
00497       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00498          Error("AMultB","this = &b");
00499          return;
00500       }
00501    }
00502 
00503    const Int_t * const pRowIndexa = a.GetRowIndexArray();
00504    const Int_t * const pColIndexa = a.GetColIndexArray();
00505    const Int_t * const pRowIndexb = b.GetRowIndexArray();
00506    const Int_t * const pColIndexb = b.GetColIndexArray();
00507 
00508    Int_t *pRowIndexc;
00509    Int_t *pColIndexc;
00510    if (constr) {
00511       // make a best guess of the sparse structure; it will guarantee
00512       // enough allocated space !
00513 
00514       Int_t nr_nonzero_rowa = 0;
00515       {
00516          for (Int_t irowa = 0; irowa < a.GetNrows(); irowa++)
00517             if (pRowIndexa[irowa] < pRowIndexa[irowa+1])
00518                nr_nonzero_rowa++;
00519       }
00520       Int_t nr_nonzero_rowb = 0;
00521       {
00522          for (Int_t irowb = 0; irowb < b.GetNrows(); irowb++)
00523             if (pRowIndexb[irowb] < pRowIndexb[irowb+1])
00524               nr_nonzero_rowb++;
00525       }
00526 
00527       const Int_t nc = nr_nonzero_rowa*nr_nonzero_rowb; // best guess
00528       Allocate(a.GetNrows(),b.GetNrows(),a.GetRowLwb(),b.GetRowLwb(),1,nc);
00529 
00530       pRowIndexc = this->GetRowIndexArray();
00531       pColIndexc = this->GetColIndexArray();
00532 
00533       pRowIndexc[0] = 0;
00534       Int_t ielem = 0;
00535       for (Int_t irowa = 0; irowa < a.GetNrows(); irowa++) {
00536          pRowIndexc[irowa+1] = pRowIndexc[irowa];
00537          if (pRowIndexa[irowa] >= pRowIndexa[irowa+1]) continue;
00538          for (Int_t irowb = 0; irowb < b.GetNrows(); irowb++) {
00539             if (pRowIndexb[irowb] >= pRowIndexb[irowb+1]) continue;
00540             pRowIndexc[irowa+1]++;
00541             pColIndexc[ielem++] = irowb;
00542          }
00543       }
00544    } else {
00545       pRowIndexc = this->GetRowIndexArray();
00546       pColIndexc = this->GetColIndexArray();
00547    }
00548 
00549    const Element * const pDataa = a.GetMatrixArray();
00550    const Element * const pDatab = b.GetMatrixArray();
00551    Element * const pDatac = this->GetMatrixArray();
00552    Int_t indexc_r = 0;
00553    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00554       const Int_t sIndexa = pRowIndexa[irowc];
00555       const Int_t eIndexa = pRowIndexa[irowc+1];
00556       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00557          const Int_t sIndexb = pRowIndexb[icolc];
00558          const Int_t eIndexb = pRowIndexb[icolc+1];
00559          Element sum = 0.0;
00560          Int_t indexb = sIndexb;
00561          for (Int_t indexa = sIndexa; indexa < eIndexa && indexb < eIndexb; indexa++) {
00562             const Int_t icola = pColIndexa[indexa];
00563             while (indexb < eIndexb && pColIndexb[indexb] <= icola) {
00564                if (icola == pColIndexb[indexb]) {
00565                  sum += pDataa[indexa]*pDatab[indexb];
00566                  break;
00567                }
00568                indexb++;
00569             }
00570          }
00571          if (sum != 0.0) {
00572             pColIndexc[indexc_r] = icolc;
00573             pDatac[indexc_r] = sum;
00574             indexc_r++;
00575          }
00576       }
00577       pRowIndexc[irowc+1] = indexc_r;
00578    }
00579 
00580    if (constr)
00581       SetSparseIndex(indexc_r);
00582 }
00583 
00584 //______________________________________________________________________________
00585 template<class Element>
00586 void TMatrixTSparse<Element>::AMultBt(const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr)
00587 {
00588   // General matrix multiplication. Create a matrix C such that C = A * B'.
00589   // Note, matrix C is allocated for constr=1.
00590 
00591    if (gMatrixCheck) {
00592       R__ASSERT(a.IsValid());
00593       R__ASSERT(b.IsValid());
00594 
00595       if (a.GetNcols() != b.GetNcols() || a.GetColLwb() != b.GetColLwb()) {
00596          Error("AMultBt","A and B columns incompatible");
00597          return;
00598       }
00599 
00600       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00601          Error("AMultB","this = &a");
00602          return;
00603       }
00604 
00605       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00606          Error("AMultB","this = &b");
00607          return;
00608       }
00609    }
00610 
00611    const Int_t * const pRowIndexa = a.GetRowIndexArray();
00612    const Int_t * const pColIndexa = a.GetColIndexArray();
00613 
00614    Int_t *pRowIndexc;
00615    Int_t *pColIndexc;
00616    if (constr) {
00617       // make a best guess of the sparse structure; it will guarantee
00618       // enough allocated space !
00619 
00620       Int_t nr_nonzero_rowa = 0;
00621       {
00622          for (Int_t irowa = 0; irowa < a.GetNrows(); irowa++)
00623             if (pRowIndexa[irowa] < pRowIndexa[irowa+1])
00624                nr_nonzero_rowa++;
00625       }
00626       Int_t nr_nonzero_rowb = b.GetNrows();
00627 
00628       const Int_t nc = nr_nonzero_rowa*nr_nonzero_rowb; // best guess
00629       Allocate(a.GetNrows(),b.GetNrows(),a.GetRowLwb(),b.GetRowLwb(),1,nc);
00630 
00631       pRowIndexc = this->GetRowIndexArray();
00632       pColIndexc = this->GetColIndexArray();
00633 
00634       pRowIndexc[0] = 0;
00635       Int_t ielem = 0;
00636       for (Int_t irowa = 0; irowa < a.GetNrows(); irowa++) {
00637          pRowIndexc[irowa+1] = pRowIndexc[irowa];
00638          for (Int_t irowb = 0; irowb < b.GetNrows(); irowb++) {
00639             pRowIndexc[irowa+1]++;
00640             pColIndexc[ielem++] = irowb;
00641          }
00642       }
00643    } else {
00644       pRowIndexc = this->GetRowIndexArray();
00645       pColIndexc = this->GetColIndexArray();
00646    }
00647 
00648    const Element * const pDataa = a.GetMatrixArray();
00649    const Element * const pDatab = b.GetMatrixArray();
00650    Element * const pDatac = this->GetMatrixArray();
00651    Int_t indexc_r = 0;
00652    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00653       const Int_t sIndexa = pRowIndexa[irowc];
00654       const Int_t eIndexa = pRowIndexa[irowc+1];
00655       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00656          const Int_t off = icolc*b.GetNcols();
00657          Element sum = 0.0;
00658          for (Int_t indexa = sIndexa; indexa < eIndexa; indexa++) {
00659             const Int_t icola = pColIndexa[indexa];
00660             sum += pDataa[indexa]*pDatab[off+icola];
00661          }
00662          if (sum != 0.0) {
00663             pColIndexc[indexc_r] = icolc;
00664             pDatac[indexc_r] = sum;
00665             indexc_r++;
00666          }
00667       }
00668       pRowIndexc[irowc+1]= indexc_r;
00669    }
00670 
00671    if (constr)
00672       SetSparseIndex(indexc_r);
00673 }
00674 
00675 //______________________________________________________________________________
00676 template<class Element>
00677 void TMatrixTSparse<Element>::AMultBt(const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr)
00678 {
00679   // General matrix multiplication. Create a matrix C such that C = A * B'.
00680   // Note, matrix C is allocated for constr=1.
00681 
00682    if (gMatrixCheck) {
00683       R__ASSERT(a.IsValid());
00684       R__ASSERT(b.IsValid());
00685 
00686       if (a.GetNcols() != b.GetNcols() || a.GetColLwb() != b.GetColLwb()) {
00687          Error("AMultBt","A and B columns incompatible");
00688          return;
00689       }
00690 
00691       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00692          Error("AMultB","this = &a");
00693          return;
00694       }
00695 
00696       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00697          Error("AMultB","this = &b");
00698          return;
00699       }
00700    }
00701 
00702    const Int_t * const pRowIndexb = b.GetRowIndexArray();
00703    const Int_t * const pColIndexb = b.GetColIndexArray();
00704 
00705    Int_t *pRowIndexc;
00706    Int_t *pColIndexc;
00707    if (constr) {
00708       // make a best guess of the sparse structure; it will guarantee
00709       // enough allocated space !
00710 
00711       Int_t nr_nonzero_rowa = a.GetNrows();
00712       Int_t nr_nonzero_rowb = 0;
00713       {
00714          for (Int_t irowb = 0; irowb < b.GetNrows(); irowb++)
00715             if (pRowIndexb[irowb] < pRowIndexb[irowb+1])
00716                nr_nonzero_rowb++;
00717       }
00718 
00719       const Int_t nc = nr_nonzero_rowa*nr_nonzero_rowb; // best guess
00720       Allocate(a.GetNrows(),b.GetNrows(),a.GetRowLwb(),b.GetRowLwb(),1,nc);
00721 
00722       pRowIndexc = this->GetRowIndexArray();
00723       pColIndexc = this->GetColIndexArray();
00724 
00725       pRowIndexc[0] = 0;
00726       Int_t ielem = 0;
00727       for (Int_t irowa = 0; irowa < a.GetNrows(); irowa++) {
00728          pRowIndexc[irowa+1] = pRowIndexc[irowa];
00729          for (Int_t irowb = 0; irowb < b.GetNrows(); irowb++) {
00730             if (pRowIndexb[irowb] >= pRowIndexb[irowb+1]) continue;
00731             pRowIndexc[irowa+1]++;
00732             pColIndexc[ielem++] = irowb;
00733          }
00734       }
00735    } else {
00736       pRowIndexc = this->GetRowIndexArray();
00737       pColIndexc = this->GetColIndexArray();
00738    }
00739 
00740    const Element * const pDataa = a.GetMatrixArray();
00741    const Element * const pDatab = b.GetMatrixArray();
00742    Element * const pDatac = this->GetMatrixArray();
00743    Int_t indexc_r = 0;
00744    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00745       const Int_t off = irowc*a.GetNcols();
00746       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00747          const Int_t sIndexb = pRowIndexb[icolc];
00748          const Int_t eIndexb = pRowIndexb[icolc+1];
00749          Element sum = 0.0;
00750          for (Int_t indexb = sIndexb; indexb < eIndexb; indexb++) {
00751             const Int_t icolb = pColIndexb[indexb];
00752             sum += pDataa[off+icolb]*pDatab[indexb];
00753          }
00754          if (sum != 0.0) {
00755             pColIndexc[indexc_r] = icolc;
00756             pDatac[indexc_r] = sum;
00757             indexc_r++;
00758          }
00759       }
00760       pRowIndexc[irowc+1] = indexc_r;
00761    }
00762 
00763    if (constr)
00764       SetSparseIndex(indexc_r);
00765 }
00766 
00767 //______________________________________________________________________________
00768 template<class Element>
00769 void TMatrixTSparse<Element>::APlusB(const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr)
00770 {
00771   // General matrix addition. Create a matrix C such that C = A + B.
00772   // Note, matrix C is allocated for constr=1.
00773 
00774    if (gMatrixCheck) {
00775       R__ASSERT(a.IsValid());
00776       R__ASSERT(b.IsValid());
00777 
00778       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
00779           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
00780          Error("APlusB(const TMatrixTSparse &,const TMatrixTSparse &","matrices not compatible");
00781          return;
00782       }
00783 
00784       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00785          Error("APlusB","this = &a");
00786          return;
00787       }
00788 
00789       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00790          Error("APlusB","this = &b");
00791          return;
00792       }
00793    }
00794 
00795    const Int_t * const pRowIndexa = a.GetRowIndexArray();
00796    const Int_t * const pRowIndexb = b.GetRowIndexArray();
00797    const Int_t * const pColIndexa = a.GetColIndexArray();
00798    const Int_t * const pColIndexb = b.GetColIndexArray();
00799 
00800    if (constr) {
00801       Allocate(a.GetNrows(),a.GetNcols(),a.GetRowLwb(),a.GetColLwb());
00802       SetSparseIndexAB(a,b);
00803    }
00804 
00805    Int_t * const pRowIndexc = this->GetRowIndexArray();
00806    Int_t * const pColIndexc = this->GetColIndexArray();
00807 
00808    const Element * const pDataa = a.GetMatrixArray();
00809    const Element * const pDatab = b.GetMatrixArray();
00810    Element * const pDatac = this->GetMatrixArray();
00811    Int_t indexc_r = 0;
00812    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00813       const Int_t sIndexa = pRowIndexa[irowc];
00814       const Int_t eIndexa = pRowIndexa[irowc+1];
00815       const Int_t sIndexb = pRowIndexb[irowc];
00816       const Int_t eIndexb = pRowIndexb[irowc+1];
00817       Int_t indexa = sIndexa;
00818       Int_t indexb = sIndexb;
00819       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00820          Element sum = 0.0;
00821          while (indexa < eIndexa && pColIndexa[indexa] <= icolc) {
00822             if (icolc == pColIndexa[indexa]) {
00823                sum += pDataa[indexa];
00824                break;
00825             }
00826             indexa++;
00827          }
00828          while (indexb < eIndexb && pColIndexb[indexb] <= icolc) {
00829             if (icolc == pColIndexb[indexb]) {
00830                sum += pDatab[indexb];
00831                break;
00832             }
00833             indexb++;
00834          }
00835 
00836          if (sum != 0.0) {
00837             pColIndexc[indexc_r] = icolc;
00838             pDatac[indexc_r] = sum;
00839             indexc_r++;
00840          }
00841       }
00842       pRowIndexc[irowc+1] = indexc_r;
00843    }
00844 
00845    if (constr)
00846       SetSparseIndex(indexc_r);
00847 }
00848 
00849 //______________________________________________________________________________
00850 template<class Element>
00851 void TMatrixTSparse<Element>::APlusB(const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr)
00852 {
00853   // General matrix addition. Create a matrix C such that C = A + B.
00854   // Note, matrix C is allocated for constr=1.
00855 
00856    if (gMatrixCheck) {
00857       R__ASSERT(a.IsValid());
00858       R__ASSERT(b.IsValid());
00859 
00860       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
00861           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
00862          Error("APlusB(const TMatrixTSparse &,const TMatrixT &","matrices not compatible");
00863          return;
00864       }
00865 
00866       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00867          Error("APlusB","this = &a");
00868          return;
00869       }
00870 
00871       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00872          Error("APlusB","this = &b");
00873          return;
00874       }
00875    }
00876 
00877    if (constr) {
00878       Allocate(a.GetNrows(),a.GetNcols(),a.GetRowLwb(),a.GetColLwb());
00879       SetSparseIndexAB(a,b);
00880    }
00881 
00882    Int_t * const pRowIndexc = this->GetRowIndexArray();
00883    Int_t * const pColIndexc = this->GetColIndexArray();
00884 
00885    const Int_t * const pRowIndexa = a.GetRowIndexArray();
00886    const Int_t * const pColIndexa = a.GetColIndexArray();
00887 
00888    const Element * const pDataa = a.GetMatrixArray();
00889    const Element * const pDatab = b.GetMatrixArray();
00890    Element * const pDatac = this->GetMatrixArray();
00891    Int_t indexc_r = 0;
00892    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00893       const Int_t sIndexa = pRowIndexa[irowc];
00894       const Int_t eIndexa = pRowIndexa[irowc+1];
00895       const Int_t off = irowc*this->GetNcols();
00896       Int_t indexa = sIndexa;
00897       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00898          Element sum = pDatab[off+icolc];
00899          while (indexa < eIndexa && pColIndexa[indexa] <= icolc) {
00900             if (icolc == pColIndexa[indexa]) {
00901                sum += pDataa[indexa];
00902                break;
00903             }
00904             indexa++;
00905          }
00906 
00907          if (sum != 0.0) {
00908             pColIndexc[indexc_r] = icolc;
00909             pDatac[indexc_r] = sum;
00910             indexc_r++;
00911          }
00912       }
00913       pRowIndexc[irowc+1] = indexc_r;
00914    }
00915 
00916    if (constr)
00917       SetSparseIndex(indexc_r);
00918 }
00919 
00920 //______________________________________________________________________________
00921 template<class Element>
00922 void TMatrixTSparse<Element>::AMinusB(const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr)
00923 {
00924   // General matrix subtraction. Create a matrix C such that C = A - B.
00925   // Note, matrix C is allocated for constr=1.
00926 
00927    if (gMatrixCheck) {
00928       R__ASSERT(a.IsValid());
00929       R__ASSERT(b.IsValid());
00930 
00931       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
00932           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
00933          Error("AMinusB(const TMatrixTSparse &,const TMatrixTSparse &","matrices not compatible");
00934          return;
00935       }
00936 
00937       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
00938          Error("AMinusB","this = &a");
00939          return;
00940       }
00941 
00942       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
00943          Error("AMinusB","this = &b");
00944          return;
00945       }
00946    }
00947 
00948    const Int_t * const pRowIndexa = a.GetRowIndexArray();
00949    const Int_t * const pRowIndexb = b.GetRowIndexArray();
00950    const Int_t * const pColIndexa = a.GetColIndexArray();
00951    const Int_t * const pColIndexb = b.GetColIndexArray();
00952 
00953    if (constr) {
00954       Allocate(a.GetNrows(),a.GetNcols(),a.GetRowLwb(),a.GetColLwb());
00955       SetSparseIndexAB(a,b);
00956    }
00957 
00958    Int_t * const pRowIndexc = this->GetRowIndexArray();
00959    Int_t * const pColIndexc = this->GetColIndexArray();
00960 
00961    const Element * const pDataa = a.GetMatrixArray();
00962    const Element * const pDatab = b.GetMatrixArray();
00963    Element * const pDatac = this->GetMatrixArray();
00964    Int_t indexc_r = 0;
00965    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
00966       const Int_t sIndexa = pRowIndexa[irowc];
00967       const Int_t eIndexa = pRowIndexa[irowc+1];
00968       const Int_t sIndexb = pRowIndexb[irowc];
00969       const Int_t eIndexb = pRowIndexb[irowc+1];
00970       Int_t indexa = sIndexa;
00971       Int_t indexb = sIndexb;
00972       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
00973          Element sum = 0.0;
00974          while (indexa < eIndexa && pColIndexa[indexa] <= icolc) {
00975             if (icolc == pColIndexa[indexa]) {
00976                sum += pDataa[indexa];
00977                break;
00978             }
00979             indexa++;
00980          }
00981          while (indexb < eIndexb && pColIndexb[indexb] <= icolc) {
00982             if (icolc == pColIndexb[indexb]) {
00983                sum -= pDatab[indexb];
00984                break;
00985             }
00986             indexb++;
00987          }
00988 
00989          if (sum != 0.0) {
00990             pColIndexc[indexc_r] = icolc;
00991             pDatac[indexc_r] = sum;
00992             indexc_r++;
00993          }
00994       }
00995       pRowIndexc[irowc+1] = indexc_r;
00996    }
00997 
00998    if (constr)
00999       SetSparseIndex(indexc_r);
01000 }
01001 
01002 //______________________________________________________________________________
01003 template<class Element>
01004 void TMatrixTSparse<Element>::AMinusB(const TMatrixTSparse<Element> &a,const TMatrixT<Element> &b,Int_t constr)
01005 {
01006   // General matrix subtraction. Create a matrix C such that C = A - B.
01007   // Note, matrix C is allocated for constr=1.
01008 
01009    if (gMatrixCheck) {
01010       R__ASSERT(a.IsValid());
01011       R__ASSERT(b.IsValid());
01012 
01013       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
01014           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
01015           Error("AMinusB(const TMatrixTSparse &,const TMatrixT &","matrices not compatible");
01016           return;
01017       }
01018 
01019       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
01020          Error("AMinusB","this = &a");
01021          return;
01022       }
01023 
01024       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
01025          Error("AMinusB","this = &b");
01026          return;
01027       }
01028    }
01029 
01030    if (constr) {
01031       Allocate(a.GetNrows(),a.GetNcols(),a.GetRowLwb(),a.GetColLwb());
01032       SetSparseIndexAB(a,b);
01033    }
01034 
01035    Int_t * const pRowIndexc = this->GetRowIndexArray();
01036    Int_t * const pColIndexc = this->GetColIndexArray();
01037 
01038    const Int_t * const pRowIndexa = a.GetRowIndexArray();
01039    const Int_t * const pColIndexa = a.GetColIndexArray();
01040 
01041    const Element * const pDataa = a.GetMatrixArray();
01042    const Element * const pDatab = b.GetMatrixArray();
01043    Element * const pDatac = this->GetMatrixArray();
01044    Int_t indexc_r = 0;
01045    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
01046       const Int_t sIndexa = pRowIndexa[irowc];
01047       const Int_t eIndexa = pRowIndexa[irowc+1];
01048       const Int_t off = irowc*this->GetNcols();
01049       Int_t indexa = sIndexa;
01050       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
01051          Element sum = -pDatab[off+icolc];
01052          while (indexa < eIndexa && pColIndexa[indexa] <= icolc) {
01053             if (icolc == pColIndexa[indexa]) {
01054                sum += pDataa[indexa];
01055                break;
01056             }
01057             indexa++;
01058          }
01059 
01060          if (sum != 0.0) {
01061             pColIndexc[indexc_r] = icolc;
01062             pDatac[indexc_r] = sum;
01063             indexc_r++;
01064          }
01065       }
01066       pRowIndexc[irowc+1] = indexc_r;
01067    }
01068 
01069    if (constr)
01070       SetSparseIndex(indexc_r);
01071 }
01072 
01073 //______________________________________________________________________________
01074 template<class Element>
01075 void TMatrixTSparse<Element>::AMinusB(const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b,Int_t constr)
01076 {
01077   // General matrix subtraction. Create a matrix C such that C = A - B.
01078   // Note, matrix C is allocated for constr=1.
01079 
01080    if (gMatrixCheck) {
01081       R__ASSERT(a.IsValid());
01082       R__ASSERT(b.IsValid());
01083 
01084       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
01085           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
01086          Error("AMinusB(const TMatrixT &,const TMatrixTSparse &","matrices not compatible");
01087          return;
01088       }
01089 
01090       if (!constr && this->GetMatrixArray() == a.GetMatrixArray()) {
01091          Error("AMinusB","this = &a");
01092          return;
01093       }
01094 
01095       if (!constr && this->GetMatrixArray() == b.GetMatrixArray()) {
01096          Error("AMinusB","this = &b");
01097          return;
01098       }
01099    }
01100 
01101    if (constr) {
01102       Allocate(a.GetNrows(),a.GetNcols(),a.GetRowLwb(),a.GetColLwb());
01103       SetSparseIndexAB(a,b);
01104    }
01105 
01106    Int_t * const pRowIndexc = this->GetRowIndexArray();
01107    Int_t * const pColIndexc = this->GetColIndexArray();
01108 
01109    const Int_t * const pRowIndexb = b.GetRowIndexArray();
01110    const Int_t * const pColIndexb = b.GetColIndexArray();
01111 
01112    const Element * const pDataa = a.GetMatrixArray();
01113    const Element * const pDatab = b.GetMatrixArray();
01114    Element * const pDatac = this->GetMatrixArray();
01115    Int_t indexc_r = 0;
01116    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
01117       const Int_t sIndexb = pRowIndexb[irowc];
01118       const Int_t eIndexb = pRowIndexb[irowc+1];
01119       const Int_t off = irowc*this->GetNcols();
01120       Int_t indexb = sIndexb;
01121       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
01122          Element sum = pDataa[off+icolc];
01123          while (indexb < eIndexb && pColIndexb[indexb] <= icolc) {
01124             if (icolc == pColIndexb[indexb]) {
01125                sum -= pDatab[indexb];
01126                break;
01127             }
01128             indexb++;
01129          }
01130 
01131          if (sum != 0.0) {
01132             pColIndexc[indexc_r] = icolc;
01133             pDatac[indexc_r] = sum;
01134             indexc_r++;
01135          }
01136       }
01137       pRowIndexc[irowc+1] = indexc_r;
01138    }
01139 
01140    if (constr)
01141       SetSparseIndex(indexc_r);
01142 }
01143 
01144 //______________________________________________________________________________
01145 template<class Element>
01146 void TMatrixTSparse<Element>::GetMatrix2Array(Element *data,Option_t * /*option*/) const
01147 {
01148   // Copy matrix data to array . It is assumed that array is of size >= fNelems
01149 
01150    R__ASSERT(this->IsValid());
01151 
01152    const Element * const elem = GetMatrixArray();
01153    memcpy(data,elem,this->fNelems*sizeof(Element));
01154 }
01155 
01156 //______________________________________________________________________________
01157 template<class Element>
01158 TMatrixTBase<Element> &TMatrixTSparse<Element>::SetMatrixArray(Int_t nr,Int_t *row,Int_t *col,Element *data)
01159 {
01160   // Copy nr elements from row/col index and data array to matrix . It is assumed
01161   // that arrays are of size >= nr
01162 
01163    R__ASSERT(this->IsValid());
01164    if (nr <= 0) {
01165       Error("SetMatrixArray(Int_t,Int_t*,Int_t*,Element*","nr <= 0");
01166       return *this;
01167    }
01168 
01169    const Int_t irowmin = TMath::LocMin(nr,row);
01170    const Int_t irowmax = TMath::LocMax(nr,row);
01171    const Int_t icolmin = TMath::LocMin(nr,col);
01172    const Int_t icolmax = TMath::LocMax(nr,col);
01173 
01174    R__ASSERT(row[irowmin] >= this->fRowLwb && row[irowmax] <= this->fRowLwb+this->fNrows-1);
01175    R__ASSERT(col[icolmin] >= this->fColLwb && col[icolmax] <= this->fColLwb+this->fNcols-1);
01176 
01177    if (row[irowmin] < this->fRowLwb || row[irowmax] > this->fRowLwb+this->fNrows-1) {
01178       Error("SetMatrixArray","Inconsistency between row index and its range");
01179       if (row[irowmin] < this->fRowLwb) {
01180          Info("SetMatrixArray","row index lower bound adjusted to %d",row[irowmin]);
01181          this->fRowLwb = row[irowmin];
01182       }
01183       if (row[irowmax] > this->fRowLwb+this->fNrows-1) {
01184          Info("SetMatrixArray","row index upper bound adjusted to %d",row[irowmax]);
01185          this->fNrows = row[irowmax]-this->fRowLwb+1;
01186       }
01187    }
01188    if (col[icolmin] < this->fColLwb || col[icolmax] > this->fColLwb+this->fNcols-1) {
01189       Error("SetMatrixArray","Inconsistency between column index and its range");
01190       if (col[icolmin] < this->fColLwb) {
01191          Info("SetMatrixArray","column index lower bound adjusted to %d",col[icolmin]);
01192          this->fColLwb = col[icolmin];
01193       }
01194       if (col[icolmax] > this->fColLwb+this->fNcols-1) {
01195          Info("SetMatrixArray","column index upper bound adjusted to %d",col[icolmax]);
01196          this->fNcols = col[icolmax]-this->fColLwb+1;
01197       }
01198    }
01199 
01200    TMatrixTBase<Element>::DoubleLexSort(nr,row,col,data);
01201 
01202    Int_t nr_nonzeros = 0;
01203    const Element *ep        = data;
01204    const Element * const fp = data+nr;
01205 
01206    while (ep < fp)
01207       if (*ep++ != 0.0) nr_nonzeros++;
01208 
01209    if (nr_nonzeros != this->fNelems) {
01210       if (fColIndex) { delete [] fColIndex; fColIndex = 0; }
01211       if (fElements) { delete [] fElements; fElements = 0; }
01212       this->fNelems = nr_nonzeros;
01213       if (this->fNelems > 0) {
01214          fColIndex = new Int_t[nr_nonzeros];
01215          fElements = new Element[nr_nonzeros];
01216       } else {
01217          fColIndex = 0;
01218          fElements = 0;
01219       }
01220    }
01221 
01222    if (this->fNelems <= 0)
01223       return *this;
01224 
01225    fRowIndex[0] = 0;
01226    Int_t ielem = 0;
01227    nr_nonzeros = 0;
01228    for (Int_t irow = 1; irow < this->fNrows+1; irow++) {
01229       if (ielem < nr && row[ielem] < irow) {
01230          while (ielem < nr) {
01231             if (data[ielem] != 0.0) {
01232                fColIndex[nr_nonzeros] = col[ielem]-this->fColLwb;
01233                fElements[nr_nonzeros] = data[ielem];
01234                nr_nonzeros++;
01235             }
01236             ielem++;
01237             if (ielem >= nr || row[ielem] != row[ielem-1])
01238                break;
01239          }
01240       }
01241       fRowIndex[irow] = nr_nonzeros;
01242    }
01243 
01244    return *this;
01245 }
01246 
01247 //______________________________________________________________________________
01248 template<class Element>
01249 TMatrixTSparse<Element> &TMatrixTSparse<Element>::SetSparseIndex(Int_t nelems_new)
01250 {
01251   // Increase/decrease the number of non-zero elements to nelems_new
01252 
01253    if (nelems_new != this->fNelems) {
01254       Int_t nr = TMath::Min(nelems_new,this->fNelems);
01255       Int_t *oIp = fColIndex;
01256       fColIndex = new Int_t[nelems_new];
01257       if (oIp) { 
01258          memmove(fColIndex,oIp,nr*sizeof(Int_t));
01259          delete [] oIp;
01260       }
01261       Element *oDp = fElements;
01262       fElements = new Element[nelems_new];
01263       if (oDp) { 
01264          memmove(fElements,oDp,nr*sizeof(Element));
01265          delete [] oDp;
01266       }
01267       this->fNelems = nelems_new;
01268       if (nelems_new > nr) {
01269          memset(fElements+nr,0,(nelems_new-nr)*sizeof(Element));
01270          memset(fColIndex+nr,0,(nelems_new-nr)*sizeof(Int_t));
01271       } else {
01272          for (Int_t irow = 0; irow < this->fNrowIndex; irow++)
01273             if (fRowIndex[irow] > nelems_new)
01274                fRowIndex[irow] = nelems_new;
01275       }
01276    }
01277 
01278    return *this;
01279 }
01280 
01281 //______________________________________________________________________________
01282 template<class Element>
01283 TMatrixTSparse<Element> &TMatrixTSparse<Element>::SetSparseIndex(const TMatrixTBase<Element> &source)
01284 {
01285   // Use non-zero data of matrix source to set the sparse structure
01286 
01287    if (gMatrixCheck) {
01288       R__ASSERT(source.IsValid());
01289       if (this->GetNrows()  != source.GetNrows()  || this->GetNcols()  != source.GetNcols() ||
01290           this->GetRowLwb() != source.GetRowLwb() || this->GetColLwb() != source.GetColLwb()) {
01291          Error("SetSparseIndex","matrices not compatible");
01292          return *this;
01293       }
01294    }
01295 
01296    const Int_t nr_nonzeros = source.NonZeros();
01297 
01298    if (nr_nonzeros != this->fNelems)
01299       SetSparseIndex(nr_nonzeros);
01300 
01301    if (source.GetRowIndexArray() && source.GetColIndexArray()) {
01302       memmove(fRowIndex,source.GetRowIndexArray(),this->fNrowIndex*sizeof(Int_t));
01303       memmove(fColIndex,source.GetColIndexArray(),this->fNelems*sizeof(Int_t));
01304    } else {
01305       const Element *ep = source.GetMatrixArray();
01306       Int_t nr = 0;
01307       for (Int_t irow = 0; irow < this->fNrows; irow++) {
01308          fRowIndex[irow] = nr;
01309          for (Int_t icol = 0; icol < this->fNcols; icol++) {
01310             if (*ep != 0.0) {
01311                fColIndex[nr] = icol;
01312                nr++;
01313             }
01314             ep++;
01315          }
01316       }
01317       fRowIndex[this->fNrows] = nr;
01318    }
01319 
01320    return *this;
01321 }
01322 
01323 //______________________________________________________________________________
01324 template<class Element>
01325 TMatrixTSparse<Element> &TMatrixTSparse<Element>::SetSparseIndexAB(const TMatrixTSparse<Element> &a,const TMatrixTSparse<Element> &b)
01326 {
01327   // Set the row/column indices to the "sum" of matrices a and b
01328   // It is checked that enough space has been allocated
01329 
01330    if (gMatrixCheck) {
01331       R__ASSERT(a.IsValid());
01332       R__ASSERT(b.IsValid());
01333 
01334       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
01335           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
01336          Error("SetSparseIndexAB","source matrices not compatible");
01337          return *this;
01338       }
01339 
01340       if (this->GetNrows()  != a.GetNrows()  || this->GetNcols()  != a.GetNcols() ||
01341           this->GetRowLwb() != a.GetRowLwb() || this->GetColLwb() != a.GetColLwb()) {
01342          Error("SetSparseIndexAB","matrix not compatible with source matrices");
01343          return *this;
01344       }
01345    }
01346 
01347    const Int_t * const pRowIndexa = a.GetRowIndexArray();
01348    const Int_t * const pRowIndexb = b.GetRowIndexArray();
01349    const Int_t * const pColIndexa = a.GetColIndexArray();
01350    const Int_t * const pColIndexb = b.GetColIndexArray();
01351 
01352    // Count first the number of non-zero slots that are needed
01353    Int_t nc = 0;
01354    for (Int_t irowc = 0; irowc < a.GetNrows(); irowc++) {
01355       const Int_t sIndexa = pRowIndexa[irowc];
01356       const Int_t eIndexa = pRowIndexa[irowc+1];
01357       const Int_t sIndexb = pRowIndexb[irowc];
01358       const Int_t eIndexb = pRowIndexb[irowc+1];
01359       nc += eIndexa-sIndexa;
01360       Int_t indexb = sIndexb;
01361       for (Int_t indexa = sIndexa; indexa < eIndexa; indexa++) {
01362          const Int_t icola = pColIndexa[indexa];
01363          for (; indexb < eIndexb; indexb++) {
01364             if (pColIndexb[indexb] >= icola) {
01365                if (pColIndexb[indexb] == icola)
01366                   indexb++;
01367                break;
01368             }
01369             nc++;
01370          }
01371       }
01372       while (indexb < eIndexb) {
01373          const Int_t icola = (eIndexa > sIndexa && eIndexa > 0) ? pColIndexa[eIndexa-1] : -1;
01374          if (pColIndexb[indexb++] > icola)
01375             nc++;
01376       }
01377    }
01378 
01379    // Allocate the necessary space in fRowIndex and fColIndex
01380    if (this->NonZeros() != nc)
01381       SetSparseIndex(nc);
01382 
01383    Int_t * const pRowIndexc = this->GetRowIndexArray();
01384    Int_t * const pColIndexc = this->GetColIndexArray();
01385 
01386    nc = 0;
01387    pRowIndexc[0] = 0;
01388    for (Int_t irowc = 0; irowc < a.GetNrows(); irowc++) {
01389       const Int_t sIndexa = pRowIndexa[irowc];
01390       const Int_t eIndexa = pRowIndexa[irowc+1];
01391       const Int_t sIndexb = pRowIndexb[irowc];
01392       const Int_t eIndexb = pRowIndexb[irowc+1];
01393       Int_t indexb = sIndexb;
01394       for (Int_t indexa = sIndexa; indexa < eIndexa; indexa++) {
01395          const Int_t icola = pColIndexa[indexa];
01396          for (; indexb < eIndexb; indexb++) {
01397             if (pColIndexb[indexb] >= icola) {
01398                if (pColIndexb[indexb] == icola)
01399                   indexb++;
01400                break;
01401             }
01402             pColIndexc[nc++] = pColIndexb[indexb];
01403          }
01404          pColIndexc[nc++] = pColIndexa[indexa];
01405       }
01406       while (indexb < eIndexb) {
01407          const Int_t icola = (eIndexa > 0) ? pColIndexa[eIndexa-1] : -1;
01408          if (pColIndexb[indexb++] > icola)
01409             pColIndexc[nc++] = pColIndexb[indexb-1];
01410       }
01411       pRowIndexc[irowc+1] = nc;
01412    }
01413 
01414    return *this;
01415 }
01416 
01417 //______________________________________________________________________________
01418 template<class Element>
01419 TMatrixTSparse<Element> &TMatrixTSparse<Element>::SetSparseIndexAB(const TMatrixT<Element> &a,const TMatrixTSparse<Element> &b)
01420 {
01421   // Set the row/column indices to the "sum" of matrices a and b
01422   // It is checked that enough space has been allocated
01423 
01424    if (gMatrixCheck) {
01425       R__ASSERT(a.IsValid());
01426       R__ASSERT(b.IsValid());
01427 
01428       if (a.GetNrows()  != b.GetNrows()  || a.GetNcols()  != b.GetNcols() ||
01429           a.GetRowLwb() != b.GetRowLwb() || a.GetColLwb() != b.GetColLwb()) {
01430          Error("SetSparseIndexAB","source matrices not compatible");
01431          return *this;
01432       }
01433 
01434       if (this->GetNrows()  != a.GetNrows()  || this->GetNcols()  != a.GetNcols() ||
01435           this->GetRowLwb() != a.GetRowLwb() || this->GetColLwb() != a.GetColLwb()) {
01436          Error("SetSparseIndexAB","matrix not compatible with source matrices");
01437          return *this;
01438       }
01439    }
01440 
01441    const Element * const pDataa = a.GetMatrixArray();
01442 
01443    const Int_t * const pRowIndexb = b.GetRowIndexArray();
01444    const Int_t * const pColIndexb = b.GetColIndexArray();
01445 
01446    // Count first the number of non-zero slots that are needed
01447    Int_t nc = a.NonZeros();
01448    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
01449       const Int_t sIndexb = pRowIndexb[irowc];
01450       const Int_t eIndexb = pRowIndexb[irowc+1];
01451       const Int_t off = irowc*this->GetNcols();
01452       Int_t indexb = sIndexb;
01453       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
01454          if (pDataa[off+icolc] != 0.0 || pColIndexb[indexb] > icolc) continue;
01455          for (; indexb < eIndexb; indexb++) {
01456             if (pColIndexb[indexb] >= icolc) {
01457                if (pColIndexb[indexb] == icolc) {
01458                   nc++;
01459                   indexb++;
01460                }
01461                break;
01462             }
01463          }
01464       }
01465    }
01466 
01467    // Allocate thre necessary space in fRowIndex and fColIndex
01468    if (this->NonZeros() != nc)
01469       SetSparseIndex(nc);
01470 
01471    Int_t * const pRowIndexc = this->GetRowIndexArray();
01472    Int_t * const pColIndexc = this->GetColIndexArray();
01473 
01474    nc = 0;
01475    pRowIndexc[0] = 0;
01476    for (Int_t irowc = 0; irowc < this->GetNrows(); irowc++) {
01477       const Int_t sIndexb = pRowIndexb[irowc];
01478       const Int_t eIndexb = pRowIndexb[irowc+1];
01479       const Int_t off = irowc*this->GetNcols();
01480       Int_t indexb = sIndexb;
01481       for (Int_t icolc = 0; icolc < this->GetNcols(); icolc++) {
01482          if (pDataa[off+icolc] != 0.0)
01483             pColIndexc[nc++] = icolc;
01484          else if (pColIndexb[indexb] <= icolc) {
01485             for (; indexb < eIndexb; indexb++) {
01486                if (pColIndexb[indexb] >= icolc) {
01487                   if (pColIndexb[indexb] == icolc)
01488                      pColIndexc[nc++] = pColIndexb[indexb++];
01489                   break;
01490                }
01491             }
01492          }
01493       }
01494       pRowIndexc[irowc+1] = nc;
01495    }
01496 
01497    return *this;
01498 }
01499 
01500 //______________________________________________________________________________
01501 template<class Element>
01502 TMatrixTBase<Element> &TMatrixTSparse<Element>::ResizeTo(Int_t nrows,Int_t ncols,Int_t nr_nonzeros)
01503 {
01504   // Set size of the matrix to nrows x ncols with nr_nonzeros non-zero entries
01505   // if nr_nonzeros > 0 .
01506   // New dynamic elements are created, the overlapping part of the old ones are
01507   // copied to the new structures, then the old elements are deleted.
01508 
01509    R__ASSERT(this->IsValid());
01510    if (!this->fIsOwner) {
01511       Error("ResizeTo(Int_t,Int_t,Int_t)","Not owner of data array,cannot resize");
01512       return *this;
01513    }
01514 
01515    if (this->fNelems > 0) {
01516       if (this->fNrows == nrows && this->fNcols == ncols &&
01517          (this->fNelems == nr_nonzeros || nr_nonzeros < 0))
01518          return *this;
01519       else if (nrows == 0 || ncols == 0 || nr_nonzeros == 0) {
01520          this->fNrows = nrows; this->fNcols = ncols;
01521          Clear();
01522          return *this;
01523       }
01524 
01525       const Element *elements_old = GetMatrixArray();
01526       const Int_t   *rowIndex_old = GetRowIndexArray();
01527       const Int_t   *colIndex_old = GetColIndexArray();
01528 
01529       Int_t nrows_old     = this->fNrows;
01530       Int_t nrowIndex_old = this->fNrowIndex;
01531 
01532       Int_t nelems_new;
01533       if (nr_nonzeros > 0)
01534          nelems_new = nr_nonzeros;
01535       else {
01536          nelems_new = 0;
01537          for (Int_t irow = 0; irow < nrows_old; irow++) {
01538             if (irow >= nrows) continue;
01539             const Int_t sIndex = rowIndex_old[irow];
01540             const Int_t eIndex = rowIndex_old[irow+1];
01541             for (Int_t index = sIndex; index < eIndex; index++) {
01542                const Int_t icol = colIndex_old[index];
01543                if (icol < ncols)
01544                   nelems_new++;
01545             }
01546          }
01547       }
01548 
01549       Allocate(nrows,ncols,0,0,1,nelems_new);
01550       R__ASSERT(this->IsValid());
01551 
01552       Element *elements_new = GetMatrixArray();
01553       Int_t   *rowIndex_new = GetRowIndexArray();
01554       Int_t   *colIndex_new = GetColIndexArray();
01555 
01556       Int_t nelems_copy = 0;
01557       rowIndex_new[0] = 0;
01558       Bool_t cont = kTRUE;
01559       for (Int_t irow = 0; irow < nrows_old && cont; irow++) {
01560          if (irow >= nrows) continue;
01561          const Int_t sIndex = rowIndex_old[irow];
01562          const Int_t eIndex = rowIndex_old[irow+1];
01563          for (Int_t index = sIndex; index < eIndex; index++) {
01564             const Int_t icol = colIndex_old[index];
01565             if (icol < ncols) {
01566                rowIndex_new[irow+1]      = nelems_copy+1;
01567                colIndex_new[nelems_copy] = icol;
01568                elements_new[nelems_copy] = elements_old[index];
01569                nelems_copy++;
01570             }
01571             if (nelems_copy >= nelems_new) {
01572                cont = kFALSE;
01573                break;
01574             }
01575          }
01576       }
01577 
01578       if (rowIndex_old) delete [] (Int_t*)   rowIndex_old;
01579       if (colIndex_old) delete [] (Int_t*)   colIndex_old;
01580       if (elements_old) delete [] (Element*) elements_old;
01581 
01582       if (nrowIndex_old < this->fNrowIndex) {
01583           for (Int_t irow = nrowIndex_old; irow < this->fNrowIndex; irow++)
01584             rowIndex_new[irow] = rowIndex_new[nrowIndex_old-1];
01585       }
01586    } else {
01587       const Int_t nelems_new = (nr_nonzeros >= 0) ? nr_nonzeros : 0;
01588       Allocate(nrows,ncols,0,0,1,nelems_new);
01589    }
01590 
01591    return *this;
01592 }
01593 
01594 //______________________________________________________________________________
01595 template<class Element>
01596 TMatrixTBase<Element> &TMatrixTSparse<Element>::ResizeTo(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
01597                                                          Int_t nr_nonzeros)
01598 {
01599   // Set size of the matrix to [row_lwb:row_upb] x [col_lwb:col_upb] with nr_nonzeros
01600   // non-zero entries if nr_nonzeros > 0 .
01601   // New dynamic elemenst are created, the overlapping part of the old ones are
01602   // copied to the new structures, then the old elements are deleted.
01603 
01604    R__ASSERT(this->IsValid());
01605    if (!this->fIsOwner) {
01606       Error("ResizeTo(Int_t,Int_t,Int_t,Int_t,Int_t)","Not owner of data array,cannot resize");
01607       return *this;
01608    }
01609 
01610    const Int_t new_nrows = row_upb-row_lwb+1;
01611    const Int_t new_ncols = col_upb-col_lwb+1;
01612 
01613    if (this->fNelems > 0) {
01614       if (this->fNrows  == new_nrows && this->fNcols  == new_ncols &&
01615           this->fRowLwb == row_lwb   && this->fColLwb == col_lwb &&
01616           (this->fNelems == nr_nonzeros || nr_nonzeros < 0))
01617           return *this;
01618       else if (new_nrows == 0 || new_ncols == 0 || nr_nonzeros == 0) {
01619          this->fNrows = new_nrows; this->fNcols = new_ncols;
01620          this->fRowLwb = row_lwb; this->fColLwb = col_lwb;
01621          Clear();
01622          return *this;
01623       }
01624 
01625       const Int_t   *rowIndex_old = GetRowIndexArray();
01626       const Int_t   *colIndex_old = GetColIndexArray();
01627       const Element *elements_old = GetMatrixArray();
01628 
01629       const Int_t  nrowIndex_old = this->fNrowIndex;
01630 
01631       const Int_t  nrows_old     = this->fNrows;
01632       const Int_t  rowLwb_old    = this->fRowLwb;
01633       const Int_t  colLwb_old    = this->fColLwb;
01634 
01635       Int_t nelems_new;
01636       if (nr_nonzeros > 0)
01637          nelems_new = nr_nonzeros;
01638       else {
01639          nelems_new = 0;
01640          for (Int_t irow = 0; irow < nrows_old; irow++) {
01641             if (irow+rowLwb_old > row_upb || irow+rowLwb_old < row_lwb) continue;
01642             const Int_t sIndex = rowIndex_old[irow];
01643             const Int_t eIndex = rowIndex_old[irow+1];
01644             for (Int_t index = sIndex; index < eIndex; index++) {
01645                const Int_t icol = colIndex_old[index];
01646                if (icol+colLwb_old <= col_upb && icol+colLwb_old >= col_lwb)
01647                   nelems_new++;
01648             }
01649          }
01650       }
01651 
01652       Allocate(new_nrows,new_ncols,row_lwb,col_lwb,1,nelems_new);
01653       R__ASSERT(this->IsValid());
01654 
01655       Int_t    *rowIndex_new = GetRowIndexArray();
01656       Int_t    *colIndex_new = GetColIndexArray();
01657       Element *elements_new = GetMatrixArray();
01658 
01659       Int_t nelems_copy = 0;
01660       rowIndex_new[0] = 0;
01661       Bool_t cont = kTRUE;
01662       const Int_t row_off = rowLwb_old-row_lwb;
01663       const Int_t col_off = colLwb_old-col_lwb;
01664       for (Int_t irow = 0; irow < nrows_old; irow++) {
01665          if (irow+rowLwb_old > row_upb || irow+rowLwb_old < row_lwb) continue;
01666          const Int_t sIndex = rowIndex_old[irow];
01667          const Int_t eIndex = rowIndex_old[irow+1];
01668          for (Int_t index = sIndex; index < eIndex; index++) {
01669             const Int_t icol = colIndex_old[index];
01670             if (icol+colLwb_old <= col_upb && icol+colLwb_old >= col_lwb) {
01671                rowIndex_new[irow+row_off+1] = nelems_copy+1;
01672                colIndex_new[nelems_copy]    = icol+col_off;
01673                elements_new[nelems_copy]    = elements_old[index];
01674                nelems_copy++;
01675             }
01676             if (nelems_copy >= nelems_new) {
01677                cont = kFALSE;
01678                break;
01679             }
01680          }
01681       }
01682 
01683       if (rowIndex_old) delete [] (Int_t*)   rowIndex_old;
01684       if (colIndex_old) delete [] (Int_t*)   colIndex_old;
01685       if (elements_old) delete [] (Element*) elements_old;
01686 
01687       if (nrowIndex_old < this->fNrowIndex) {
01688          for (Int_t irow = nrowIndex_old; irow < this->fNrowIndex; irow++)
01689             rowIndex_new[irow] = rowIndex_new[nrowIndex_old-1];
01690       }
01691    } else {
01692       const Int_t nelems_new = (nr_nonzeros >= 0) ? nr_nonzeros : 0;
01693       Allocate(new_nrows,new_ncols,row_lwb,col_lwb,1,nelems_new);
01694    }
01695 
01696    return *this;
01697 }
01698 
01699 //______________________________________________________________________________
01700 template<class Element>
01701 TMatrixTSparse<Element> &TMatrixTSparse<Element>::Use(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
01702                                                       Int_t nr_nonzeros,Int_t *pRowIndex,Int_t *pColIndex,Element *pData)
01703 {
01704    if (gMatrixCheck) {
01705       if (row_upb < row_lwb) {
01706          Error("Use","row_upb=%d < row_lwb=%d",row_upb,row_lwb);
01707          return *this;
01708       }
01709       if (col_upb < col_lwb) {
01710          Error("Use","col_upb=%d < col_lwb=%d",col_upb,col_lwb);
01711          return *this;
01712       }
01713    }
01714 
01715    Clear();
01716 
01717    this->fNrows     = row_upb-row_lwb+1;
01718    this->fNcols     = col_upb-col_lwb+1;
01719    this->fRowLwb    = row_lwb;
01720    this->fColLwb    = col_lwb;
01721    this->fNrowIndex = this->fNrows+1;
01722    this->fNelems    = nr_nonzeros;
01723    this->fIsOwner   = kFALSE;
01724    this->fTol       = std::numeric_limits<Element>::epsilon();
01725 
01726    fElements  = pData;
01727    fRowIndex  = pRowIndex;
01728    fColIndex  = pColIndex;
01729 
01730    return *this;
01731 }
01732 
01733 //______________________________________________________________________________
01734 template<class Element>
01735 TMatrixTBase<Element> &TMatrixTSparse<Element>::GetSub(Int_t row_lwb,Int_t row_upb,Int_t col_lwb,Int_t col_upb,
01736                                                        TMatrixTBase<Element> &target,Option_t *option) const
01737 {
01738   // Get submatrix [row_lwb..row_upb][col_lwb..col_upb]; The indexing range of the
01739   // returned matrix depends on the argument option:
01740   //
01741   // option == "S" : return [0..row_upb-row_lwb+1][0..col_upb-col_lwb+1] (default)
01742   // else          : return [row_lwb..row_upb][col_lwb..col_upb]
01743 
01744    if (gMatrixCheck) {
01745       R__ASSERT(this->IsValid());
01746       if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
01747          Error("GetSub","row_lwb out-of-bounds");
01748          return target;
01749       }
01750       if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
01751          Error("GetSub","col_lwb out-of-bounds");
01752          return target;
01753       }
01754       if (row_upb < this->fRowLwb || row_upb > this->fRowLwb+this->fNrows-1) {
01755          Error("GetSub","row_upb out-of-bounds");
01756          return target;
01757       }
01758       if (col_upb < this->fColLwb || col_upb > this->fColLwb+this->fNcols-1) {
01759          Error("GetSub","col_upb out-of-bounds");
01760          return target;
01761       }
01762       if (row_upb < row_lwb || col_upb < col_lwb) {
01763          Error("GetSub","row_upb < row_lwb || col_upb < col_lwb");
01764          return target;
01765       }
01766    }
01767 
01768    TString opt(option);
01769    opt.ToUpper();
01770    const Int_t shift = (opt.Contains("S")) ? 1 : 0;
01771 
01772    const Int_t row_lwb_sub = (shift) ? 0               : row_lwb;
01773    const Int_t row_upb_sub = (shift) ? row_upb-row_lwb : row_upb;
01774    const Int_t col_lwb_sub = (shift) ? 0               : col_lwb;
01775    const Int_t col_upb_sub = (shift) ? col_upb-col_lwb : col_upb;
01776 
01777    Int_t nr_nonzeros = 0;
01778    for (Int_t irow = 0; irow < this->fNrows; irow++) {
01779       if (irow+this->fRowLwb > row_upb || irow+this->fRowLwb < row_lwb) continue;
01780       const Int_t sIndex = fRowIndex[irow];
01781       const Int_t eIndex = fRowIndex[irow+1];
01782       for (Int_t index = sIndex; index < eIndex; index++) {
01783          const Int_t icol = fColIndex[index];
01784          if (icol+this->fColLwb <= col_upb && icol+this->fColLwb >= col_lwb)
01785             nr_nonzeros++;
01786       }
01787    }
01788 
01789    target.ResizeTo(row_lwb_sub,row_upb_sub,col_lwb_sub,col_upb_sub,nr_nonzeros);
01790 
01791    const Element *ep = this->GetMatrixArray();
01792 
01793    Int_t   *rowIndex_sub = target.GetRowIndexArray();
01794    Int_t   *colIndex_sub = target.GetColIndexArray();
01795    Element *ep_sub       = target.GetMatrixArray();
01796 
01797    if (target.GetRowIndexArray() && target.GetColIndexArray()) {
01798       Int_t nelems_copy = 0;
01799       rowIndex_sub[0] = 0;
01800       const Int_t row_off = this->fRowLwb-row_lwb;
01801       const Int_t col_off = this->fColLwb-col_lwb;
01802       for (Int_t irow = 0; irow < this->fNrows; irow++) {
01803          if (irow+this->fRowLwb > row_upb || irow+this->fRowLwb < row_lwb) continue;
01804          const Int_t sIndex = fRowIndex[irow];
01805          const Int_t eIndex = fRowIndex[irow+1];
01806          for (Int_t index = sIndex; index < eIndex; index++) {
01807             const Int_t icol = fColIndex[index];
01808             if (icol+this->fColLwb <= col_upb && icol+this->fColLwb >= col_lwb) {
01809                rowIndex_sub[irow+row_off+1] = nelems_copy+1;
01810                colIndex_sub[nelems_copy]    = icol+col_off;
01811                ep_sub[nelems_copy]          = ep[index];
01812                nelems_copy++;
01813             }
01814          }
01815       }
01816    } else {
01817       const Int_t row_off = this->fRowLwb-row_lwb;
01818       const Int_t col_off = this->fColLwb-col_lwb;
01819       const Int_t ncols_sub = col_upb_sub-col_lwb_sub+1;
01820       for (Int_t irow = 0; irow < this->fNrows; irow++) {
01821          if (irow+this->fRowLwb > row_upb || irow+this->fRowLwb < row_lwb) continue;
01822          const Int_t sIndex = fRowIndex[irow];
01823          const Int_t eIndex = fRowIndex[irow+1];
01824          const Int_t off = (irow+row_off)*ncols_sub;
01825          for (Int_t index = sIndex; index < eIndex; index++) {
01826             const Int_t icol = fColIndex[index];
01827             if (icol+this->fColLwb <= col_upb && icol+this->fColLwb >= col_lwb)
01828                ep_sub[off+icol+col_off] = ep[index];
01829          }
01830       }
01831     }
01832 
01833    return target;
01834 }
01835 
01836 //______________________________________________________________________________
01837 template<class Element>
01838 TMatrixTBase<Element> &TMatrixTSparse<Element>::SetSub(Int_t row_lwb,Int_t col_lwb,const TMatrixTBase<Element> &source)
01839 {
01840   // Insert matrix source starting at [row_lwb][col_lwb], thereby overwriting the part
01841   // [row_lwb..row_lwb+nrows_source-1][col_lwb..col_lwb+ncols_source-1];
01842 
01843    if (gMatrixCheck) {
01844       R__ASSERT(this->IsValid());
01845       R__ASSERT(source.IsValid());
01846 
01847       if (row_lwb < this->fRowLwb || row_lwb > this->fRowLwb+this->fNrows-1) {
01848          Error("SetSub","row_lwb out-of-bounds");
01849          return *this;
01850       }
01851       if (col_lwb < this->fColLwb || col_lwb > this->fColLwb+this->fNcols-1) {
01852          Error("SetSub","col_lwb out-of-bounds");
01853          return *this;
01854       }
01855       if (row_lwb+source.GetNrows() > this->fRowLwb+this->fNrows || col_lwb+source.GetNcols() > this->fColLwb+this->fNcols) {
01856          Error("SetSub","source matrix too large");
01857          return *this;
01858       }
01859    }
01860 
01861    const Int_t nRows_source = source.GetNrows();
01862    const Int_t nCols_source = source.GetNcols();
01863 
01864    // Determine how many non-zero's are already available in
01865    // [row_lwb..row_lwb+nrows_source-1][col_lwb..col_lwb+ncols_source-1]
01866    Int_t nr_nonzeros = 0;
01867    for (Int_t irow = 0; irow < this->fNrows; irow++) {
01868       if (irow+this->fRowLwb >= row_lwb+nRows_source || irow+this->fRowLwb < row_lwb) continue;
01869       const Int_t sIndex = fRowIndex[irow];
01870       const Int_t eIndex = fRowIndex[irow+1];
01871       for (Int_t index = sIndex; index < eIndex; index++) {
01872          const Int_t icol = fColIndex[index];
01873          if (icol+this->fColLwb < col_lwb+nCols_source && icol+this->fColLwb >= col_lwb)
01874             nr_nonzeros++;
01875       }
01876    }
01877 
01878    const Int_t   *rowIndex_s = source.GetRowIndexArray();
01879    const Int_t   *colIndex_s = source.GetColIndexArray();
01880    const Element *elements_s = source.GetMatrixArray();
01881 
01882    const Int_t nelems_old = this->fNelems;
01883    const Int_t   *rowIndex_old = GetRowIndexArray();
01884    const Int_t   *colIndex_old = GetColIndexArray();
01885    const Element *elements_old = GetMatrixArray();
01886 
01887    const Int_t nelems_new = nelems_old+source.NonZeros()-nr_nonzeros;
01888    fRowIndex = new Int_t[this->fNrowIndex];
01889    fColIndex = new Int_t[nelems_new];
01890    fElements = new Element[nelems_new];
01891    this->fNelems   = nelems_new;
01892 
01893    Int_t   *rowIndex_new = GetRowIndexArray();
01894    Int_t   *colIndex_new = GetColIndexArray();
01895    Element *elements_new = GetMatrixArray();
01896 
01897    const Int_t row_off = row_lwb-this->fRowLwb;
01898    const Int_t col_off = col_lwb-this->fColLwb;
01899 
01900    Int_t nr = 0;
01901    rowIndex_new[0] = 0;
01902    for (Int_t irow = 0; irow < this->fNrows; irow++) {
01903       rowIndex_new[irow+1] = rowIndex_new[irow];
01904       Bool_t flagRow = kFALSE;
01905       if (irow+this->fRowLwb < row_lwb+nRows_source && irow+this->fRowLwb >= row_lwb)
01906          flagRow = kTRUE;
01907 
01908       const Int_t sIndex_o = rowIndex_old[irow];
01909       const Int_t eIndex_o = rowIndex_old[irow+1];
01910 
01911       if (flagRow) {
01912          const Int_t icol_left = col_off-1;
01913          const Int_t left = TMath::BinarySearch(eIndex_o-sIndex_o,colIndex_old+sIndex_o,icol_left)+sIndex_o;
01914          for (Int_t index = sIndex_o; index <= left; index++) {
01915             rowIndex_new[irow+1]++;
01916             colIndex_new[nr] = colIndex_old[index];
01917             elements_new[nr] = elements_old[index];
01918             nr++;
01919          }
01920 
01921          if (rowIndex_s && colIndex_s) {
01922             const Int_t sIndex_s = rowIndex_s[irow-row_off];
01923             const Int_t eIndex_s = rowIndex_s[irow-row_off+1];
01924             for (Int_t index = sIndex_s; index < eIndex_s; index++) {
01925                rowIndex_new[irow+1]++;
01926                colIndex_new[nr] = colIndex_s[index]+col_off;
01927                elements_new[nr] = elements_s[index];
01928                nr++;
01929             }
01930          } else {
01931             const Int_t off = (irow-row_off)*nCols_source;
01932             for (Int_t icol = 0; icol < nCols_source; icol++) {
01933                 rowIndex_new[irow+1]++;
01934                 colIndex_new[nr] = icol+col_off;
01935                 elements_new[nr] = elements_s[off+icol];
01936                 nr++;
01937             }
01938          }
01939 
01940          const Int_t icol_right = col_off+nCols_source-1;
01941          if (colIndex_old) {
01942             Int_t right = TMath::BinarySearch(eIndex_o-sIndex_o,colIndex_old+sIndex_o,icol_right)+sIndex_o;
01943             while (right < eIndex_o-1 && colIndex_old[right+1] <= icol_right)
01944                right++;
01945             right++;
01946 
01947             for (Int_t index = right; index < eIndex_o; index++) {
01948                rowIndex_new[irow+1]++;
01949                colIndex_new[nr] = colIndex_old[index];
01950                elements_new[nr] = elements_old[index];
01951                nr++;
01952             }
01953          }
01954       } else {
01955          for (Int_t index = sIndex_o; index < eIndex_o; index++) {
01956             rowIndex_new[irow+1]++;
01957             colIndex_new[nr] = colIndex_old[index];
01958             elements_new[nr] = elements_old[index];
01959             nr++;
01960          }
01961       }
01962    }
01963 
01964    R__ASSERT(this->fNelems == fRowIndex[this->fNrowIndex-1]);
01965 
01966    if (rowIndex_old) delete [] (Int_t*)    rowIndex_old;
01967    if (colIndex_old) delete [] (Int_t*)    colIndex_old;
01968    if (elements_old) delete [] (Element*) elements_old;
01969 
01970    return *this;
01971 }
01972 
01973 //______________________________________________________________________________
01974 template<class Element>
01975 TMatrixTSparse<Element> &TMatrixTSparse<Element>::Transpose(const TMatrixTSparse<Element> &source)
01976 {
01977   // Transpose a matrix.
01978 
01979    if (gMatrixCheck) {
01980       R__ASSERT(this->IsValid());
01981       R__ASSERT(source.IsValid());
01982 
01983       if (this->fNrows  != source.GetNcols()  || this->fNcols  != source.GetNrows() ||
01984           this->fRowLwb != source.GetColLwb() || this->fColLwb != source.GetRowLwb()) {
01985          Error("Transpose","matrix has wrong shape");
01986          return *this;
01987       }
01988 
01989       if (source.NonZeros() <= 0)
01990          return *this;
01991    }
01992 
01993    const Int_t nr_nonzeros = source.NonZeros();
01994 
01995    const Int_t   * const pRowIndex_s = source.GetRowIndexArray();
01996    const Int_t   * const pColIndex_s = source.GetColIndexArray();
01997    const Element * const pData_s     = source.GetMatrixArray();
01998 
01999    Int_t   *rownr   = new Int_t  [nr_nonzeros];
02000    Int_t   *colnr   = new Int_t  [nr_nonzeros];
02001    Element *pData_t = new Element[nr_nonzeros];
02002 
02003    Int_t ielem = 0;
02004    for (Int_t irow_s = 0; irow_s < source.GetNrows(); irow_s++) {
02005       const Int_t sIndex = pRowIndex_s[irow_s];
02006       const Int_t eIndex = pRowIndex_s[irow_s+1];
02007       for (Int_t index = sIndex; index < eIndex; index++) {
02008          rownr[ielem]   = pColIndex_s[index]+this->fRowLwb;
02009          colnr[ielem]   = irow_s+this->fColLwb;
02010          pData_t[ielem] = pData_s[index];
02011          ielem++;
02012       }
02013    }
02014 
02015    R__ASSERT(nr_nonzeros >= ielem);
02016 
02017    TMatrixTBase<Element>::DoubleLexSort(nr_nonzeros,rownr,colnr,pData_t);
02018    SetMatrixArray(nr_nonzeros,rownr,colnr,pData_t);
02019 
02020    R__ASSERT(this->fNelems == fRowIndex[this->fNrowIndex-1]);
02021 
02022    if (pData_t) delete [] pData_t;
02023    if (rownr)   delete [] rownr;
02024    if (colnr)   delete [] colnr;
02025 
02026    return *this;
02027 }
02028 
02029 //______________________________________________________________________________
02030 template<class Element>
02031 TMatrixTBase<Element> &TMatrixTSparse<Element>::Zero()
02032 {
02033    R__ASSERT(this->IsValid());
02034 
02035    if (fElements) { delete [] fElements; fElements = 0; }
02036    if (fColIndex) { delete [] fColIndex; fColIndex = 0; }
02037    this->fNelems = 0;
02038    memset(this->GetRowIndexArray(),0,this->fNrowIndex*sizeof(Int_t));
02039 
02040    return *this;
02041 }
02042 
02043 //______________________________________________________________________________
02044 template<class Element>
02045 TMatrixTBase<Element> &TMatrixTSparse<Element>::UnitMatrix()
02046 {
02047   // Make a unit matrix (matrix need not be a square one).
02048 
02049    R__ASSERT(this->IsValid());
02050 
02051    Int_t i;
02052 
02053    Int_t nr_nonzeros = 0;
02054    for (i = this->fRowLwb; i <= this->fRowLwb+this->fNrows-1; i++)
02055       for (Int_t j = this->fColLwb; j <= this->fColLwb+this->fNcols-1; j++)
02056          if (i == j) nr_nonzeros++;
02057 
02058    if (nr_nonzeros != this->fNelems) {
02059       this->fNelems = nr_nonzeros;
02060       Int_t *oIp = fColIndex;
02061       fColIndex = new Int_t[nr_nonzeros];
02062       if (oIp) delete [] oIp;
02063       Element *oDp = fElements;
02064       fElements = new Element[nr_nonzeros];
02065       if (oDp) delete [] oDp;
02066    }
02067 
02068    Int_t ielem = 0;
02069    fRowIndex[0] = 0;
02070    for (i = this->fRowLwb; i <= this->fRowLwb+this->fNrows-1; i++) {
02071       for (Int_t j = this->fColLwb; j <= this->fColLwb+this->fNcols-1; j++) {
02072          if (i == j) {
02073             const Int_t irow = i-this->fRowLwb;
02074             fRowIndex[irow+1]  = ielem+1;
02075             fElements[ielem]   = 1.0;
02076             fColIndex[ielem++] = j-this->fColLwb;
02077          }
02078       }
02079    }
02080 
02081    return *this;
02082 }
02083 
02084 //______________________________________________________________________________
02085 template<class Element>
02086 Element TMatrixTSparse<Element>::RowNorm() const
02087 {
02088   // Row matrix norm, MAX{ SUM{ |M(i,j)|, over j}, over i}.
02089   // The norm is induced by the infinity vector norm.
02090 
02091    R__ASSERT(this->IsValid());
02092 
02093    const Element *       ep = GetMatrixArray();
02094    const Element * const fp = ep+this->fNelems;
02095    const Int_t   * const pR = GetRowIndexArray();
02096          Element norm = 0;
02097 
02098    // Scan the matrix row-after-row
02099    for (Int_t irow = 0; irow < this->fNrows; irow++) {
02100       const Int_t sIndex = pR[irow];
02101       const Int_t eIndex = pR[irow+1];
02102       Element sum = 0;
02103       for (Int_t index = sIndex; index < eIndex; index++)
02104           sum += TMath::Abs(*ep++);
02105       norm = TMath::Max(norm,sum);
02106    }
02107 
02108    R__ASSERT(ep == fp);
02109 
02110    return norm;
02111 }
02112 
02113 //______________________________________________________________________________
02114 template<class Element>
02115 Element TMatrixTSparse<Element>::ColNorm() const
02116 {
02117   // Column matrix norm, MAX{ SUM{ |M(i,j)|, over i}, over j}.
02118   // The norm is induced by the 1 vector norm.
02119 
02120    R__ASSERT(this->IsValid());
02121 
02122    const TMatrixTSparse<Element> mt(kTransposed,*this);
02123 
02124    const Element *       ep = mt.GetMatrixArray();
02125    const Element * const fp = ep+this->fNcols;
02126          Element norm = 0;
02127 
02128    // Scan the matrix col-after-col
02129    while (ep < fp) {
02130       Element sum = 0;
02131       // Scan a col to compute the sum
02132       for (Int_t i = 0; i < this->fNrows; i++,ep += this->fNcols)
02133          sum += TMath::Abs(*ep);
02134       ep -= this->fNelems-1;         // Point ep to the beginning of the next col
02135       norm = TMath::Max(norm,sum);
02136    }
02137 
02138    R__ASSERT(ep == fp);
02139 
02140    return norm;
02141 }
02142 
02143 //______________________________________________________________________________
02144 template<class Element>
02145 Element &TMatrixTSparse<Element>::operator()(Int_t rown,Int_t coln)
02146 {
02147    R__ASSERT(this->IsValid());
02148 
02149    const Int_t arown = rown-this->fRowLwb;
02150    const Int_t acoln = coln-this->fColLwb;
02151    if (arown >= this->fNrows || arown < 0) {
02152       Error("operator()","Request row(%d) outside matrix range of %d - %d",rown,this->fRowLwb,this->fRowLwb+this->fNrows);
02153       return fElements[0];
02154    }
02155    if (acoln >= this->fNcols || acoln < 0) {
02156       Error("operator()","Request column(%d) outside matrix range of %d - %d",coln,this->fColLwb,this->fColLwb+this->fNcols);
02157       return fElements[0];
02158    }
02159 
02160    Int_t index = -1;
02161    Int_t sIndex = 0;
02162    Int_t eIndex = 0;
02163    if (this->fNrowIndex > 0 && fRowIndex[this->fNrowIndex-1] != 0) {
02164       sIndex = fRowIndex[arown];
02165       eIndex = fRowIndex[arown+1];
02166       index = TMath::BinarySearch(eIndex-sIndex,fColIndex+sIndex,acoln)+sIndex;
02167    }
02168 
02169    if (index >= sIndex && fColIndex[index] == acoln)
02170       return fElements[index];
02171    else {
02172       Element val = 0.;
02173       InsertRow(rown,coln,&val,1);
02174       sIndex = fRowIndex[arown];
02175       eIndex = fRowIndex[arown+1];
02176       index = TMath::BinarySearch(eIndex-sIndex,fColIndex+sIndex,acoln)+sIndex;
02177       if (index >= sIndex && fColIndex[index] == acoln)
02178          return fElements[index];
02179       else {
02180          Error("operator()(Int_t,Int_t","Insert row failed");
02181          return fElements[0];
02182       }
02183    }
02184 }
02185 
02186 //______________________________________________________________________________
02187 template <class Element>
02188 Element TMatrixTSparse<Element>::operator()(Int_t rown,Int_t coln) const
02189 {
02190    R__ASSERT(this->IsValid());
02191    if (this->fNrowIndex > 0 && this->fRowIndex[this->fNrowIndex-1] == 0) {
02192       Error("operator()(Int_t,Int_t) const","row/col indices are not set");
02193       Info("operator()","fNrowIndex = %d fRowIndex[fNrowIndex-1] = %d\n",this->fNrowIndex,this->fRowIndex[this->fNrowIndex-1]);
02194       return 0.0;
02195    }
02196    const Int_t arown = rown-this->fRowLwb;
02197    const Int_t acoln = coln-this->fColLwb;
02198 
02199    if (arown >= this->fNrows || arown < 0) {
02200       Error("operator()","Request row(%d) outside matrix range of %d - %d",rown,this->fRowLwb,this->fRowLwb+this->fNrows);
02201       return 0.0;
02202    }
02203    if (acoln >= this->fNcols || acoln < 0) {
02204       Error("operator()","Request column(%d) outside matrix range of %d - %d",coln,this->fColLwb,this->fColLwb+this->fNcols);
02205       return 0.0;
02206    }
02207 
02208    const Int_t sIndex = fRowIndex[arown];
02209    const Int_t eIndex = fRowIndex[arown+1];
02210    const Int_t index = (Int_t)TMath::BinarySearch(eIndex-sIndex,fColIndex+sIndex,acoln)+sIndex;
02211    if (index < sIndex || fColIndex[index] != acoln) return 0.0;
02212    else                                             return fElements[index];
02213 }
02214 
02215 //______________________________________________________________________________
02216 template<class Element>
02217 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator=(const TMatrixTSparse<Element> &source)
02218 {
02219   // Notice that the sparsity of the matrix is NOT changed : its fRowIndex/fColIndex
02220   // are used !
02221 
02222    if (gMatrixCheck && !AreCompatible(*this,source)) {
02223       Error("operator=(const TMatrixTSparse &)","matrices not compatible");
02224       return *this;
02225    }
02226 
02227    if (this->GetMatrixArray() != source.GetMatrixArray()) {
02228       TObject::operator=(source);
02229 
02230       const Element * const sp = source.GetMatrixArray();
02231             Element * const tp = this->GetMatrixArray();
02232       memcpy(tp,sp,this->fNelems*sizeof(Element));
02233       this->fTol = source.GetTol();
02234    }
02235    return *this;
02236 }
02237 
02238 //______________________________________________________________________________
02239 template<class Element>
02240 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator=(const TMatrixT<Element> &source)
02241 {
02242   // Notice that the sparsity of the matrix is NOT changed : its fRowIndex/fColIndex
02243   // are used !
02244 
02245    if (gMatrixCheck && !AreCompatible(*this,(TMatrixTBase<Element> &)source)) {
02246       Error("operator=(const TMatrixT &)","matrices not compatible");
02247       return *this;
02248    }
02249 
02250    if (this->GetMatrixArray() != source.GetMatrixArray()) {
02251       TObject::operator=(source);
02252 
02253       const Element * const sp = source.GetMatrixArray();
02254             Element * const tp = this->GetMatrixArray();
02255       for (Int_t irow = 0; irow < this->fNrows; irow++) {
02256          const Int_t sIndex = fRowIndex[irow];
02257          const Int_t eIndex = fRowIndex[irow+1];
02258          const Int_t off = irow*this->fNcols;
02259          for (Int_t index = sIndex; index < eIndex; index++) {
02260             const Int_t icol = fColIndex[index];
02261             tp[index] = sp[off+icol];
02262          }
02263       }
02264       this->fTol = source.GetTol();
02265    }
02266    return *this;
02267 }
02268 
02269 //______________________________________________________________________________
02270 template<class Element>
02271 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator=(Element val)
02272 {
02273   // Assign val to every element of the matrix. Check that the row/col
02274   // indices are set !
02275 
02276    R__ASSERT(this->IsValid());
02277 
02278    if (fRowIndex[this->fNrowIndex-1] == 0) {
02279       Error("operator=(Element","row/col indices are not set");
02280       return *this;
02281    }
02282 
02283    Element *ep = this->GetMatrixArray();
02284    const Element * const ep_last = ep+this->fNelems;
02285    while (ep < ep_last)
02286       *ep++ = val;
02287 
02288    return *this;
02289 }
02290 
02291 //______________________________________________________________________________
02292 template<class Element>
02293 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator+=(Element val)
02294 {
02295   // Add val to every element of the matrix.
02296 
02297    R__ASSERT(this->IsValid());
02298 
02299    Element *ep = this->GetMatrixArray();
02300    const Element * const ep_last = ep+this->fNelems;
02301    while (ep < ep_last)
02302       *ep++ += val;
02303 
02304    return *this;
02305 }
02306 
02307 //______________________________________________________________________________
02308 template<class Element>
02309 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator-=(Element val)
02310 {
02311   // Subtract val from every element of the matrix.
02312 
02313    R__ASSERT(this->IsValid());
02314 
02315    Element *ep = this->GetMatrixArray();
02316    const Element * const ep_last = ep+this->fNelems;
02317    while (ep < ep_last)
02318       *ep++ -= val;
02319 
02320    return *this;
02321 }
02322 
02323 //______________________________________________________________________________
02324 template<class Element>
02325 TMatrixTSparse<Element> &TMatrixTSparse<Element>::operator*=(Element val)
02326 {
02327   // Multiply every element of the matrix with val.
02328 
02329    R__ASSERT(this->IsValid());
02330 
02331    Element *ep = this->GetMatrixArray();
02332    const Element * const ep_last = ep+this->fNelems;
02333    while (ep < ep_last)
02334       *ep++ *= val;
02335 
02336    return *this;
02337 }
02338 
02339 //______________________________________________________________________________
02340 template<class Element>
02341 TMatrixTBase<Element> &TMatrixTSparse<Element>::Randomize(Element alpha,Element beta,Double_t &seed)
02342 {
02343   // randomize matrix element values
02344 
02345    R__ASSERT(this->IsValid());
02346 
02347    const Element scale = beta-alpha;
02348    const Element shift = alpha/scale;
02349 
02350    Int_t   * const pRowIndex = GetRowIndexArray();
02351    Int_t   * const pColIndex = GetColIndexArray();
02352    Element * const ep        = GetMatrixArray();
02353 
02354    const Int_t m = this->GetNrows();
02355    const Int_t n = this->GetNcols();
02356 
02357    // Knuth's algorithm for choosing "length" elements out of nn .
02358    const Int_t nn     = this->GetNrows()*this->GetNcols();
02359    const Int_t length = (this->GetNoElements() <= nn) ? this->GetNoElements() : nn;
02360    Int_t chosen   = 0;
02361    Int_t icurrent = 0;
02362    pRowIndex[0] = 0;
02363    for (Int_t k = 0; k < nn; k++) {
02364       const Element r = Drand(seed);
02365 
02366       if ((nn-k)*r < length-chosen) {
02367          pColIndex[chosen] = k%n;
02368          const Int_t irow  = k/n;
02369 
02370          if (irow > icurrent) {
02371             for ( ; icurrent < irow; icurrent++)
02372                pRowIndex[icurrent+1] = chosen;
02373          }
02374          ep[chosen] = scale*(Drand(seed)+shift);
02375          chosen++;
02376       }
02377    }
02378    for ( ; icurrent < m; icurrent++)
02379       pRowIndex[icurrent+1] = length;
02380 
02381    R__ASSERT(chosen == length);
02382 
02383    return *this;
02384 }
02385 
02386 //______________________________________________________________________________
02387 template<class Element>
02388 TMatrixTSparse<Element> &TMatrixTSparse<Element>::RandomizePD(Element alpha,Element beta,Double_t &seed)
02389 {
02390   // randomize matrix element values but keep matrix symmetric positive definite
02391 
02392    const Element scale = beta-alpha;
02393    const Element shift = alpha/scale;
02394 
02395    if (gMatrixCheck) {
02396       R__ASSERT(this->IsValid());
02397 
02398       if (this->fNrows != this->fNcols || this->fRowLwb != this->fColLwb) {
02399          Error("RandomizePD(Element &","matrix should be square");
02400          return *this;
02401       }
02402    }
02403 
02404    const Int_t n = this->fNcols;
02405 
02406    Int_t   * const pRowIndex = this->GetRowIndexArray();
02407    Int_t   * const pColIndex = this->GetColIndexArray();
02408    Element * const ep        = GetMatrixArray();
02409 
02410    // We will always have non-zeros on the diagonal, so there
02411    // is no randomness there. In fact, choose the (0,0) element now
02412    pRowIndex[0] = 0;
02413    pColIndex[0] = 0;
02414    pRowIndex[1] = 1;
02415    ep[0]        = 1e-8+scale*(Drand(seed)+shift);
02416 
02417    // Knuth's algorithm for choosing length elements out of nn .
02418    // nn here is the number of elements in the strict lower triangle.
02419    const Int_t nn = n*(n-1)/2;
02420 
02421    // length is the number of elements that can be stored, minus the number
02422    // of elements in the diagonal, which will always be in the matrix.
02423    //  Int_t length = (this->fNelems-n)/2;
02424    Int_t length = this->fNelems-n;
02425    length = (length <= nn) ? length : nn;
02426 
02427    // chosen   : the number of elements that have already been chosen (now 0)
02428    // nnz      : the number of non-zeros in the matrix (now 1, because the
02429    //            (0,0) element is already in the matrix.
02430    // icurrent : the index of the last row whose start has been stored in pRowIndex;
02431 
02432    Int_t chosen   = 0;
02433    Int_t icurrent = 1;
02434    Int_t nnz      = 1;
02435    for (Int_t k = 0; k < nn; k++ ) {
02436       const Element r = Drand(seed);
02437 
02438       if( (nn-k)*r < length-chosen) {
02439          // Element k is chosen. What row is it in?
02440          // In a lower triangular matrix (including a diagonal), it will be in
02441          // the largest row such that row*(row+1)/2 < k. In other words
02442 
02443          Int_t row = (int) TMath::Floor((-1+TMath::Sqrt(1.0+8.0*k))/2);
02444          // and its column will be the remainder
02445          Int_t col = k-row*(row+1)/2;
02446          // but since we are only filling in the *strict* lower triangle of
02447          // the matrix, we shift the row by 1
02448          row++;
02449 
02450          if (row > icurrent) {
02451             // We have chosen a row beyond the current row.
02452             // Choose a diagonal element for each intermediate row and fix the
02453             // data structure.
02454             for ( ; icurrent < row; icurrent++) {
02455                // Choose the diagonal
02456                ep[nnz] = 0.0;
02457                for (Int_t ll = pRowIndex[icurrent]; ll < nnz; ll++)
02458                   ep[nnz] += TMath::Abs(ep[ll]);
02459                ep[nnz] +=  1e-8+scale*(Drand(seed)+shift);
02460                pColIndex[nnz] = icurrent;
02461 
02462                nnz++;
02463                pRowIndex[icurrent+1] = nnz;
02464             }
02465          } // end if we have chosen a row beyond the current row;
02466          ep[nnz] = scale*(Drand(seed)+shift);
02467          pColIndex[nnz] = col;
02468          // add the value of this element (which occurs symmetrically in the
02469          // upper triangle) to the appropriate diagonal element
02470          ep[pRowIndex[col+1]-1] += TMath::Abs(ep[nnz]);
02471 
02472          nnz++; // We have added another element to the matrix
02473          chosen++; // And finished choosing another element.
02474       }
02475    }
02476 
02477    R__ASSERT(chosen == length);
02478 
02479    // and of course, we must choose all remaining diagonal elements .
02480    for ( ; icurrent < n; icurrent++) {
02481       // Choose the diagonal
02482       ep[nnz] = 0.0;
02483       for(Int_t ll = pRowIndex[icurrent]; ll < nnz; ll++)
02484          ep[nnz] += TMath::Abs(ep[ll]);
02485       ep[nnz] += 1e-8+scale*(Drand(seed)+shift);
02486       pColIndex[nnz] = icurrent;
02487 
02488       nnz++;
02489       pRowIndex[icurrent+1] = nnz;
02490    }
02491    this->fNelems = nnz;
02492 
02493    TMatrixTSparse<Element> tmp(TMatrixTSparse<Element>::kTransposed,*this);
02494    *this += tmp;
02495 
02496    // make sure to divide the diagonal by 2 becuase the operation
02497    // *this += tmp; adds the diagonal again
02498    {
02499       const Int_t    * const pR = this->GetRowIndexArray();
02500       const Int_t    * const pC = this->GetColIndexArray();
02501              Element * const pD = GetMatrixArray();
02502       for (Int_t irow = 0; irow < this->fNrows+1; irow++) {
02503          const Int_t sIndex = pR[irow];
02504          const Int_t eIndex = pR[irow+1];
02505          for (Int_t index = sIndex; index < eIndex; index++) {
02506             const Int_t icol = pC[index];
02507             if (irow == icol)
02508                pD[index] /= 2.;
02509          }
02510       }
02511    }
02512 
02513    return *this;
02514 }
02515 
02516 //______________________________________________________________________________
02517 template<class Element>
02518 TMatrixTSparse<Element> operator+(const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2)
02519 {
02520   TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kPlus,source2);
02521   return target;
02522 }
02523 
02524 //______________________________________________________________________________
02525 template<class Element>
02526 TMatrixTSparse<Element> operator+(const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2)
02527 {
02528    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kPlus,source2);
02529    return target;
02530 }
02531 
02532 //______________________________________________________________________________
02533 template<class Element>
02534 TMatrixTSparse<Element> operator+(const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2)
02535 {
02536    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kPlus,source2);
02537    return target;
02538 }
02539 
02540 //______________________________________________________________________________
02541 template<class Element>
02542 TMatrixTSparse<Element> operator+(const TMatrixTSparse<Element> &source,Element val)
02543 {
02544    TMatrixTSparse<Element> target(source);
02545    target += val;
02546    return target;
02547 }
02548 
02549 //______________________________________________________________________________
02550 template<class Element>
02551 TMatrixTSparse<Element> operator+(Element val,const TMatrixTSparse<Element> &source)
02552 {
02553    TMatrixTSparse<Element> target(source);
02554    target += val;
02555    return target;
02556 }
02557 
02558 //______________________________________________________________________________
02559 template<class Element>
02560 TMatrixTSparse<Element> operator-(const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2)
02561 {
02562    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMinus,source2);
02563    return target;
02564 }
02565 
02566 //______________________________________________________________________________
02567 template<class Element>
02568 TMatrixTSparse<Element> operator-(const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2)
02569 {
02570    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMinus,source2);
02571    return target;
02572 }
02573 
02574 //______________________________________________________________________________
02575 template<class Element>
02576 TMatrixTSparse<Element> operator-(const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2)
02577 {
02578    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMinus,source2);
02579    return target;
02580 }
02581 
02582 //______________________________________________________________________________
02583 template<class Element>
02584 TMatrixTSparse<Element> operator-(const TMatrixTSparse<Element> &source,Element val)
02585 {
02586    TMatrixTSparse<Element> target(source);
02587    target -= val;
02588    return target;
02589 }
02590 
02591 //______________________________________________________________________________
02592 template<class Element>
02593 TMatrixTSparse<Element> operator-(Element val,const TMatrixTSparse<Element> &source)
02594 {
02595    TMatrixTSparse<Element> target(source);
02596    target -= val;
02597    return target;
02598 }
02599 
02600 //______________________________________________________________________________
02601 template<class Element>
02602 TMatrixTSparse<Element> operator*(const TMatrixTSparse<Element> &source1,const TMatrixTSparse<Element> &source2)
02603 {
02604    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMult,source2);
02605    return target;
02606 }
02607 
02608 //______________________________________________________________________________
02609 template<class Element>
02610 TMatrixTSparse<Element> operator*(const TMatrixTSparse<Element> &source1,const TMatrixT<Element> &source2)
02611 {
02612    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMult,source2);
02613    return target;
02614 }
02615 
02616 //______________________________________________________________________________
02617 template<class Element>
02618 TMatrixTSparse<Element> operator*(const TMatrixT<Element> &source1,const TMatrixTSparse<Element> &source2)
02619 {
02620    TMatrixTSparse<Element> target(source1,TMatrixTSparse<Element>::kMult,source2);
02621    return target;
02622 }
02623 
02624 //______________________________________________________________________________
02625 template<class Element>
02626 TMatrixTSparse<Element> operator*(Element val,const TMatrixTSparse<Element> &source)
02627 {
02628    TMatrixTSparse<Element> target(source);
02629    target *= val;
02630    return target;
02631 }
02632 
02633 //______________________________________________________________________________
02634 template<class Element>
02635 TMatrixTSparse<Element> operator*(const TMatrixTSparse<Element> &source,Element val)
02636 {
02637    TMatrixTSparse<Element> target(source);
02638    target *= val;
02639    return target;
02640 }
02641 
02642 //______________________________________________________________________________
02643 template<class Element>
02644 TMatrixTSparse<Element> &Add(TMatrixTSparse<Element> &target,Element scalar,const TMatrixTSparse<Element> &source)
02645 {
02646   // Modify addition: target += scalar * source.
02647 
02648    target += scalar * source;
02649    return target;
02650 }
02651 
02652 //______________________________________________________________________________
02653 template<class Element>
02654 TMatrixTSparse<Element> &ElementMult(TMatrixTSparse<Element> &target,const TMatrixTSparse<Element> &source)
02655 {
02656   // Multiply target by the source, element-by-element.
02657 
02658    if (gMatrixCheck && !AreCompatible(target,source)) {
02659       ::Error("ElementMult(TMatrixTSparse &,const TMatrixTSparse &)","matrices not compatible");
02660       return target;
02661    }
02662 
02663    const Element *sp  = source.GetMatrixArray();
02664          Element *tp  = target.GetMatrixArray();
02665    const Element *ftp = tp+target.GetNoElements();
02666    while ( tp < ftp )
02667       *tp++ *= *sp++;
02668 
02669    return target;
02670 }
02671 
02672 //______________________________________________________________________________
02673 template<class Element>
02674 TMatrixTSparse<Element> &ElementDiv(TMatrixTSparse<Element> &target,const TMatrixTSparse<Element> &source)
02675 {
02676    // Divide target by the source, element-by-element.
02677 
02678    if (gMatrixCheck && !AreCompatible(target,source)) {
02679       ::Error("ElementDiv(TMatrixT &,const TMatrixT &)","matrices not compatible");
02680       return target;
02681    }
02682 
02683    const Element *sp  = source.GetMatrixArray();
02684          Element *tp  = target.GetMatrixArray();
02685    const Element *ftp = tp+target.GetNoElements();
02686    while ( tp < ftp ) {
02687       if (*sp != 0.0)
02688          *tp++ /= *sp++;
02689       else {
02690          Error("ElementDiv","source element is zero");
02691          tp++;
02692       }
02693    }
02694 
02695    return target;
02696 }
02697 
02698 //______________________________________________________________________________
02699 template<class Element>
02700 Bool_t AreCompatible(const TMatrixTSparse<Element> &m1,const TMatrixTSparse<Element> &m2,Int_t verbose)
02701 {
02702    if (!m1.IsValid()) {
02703       if (verbose)
02704          ::Error("AreCompatible", "matrix 1 not valid");
02705       return kFALSE;
02706    }
02707    if (!m2.IsValid()) {
02708       if (verbose)
02709          ::Error("AreCompatible", "matrix 2 not valid");
02710       return kFALSE;
02711    }
02712 
02713    if (m1.GetNrows()  != m2.GetNrows()  || m1.GetNcols()  != m2.GetNcols() ||
02714        m1.GetRowLwb() != m2.GetRowLwb() || m1.GetColLwb() != m2.GetColLwb()) {
02715       if (verbose)
02716          ::Error("AreCompatible", "matrices 1 and 2 not compatible");
02717       return kFALSE;
02718    }
02719 
02720    const Int_t *pR1 = m1.GetRowIndexArray();
02721    const Int_t *pR2 = m2.GetRowIndexArray();
02722    const Int_t nRows = m1.GetNrows();
02723    if (memcmp(pR1,pR2,(nRows+1)*sizeof(Int_t))) {
02724       if (verbose)
02725          ::Error("AreCompatible", "matrices 1 and 2 have different rowIndex");
02726       for (Int_t i = 0; i < nRows+1; i++)
02727          printf("%d: %d %d\n",i,pR1[i],pR2[i]);
02728       return kFALSE;
02729    }
02730    const Int_t *pD1 = m1.GetColIndexArray();
02731    const Int_t *pD2 = m2.GetColIndexArray();
02732    const Int_t nData = m1.GetNoElements();
02733    if (memcmp(pD1,pD2,nData*sizeof(Int_t))) {
02734       if (verbose)
02735          ::Error("AreCompatible", "matrices 1 and 2 have different colIndex");
02736       for (Int_t i = 0; i < nData; i++)
02737          printf("%d: %d %d\n",i,pD1[i],pD2[i]);
02738       return kFALSE;
02739    }
02740 
02741    return kTRUE;
02742 }
02743 
02744 //______________________________________________________________________________
02745 template<class Element>
02746 void TMatrixTSparse<Element>::Streamer(TBuffer &R__b)
02747 {
02748 // Stream an object of class TMatrixTSparse.
02749 
02750    if (R__b.IsReading()) {
02751       UInt_t R__s, R__c;
02752       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
02753       Clear();
02754       R__b.ReadClassBuffer(TMatrixTSparse<Element>::Class(),this,R__v,R__s,R__c);
02755       if (this->fNelems < 0)
02756          this->Invalidate();
02757       } else {
02758          R__b.WriteClassBuffer(TMatrixTSparse<Element>::Class(),this);
02759    }
02760 }
02761 
02762 template class TMatrixTSparse<Float_t>;
02763 
02764 #ifndef ROOT_TMatrixFSparsefwd
02765 #include "TMatrixFSparsefwd.h"
02766 #endif
02767 #ifndef ROOT_TMatrixFfwd
02768 #include "TMatrixFfwd.h"
02769 #endif
02770 
02771 template TMatrixFSparse  operator+    <Float_t>(const TMatrixFSparse &source1,const TMatrixFSparse &source2);
02772 template TMatrixFSparse  operator+    <Float_t>(const TMatrixFSparse &source1,const TMatrixF       &source2);
02773 template TMatrixFSparse  operator+    <Float_t>(const TMatrixF       &source1,const TMatrixFSparse &source2);
02774 template TMatrixFSparse  operator+    <Float_t>(const TMatrixFSparse &source ,      Float_t         val    );
02775 template TMatrixFSparse  operator+    <Float_t>(      Float_t         val    ,const TMatrixFSparse &source );
02776 template TMatrixFSparse  operator-    <Float_t>(const TMatrixFSparse &source1,const TMatrixFSparse &source2);
02777 template TMatrixFSparse  operator-    <Float_t>(const TMatrixFSparse &source1,const TMatrixF       &source2);
02778 template TMatrixFSparse  operator-    <Float_t>(const TMatrixF       &source1,const TMatrixFSparse &source2);
02779 template TMatrixFSparse  operator-    <Float_t>(const TMatrixFSparse &source ,      Float_t         val    );
02780 template TMatrixFSparse  operator-    <Float_t>(      Float_t         val    ,const TMatrixFSparse &source );
02781 template TMatrixFSparse  operator*    <Float_t>(const TMatrixFSparse &source1,const TMatrixFSparse &source2);
02782 template TMatrixFSparse  operator*    <Float_t>(const TMatrixFSparse &source1,const TMatrixF       &source2);
02783 template TMatrixFSparse  operator*    <Float_t>(const TMatrixF       &source1,const TMatrixFSparse &source2);
02784 template TMatrixFSparse  operator*    <Float_t>(      Float_t         val    ,const TMatrixFSparse &source );
02785 template TMatrixFSparse  operator*    <Float_t>(const TMatrixFSparse &source,       Float_t         val    );
02786 
02787 template TMatrixFSparse &Add          <Float_t>(TMatrixFSparse &target,      Float_t         scalar,const TMatrixFSparse &source);
02788 template TMatrixFSparse &ElementMult  <Float_t>(TMatrixFSparse &target,const TMatrixFSparse &source);
02789 template TMatrixFSparse &ElementDiv   <Float_t>(TMatrixFSparse &target,const TMatrixFSparse &source);
02790 
02791 template Bool_t          AreCompatible<Float_t>(const TMatrixFSparse &m1,const TMatrixFSparse &m2,Int_t verbose);
02792 
02793 #ifndef ROOT_TMatrixDSparsefwd
02794 #include "TMatrixDSparsefwd.h"
02795 #endif
02796 #ifndef ROOT_TMatrixDfwd
02797 #include "TMatrixDfwd.h"
02798 #endif
02799 
02800 template class TMatrixTSparse<Double_t>;
02801 
02802 template TMatrixDSparse  operator+    <Double_t>(const TMatrixDSparse &source1,const TMatrixDSparse &source2);
02803 template TMatrixDSparse  operator+    <Double_t>(const TMatrixDSparse &source1,const TMatrixD       &source2);
02804 template TMatrixDSparse  operator+    <Double_t>(const TMatrixD       &source1,const TMatrixDSparse &source2);
02805 template TMatrixDSparse  operator+    <Double_t>(const TMatrixDSparse &source ,      Double_t        val    );
02806 template TMatrixDSparse  operator+    <Double_t>(      Double_t        val    ,const TMatrixDSparse &source );
02807 template TMatrixDSparse  operator-    <Double_t>(const TMatrixDSparse &source1,const TMatrixDSparse &source2);
02808 template TMatrixDSparse  operator-    <Double_t>(const TMatrixDSparse &source1,const TMatrixD       &source2);
02809 template TMatrixDSparse  operator-    <Double_t>(const TMatrixD       &source1,const TMatrixDSparse &source2);
02810 template TMatrixDSparse  operator-    <Double_t>(const TMatrixDSparse &source ,      Double_t        val    );
02811 template TMatrixDSparse  operator-    <Double_t>(      Double_t        val    ,const TMatrixDSparse &source );
02812 template TMatrixDSparse  operator*    <Double_t>(const TMatrixDSparse &source1,const TMatrixDSparse &source2);
02813 template TMatrixDSparse  operator*    <Double_t>(const TMatrixDSparse &source1,const TMatrixD       &source2);
02814 template TMatrixDSparse  operator*    <Double_t>(const TMatrixD       &source1,const TMatrixDSparse &source2);
02815 template TMatrixDSparse  operator*    <Double_t>(      Double_t        val    ,const TMatrixDSparse &source );
02816 template TMatrixDSparse  operator*    <Double_t>(const TMatrixDSparse &source,       Double_t        val    );
02817 
02818 template TMatrixDSparse &Add          <Double_t>(TMatrixDSparse &target,      Double_t        scalar,const TMatrixDSparse &source);
02819 template TMatrixDSparse &ElementMult  <Double_t>(TMatrixDSparse &target,const TMatrixDSparse &source);
02820 template TMatrixDSparse &ElementDiv   <Double_t>(TMatrixDSparse &target,const TMatrixDSparse &source);
02821 
02822 template Bool_t          AreCompatible<Double_t>(const TMatrixDSparse &m1,const TMatrixDSparse &m2,Int_t verbose);

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