TDecompSparse.cxx

Go to the documentation of this file.
00001 // @(#)root/matrix:$Id: TDecompSparse.cxx 34976 2010-08-25 04:11:07Z pcanal $
00002 // Authors: Fons Rademakers, Eddy Offermann  Apr 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 #include "TDecompSparse.h"
00013 #include "TMath.h"
00014 
00015 ClassImp(TDecompSparse)
00016 
00017 ///////////////////////////////////////////////////////////////////////////
00018 //                                                                       //
00019 // Sparse Symmetric Decomposition class                                  //
00020 //                                                                       //
00021 // Solve a sparse symmetric system of linear equations using a method    //
00022 // based on Gaussian elimination as discussed in Duff and Reid,          //
00023 // ACM Trans. Math. Software 9 (1983), 302-325.                          //
00024 //                                                                       //
00025 ///////////////////////////////////////////////////////////////////////////
00026 
00027 //______________________________________________________________________________
00028 TDecompSparse::TDecompSparse()
00029 {
00030 // Default constructor
00031 
00032    fVerbose = 0;
00033    InitParam();
00034    memset(fInfo,0,21*sizeof(Int_t));
00035 }
00036 
00037 //______________________________________________________________________________
00038 TDecompSparse::TDecompSparse(Int_t nRows,Int_t nr_nonZeros,Int_t verbose)
00039 {
00040 // Constructor for a matrix with nrows and unspecified number of columns .
00041 // nr_nonZeros is the total number of non-zero entries in the matrix .
00042 
00043    fVerbose = verbose;
00044    InitParam();
00045 
00046    fNrows     = nRows;
00047    fNnonZeros = nr_nonZeros;
00048 
00049    fRowFact.Set(nr_nonZeros+1);
00050    fColFact.Set(nr_nonZeros+1);
00051    fW      .Set(fNrows+1);
00052    fIkeep  .Set(3*(fNrows+1));
00053    fIw     .Set((Int_t)(1.3 * (2*fNnonZeros+3*fNrows+1)+1));
00054    fIw1    .Set(2*(fNrows+1));
00055 
00056    memset(fInfo,0,21*sizeof(Int_t));
00057 
00058    // These parameters can only be set after sparsity/pivoting pattern is known
00059    fNsteps = 0;
00060    fMaxfrt = 0;
00061 }
00062 
00063 //______________________________________________________________________________
00064 TDecompSparse::TDecompSparse(Int_t row_lwb,Int_t row_upb,Int_t nr_nonZeros,Int_t verbose)
00065 {
00066 // Constructor for a matrix with row range, [row_lwb..row_upb] and unspecified column
00067 // range . nr_nonZeros is the total number of non-zero entries in the matrix .
00068 
00069    fVerbose = verbose;
00070    InitParam();
00071 
00072    fRowLwb    = row_lwb;
00073    fColLwb    = row_lwb;
00074    fNrows     = row_upb-row_lwb+1;
00075    fNnonZeros = nr_nonZeros;
00076 
00077    fRowFact.Set(nr_nonZeros+1);
00078    fColFact.Set(nr_nonZeros+1);
00079    fW      .Set(fNrows+1);
00080    fIkeep  .Set(3*(fNrows+1));
00081    fIw     .Set((Int_t)(1.3 * (2*fNnonZeros+3*fNrows+1)+1));
00082    fIw1    .Set(2*(fNrows+1));
00083 
00084    memset(fInfo,0,21*sizeof(Int_t));
00085 
00086    // These parameters can only be set after sparsity/pivoting pattern is known
00087    fNsteps = 0;
00088    fMaxfrt = 0;
00089 }
00090 
00091 //______________________________________________________________________________
00092 TDecompSparse::TDecompSparse(const TMatrixDSparse &a,Int_t verbose)
00093 {
00094 // Constructor for matrix A .
00095 
00096    fVerbose = verbose;
00097 
00098    InitParam();
00099    SetMatrix(a);
00100 
00101    memset(fInfo,0,21*sizeof(Int_t));
00102 }
00103 
00104 //______________________________________________________________________________
00105 TDecompSparse::TDecompSparse(const TDecompSparse &another) : TDecompBase(another)
00106 {
00107 // Copy constructor
00108 
00109    *this = another;
00110 }
00111 
00112 //______________________________________________________________________________
00113 Int_t TDecompSparse::NonZerosUpperTriang(const TMatrixDSparse &a)
00114 {
00115 // Static function, returning the number of non-zero entries in the upper triangular matrix .
00116 
00117    const Int_t  rowLwb   = a.GetRowLwb();
00118    const Int_t  colLwb   = a.GetColLwb();
00119    const Int_t  nrows    = a.GetNrows();;
00120    const Int_t *pRowIndex = a.GetRowIndexArray();
00121    const Int_t *pColIndex = a.GetColIndexArray();
00122 
00123    Int_t nr_nonzeros = 0;
00124    for (Int_t irow = 0; irow < nrows; irow++ ) {
00125       const Int_t rown = irow+rowLwb;
00126       for (Int_t index = pRowIndex[irow]; index < pRowIndex[irow+1]; index++ ) {
00127          const Int_t coln = pColIndex[index]+colLwb;
00128          if (coln >= rown) nr_nonzeros++;
00129       }
00130    }
00131 
00132    return nr_nonzeros;
00133 }
00134 
00135 //______________________________________________________________________________
00136 void TDecompSparse::CopyUpperTriang(const TMatrixDSparse &a,Double_t *b)
00137 {
00138 // Static function, copying the non-zero entries in the upper triangle to
00139 // array b . User should allocate enough memory for array b .
00140 
00141    const Int_t     rowLwb    = a.GetRowLwb();
00142    const Int_t     colLwb    = a.GetColLwb();
00143    const Int_t     nrows     = a.GetNrows();;
00144    const Int_t    *pRowIndex = a.GetRowIndexArray();
00145    const Int_t    *pColIndex = a.GetColIndexArray();
00146    const Double_t *pData     = a.GetMatrixArray();
00147 
00148    Int_t nr = 0;
00149    for (Int_t irow = 0; irow < nrows; irow++ ) {
00150       const Int_t rown = irow+rowLwb;
00151       for (Int_t index = pRowIndex[irow]; index < pRowIndex[irow+1]; index++ ) {
00152          const Int_t coln = pColIndex[index]+colLwb;
00153          if (coln >= rown) b[nr++] = pData[index];
00154       }
00155    }
00156 }
00157 
00158 //______________________________________________________________________________
00159 void TDecompSparse::SetMatrix(const TMatrixDSparse &a)
00160 {
00161 // Set matrix to be decomposed .
00162 
00163    ResetStatus();
00164 
00165    fA.Use(*const_cast<TMatrixDSparse *>(&a));
00166    fRowLwb    = fA.GetRowLwb();
00167    fColLwb    = fA.GetColLwb();
00168    fNrows     = fA.GetNrows();
00169    fNnonZeros = NonZerosUpperTriang(a);
00170 
00171    fRowFact.Set(fNnonZeros+1);
00172    fColFact.Set(fNnonZeros+1);
00173 
00174    const Int_t *rowIndex = a.GetRowIndexArray();
00175    const Int_t *colIndex = a.GetColIndexArray();
00176 
00177    Int_t nr = 0;
00178    for (Int_t irow = 0; irow < fNrows; irow++ ) {
00179       const Int_t rown = irow+fRowLwb;
00180       for (Int_t index = rowIndex[irow]; index < rowIndex[irow+1]; index++ ) {
00181          const Int_t coln = colIndex[index]+fColLwb;
00182          if (coln >= rown) {
00183             fRowFact[nr+1] = irow+1;
00184             fColFact[nr+1] = colIndex[index]+1;
00185             nr++;
00186          }
00187       }
00188    }
00189 
00190    fW    .Set(fNrows+1);
00191    fIkeep.Set(3*(fNrows+1));
00192    fIw   .Set((Int_t)(1.3 * (2*fNnonZeros+3*fNrows+1)+1));
00193    fIw1  .Set(2*(fNrows+1));
00194 
00195    // Determine pivot sequence, set iflag = 0 in order to make InitPivot choose the order.
00196    Int_t iflag = 0;
00197    Double_t ops;
00198    InitPivot(fNrows,fNnonZeros,fRowFact,fColFact,fIw,fIkeep,fIw1,fNsteps,iflag,
00199              fIcntl,fCntl,fInfo,ops);
00200 
00201    switch ( this->ErrorFlag() ) {
00202       case -1 :
00203          Error("SetMatrix(const TMatrixDSparse &","nRows  = %d out of range",fNrows);
00204          return;
00205       case -2 :
00206          Error("SetMatrix(const TMatrixDSparse &","nr_nonzeros  = %d out of range",fNnonZeros);
00207          return;
00208       case -3 :
00209          Error("SetMatrix(const TMatrixDSparse &",
00210                "insufficient space in fIw of %d suggest reset to %d",fIw.GetSize(),this->IError());
00211          return;
00212       case 1 :
00213          Error("SetMatrix(const TMatrixDSparse &",
00214                "detected %d entries out of rage in row/col indices; ignored",this->IError());
00215          return;
00216    }
00217 
00218    // set fIw and fIw1 in prep for calls to Factor and Solve
00219 
00220 //   fIw  .Set((Int_t) 1.2*this->MinRealWorkspace()+1);
00221    fIw  .Set((Int_t) 3*this->MinRealWorkspace()+1);
00222    fIw1 .Set(fNrows+1);
00223    fIw2 .Set(fNsteps+1);
00224 //   fFact.Set((Int_t) 1.2*this->MinRealWorkspace()+1);
00225    fFact.Set((Int_t) 3*this->MinRealWorkspace()+1);
00226 
00227    SetBit(kMatrixSet);
00228 }
00229 
00230 //______________________________________________________________________________
00231 Bool_t TDecompSparse::Decompose()
00232 {
00233 // Decomposition engine .
00234 // If the decomposition succeeds, bit kDecomposed is set .
00235 
00236    if (TestBit(kDecomposed)) return kTRUE;
00237 
00238    if ( !TestBit(kMatrixSet) ) {
00239       Error("Decompose()","Matrix has not been set");
00240       return kFALSE;
00241    }
00242 
00243    Int_t done = 0; Int_t tries = 0;
00244    do {
00245       fFact[0] = 0.;
00246       CopyUpperTriang(fA,fFact.GetArray()+1);
00247 
00248       Factor(fNrows,fNnonZeros,fRowFact,fColFact,fFact,fIw,fIkeep,
00249              fNsteps,fMaxfrt,fIw1,fIcntl,fCntl,fInfo);
00250 
00251       switch ( this->ErrorFlag() ) {
00252          case 0 :
00253             done = 1;
00254             break;
00255          case -1 :
00256             Error("Decompose()","nRows  = %d out of range",fNrows);
00257             return kFALSE;
00258          case -2 :
00259             Error("Decompose()","nr_nonzeros  = %d out of range",fNnonZeros);
00260             return kFALSE;
00261          case -3 :
00262             {
00263                if (fVerbose)
00264                   Info("Decompose()","insufficient space of fIw: %d",fIw.GetSize());
00265                const Int_t nIw_old = fIw.GetSize();
00266                const Int_t nIw = (this->IError() > fIPessimism*nIw_old) ? this->IError() :
00267                                                                          (Int_t)(fIPessimism*nIw_old);
00268                fIw.Set(nIw);
00269                if (fVerbose)
00270                   Info("Decompose()","resetting to fIw: %d",nIw);
00271                fIPessimism *= 1.1;
00272                break;
00273             }
00274          case -4 :
00275             {
00276                if (fVerbose)
00277                   Info("Decompose()","insufficient factorization space: %d",fFact.GetSize());
00278                const Int_t nFact_old = fFact.GetSize();
00279                const Int_t nFact = (this->IError() > fRPessimism*nFact_old) ? this->IError() :
00280                                                                              (Int_t) (fRPessimism*nFact_old);
00281                fFact.Set(nFact); fFact.Reset(0.0);
00282                CopyUpperTriang(fA,fFact.GetArray()+1);
00283                if (fVerbose)
00284                   Info("Decompose()","reseting to: %d",nFact);
00285                fRPessimism *= 1.1;
00286                break;
00287             }
00288          case -5 :
00289             if (fVerbose) {
00290                Info("Decompose()","matrix apparently numerically singular");
00291                Info("Decompose()","detected at stage %d",this->IError());
00292                Info("Decompose()","accept this factorization and hope for the best..");
00293             }
00294             done = 1;
00295             break;
00296          case -6 :
00297             if (fVerbose) {
00298                Info("Decompose()","change of sign of pivots detected at stage %d",this->IError());
00299                Info("Decompose()","but who cares ");
00300             }
00301             done = 1;
00302             break;
00303          case -7 :
00304             Error("Decompose()","value of fNsteps out of range: %d",fNsteps);
00305             return kFALSE;
00306          case 1 :
00307             if (fVerbose) {
00308                Info("Decompose()","detected %d entries out of range in row/column index",this->IError());
00309                Info("Decompose()","they are ignored");
00310             }
00311             done = 1;
00312             break;
00313          case 3 :
00314             if (fVerbose)
00315                Info("Decompose()","rank deficient matrix detected; apparent rank = %d",this->IError());
00316             done = 1;
00317             break;
00318          default:
00319             break;
00320       }
00321 
00322       tries++;
00323    } while (!done && tries < 10);
00324 
00325    Int_t ok;
00326    if ( !done && tries >= 10) {
00327       ok = kFALSE;
00328       if (fVerbose)
00329          Error("Decompose()","did not get a factorization after 10 tries");
00330    } else {
00331       ok = kTRUE;
00332       SetBit(kDecomposed);
00333    }
00334 
00335    return ok;
00336 }
00337 
00338 //______________________________________________________________________________
00339 Bool_t TDecompSparse::Solve(TVectorD &b)
00340 {
00341 // Solve Ax=b . Solution returned in b.
00342 
00343    R__ASSERT(b.IsValid());
00344    if (TestBit(kSingular)) {
00345       Error("Solve()","Matrix is singular");
00346       return kFALSE;
00347    }
00348    if ( !TestBit(kDecomposed) ) {
00349       if (!Decompose()) {
00350          Error("Solve()","Decomposition failed");
00351          return kFALSE;
00352       }
00353    }
00354 
00355    if (fNrows != b.GetNrows() || fRowLwb != b.GetLwb())
00356    {
00357       Error("Solve(TVectorD &","vector and matrix incompatible");
00358       return kFALSE;
00359    }
00360    b.Shift(-fRowLwb); // make sure rowlwb = 0
00361 
00362    // save bs and store residuals
00363    TVectorD resid = b;
00364    TVectorD bSave = b;
00365 
00366    Double_t bnorm = b.NormInf();
00367    Double_t rnorm = 0.0;
00368 
00369    Int_t done = 0;
00370    Int_t refactorizations = 0;
00371 
00372    while (!done && refactorizations < 10) {
00373 
00374       Solve(fNrows,fFact,fIw,fW,fMaxfrt,b,fIw1,fNsteps,fIcntl,fInfo);
00375 
00376       // compute residuals
00377       resid = fA*b-resid;
00378       rnorm = resid.NormInf();
00379 
00380       if (rnorm < fPrecision*(1.+bnorm)) {
00381          // residuals are small enough, use this solution
00382          done = 1;
00383       } else if (this->GetThresholdPivoting() >= kThresholdPivotingMax
00384                   || refactorizations > 10)  {
00385          // ThresholdPivoting parameter is already too high; give up and
00386          // use this solution, whatever it is (the algorithm may bomb in
00387          // an outer loop).
00388          done = 1;
00389       } else {
00390          // refactor with a higher Threshold Pivoting parameter
00391          Double_t tp = this->GetThresholdPivoting();
00392          tp *= kThresholdPivotingFactor;
00393          if (tp > kThresholdPivotingMax) tp = kThresholdPivotingMax;
00394          this->SetThresholdPivoting(tp);
00395          if (fVerbose)
00396             Info("Solve","Setting ThresholdPivoting parameter to %.4e for future factorizations",
00397                   this->GetThresholdPivoting());
00398 
00399          SetMatrix(fA);
00400          refactorizations++;
00401          resid = bSave;
00402          b     = bSave;
00403       }
00404    }
00405 
00406    b.Shift(fRowLwb);
00407    return kTRUE;
00408 }
00409 
00410 //______________________________________________________________________________
00411 void TDecompSparse::InitParam()
00412 {
00413 // initializing control parameters
00414 
00415    fPrecision  = kInitPrecision;
00416    fIPessimism = 1.2;
00417    fRPessimism = 1.2;
00418 
00419    const Int_t ifrlvl = 5;
00420 
00421    SetVerbose(fVerbose);
00422    fIcntl[4] = 2139062143;
00423    fIcntl[5] = 1;
00424    fIcntl[ifrlvl+1]  = 32639;
00425    fIcntl[ifrlvl+2]  = 32639;
00426    fIcntl[ifrlvl+3]  = 32639;
00427    fIcntl[ifrlvl+4]  = 32639;
00428    fIcntl[ifrlvl+5]  = 14;
00429    fIcntl[ifrlvl+6]  = 9;
00430    fIcntl[ifrlvl+7]  = 8;
00431    fIcntl[ifrlvl+8]  = 8;
00432    fIcntl[ifrlvl+9]  = 9;
00433    fIcntl[ifrlvl+10] = 10;
00434    fIcntl[ifrlvl+11] = 32639;
00435    fIcntl[ifrlvl+12] = 32639;
00436    fIcntl[ifrlvl+13] = 32639;
00437    fIcntl[ifrlvl+14] = 32689;
00438    fIcntl[ifrlvl+15] = 24;
00439    fIcntl[ifrlvl+16] = 11;
00440    fIcntl[ifrlvl+17] = 9;
00441    fIcntl[ifrlvl+18] = 8;
00442    fIcntl[ifrlvl+19] = 9;
00443    fIcntl[ifrlvl+20] = 10;
00444    fIcntl[26] = 0;
00445    fIcntl[27] = 0;
00446    fIcntl[28] = 0;
00447    fIcntl[29] = 0;
00448    fIcntl[30] = 0;
00449    fCntl[1] = 0.10;
00450    fCntl[2] = 1.00;
00451    fCntl[3] = 0.00;
00452    fCntl[4] = 0.0;
00453    fCntl[5] = 0.0;
00454 
00455    // set initial value of "Treat As Zero" parameter
00456    this->SetTreatAsZero(kInitTreatAsZero);
00457 
00458    // set initial value of Threshold parameter
00459    this->SetThresholdPivoting(kInitThresholdPivoting);
00460 
00461    fNsteps    = 0;
00462    fMaxfrt    = 0;
00463    fNrows     = 0;
00464    fNnonZeros = 0;
00465 }
00466 
00467 //______________________________________________________________________________
00468 void TDecompSparse::InitPivot(const Int_t n,const Int_t nz,TArrayI &Airn,TArrayI &Aicn,
00469                               TArrayI &Aiw,TArrayI &Aikeep,TArrayI &Aiw1,Int_t &nsteps,
00470                               const Int_t iflag,Int_t *icntl,Double_t *cntl,Int_t *info,
00471                               Double_t &ops)
00472 {
00473 // Setup Pivoting variables
00474 
00475    Int_t i,iwfr,k,l1,l2,lliw;
00476 
00477    Int_t *irn      = Airn.GetArray();
00478    Int_t *icn      = Aicn.GetArray();
00479    Int_t *iw       = Aiw.GetArray();
00480    Int_t *ikeep    = Aikeep.GetArray();
00481    Int_t *iw1      = Aiw1.GetArray();
00482    const Int_t liw = Aiw.GetSize()-1;
00483 
00484    for (i = 1; i < 16; i++)
00485       info[i] = 0;
00486 
00487    if (icntl[3] > 0 && icntl[2] > 0) {
00488       ::Info("TDecompSparse::InitPivot","Start with n = %d  nz = %d  liw = %d  iflag = %d",n,nz,liw,iflag);
00489       nsteps = 0;
00490       k = TMath::Min(8,nz);
00491       if (icntl[3] > 1) k = nz;
00492       if (k > 0) {
00493          printf("matrix non-zeros:\n");
00494          for (i = 1; i < k+1; i++) {
00495             printf("%d %d ",irn[i],icn[i]);
00496             if (i%5 == 0 || i == k) printf("\n");
00497          }
00498       }
00499 
00500       k = TMath::Min(10,n);
00501       if (icntl[3] > 1) k = n;
00502       if (iflag == 1 && k > 0) {
00503          for (i = 1; i < k+1; i++) {
00504             printf("%d ",ikeep[i]);
00505             if (i%10 == 0 || i == k) printf("\n");
00506          }
00507       }
00508    }
00509 
00510    if (n >= 1 && n <= icntl[4]) {
00511       if (nz < 0) {
00512          info[1] = -2;
00513          if (icntl[1] > 0)
00514             ::Error("TDecompSparse::InitPivot","info[1]= %d; value of nz out of range .. = %d",info[1],nz);
00515          return;
00516       }
00517       lliw = liw-2*n;
00518       l1 = lliw+1;
00519       l2 = l1+n;
00520       if (iflag != 1) {
00521          if (liw < 2*nz+3*n+1) {
00522             info[1] = -3;
00523             info[2] = 2*nz+3*n+1;
00524             if (icntl[1] > 0)
00525                ::Error("TDecompSparse::InitPivot","info[1]= %d; liw too small, must be increased from %d to at least %d",info[1],liw,info[2]);
00526             return;
00527          }
00528          InitPivot_sub1(n,nz,irn,icn,iw,iw1,iw1+n+1,iw+l1-1,iwfr,icntl,info);
00529          InitPivot_sub2(n,iw1,iw,lliw,iwfr,iw+l1-1,iw+l2-1,ikeep+n+1,
00530                         ikeep+2*(n+1),ikeep,icntl[4],info[11],cntl[2]);
00531       } else {
00532          if (liw < nz+3*n+1) {
00533             info[1] = -3;
00534             info[2] = nz+3*n+1;
00535             if (icntl[1] > 0)
00536                ::Error("TDecompSparse::InitPivot","info[1]= %d; liw too small, must be increased from %d to at least %d",info[1],liw,info[2]);
00537             return;
00538          }
00539          InitPivot_sub3(n,nz,irn,icn,ikeep,iw,iw1,iw1+n+1,iw+l1-1,iwfr,icntl,info);
00540          InitPivot_sub4(n,iw1,iw,lliw,iwfr,ikeep,ikeep+n+1,iw+l1-1,iw+l2-1,info[11]);
00541       }
00542       InitPivot_sub5(n,iw1,iw+l1-1,ikeep,ikeep+n+1,ikeep+2*(n+1),iw+l2-1,nsteps,icntl[5]);
00543       if (nz >= 1) iw[1] = irn[1]+1;
00544       InitPivot_sub6(n,nz,irn,icn,ikeep,ikeep+2*(n+1),ikeep+n+1,iw+l2-1,
00545                      nsteps,iw1,iw1+n+1,iw,info,ops);
00546    } else {
00547       info[1] = -1;
00548       if (icntl[1] > 0)
00549          ::Error("TDecompSparse::InitPivot","info[1]= %d; value of n out of range ... = %d",info[1],n);
00550       return;
00551    }
00552 
00553    if (icntl[3] <= 0 || icntl[2] <= 0) return;
00554 
00555    printf("Leaving with nsteps =%d info(1)=%d ops=%14.5e ierror=%d\n",nsteps,info[1],ops,info[2]);
00556    printf("nrltot=%d nirtot=%d nrlnec=%d nirnec=%d nrladu=%d niradu=%d ncmpa=%d\n",
00557            info[3],info[4],info[5],info[6],info[7],info[8],info[11]);
00558 
00559    k = TMath::Min(9,n);
00560    if (icntl[3] > 1) k = n;
00561    if (k > 0) {
00562       printf("ikeep[0][.]=\n");
00563       for (i = 1; i < k+1; i++) {
00564          printf("%d ",ikeep[i]);
00565          if (k%10 == 0 || i == k) printf("\n");
00566       }
00567    }
00568    k = TMath::Min(k,nsteps);
00569    if (k > 0) {
00570       printf("ikeep[2][.]=\n");
00571       for (i = 1; i < k+1; i++) {
00572          printf("%d ",ikeep[2*(n+1)+i]);
00573          if (k%10 == 0 || i == k) printf("\n");
00574       }
00575    }
00576 }
00577 
00578 //______________________________________________________________________________
00579 void TDecompSparse::Factor(const Int_t n,const Int_t nz,TArrayI &Airn,TArrayI &Aicn,TArrayD &Aa,
00580                            TArrayI &Aiw,TArrayI &Aikeep,const Int_t nsteps,Int_t &maxfrt,
00581                            TArrayI &Aiw1,Int_t *icntl,Double_t *cntl,Int_t *info)
00582 {
00583 // Factorization routine, the workhorse for the decompostion step
00584 
00585    Int_t i,iapos,iblk,ipos,irows,j1,j2,jj,k,kblk,kz,len,ncols,nrows,nz1;
00586 
00587    Int_t    *irn   = Airn.GetArray();
00588    Int_t    *icn   = Aicn.GetArray();
00589    Int_t    *iw    = Aiw.GetArray();
00590    Int_t    *ikeep = Aikeep.GetArray();
00591    Int_t    *iw1   = Aiw1.GetArray();
00592    Double_t *a     = Aa.GetArray();
00593 
00594    const Int_t la = Aa.GetSize()-1;
00595    const Int_t liw = Aiw.GetSize()-1;
00596 
00597    info[1] = 0;
00598    if (icntl[3] > 0 && icntl[2] > 0) {
00599       printf("entering Factor with n=%d nz=%d la=%d liw=%d nsteps=%d u=%10.2e\n",
00600                n,nz,la,liw,nsteps,cntl[1]);
00601       kz = TMath::Min(6,nz);
00602       if (icntl[3] > 1) kz = nz;
00603       if (nz > 0) {
00604          printf("matrix non-zeros:\n");
00605          for (i = 1; i < kz+1; i++) {
00606             printf("%16.3e %d %d ",a[i],irn[i],icn[i]);
00607             if (i%2 == 0 || i==kz) printf("\n");
00608          }
00609       }
00610       k = TMath::Min(9,n);
00611       if (icntl[3] > 1) k = n;
00612       if (k > 0) {
00613          printf("ikeep(0,.)=\n");
00614          for (i = 1; i < k+1; i++) {
00615             printf("%d ",ikeep[i]);
00616             if (i%10 == 0 || i == k) printf("\n");
00617          }
00618       }
00619       k = TMath::Min(k,nsteps);
00620       if (k > 0) {
00621          printf("ikeep(1,.)=\n");
00622          for (i = 1; i < k+1; i++) {
00623             printf("%d ",ikeep[n+1+i]);
00624             if (i%10 == 0 || i == k) printf("\n");
00625          }
00626          printf("ikeep(2,.)=\n");
00627          for (i = 1; i < k+1; i++) {
00628             printf("%d ",ikeep[2*(n+1)+i]);
00629             if (i%10 == 0 || i == k) printf("\n");
00630          }
00631       }
00632    }
00633 
00634    if (n < 1 || n > icntl[4])
00635       info[1] = -1;
00636    else if (nz < 0)
00637       info[1] = -2;
00638    else if (liw < nz) {
00639       info[1] = -3;
00640       info[2] = nz;
00641    } else if (la < nz+n) {
00642       info[1] = -4;
00643       info[2] = nz+n;
00644    } else if (nsteps < 1 || nsteps > n)
00645       info[1] = -7;
00646    else {
00647       Factor_sub1(n,nz,nz1,a,la,irn,icn,iw,liw,ikeep,iw1,icntl,info);
00648       if (info[1] != -3 && info[1] != -4) {
00649          Factor_sub2(n,nz1,a,la,iw,liw,ikeep,ikeep+2*(n+1),nsteps,maxfrt,ikeep+n+1,iw1,icntl,cntl,info);
00650          if (info[1] == 3 && icntl[2] > 0)
00651             ::Warning("TDecompSparse::Factor","info[1]= %d; matrix is singular. rank=%d",info[1],info[2]);
00652       }
00653    }
00654 
00655    if (icntl[1] > 0) {
00656       switch(info[1]) {
00657          case -1:
00658             ::Error("TDecompSparse::Factor","info[1]= %d; value of n out of range ... =%d",info[1],n);
00659             break;
00660 
00661          case -2:
00662             ::Error("TDecompSparse::Factor","info[1]= %d; value of nz out of range ... =%d",info[1],nz);
00663             break;
00664 
00665          case -3:
00666             ::Error("TDecompSparse::Factor","info[1]= %d; liw too small, must be increased from %d to at least %d",info[1],liw,info[2]);
00667             break;
00668 
00669          case -4:
00670             ::Error("TDecompSparse::Factor","info[1]= %d; la too small, must be increased from %d to at least %d",info[1],la,info[2]);
00671             break;
00672 
00673          case -5:
00674             ::Error("TDecompSparse::Factor","info[1]= %d; zero pivot at stage %d zero pivot at stage",info[1],info[2]);
00675             break;
00676 
00677          case -6:
00678             ::Error("TDecompSparse::Factor","info[1]= %d; change in sign of pivot encountered when factoring allegedly definite matrix",info[1]);
00679             break;
00680 
00681          case -7:
00682             ::Error("TDecompSparse::Factor","info[1]= %d; nsteps is out of range",info[1]);
00683             break;
00684       }
00685    }
00686 
00687    if (icntl[3] <= 0 || icntl[2] <= 0 || info[1] < 0)
00688       return;
00689 
00690    ::Info("TDecompSparse::Factor","leaving Factor with maxfrt=%d info[1]=%d nrlbdu=%d nirbdu=%d ncmpbr=%d ncmpbi=%d ntwo=%d ierror=%d",
00691           maxfrt,info[1],info[9],info[10],info[12],info[13],info[14],info[2]);
00692 
00693    if (info[1] < 0) return;
00694 
00695    kblk = TMath::Abs(iw[1]+0);
00696    if (kblk == 0) return;
00697    if (icntl[3] == 1) kblk = 1;
00698    ipos = 2;
00699    iapos = 1;
00700 
00701    for (iblk = 1; iblk < kblk+1; iblk++) {
00702       ncols = iw[ipos];
00703       nrows = iw[ipos+1];
00704       j1 = ipos+2;
00705       if (ncols <= 0) {
00706          ncols = -ncols;
00707          nrows = 1;
00708          j1 = j1-1;
00709       }
00710       ::Info("TDecompSparse::Factor","block pivot =%d nrows =%d ncols =%d",iblk,nrows,ncols);
00711       j2 = j1+ncols-1;
00712       ipos = j2+1;
00713 
00714       printf(" column indices =\n");
00715       for (jj = j1; jj < j2+1; jj++) {
00716          printf("%d ",iw[jj]);
00717          if (jj%10 == 0 || jj == j2) printf("\n");
00718       }
00719 
00720       printf(" real entries .. each row starts on a new line\n");
00721       len = ncols;
00722       for (irows = 1; irows < nrows+1; irows++) {
00723          j1 = iapos;
00724          j2 = iapos+len-1;
00725          for (jj = j1; jj < j2+1; jj++) {
00726             printf("%13.4e ",a[jj]);
00727             if (jj%5 == 0 || jj == j2) printf("\n");
00728          }
00729          len = len-1;
00730          iapos = j2+1;
00731       }
00732    }
00733 }
00734 
00735 //______________________________________________________________________________
00736 void TDecompSparse::Solve(const Int_t n,TArrayD &Aa,TArrayI &Aiw,
00737                           TArrayD &Aw,const Int_t maxfrt,TVectorD &b,TArrayI &Aiw1,
00738                           const Int_t nsteps,Int_t *icntl,Int_t *info)
00739 {
00740 // Main routine for solving Ax=b
00741 
00742    Int_t i,iapos,iblk,ipos,irows,j1,j2,jj,k,kblk,latop,len,nblk,ncols,nrows;
00743 
00744    Double_t *a   = Aa.GetArray();
00745    Double_t *w   = Aw.GetArray();
00746    Int_t    *iw  = Aiw.GetArray();
00747    Int_t    *iw1 = Aiw1.GetArray();
00748    Double_t *rhs = new Double_t[n+1];
00749    rhs[0] = 0.;
00750    memcpy(rhs+1,b.GetMatrixArray(),n*sizeof(Double_t));
00751    const Int_t la  = Aa.GetSize()-1;
00752    const Int_t liw = Aiw.GetSize()-1;
00753 
00754    info[1] = 0;
00755    k = 0;
00756    if (icntl[3] > 0 && icntl[2] > 0) {
00757       printf("nentering Solve with n=%d la=%d liw=%d maxfrt=%d nsteps=%d",n,la,liw,maxfrt,nsteps);
00758 
00759       kblk = TMath::Abs(iw[1]+0);
00760       if (kblk != 0) {
00761          if (icntl[3] == 1) kblk = 1;
00762          ipos = 2;
00763          iapos = 1;
00764          for (iblk = 1; iblk < kblk+1; iblk++) {
00765             ncols = iw[ipos];
00766             nrows = iw[ipos+1];
00767             j1 = ipos+2;
00768             if (ncols <= 0) {
00769                ncols = -ncols;
00770                nrows = 1;
00771                j1 = j1-1;
00772             }
00773             printf("block pivot=%d nrows=%d ncols=%d\n",iblk,nrows,ncols);
00774             j2 = j1+ncols-1;
00775             ipos = j2+1;
00776             printf("column indices =\n");
00777             for (jj = j1; jj < j2+1; jj++) {
00778                printf("%d ",iw[jj]);
00779                if (jj%10 == 0 || jj == j2) printf("\n");
00780             }
00781             printf("real entries .. each row starts on a new line\n");
00782             len = ncols;
00783             for (irows = 1; irows < nrows+1; irows++) {
00784                j1 = iapos;
00785                j2 = iapos+len-1;
00786                for (jj = j1; jj < j2+1; jj++) {
00787                   printf("%13.3e ",a[jj]);
00788                   if (jj%5 == 0 || jj == j2) printf("\n");
00789                }
00790                len = len-1;
00791                iapos = j2+1;
00792             }
00793          }
00794       }
00795 
00796       k = TMath::Min(10,n);
00797       if (icntl[3] > 1) k = n;
00798       if (n > 0) {
00799          printf("rhs =\n");
00800          for (i = 1; i < k+1; i++) {
00801             printf("%13.3e ",rhs[i]);
00802             if (i%5 == 0 || i == k) printf("\n");
00803          }
00804       }
00805    }
00806 
00807    nblk = 0;
00808    if (iw[1] == 0) {
00809       nblk = 0;
00810       for (i = 1; i < n+1; i++)
00811          rhs[i] = 0.0;
00812    } else {
00813       nblk = (iw[1] <= 0) ? -iw[1] : iw[1];
00814       Solve_sub1(n,a,iw+1,w,rhs,iw1,nblk,latop,icntl);
00815       Solve_sub2(n,a,iw+1,w,rhs,iw1,nblk,latop,icntl);
00816    }
00817 
00818    if (icntl[3] > 0 && icntl[2] > 0) {
00819       printf("leaving Solve with:\n");
00820       if (n > 0) {
00821          printf("rhs =\n");
00822          for (i = 1; i < k+1; i++) {
00823             printf("%13.3e ",rhs[i]);
00824             if (i%5 == 0 || i == k) printf("\n");
00825          }
00826       }
00827    }
00828 
00829    memcpy(b.GetMatrixArray(),rhs+1,n*sizeof(Double_t));
00830    delete [] rhs;
00831 }
00832 
00833 //______________________________________________________________________________
00834 void TDecompSparse::InitPivot_sub1(const Int_t n,const Int_t nz,Int_t *irn,Int_t *icn,
00835                                    Int_t *iw,Int_t *ipe,Int_t *iq,Int_t *flag,
00836                                    Int_t &iwfr,Int_t *icntl,Int_t *info)
00837 {
00838 // Help routine for pivoting setup
00839 
00840    Int_t i,id,j,jn,k,k1,k2,l,last,lr,n1,ndup;
00841 
00842    info[2] = 0;
00843    for (i = 1; i < n+1; i++)
00844       ipe[i] = 0;
00845    lr = nz;
00846 
00847    if (nz != 0) {
00848       for (k = 1; k < nz+1; k++) {
00849          i = irn[k];
00850          j = icn[k];
00851 
00852          Bool_t outRange = (i < 1 || i > n || j < 1 || j > n);
00853          if (outRange) {
00854             info[2] = info[2]+1;
00855             info[1] = 1;
00856             if (info[2] <= 1 && icntl[2]> 0)
00857                ::Warning("TDecompSparse::InitPivot_sub1","info[1]= %d; %d th non-zero (in row=%d and column=%d) ignored",info[1],k,i,j);
00858          }
00859 
00860          if (outRange || i == j) {
00861             i = 0;
00862             j = 0;
00863          } else {
00864             ipe[i] = ipe[i]+1;
00865             ipe[j] = ipe[j]+1;
00866          }
00867          iw[k] = j;
00868          lr = lr+1;
00869          iw[lr] = i;
00870       }
00871    }
00872 
00873    iq[1] = 1;
00874    n1 = n-1;
00875    if (n1 > 0) {
00876       for (i = 1; i < n1+1; i++) {
00877          flag[i] = 0;
00878          if (ipe[i] == 0) ipe[i] = -1;
00879          iq[i+1] = ipe[i]+iq[i]+1;
00880          ipe[i] = iq[i];
00881       }
00882    }
00883 
00884    last = ipe[n]+iq[n];
00885    flag[n] = 0;
00886    if (lr < last) {
00887       k1 = lr+1;
00888       for (k = k1; k < last+1; k++)
00889          iw[k] = 0;
00890    }
00891    ipe[n] = iq[n];
00892    iwfr = last+1;
00893    if (nz != 0) {
00894       for (k = 1; k < nz+1; k++) {
00895          j = iw[k];
00896          if (j <= 0)  continue;
00897          l = k;
00898          iw[k] = 0;
00899          for (id = 1; id < nz+1; id++) {
00900             if (l <= nz)
00901                l = l+nz;
00902             else
00903                l = l-nz;
00904             i = iw[l];
00905             iw[l] = 0;
00906             if (i >= j) {
00907                l = iq[j]+1;
00908                iq[j] = l;
00909                jn = iw[l];
00910                iw[l] = -i;
00911             } else {
00912                l = iq[i]+1;
00913                iq[i] = l;
00914                jn = iw[l];
00915                iw[l] = -j;
00916             }
00917             j = jn;
00918             if (j <= 0) break;
00919          }
00920       }
00921    }
00922 
00923    ndup = 0;
00924 
00925    for (i = 1; i < n+1; i++) {
00926       k1 = ipe[i]+1;
00927       k2 = iq[i];
00928       if (k1 > k2) {
00929          ipe[i] = 0;
00930          iq[i] = 0;
00931       } else {
00932          for (k = k1; k < k2+1; k++) {
00933             j = -iw[k];
00934             if (j <= 0) break;
00935             l = iq[j]+1;
00936             iq[j] = l;
00937             iw[l] = i;
00938             iw[k] = j;
00939             if (flag[j] == i) {
00940                ndup = ndup + 1;
00941                iw[l] = 0;
00942                iw[k] = 0;
00943             }
00944             flag[j] = i;
00945          }
00946 
00947          iq[i] = iq[i]-ipe[i];
00948          if (ndup == 0) iw[k1-1] = iq[i];
00949       }
00950    }
00951 
00952    if (ndup != 0) {
00953       iwfr = 1;
00954       for (i = 1; i < n+1; i++) {
00955          k1 = ipe[i]+1;
00956          if (k1 == 1) continue;
00957          k2 = iq[i]+ipe[i];
00958          l = iwfr;
00959          ipe[i] = iwfr;
00960          iwfr = iwfr+1;
00961          for (k = k1; k < k2+1; k++) {
00962             if (iw[k] == 0) continue;
00963             iw[iwfr] = iw[k];
00964             iwfr = iwfr+1;
00965          }
00966          iw[l] = iwfr-l-1;
00967       }
00968    }
00969 
00970 }
00971 
00972 //______________________________________________________________________________
00973 void TDecompSparse::InitPivot_sub2(const Int_t n,Int_t *ipe,Int_t *iw,const Int_t lw,
00974                                    Int_t &iwfr,Int_t *nv,Int_t *nxt,Int_t *lst,Int_t *ipd,
00975                                    Int_t *flag,const Int_t iovflo,Int_t &ncmpa,
00976                                    const Double_t fratio)
00977 {
00978 // Help routine for pivoting setup
00979 
00980    Int_t i,id,idl,idn,ie,ip,is,jp,jp1,jp2,js,k,k1,k2,ke,kp,kp0,kp1,
00981          kp2,ks,l,len,limit,ln,ls,lwfr,md,me,ml,ms,nel,nflg,np,
00982          np0,ns,nvpiv,nvroot,root;
00983 
00984    for (i = 1; i < n+1; i++) {
00985       ipd[i]  = 0;
00986       nv[i]   = 1;
00987       flag[i] = iovflo;
00988    }
00989 
00990    js = 0;
00991    ms = 0;
00992    ncmpa = 0;
00993    md     = 1;
00994    nflg   = iovflo;
00995    nel    = 0;
00996    root   = n+1;
00997    nvroot = 0;
00998    for (is = 1; is < n+1; is++) {
00999       k = ipe[is];
01000       if (k > 0) {
01001          id = iw[k]+1;
01002          ns = ipd[id];
01003          if (ns > 0) lst[ns] = is;
01004          nxt[is] = ns;
01005          ipd[id] = is;
01006          lst[is] = -id;
01007       } else {
01008          nel = nel+1;
01009          flag[is] = -1;
01010          nxt[is]  = 0;
01011          lst[is]  = 0;
01012       }
01013    }
01014 
01015    for (ml = 1; ml < n+1; ml++) {
01016       if (nel+nvroot+1 >= n) break;
01017       for (id = md; id < n+1; id++) {
01018          ms = ipd[id];
01019          if (ms > 0) break;
01020       }
01021 
01022       md = id;
01023       nvpiv = nv[ms];
01024       ns = nxt[ms];
01025       nxt[ms] = 0;
01026       lst[ms] = 0;
01027       if (ns > 0) lst[ns] = -id;
01028       ipd[id] = ns;
01029       me = ms;
01030       nel = nel+nvpiv;
01031       idn = 0;
01032       kp = ipe[me];
01033       flag[ms] = -1;
01034       ip = iwfr;
01035       len = iw[kp];
01036       jp = 0;
01037       for (kp1 = 1; kp1 < len+1; kp1++) {
01038          kp = kp+1;
01039          ke = iw[kp];
01040          if (flag[ke] > -2) {
01041             if (flag[ke] <= 0) {
01042                if (ipe[ke] != -root) continue;
01043                ke = root;
01044                if (flag[ke] <= 0) continue;
01045             }
01046             jp = kp-1;
01047             ln = len-kp1+1;
01048             ie = ms;
01049          } else {
01050             ie = ke;
01051             jp = ipe[ie];
01052             ln = iw[jp];
01053          }
01054 
01055          for (jp1 = 1; jp1 < ln+1; jp1++) {
01056             jp = jp+1;
01057             is = iw[jp];
01058             if (flag[is] <= 0) {
01059                if (ipe[is] == -root) {
01060                   is = root;
01061                   iw[jp] = root;
01062                   if (flag[is] <= 0) continue;
01063                } else
01064                continue;
01065             }
01066             flag[is] = 0;
01067             if (iwfr >= lw) {
01068                ipe[ms] = kp;
01069                iw[kp] = len-kp1;
01070                ipe[ie] = jp;
01071                iw[jp] = ln-jp1;
01072                InitPivot_sub2a(n,ipe,iw,ip-1,lwfr,ncmpa);
01073                jp2 = iwfr-1;
01074                iwfr = lwfr;
01075                if (ip <= jp2) {
01076                   for (jp = ip; jp < jp2+1; jp++) {
01077                      iw[iwfr] = iw[jp];
01078                      iwfr = iwfr+1;
01079                   }
01080                }
01081                ip = lwfr;
01082                jp = ipe[ie];
01083                kp = ipe[me];
01084             }
01085             iw[iwfr] = is;
01086             idn = idn+nv[is];
01087             iwfr = iwfr+1;
01088             ls = lst[is];
01089             lst[is] = 0;
01090             ns = nxt[is];
01091             nxt[is] = 0;
01092             if (ns > 0) lst[ns] = ls;
01093             if (ls < 0) {
01094                ls = -ls;
01095                ipd[ls] = ns;
01096             } else if (ls > 0)
01097                nxt[ls] = ns;
01098          }
01099 
01100          if (ie == ms)
01101             break;
01102          ipe[ie] = -me;
01103          flag[ie] = -1;
01104       }
01105 
01106       nv[ms] = idn+nvpiv;
01107       if (iwfr != ip) {
01108          k1 = ip;
01109          k2 = iwfr-1;
01110          limit = TMath::Nint(fratio*(n-nel));
01111 
01112          for (k = k1; k < k2+1; k++) {
01113             is = iw[k];
01114             if (is == root) continue;
01115             if (nflg <= 2)  {
01116                for (i = 1; i < n+1; i++) {
01117                   if (flag[i] > 0)   flag[i] =  iovflo;
01118                   if (flag[i] <= -2) flag[i] = -iovflo;
01119                }
01120                nflg = iovflo;
01121             }
01122             nflg = nflg-1;
01123             id = idn;
01124             kp1 = ipe[is]+1;
01125             np  = kp1;
01126             kp2 = iw[kp1-1]+kp1-1;
01127 
01128             Int_t skip = 0;
01129             for (kp = kp1; kp < kp2+1; kp++) {
01130                ke = iw[kp];
01131                if (flag[ke] == -1) {
01132                   if (ipe[ke] != -root) continue;
01133                   ke = root;
01134                   iw[kp] = root;
01135                   if (flag[ke] == -1) continue;
01136                }
01137                if (flag[ke] >= 0) {
01138                   skip = 1;
01139                   break;
01140                }
01141                jp1 = ipe[ke]+1;
01142                jp2 = iw[jp1-1]+jp1-1;
01143                idl = id;
01144                for (jp = jp1; jp < jp2+1; jp++) {
01145                   js = iw[jp];
01146                   if (flag[js] <= nflg) continue;
01147                   id = id+nv[js];
01148                   flag[js] = nflg;
01149                }
01150                if (id <= idl) {
01151                   Int_t skip2 = 0;
01152                   for (jp = jp1; jp < jp2+1; jp++) {
01153                      js = iw[jp];
01154                      if (flag[js] != 0) {
01155                         skip2 = 1;
01156                         break;
01157                      }
01158                   }
01159                   if (skip2) {
01160                      iw[np] = ke;
01161                      flag[ke] = -nflg;
01162                      np = np+1;
01163                   } else {
01164                      ipe[ke] = -me;
01165                      flag[ke] = -1;
01166                   }
01167                } else {
01168                   iw[np] = ke;
01169                   flag[ke] = -nflg;
01170                   np = np+1;
01171                }
01172             }
01173 
01174             if (!skip)
01175                np0 = np;
01176             else {
01177                np0 = np;
01178                kp0 = kp;
01179                for (kp = kp0; kp < kp2+1; kp++) {
01180                   ks = iw[kp];
01181                   if (flag[ks] <= nflg) {
01182                      if (ipe[ks] == -root) {
01183                         ks = root;
01184                         iw[kp] = root;
01185                         if (flag[ks] <= nflg) continue;
01186                      } else
01187                         continue;
01188                   }
01189                   id = id+nv[ks];
01190                   flag[ks] = nflg;
01191                   iw[np] = ks;
01192                   np = np+1;
01193                }
01194             }
01195  
01196             Int_t doit = 2;
01197             if (id < limit) {
01198                iw[np] = iw[np0];
01199                iw[np0] = iw[kp1];
01200                iw[kp1] = me;
01201                iw[kp1-1] = np-kp1+1;
01202                js = ipd[id];
01203                for (l = 1; l < n+1; l++) {
01204                   if (js <= 0) {
01205                      doit = 3;
01206                      break;
01207                   }
01208                   kp1 = ipe[js]+1;
01209                   if (iw[kp1] != me) {
01210                      doit = 3;
01211                      break;
01212                   }
01213                   kp2 = kp1-1+iw[kp1-1];
01214                   Int_t stayInLoop = 0;
01215                   for (kp = kp1; kp < kp2+1; kp++) {
01216                      ie = iw[kp];
01217                      if (TMath::Abs(flag[ie]+0) > nflg) {
01218                         stayInLoop = 1;
01219                         break;
01220                      }
01221                   }
01222                   if (!stayInLoop) {
01223                      doit = 1;
01224                      break;
01225                   }
01226                   js = nxt[js];
01227                }
01228             }
01229  
01230             if (doit == 1) {
01231                ipe[js] = -is;
01232                nv[is] = nv[is]+nv[js];
01233                nv[js] = 0;
01234                flag[js] = -1;
01235                ns = nxt[js];
01236                ls = lst[js];
01237                if (ns > 0) lst[ns] = is;
01238                if (ls > 0) nxt[ls] = is;
01239                lst[is] = ls;
01240                nxt[is] = ns;
01241                lst[js] = 0;
01242                nxt[js] = 0;
01243                if (ipd[id] == js) ipd[id] = is;
01244             } else if (doit == 2) {
01245                if (nvroot == 0) {
01246                   root = is;
01247                   ipe[is] = 0;
01248                } else {
01249                   iw[k] = root;
01250                   ipe[is] = -root;
01251                   nv[root] = nv[root]+nv[is];
01252                   nv[is] = 0;
01253                   flag[is] = -1;
01254                }
01255                nvroot = nv[root];
01256             } else if (doit == 3) {
01257                ns = ipd[id];
01258                if (ns > 0) lst[ns] = is;
01259                nxt[is] = ns;
01260                ipd[id] = is;
01261                lst[is] = -id;
01262                md = TMath::Min(md,id);
01263             }
01264          }
01265 
01266          for (k = k1; k < k2+1; k++) {
01267             is = iw[k];
01268             if (nv[is] == 0) continue;
01269             flag[is] = nflg;
01270             iw[ip] = is;
01271             ip = ip+1;
01272          }
01273          iwfr = k1;
01274          flag[me] = -nflg;
01275          iw[ip] = iw[k1];
01276          iw[k1] = ip-k1;
01277          ipe[me] = k1;
01278          iwfr = ip+1;
01279       } else
01280          ipe[me] = 0;
01281    }
01282 
01283    for (is = 1; is < n+1; is++) {
01284       if (nxt[is] != 0 || lst[is] != 0) {
01285          if (nvroot == 0) {
01286             root = is;
01287             ipe[is] = 0;
01288          } else {
01289             ipe[is] = -root;
01290          }
01291          nvroot = nvroot+nv[is];
01292          nv[is] = 0;
01293       }
01294    }
01295 
01296    for (ie = 1; ie < n+1; ie++)
01297       if (ipe[ie] > 0) ipe[ie] = -root;
01298 
01299    if (nvroot> 0) nv[root] = nvroot;
01300 }
01301 
01302 //______________________________________________________________________________
01303 void TDecompSparse::InitPivot_sub2a(const Int_t n,Int_t *ipe,Int_t *iw,const Int_t lw,
01304                                     Int_t &iwfr,Int_t &ncmpa)
01305 {
01306 // Help routine for pivoting setup
01307 
01308    Int_t i,ir,k,k1,k2,lwfr;
01309 
01310    ncmpa = ncmpa+1;
01311    for (i = 1; i < n+1; i++) {
01312       k1 = ipe[i];
01313       if (k1 <= 0)  continue;
01314       ipe[i] = iw[k1];
01315       iw[k1] = -i;
01316    }
01317 
01318    iwfr = 1;
01319    lwfr = iwfr;
01320    for (ir = 1; ir < n+1; ir++) {
01321       if (lwfr > lw) break;
01322       Int_t skip = 1;
01323       for (k = lwfr; k < lw+1; k++) {
01324          if (iw[k] < 0) {
01325             skip = 0;
01326             break;
01327          }
01328       }
01329       if (skip) break;
01330       i = -iw[k];
01331       iw[iwfr] = ipe[i];
01332       ipe[i] = iwfr;
01333       k1 = k+1;
01334       k2 = k+iw[iwfr];
01335       iwfr = iwfr+1;
01336       if (k1 <= k2)  {
01337          for (k = k1; k < k2+1; k++) {
01338             iw[iwfr] = iw[k];
01339             iwfr = iwfr+1;
01340          }
01341       }
01342       lwfr = k2+1;
01343    }
01344 }
01345 
01346 //______________________________________________________________________________
01347 void TDecompSparse::InitPivot_sub3(const Int_t n,const Int_t nz,Int_t *irn,Int_t *icn,
01348                                    Int_t *perm,Int_t *iw,Int_t *ipe,Int_t *iq,
01349                                    Int_t *flag,Int_t &iwfr,Int_t *icntl,Int_t *info)
01350 {
01351 // Help routine for pivoting setup
01352 
01353    Int_t i,id,in,j,jdummy,k,k1,k2,l,lbig,len;
01354 
01355    info[1] = 0;
01356    info[2] = 0;
01357    for (i = 1; i < n+1; i++)
01358       iq[i] = 0;
01359 
01360    if (nz != 0) {
01361       for (k = 1; k < nz+1; k++) {
01362          i = irn[k];
01363          j = icn[k];
01364          iw[k] = -i;
01365 
01366          Bool_t outRange = (i < 1 || i > n || j < 1 || j > n);
01367          if (outRange) {
01368             info[2] = info[2]+1;
01369             info[1] = 1;
01370             if (info[2] <= 1 && icntl[2] > 0)
01371                ::Warning("TDecompSparse::InitPivot_sub3","info[1]= %d; %d 'th non-zero (in row %d and column %d) ignored",info[1],k,i,j);
01372          }
01373 
01374          if (outRange || i==j) {
01375             iw[k] = 0;
01376          } else {
01377             if (perm[j] <= perm[i])
01378               iq[j] = iq[j]+1;
01379             else
01380                iq[i] = iq[i]+1;
01381          }
01382       }
01383    }
01384 
01385    iwfr = 1;
01386    lbig = 0;
01387    for (i = 1; i < n+1; i++) {
01388       l = iq[i];
01389       lbig = TMath::Max(l,lbig);
01390       iwfr = iwfr+l;
01391       ipe[i] = iwfr-1;
01392    }
01393 
01394    if (nz != 0) {
01395       for (k = 1; k < nz+1; k++) {
01396          i = -iw[k];
01397          if (i <= 0) continue;
01398          l = k;
01399          iw[k] = 0;
01400          for (id = 1; id < nz+1; id++) {
01401             j = icn[l];
01402             if (perm[i] >= perm[j]) {
01403                l = ipe[j];
01404                ipe[j] = l-1;
01405                in = iw[l];
01406                iw[l] = i;
01407             } else {
01408                l = ipe[i];
01409                ipe[i] = l-1;
01410                in = iw[l];
01411                iw[l] = j;
01412             }
01413             i = -in;
01414             if (i <= 0) continue;
01415          }
01416       }
01417 
01418       k = iwfr-1;
01419       l = k+n;
01420       iwfr = l+1;
01421       for (i = 1; i < n+1; i++) {
01422          flag[i] = 0;
01423          j = n+1-i;
01424          len = iq[j];
01425          if (len > 0)  {
01426             for (jdummy = 1; jdummy < len+1; jdummy++) {
01427                iw[l] = iw[k];
01428                k = k-1;
01429                l = l-1;
01430             }
01431          }
01432          ipe[j] = l;
01433          l = l-1;
01434       }
01435 
01436       if (lbig < icntl[4]) {
01437          for (i = 1; i < n+1; i++) {
01438             k = ipe[i];
01439             iw[k] = iq[i];
01440             if (iq[i] == 0) ipe[i] = 0;
01441          }
01442       } else {
01443          iwfr = 1;
01444          for (i = 1; i < n+1; i++) {
01445             k1 = ipe[i]+1;
01446             k2 = ipe[i]+iq[i];
01447             if (k1 > k2) {
01448                ipe[i] = 0;
01449             } else {
01450                ipe[i] = iwfr;
01451                iwfr = iwfr+1;
01452                for (k = k1; k < k2+1; k++) {
01453                   j = iw[k];
01454                   if (flag[j] == i) continue;
01455                   iw[iwfr] = j;
01456                   iwfr = iwfr+1;
01457                   flag[j] = i;
01458                }
01459                k = ipe[i];
01460                iw[k] = iwfr-k-1;
01461             }
01462          }
01463       }
01464    }
01465 
01466 }
01467 
01468 //______________________________________________________________________________
01469 void TDecompSparse::InitPivot_sub4(const Int_t n,Int_t *ipe,Int_t *iw,const Int_t lw,
01470                                    Int_t &iwfr,Int_t *ips,Int_t *ipv,Int_t *nv,Int_t *flag,
01471                                    Int_t &ncmpa)
01472 {
01473 // Help routine for pivoting setup
01474 
01475    Int_t i,ie,ip,j,je,jp,jp1,jp2,js,kdummy,ln,lwfr,me,minjs,ml,ms;
01476 
01477    for (i = 1; i < n+1; i++) {
01478       flag[i] = 0;
01479       nv[i] = 0;
01480       j = ips[i];
01481       ipv[j] = i;
01482    }
01483 
01484    ncmpa = 0;
01485    for (ml = 1; ml < n+1; ml++) {
01486       ms = ipv[ml];
01487       me = ms;
01488       flag[ms] = me;
01489       ip = iwfr;
01490       minjs = n;
01491       ie = me;
01492 
01493       for (kdummy = 1; kdummy < n+1; kdummy++) {
01494          jp = ipe[ie];
01495          ln = 0;
01496          if (jp > 0) {
01497             ln = iw[jp];
01498             for (jp1 = 1; jp1 < ln+1; jp1++) {
01499                jp = jp+1;
01500                js = iw[jp];
01501                if (flag[js] == me) continue;
01502                flag[js] = me;
01503                if (iwfr >= lw) {
01504                   ipe[ie] = jp;
01505                   iw[jp] = ln-jp1;
01506                   InitPivot_sub2a(n,ipe,iw,ip-1,lwfr,ncmpa);
01507                   jp2 = iwfr-1;
01508                   iwfr = lwfr;
01509                   if (ip <= jp2)  {
01510                      for (jp = ip; jp < jp2+1; jp++) {
01511                         iw[iwfr] = iw[jp];
01512                         iwfr = iwfr+1;
01513                      }
01514                   }
01515                   ip = lwfr;
01516                   jp = ipe[ie];
01517                }
01518                iw[iwfr] = js;
01519                minjs = TMath::Min(minjs,ips[js]+0);
01520                iwfr = iwfr+1;
01521             }
01522          }
01523          ipe[ie] = -me;
01524          je = nv[ie];
01525          nv[ie] = ln+1;
01526          ie = je;
01527          if (ie == 0) break;
01528       }
01529 
01530       if (iwfr <= ip) {
01531          ipe[me] = 0;
01532          nv[me] = 1;
01533       } else {
01534          minjs = ipv[minjs];
01535          nv[me] = nv[minjs];
01536          nv[minjs] = me;
01537          iw[iwfr] = iw[ip];
01538          iw[ip] = iwfr-ip;
01539          ipe[me] = ip;
01540          iwfr = iwfr+1;
01541       }
01542    }
01543 }
01544 
01545 //______________________________________________________________________________
01546 void TDecompSparse::InitPivot_sub5(const Int_t n,Int_t *ipe,Int_t *nv,Int_t *ips,Int_t *ne,
01547                                    Int_t *na,Int_t *nd,Int_t &nsteps,const Int_t nemin)
01548 {
01549 // Help routine for pivoting setup
01550 
01551    Int_t i,ib,iff,il,is,ison,k,l,nr;
01552 
01553    il = 0;
01554    for (i = 1; i < n+1; i++) {
01555       ips[i] = 0;
01556       ne[i] = 0;
01557    }
01558    for (i = 1; i < n+1; i++) {
01559       if (nv[i] > 0) continue;
01560       iff = -ipe[i];
01561       is = -ips[iff];
01562       if (is > 0) ipe[i] = is;
01563       ips[iff] = -i;
01564    }
01565 
01566    nr = n+1;
01567    for (i = 1; i < n+1; i++) {
01568       if (nv[i] <= 0) continue;
01569       iff = -ipe[i];
01570       if (iff != 0) {
01571          is = -ips[iff];
01572          if (is > 0)
01573             ipe[i] = is;
01574          ips[iff] = -i;
01575       } else {
01576          nr = nr-1;
01577          ne[nr] = i;
01578       }
01579    }
01580 
01581    is = 1;
01582    i = 0;
01583    for (k = 1; k < n+1; k++) {
01584       if (i <= 0) {
01585          i = ne[nr];
01586          ne[nr] = 0;
01587          nr = nr+1;
01588          il = n;
01589          na[n] = 0;
01590       }
01591       for (l = 1; l < n+1; l++) {
01592          if (ips[i] >= 0) break;
01593          ison = -ips[i];
01594          ips[i] = 0;
01595          i = ison;
01596          il = il-1;
01597          na[il] = 0;
01598       }
01599 
01600       ips[i] = k;
01601       ne[is] = ne[is]+1;
01602       if (nv[i] > 0) {
01603          if (il < n) na[il+1] = na[il+1]+1;
01604          na[is] = na[il];
01605          nd[is] = nv[i];
01606 
01607          Bool_t doit = (na[is] == 1 && (nd[is-1]-ne[is-1] == nd[is])) ||
01608                        (na[is] != 1 && ne[is] < nemin && na[is] != 0 && ne[is-1] < nemin);
01609 
01610          if (doit) {
01611             na[is-1] = na[is-1]+na[is]-1;
01612             nd[is-1] = nd[is]+ne[is-1];
01613             ne[is-1] = ne[is]+ne[is-1];
01614             ne[is] = 0;
01615          } else {
01616             is = is+1;
01617          }
01618       }
01619 
01620       ib = ipe[i];
01621       if (ib >= 0) {
01622          if (ib > 0)
01623             na[il] = 0;
01624          i = ib;
01625       } else {
01626          i = -ib;
01627          il = il+1;
01628       }
01629    }
01630 
01631    nsteps = is-1;
01632 }
01633 
01634 //______________________________________________________________________________
01635 void TDecompSparse::InitPivot_sub6(const Int_t n,const Int_t nz,Int_t *irn,Int_t *icn,
01636                                    Int_t *perm,Int_t *na,Int_t *ne,Int_t *nd,const Int_t nsteps,
01637                                    Int_t *lstki,Int_t *lstkr,Int_t *iw,Int_t *info,Double_t &ops)
01638 {
01639 // Help routine for pivoting setup
01640 
01641    Int_t i,inew,iold,iorg,irow,istki,istkr,itop,itree,jold,jorg,k,lstk,nassr,nelim,nfr,nstk,
01642          numorg,nz1,nz2,nrladu,niradu,nirtot,nrltot,nirnec,nrlnec;
01643    Double_t delim;
01644 
01645    if (nz != 0 && irn[1] == iw[1]) {
01646       irn[1] = iw[1]-1;
01647       nz2 = 0;
01648       for (iold = 1; iold < n+1; iold++) {
01649          inew = perm[iold];
01650          lstki[inew] = lstkr[iold]+1;
01651          nz2 = nz2+lstkr[iold];
01652       }
01653       nz1 = nz2/2+n;
01654       nz2 = nz2+n;
01655    } else {
01656       for (i = 1; i < n+1; i++)
01657          lstki[i] = 1;
01658       nz1 = n;
01659       if (nz != 0) {
01660          for (i = 1; i < nz+1; i++) {
01661             iold = irn[i];
01662             jold = icn[i];
01663             if (iold < 1 || iold > n) continue;
01664             if (jold < 1 || jold > n) continue;
01665             if (iold == jold) continue;
01666             nz1 = nz1+1;
01667             irow = TMath::Min(perm[iold]+0,perm[jold]+0);
01668             lstki[irow] = lstki[irow]+1;
01669          }
01670       }
01671       nz2 = nz1;
01672    }
01673 
01674    ops = 0.0;
01675    istki = 0;
01676    istkr = 0;
01677    nrladu = 0;
01678    niradu = 1;
01679    nirtot = nz1;
01680    nrltot = nz1;
01681    nirnec = nz2;
01682    nrlnec = nz2;
01683    numorg = 0;
01684    itop = 0;
01685    for (itree = 1; itree < nsteps+1; itree++) {
01686       nelim = ne[itree];
01687       delim = Double_t(nelim);
01688       nfr = nd[itree];
01689       nstk = na[itree];
01690       nassr = nfr*(nfr+1)/2;
01691       if (nstk != 0) nassr = nassr-lstkr[itop]+1;
01692       nrltot = TMath::Max(nrltot,nrladu+nassr+istkr+nz1);
01693       nirtot = TMath::Max(nirtot,niradu+nfr+2+istki+nz1);
01694       nrlnec = TMath::Max(nrlnec,nrladu+nassr+istkr+nz2);
01695       nirnec = TMath::Max(nirnec,niradu+nfr+2+istki+nz2);
01696       for (iorg = 1; iorg < nelim+1; iorg++) {
01697          jorg = numorg+iorg;
01698          nz2 = nz2-lstki[jorg];
01699       }
01700       numorg = numorg+nelim;
01701       if (nstk > 0) {
01702          for (k = 1; k < nstk+1; k++) {
01703             lstk = lstkr[itop];
01704             istkr = istkr-lstk;
01705             lstk = lstki[itop];
01706             istki = istki-lstk;
01707             itop = itop-1;
01708          }
01709       }
01710       nrladu = nrladu+(nelim* (2*nfr-nelim+1))/2;
01711       niradu = niradu+2+nfr;
01712       if (nelim == 1) niradu = niradu-1;
01713       ops = ops+((nfr*delim*(nfr+1)-(2*nfr+1)*delim*(delim+1)/2+delim*(delim+1)*(2*delim+1)/6)/2);
01714       if (itree == nsteps || nfr == nelim) continue;
01715       itop = itop+1;
01716       lstkr[itop] = (nfr-nelim)* (nfr-nelim+1)/2;
01717       lstki[itop] = nfr-nelim+1;
01718       istki = istki+lstki[itop];
01719       istkr = istkr+lstkr[itop];
01720       nirtot = TMath::Max(nirtot,niradu+istki+nz1);
01721       nirnec = TMath::Max(nirnec,niradu+istki+nz2);
01722    }
01723 
01724    nrlnec = TMath::Max(nrlnec,n+TMath::Max(nz,nz1));
01725    nrltot = TMath::Max(nrltot,n+TMath::Max(nz,nz1));
01726    nrlnec = TMath::Min(nrlnec,nrltot);
01727    nirnec = TMath::Max(nz,nirnec);
01728    nirtot = TMath::Max(nz,nirtot);
01729    nirnec = TMath::Min(nirnec,nirtot);
01730    info[3] = nrltot;
01731    info[4] = nirtot;
01732    info[5] = nrlnec;
01733    info[6] = nirnec;
01734    info[7] = nrladu;
01735    info[8] = niradu;
01736 }
01737 
01738 //______________________________________________________________________________
01739 void TDecompSparse::Factor_sub1(const Int_t n,const Int_t nz,Int_t &nz1,Double_t *a,
01740                                 const Int_t la,Int_t *irn,Int_t *icn,Int_t *iw,const Int_t liw,
01741                                 Int_t *perm,Int_t *iw2,Int_t *icntl,Int_t *info)
01742 {
01743 // Help routine for factorization
01744 
01745    Int_t i,ia,ich,ii,iiw,inew,iold,ipos,j1,j2,jj,jnew,jold,jpos,k;
01746    Double_t anext,anow;
01747 
01748    const Double_t zero = 0.0;
01749    info[1] = 0;
01750    ia = la;
01751    for (iold = 1; iold < n+1; iold++) {
01752       iw2[iold] = 1;
01753       a[ia] = zero;
01754       ia = ia-1;
01755    }
01756 
01757    info[2] = 0;
01758    nz1 = n;
01759    if (nz != 0) {
01760       for (k = 1; k < nz+1; k++) {
01761          iold = irn[k];
01762          jold = icn[k];
01763          Bool_t outRange = (iold < 1 || iold > n || jold < 1 || jold > n);
01764 
01765          inew = perm[iold];
01766          jnew = perm[jold];
01767 
01768          if (!outRange && inew == jnew) {
01769             ia = la-n+iold;
01770             a[ia] = a[ia]+a[k];
01771             iw[k] = 0;
01772          } else {
01773             if (!outRange) {
01774                inew = TMath::Min(inew,jnew);
01775                iw2[inew] = iw2[inew]+1;
01776                iw[k] = -iold;
01777                nz1 = nz1+1;
01778             } else {
01779                info[1] = 1;
01780                info[2] = info[2]+1;
01781                if (info[2] <= 1 && icntl[2] > 0)
01782                   ::Warning("TDecompSparse::Factor_sub1","info[1]= %d; %d 'th non-zero (in row %d and column %d) ignored",
01783                             info[1],k,irn[k],icn[k]);
01784                iw[k] = 0;
01785             }
01786          }
01787       }
01788    }
01789 
01790    if (nz >= nz1 || nz1 == n) {
01791       k = 1;
01792       for (i = 1; i < n+1; i++) {
01793          k = k+iw2[i];
01794          iw2[i] = k;
01795       }
01796    } else {
01797       k = 1;
01798       for (i = 1; i < n+1; i++) {
01799          k = k+iw2[i]-1;
01800          iw2[i] = k;
01801       }
01802    }
01803 
01804    if (nz1 > liw) {
01805       info[1] = -3;
01806       info[2] = nz1;
01807       return;
01808    }
01809 
01810    if (nz1+n > la) {
01811       info[1] = -4;
01812       info[2] = nz1+n;
01813       return;
01814    }
01815 
01816    if (nz1 != n) {
01817       for (k = 1; k < nz+1; k++) {
01818          iold = -iw[k];
01819          if (iold <= 0) continue;
01820          jold = icn[k];
01821          anow = a[k];
01822          iw[k] = 0;
01823          for (ich = 1; ich < nz+1; ich++) {
01824             inew = perm[iold];
01825             jnew = perm[jold];
01826             inew = TMath::Min(inew,jnew);
01827             if (inew == perm[jold]) jold = iold;
01828             jpos = iw2[inew]-1;
01829             iold = -iw[jpos];
01830             anext = a[jpos];
01831             a[jpos] = anow;
01832             iw[jpos] = jold;
01833             iw2[inew] = jpos;
01834             if (iold == 0) break;
01835             anow = anext;
01836             jold = icn[jpos];
01837          }
01838       }
01839 
01840       if (nz < nz1) {
01841          ipos = nz1;
01842          jpos = nz1-n;
01843          for (ii = 1; ii < n+1; ii++) {
01844             i = n-ii+1;
01845             j1 = iw2[i];
01846             j2 = jpos;
01847             if (j1 <= jpos) {
01848                for (jj = j1; jj < j2+1; jj++) {
01849                   iw[ipos] = iw[jpos];
01850                   a[ipos] = a[jpos];
01851                   ipos = ipos-1;
01852                   jpos = jpos-1;
01853                }
01854             }
01855             iw2[i] = ipos+1;
01856             ipos = ipos-1;
01857          }
01858       }
01859    }
01860 
01861    for (iold = 1; iold < n+1; iold++) {
01862       inew = perm[iold];
01863       jpos = iw2[inew]-1;
01864       ia = la-n+iold;
01865       a[jpos] = a[ia];
01866       iw[jpos] = -iold;
01867    }
01868    ipos = nz1;
01869    ia = la;
01870    iiw = liw;
01871    for (i = 1; i < nz1+1; i++) {
01872       a[ia] = a[ipos];
01873       iw[iiw] = iw[ipos];
01874       ipos = ipos-1;
01875       ia = ia-1;
01876       iiw = iiw-1;
01877    }
01878 }
01879 
01880 //______________________________________________________________________________
01881 void TDecompSparse::Factor_sub2(const Int_t n,const Int_t nz,Double_t *a,const Int_t la,
01882                                 Int_t *iw,const Int_t liw,Int_t *perm,Int_t *nstk,
01883                                 const Int_t nsteps,Int_t &maxfrt,Int_t *nelim,Int_t *iw2,
01884                                 Int_t *icntl,Double_t *cntl,Int_t *info)
01885 {
01886 // Help routine for factorization
01887 
01888    Double_t amax,amult,amult1,amult2,detpiv,rmax,swop,thresh,tmax,uu;
01889    Int_t ainput,apos,apos1,apos2,apos3,astk,astk2,azero,i,iass;
01890    Int_t ibeg,idummy,iell,iend,iexch,ifr,iinput,ioldps,iorg,ipiv;
01891    Int_t ipmnp,ipos,irow,isnpiv,istk,istk2,iswop,iwpos,j,j1;
01892    Int_t j2,jcol,jdummy,jfirst,jj,jj1,jjj,jlast,jmax,jmxmip,jnew;
01893    Int_t jnext,jpiv,jpos,k,k1,k2,kdummy,kk,kmax,krow,laell,lapos2;
01894    Int_t liell,lnass,lnpiv,lt,ltopst,nass,nblk,newel,nfront,npiv;
01895    Int_t npivp1,ntotpv,numass,numorg,numstk,pivsiz,posfac,pospv1,pospv2;
01896    Int_t ntwo,neig,ncmpbi,ncmpbr,nrlbdu,nirbdu;
01897 
01898    const Double_t zero = 0.0;
01899    const Double_t half = 0.5;
01900    const Double_t one  = 1.0;
01901 
01902    detpiv = 0.0;
01903    apos3  = 0;
01904    isnpiv = 0;
01905    jmax   = 0;
01906    jpos   = 0;
01907 
01908    nblk = 0;
01909    ntwo = 0;
01910    neig = 0;
01911    ncmpbi = 0;
01912    ncmpbr = 0;
01913    maxfrt = 0;
01914    nrlbdu = 0;
01915    nirbdu = 0;
01916    uu = TMath::Min(cntl[1],half);
01917    uu = TMath::Max(uu,-half);
01918    for (i = 1; i < n+1; i++)
01919       iw2[i] = 0;
01920    iwpos = 2;
01921    posfac = 1;
01922    istk = liw-nz+1;
01923    istk2 = istk-1;
01924    astk = la-nz+1;
01925    astk2 = astk-1;
01926    iinput = istk;
01927    ainput = astk;
01928    azero = 0;
01929    ntotpv = 0;
01930    numass = 0;
01931 
01932    for (iass = 1; iass < nsteps+1; iass++) {
01933       nass = nelim[iass];
01934       newel = iwpos+1;
01935       jfirst = n+1;
01936       nfront = 0;
01937       numstk = nstk[iass];
01938       ltopst = 1;
01939       lnass = 0;
01940       if (numstk != 0) {
01941          j2 = istk-1;
01942          lnass = nass;
01943          ltopst = ((iw[istk]+1)*iw[istk])/2;
01944          for (iell = 1; iell < numstk+1; iell++) {
01945             jnext = jfirst;
01946             jlast = n+1;
01947             j1 = j2+2;
01948             j2 = j1-1+iw[j1-1];
01949             for (jj = j1; jj < j2+1; jj++) {
01950                j = iw[jj];
01951                if (iw2[j] > 0) continue;
01952                jnew = perm[j];
01953                if (jnew <= numass) nass = nass+1;
01954                for (idummy = 1; idummy < n+1; idummy++) {
01955                   if (jnext == n+1) break;
01956                   if (perm[jnext] > jnew) break;
01957                   jlast = jnext;
01958                   jnext = iw2[jlast];
01959                }
01960                if (jlast == n+1)
01961                   jfirst = j;
01962                else
01963                   iw2[jlast] = j;
01964                iw2[j] = jnext;
01965                jlast = j;
01966                nfront = nfront+1;
01967             }
01968          }
01969          lnass = nass-lnass;
01970       }
01971 
01972       numorg = nelim[iass];
01973       j1 = iinput;
01974       for (iorg = 1; iorg < numorg+1; iorg++) {
01975          j = -iw[j1];
01976          for (idummy = 1; idummy < liw+1; idummy++) {
01977             jnew = perm[j];
01978             if (iw2[j] <= 0) {
01979                jlast = n+1;
01980                jnext = jfirst;
01981                for (jdummy = 1; jdummy < n+1; jdummy++) {
01982                   if (jnext == n+1) break;
01983                   if (perm[jnext] > jnew) break;
01984                   jlast = jnext;
01985                   jnext = iw2[jlast];
01986                }
01987                if (jlast == n+1)
01988                   jfirst = j;
01989                else
01990                   iw2[jlast] = j;
01991                iw2[j] = jnext;
01992                nfront = nfront+1;
01993             }
01994             j1 = j1+1;
01995             if (j1 > liw) break;
01996             j = iw[j1];
01997             if (j < 0) break;
01998          }
01999       }
02000 
02001       if (newel+nfront >= istk)
02002          Factor_sub3(a,iw,istk,istk2,iinput,2,ncmpbr,ncmpbi);
02003       if (newel+nfront >= istk) {
02004          info[1] = -3;
02005          info[2] = liw+1+newel+nfront-istk;
02006          goto finish;
02007       }
02008 
02009       j = jfirst;
02010       for (ifr = 1; ifr < nfront+1; ifr++) {
02011          newel = newel+1;
02012          iw[newel] = j;
02013          jnext = iw2[j];
02014          iw2[j] = newel-(iwpos+1);
02015          j = jnext;
02016       }
02017 
02018       maxfrt = TMath::Max(maxfrt,nfront);
02019       iw[iwpos] = nfront;
02020       laell = ((nfront+1)*nfront)/2;
02021       apos2 = posfac+laell-1;
02022       if (numstk != 0) lnass = lnass*(2*nfront-lnass+1)/2;
02023 
02024       if (posfac+lnass-1 >= astk || apos2 >= astk+ltopst-1) {
02025          Factor_sub3(a,iw,astk,astk2,ainput,1,ncmpbr,ncmpbi);
02026          if (posfac+lnass-1 >= astk || apos2 >= astk+ltopst-1) {
02027             info[1] = -4;
02028             info[2] = la+TMath::Max(posfac+lnass,apos2-ltopst+2)-astk;
02029             goto finish;
02030          }
02031       }
02032 
02033       if (apos2 > azero) {
02034          apos = azero+1;
02035          lapos2 = TMath::Min(apos2,astk-1);
02036          if (lapos2 >= apos) {
02037             for (k= apos; k< lapos2+1; k++)
02038             a[k] = zero;
02039          }
02040          azero = apos2;
02041       }
02042 
02043       if (numstk != 0) {
02044          for (iell = 1; iell < numstk+1; iell++) {
02045             j1 = istk+1;
02046             j2 = istk+iw[istk];
02047             for (jj = j1; jj < j2+1; jj++) {
02048                irow = iw[jj];
02049                irow = iw2[irow];
02050                apos = posfac+IDiag(nfront,irow);
02051                for (jjj = jj; jjj < j2+1; jjj++) {
02052                   j = iw[jjj];
02053                   apos2 = apos+iw2[j]-irow;
02054                   a[apos2] = a[apos2]+a[astk];
02055                   a[astk] = zero;
02056                   astk = astk+1;
02057                }
02058             }
02059             istk = j2+1;
02060          }
02061       }
02062 
02063       for (iorg = 1; iorg < numorg+1; iorg++) {
02064          j = -iw[iinput];
02065          irow = iw2[j];
02066          apos = posfac+IDiag(nfront,irow);
02067          for (idummy = 1; idummy < nz+1; idummy++) {
02068             apos2 = apos+iw2[j]-irow;
02069             a[apos2] = a[apos2]+a[ainput];
02070             ainput = ainput+1;
02071             iinput = iinput+1;
02072             if (iinput > liw) break;
02073             j = iw[iinput];
02074             if (j < 0) break;
02075          }
02076       }
02077       numass = numass+numorg;
02078       j1 = iwpos+2;
02079       j2 = iwpos+nfront+1;
02080       for (k = j1; k < j2+1; k++) {
02081          j = iw[k];
02082          iw2[j] = 0;
02083       }
02084 
02085       lnpiv = -1;
02086       npiv = 0;
02087       for (kdummy = 1; kdummy < nass+1; kdummy++) {
02088          if (npiv == nass) break;
02089          if (npiv == lnpiv) break;
02090          lnpiv = npiv;
02091          npivp1 = npiv+1;
02092          jpiv = 1;
02093          for (ipiv = npivp1; ipiv < nass+1; ipiv++) {
02094             jpiv = jpiv-1;
02095             if (jpiv == 1) continue;
02096             apos = posfac+IDiag(nfront-npiv,ipiv-npiv);
02097 
02098             if (uu <= zero) {
02099                if (TMath::Abs(a[apos]) <= cntl[3]) {
02100                   info[1] = -5;
02101                   info[2] = ntotpv+1;
02102                   goto finish;
02103                }
02104                if (ntotpv <=  0) {
02105                   if (a[apos] > zero) isnpiv = 1;
02106                   if (a[apos] < zero) isnpiv = -1;
02107                }
02108                if ((a[apos] <= zero || isnpiv !=  1) && (a[apos] >= zero || isnpiv != -1)) {
02109                   if (info[1] != 2) info[2] = 0;
02110                   info[2] = info[2]+1;
02111                   info[1] = 2;
02112                   i = ntotpv+1;
02113                   if (icntl[2] > 0 && info[2] <= 10)
02114                      ::Warning("TDecompSparse::Factor_sub2","info[1]= %d; pivot %d has different sign from the previous one",
02115                                info[1],i);
02116                   isnpiv = -isnpiv;
02117                }
02118                if ((a[apos] > zero && isnpiv ==  1) || (a[apos] < zero && isnpiv == -1) || (uu == zero)) goto hack;
02119                info[1] = -6;
02120                info[2] = ntotpv+1;
02121                goto finish;
02122             }
02123 
02124             amax = zero;
02125             tmax = amax;
02126             j1 = apos+1;
02127             j2 = apos+nass-ipiv;
02128             if (j2 >= j1) {
02129                for (jj = j1; jj < j2+1; jj++) {
02130                   if (TMath::Abs(a[jj]) <= amax) continue;
02131                   jmax = ipiv+jj-j1+1;
02132                   amax = TMath::Abs(a[jj]);
02133                }
02134             }
02135             j1 = j2+1;
02136             j2 = apos+nfront-ipiv;
02137             if (j2 >= j1) {
02138                for (jj = j1; jj < j2+1; jj++)
02139                   tmax = TMath::Max(TMath::Abs(a[jj]),tmax);
02140             }
02141             rmax = TMath::Max(tmax,amax);
02142             apos1 = apos;
02143             kk = nfront-ipiv;
02144             lt = ipiv-(npiv+1);
02145             if (lt != 0) {
02146                for (k = 1; k < lt+1; k++) {
02147                   kk = kk+1;
02148                   apos1 = apos1-kk;
02149                   rmax = TMath::Max(rmax,TMath::Abs(a[apos1]));
02150                }
02151             }
02152             if (TMath::Abs(a[apos]) <= TMath::Max(cntl[3],uu*rmax)) {
02153                if (TMath::Abs(amax) <= cntl[3]) continue;
02154                apos2 = posfac+IDiag(nfront-npiv,jmax-npiv);
02155                detpiv = a[apos]*a[apos2]-amax*amax;
02156                thresh = TMath::Abs(detpiv);
02157                thresh = thresh/(uu*TMath::Max(TMath::Abs(a[apos])+amax,TMath::Abs(a[apos2])+amax));
02158                if (thresh <= rmax) continue;
02159                rmax = zero;
02160                j1 = apos2+1;
02161                j2 = apos2+nfront-jmax;
02162                if (j2 >= j1) {
02163                   for (jj = j1; jj < j2+1; jj++)
02164                      rmax = TMath::Max(rmax,TMath::Abs(a[jj]));
02165                }
02166                kk = nfront-jmax+1;
02167                apos3 = apos2;
02168                jmxmip = jmax-ipiv-1;
02169                if (jmxmip != 0) {
02170                   for (k = 1; k < jmxmip+1; k++) {
02171                      apos2 = apos2-kk;
02172                      kk = kk+1;
02173                      rmax = TMath::Max(rmax,TMath::Abs(a[apos2]));
02174                   }
02175                }
02176                ipmnp = ipiv-npiv-1;
02177                if (ipmnp != 0) {
02178                   apos2 = apos2-kk;
02179                   kk = kk+1;
02180                   for (k = 1; k < ipmnp+1; k++) {
02181                      apos2 = apos2-kk;
02182                      kk = kk+1;
02183                      rmax = TMath::Max(rmax,TMath::Abs(a[apos2]));
02184                   }
02185                }
02186                if (thresh <= rmax) continue;
02187                pivsiz = 2;
02188             } else {
02189                pivsiz = 1;
02190             }
02191 
02192             irow = ipiv-npiv;
02193             for (krow = 1; krow < pivsiz+1; krow++) {
02194                if (irow != 1) {
02195                   j1 = posfac+irow;
02196                   j2 = posfac+nfront-(npiv+1);
02197                   if (j2 >= j1) {
02198                      apos2 = apos+1;
02199                      for (jj = j1; jj < j2+1; jj++) {
02200                         swop = a[apos2];
02201                         a[apos2] = a[jj];
02202                         a[jj] = swop;
02203                         apos2 = apos2+1;
02204                      }
02205                   }
02206                   j1 = posfac+1;
02207                   j2 = posfac+irow-2;
02208                   apos2 = apos;
02209                   kk = nfront-(irow+npiv);
02210                   if (j2 >= j1) {
02211                      for (jjj = j1; jjj < j2+1; jjj++) {
02212                         jj = j2-jjj+j1;
02213                         kk = kk+1;
02214                         apos2 = apos2-kk;
02215                         swop = a[apos2];
02216                         a[apos2] = a[jj];
02217                         a[jj] = swop;
02218                      }
02219                   }
02220                   if (npiv != 0) {
02221                      apos1 = posfac;
02222                      kk = kk+1;
02223                      apos2 = apos2-kk;
02224                      for (jj = 1; jj < npiv+1; jj++) {
02225                         kk = kk+1;
02226                         apos1 = apos1-kk;
02227                         apos2 = apos2-kk;
02228                         swop = a[apos2];
02229                         a[apos2] = a[apos1];
02230                         a[apos1] = swop;
02231                      }
02232                   }
02233                   swop = a[apos];
02234                   a[apos] = a[posfac];
02235                   a[posfac] = swop;
02236                   ipos = iwpos+npiv+2;
02237                   iexch = iwpos+irow+npiv+1;
02238                   iswop = iw[ipos];
02239                   iw[ipos] = iw[iexch];
02240                   iw[iexch] = iswop;
02241                }
02242                if (pivsiz == 1) continue;
02243                if (krow != 2)  {
02244                   irow = jmax-(npiv+1);
02245                   jpos = posfac;
02246                   posfac = posfac+nfront-npiv;
02247                   npiv = npiv+1;
02248                   apos = apos3;
02249                } else {
02250                   npiv = npiv-1;
02251                   posfac = jpos;
02252                }
02253             }
02254 
02255             if (pivsiz != 2) {
02256 hack:
02257                a[posfac] = one/a[posfac];
02258                if (a[posfac] < zero) neig = neig+1;
02259                j1 = posfac+1;
02260                j2 = posfac+nfront-(npiv+1);
02261                if (j2 >= j1) {
02262                   ibeg = j2+1;
02263                   for (jj = j1; jj < j2+1; jj++) {
02264                      amult = -a[jj]*a[posfac];
02265                      iend = ibeg+nfront-(npiv+jj-j1+2);
02266                      for (irow = ibeg; irow < iend+1; irow++) {
02267                         jcol = jj+irow-ibeg;
02268                         a[irow] = a[irow]+amult*a[jcol];
02269                      }
02270                      ibeg = iend+1;
02271                      a[jj] = amult;
02272                   }
02273                }
02274                npiv = npiv+1;
02275                ntotpv = ntotpv+1;
02276                jpiv = 1;
02277                posfac = posfac+nfront-npiv+1;
02278             } else {
02279                ipos = iwpos+npiv+2;
02280                ntwo = ntwo+1;
02281                iw[ipos] = -iw[ipos];
02282                pospv1 = posfac;
02283                pospv2 = posfac+nfront-npiv;
02284                swop = a[pospv2];
02285                if (detpiv < zero) neig = neig+1;
02286                if (detpiv > zero && swop < zero) neig = neig+2;
02287                a[pospv2] = a[pospv1]/detpiv;
02288                a[pospv1] = swop/detpiv;
02289                a[pospv1+1] = -a[pospv1+1]/detpiv;
02290                j1 = pospv1+2;
02291                j2 = pospv1+nfront-(npiv+1);
02292                if (j2 >= j1) {
02293                   jj1 = pospv2;
02294                   ibeg = pospv2+nfront-(npiv+1);
02295                   for (jj = j1; jj < j2+1; jj++) {
02296                      jj1 = jj1+1;
02297                      amult1 =-(a[pospv1]*a[jj]+a[pospv1+1]*a[jj1]);
02298                      amult2 =-(a[pospv1+1]*a[jj]+a[pospv2]*a[jj1]);
02299                      iend = ibeg+nfront-(npiv+jj-j1+3);
02300                      for (irow = ibeg; irow < iend+1; irow++) {
02301                         k1 = jj+irow-ibeg;
02302                         k2 = jj1+irow-ibeg;
02303                         a[irow] = a[irow]+amult1*a[k1]+amult2*a[k2];
02304                      }
02305                      ibeg = iend+1;
02306                      a[jj] = amult1;
02307                      a[jj1] = amult2;
02308                   }
02309                }
02310                npiv = npiv+2;
02311                ntotpv = ntotpv+2;
02312                jpiv = 2;
02313                posfac = pospv2+nfront-npiv+1;
02314             }
02315          }
02316       }
02317 
02318       if (npiv != 0) nblk = nblk+1;
02319       ioldps = iwpos;
02320       iwpos = iwpos+nfront+2;
02321       if (npiv != 0) {
02322          if (npiv <= 1) {
02323             iw[ioldps] = -iw[ioldps];
02324             for (k = 1; k < nfront+1; k++) {
02325                j1 = ioldps+k;
02326                iw[j1] = iw[j1+1];
02327             }
02328             iwpos = iwpos-1;
02329          } else {
02330             iw[ioldps+1] = npiv;
02331          }
02332       }
02333       liell = nfront-npiv;
02334 
02335       if (liell != 0 && iass != nsteps) {
02336          if (iwpos+liell >= istk)
02337             Factor_sub3(a,iw,istk,istk2,iinput,2,ncmpbr,ncmpbi);
02338          istk = istk-liell-1;
02339          iw[istk] = liell;
02340          j1 = istk;
02341          kk = iwpos-liell-1;
02342          for (k = 1; k < liell+1; k++) {
02343             j1 = j1+1;
02344             kk = kk+1;
02345             iw[j1] = iw[kk];
02346          }
02347          laell = ((liell+1)*liell)/2;
02348          kk = posfac+laell;
02349          if (kk == astk) {
02350             astk = astk-laell;
02351          } else {
02352             kmax = kk-1;
02353             for (k = 1; k < laell+1; k++) {
02354                kk = kk-1;
02355                astk = astk-1;
02356                a[astk] = a[kk];
02357             }
02358             kmax = TMath::Min(kmax,astk-1);
02359             for ( k = kk; k < kmax+1; k++)
02360                a[k] = zero;
02361          }
02362          azero = TMath::Min(azero,astk-1);
02363       }
02364       if (npiv == 0) iwpos = ioldps;
02365    }
02366 
02367    iw[1] = nblk;
02368    if (ntwo > 0) iw[1] = -nblk;
02369    nrlbdu = posfac-1;
02370    nirbdu = iwpos-1;
02371 
02372    if (ntotpv != n) {
02373       info[1] = 3;
02374       info[2] = ntotpv;
02375    }
02376 
02377 finish:
02378    info[9]  = nrlbdu;
02379    info[10] = nirbdu;
02380    info[12] = ncmpbr;
02381    info[13] = ncmpbi;
02382    info[14] = ntwo;
02383    info[15] = neig;
02384 }
02385 
02386 //______________________________________________________________________________
02387 void TDecompSparse::Factor_sub3(Double_t *a,Int_t *iw,Int_t &j1,Int_t &j2,const Int_t itop,
02388                                 const Int_t ireal,Int_t &ncmpbr,Int_t &ncmpbi)
02389 {
02390 // Help routine for factorization
02391 
02392    Int_t ipos,jj,jjj;
02393 
02394    ipos = itop-1;
02395    if (j2 != ipos) {
02396       if (ireal != 2) {
02397          ncmpbr = ncmpbr+1;
02398          if (j1 <= j2) {
02399             for (jjj = j1; jjj < j2+1; jjj++) {
02400                jj = j2-jjj+j1;
02401                a[ipos] = a[jj];
02402                ipos = ipos-1;
02403             }
02404          }
02405       } else {
02406          ncmpbi = ncmpbi+1;
02407          if (j1 <= j2) {
02408             for (jjj = j1; jjj < j2+1; jjj++) {
02409                jj = j2-jjj+j1;
02410                iw[ipos] = iw[jj];
02411                ipos = ipos-1;
02412             }
02413          }
02414       }
02415       j2 = itop-1;
02416       j1 = ipos+1;
02417    }
02418 }
02419 
02420 //______________________________________________________________________________
02421 void TDecompSparse::Solve_sub1(const Int_t n,Double_t *a,Int_t *iw,Double_t *w,
02422                                Double_t *rhs,Int_t *iw2,const Int_t nblk,Int_t &latop,
02423                                Int_t *icntl)
02424 {
02425 // Help routine for solving
02426 
02427    Int_t apos,iblk,ifr,ilvl,ipiv,ipos,irhs,irow,ist,j,j1=0,j2,j3,jj,jpiv,k,k1,k2,k3,liell,npiv;
02428    Double_t w1,w2;
02429 
02430    const Int_t ifrlvl = 5;
02431 
02432    apos = 1;
02433    ipos = 1;
02434    j2 = 0;
02435    iblk = 0;
02436    npiv = 0;
02437    for (irow = 1; irow < n+1; irow++) {
02438       if (npiv <= 0) {
02439          iblk = iblk+1;
02440          if (iblk > nblk) break;
02441          ipos = j2+1;
02442          iw2[iblk] = ipos;
02443          liell = -iw[ipos];
02444          npiv = 1;
02445          if (liell <= 0)  {
02446             liell = -liell;
02447             ipos = ipos+1;
02448             npiv = iw[ipos];
02449          }
02450          j1 = ipos+1;
02451          j2 = ipos+liell;
02452          ilvl = TMath::Min(npiv,10);
02453          if (liell < icntl[ifrlvl+ilvl]) goto hack;
02454          ifr = 0;
02455          for (jj = j1; jj < j2+1; jj++) {
02456             j = TMath::Abs(iw[jj]+0);
02457             ifr = ifr+1;
02458             w[ifr] = rhs[j];
02459          }
02460          jpiv = 1;
02461          j3 = j1;
02462 
02463          for (ipiv = 1; ipiv < npiv+1; ipiv++) {
02464             jpiv = jpiv-1;
02465             if (jpiv == 1) continue;
02466 
02467             if (iw[j3] >= 0) {
02468                jpiv = 1;
02469                j3 = j3+1;
02470                apos = apos+1;
02471                ist = ipiv+1;
02472                if (liell< ist) continue;
02473                w1 = w[ipiv];
02474                k = apos;
02475                for (j = ist; j < liell+1; j++) {
02476                   w[j] = w[j]+a[k]*w1;
02477                   k = k+1;
02478                }
02479                apos = apos+liell-ist+1;
02480             } else {
02481                jpiv = 2;
02482                j3 = j3+2;
02483                apos = apos+2;
02484                ist = ipiv+2;
02485                if (liell >= ist) {
02486                   w1 = w[ipiv];
02487                   w2 = w[ipiv+1];
02488                   k1 = apos;
02489                   k2 = apos+liell-ipiv;
02490                   for (j = ist; j < liell+1; j++) {
02491                      w[j] = w[j]+w1*a[k1]+w2*a[k2];
02492                      k1 = k1+1;
02493                      k2 = k2+1;
02494                   }
02495                }
02496                apos = apos+2*(liell-ist+1)+1;
02497             }
02498          }
02499 
02500          ifr = 0;
02501          for (jj = j1; jj < j2+1; jj++) {
02502             j = TMath::Abs(iw[jj]+0);
02503             ifr = ifr+1;
02504             rhs[j] = w[ifr];
02505          }
02506          npiv = 0;
02507       } else {
02508 hack:
02509          if (iw[j1] >= 0) {
02510             npiv = npiv-1;
02511             apos = apos+1;
02512             j1 = j1+1;
02513             if (j1 <= j2) {
02514                irhs = iw[j1-1];
02515                w1 = rhs[irhs];
02516                k = apos;
02517                for (j = j1; j < j2+1; j++) {
02518                   irhs = TMath::Abs(iw[j]+0);
02519                   rhs[irhs] = rhs[irhs]+a[k]*w1;
02520                   k = k+1;
02521                }
02522             }
02523             apos = apos+j2-j1+1;
02524          } else {
02525             npiv = npiv-2;
02526             j1 = j1+2;
02527             apos = apos+2;
02528             if (j1 <= j2) {
02529                irhs = -iw[j1-2];
02530                w1 = rhs[irhs];
02531                irhs = iw[j1-1];
02532                w2 = rhs[irhs];
02533                k1 = apos;
02534                k3 = apos+j2-j1+2;
02535                for (j = j1; j < j2+1; j++) {
02536                   irhs = TMath::Abs(iw[j]+0);
02537                   rhs[irhs] = rhs[irhs]+w1*a[k1]+w2*a[k3];
02538                   k1 = k1+1;
02539                   k3 = k3+1;
02540                }
02541             }
02542             apos = apos+2*(j2-j1+1)+1;
02543          }
02544       }
02545    }
02546 
02547    latop = apos-1;
02548 }
02549 
02550 //______________________________________________________________________________
02551 void TDecompSparse::Solve_sub2(const Int_t n,Double_t *a,Int_t *iw,Double_t *w,
02552                                Double_t *rhs,Int_t *iw2,const Int_t nblk,
02553                                const Int_t latop,Int_t *icntl)
02554 {
02555 // Help routine for solving
02556 
02557    Int_t apos,apos2,i1rhs,i2rhs,iblk,ifr,iipiv,iirhs,ilvl,ipiv,ipos,irhs,ist,
02558          j,j1=0,j2=0,jj,jj1,jj2,jpiv,jpos=0,k,liell,loop,npiv;
02559    Double_t w1,w2;
02560 
02561    const Int_t ifrlvl = 5;
02562 
02563    apos = latop+1;
02564    npiv = 0;
02565    iblk = nblk+1;
02566    for (loop = 1; loop < n+1; loop++) {
02567       if (npiv <= 0) {
02568          iblk = iblk-1;
02569          if (iblk < 1) break;
02570          ipos = iw2[iblk];
02571          liell = -iw[ipos];
02572          npiv = 1;
02573          if (liell <= 0) {
02574             liell = -liell;
02575             ipos = ipos+1;
02576             npiv = iw[ipos];
02577          }
02578          jpos = ipos+npiv;
02579          j2 = ipos+liell;
02580          ilvl = TMath::Min(10,npiv)+10;
02581          if (liell < icntl[ifrlvl+ilvl]) goto hack;
02582          j1 = ipos+1;
02583          ifr = 0;
02584          for (jj = j1; jj < j2+1; jj++) {
02585             j = TMath::Abs(iw[jj]+0);
02586             ifr = ifr+1;
02587             w[ifr] = rhs[j];
02588          }
02589          jpiv = 1;
02590          for (iipiv = 1; iipiv < npiv+1; iipiv++) {
02591             jpiv = jpiv-1;
02592             if (jpiv == 1) continue;
02593             ipiv = npiv-iipiv+1;
02594             if (ipiv == 1 || iw[jpos-1] >= 0) {
02595                jpiv = 1;
02596                apos = apos-(liell+1-ipiv);
02597                ist = ipiv+1;
02598                w1 = w[ipiv]*a[apos];
02599                if (liell >= ist) {
02600                   jj1 = apos+1;
02601                   for (j = ist; j < liell+1; j++) {
02602                      w1 = w1+a[jj1]*w[j];
02603                      jj1 = jj1+1;
02604                   }
02605                }
02606                w[ipiv] = w1;
02607                jpos = jpos-1;
02608             } else {
02609                jpiv = 2;
02610                apos2 = apos-(liell+1-ipiv);
02611                apos = apos2-(liell+2-ipiv);
02612                ist = ipiv+1;
02613                w1 = w[ipiv-1]*a[apos]+w[ipiv]*a[apos+1];
02614                w2 = w[ipiv-1]*a[apos+1]+w[ipiv]*a[apos2];
02615                if (liell >= ist) {
02616                   jj1 = apos+2;
02617                   jj2 = apos2+1;
02618                   for (j = ist; j < liell+1; j++) {
02619                      w1 = w1+w[j]*a[jj1];
02620                      w2 = w2+w[j]*a[jj2];
02621                      jj1 = jj1+1;
02622                      jj2 = jj2+1;
02623                   }
02624                }
02625                w[ipiv-1] = w1;
02626                w[ipiv] = w2;
02627                jpos = jpos-2;
02628             }
02629          }
02630          ifr = 0;
02631          for (jj = j1; jj < j2+1; jj++) {
02632             j = TMath::Abs(iw[jj]+0);
02633             ifr = ifr+1;
02634             rhs[j] = w[ifr];
02635          }
02636          npiv = 0;
02637       } else {
02638 hack:
02639          if (npiv == 1 || iw[jpos-1] >= 0) {
02640             npiv = npiv-1;
02641             apos = apos-(j2-jpos+1);
02642             iirhs = iw[jpos];
02643             w1 = rhs[iirhs]*a[apos];
02644             j1 = jpos+1;
02645             if (j1 <= j2) {
02646                k = apos+1;
02647                for (j = j1; j < j2+1; j++) {
02648                   irhs = TMath::Abs(iw[j]+0);
02649                   w1 = w1+a[k]*rhs[irhs];
02650                   k = k+1;
02651                }
02652             }
02653             rhs[iirhs] = w1;
02654             jpos = jpos-1;
02655          } else {
02656             npiv = npiv-2;
02657             apos2 = apos-(j2-jpos+1);
02658             apos = apos2-(j2-jpos+2);
02659             i1rhs = -iw[jpos-1];
02660             i2rhs = iw[jpos];
02661             w1 = rhs[i1rhs]*a[apos]+rhs[i2rhs]*a[apos+1];
02662             w2 = rhs[i1rhs]*a[apos+1]+rhs[i2rhs]*a[apos2];
02663             j1 = jpos+1;
02664             if (j1 <= j2) {
02665                jj1 = apos+2;
02666                jj2 = apos2+1;
02667                for (j = j1; j < j2+1; j++) {
02668                   irhs = TMath::Abs(iw[j]+0);
02669                   w1 = w1+rhs[irhs]*a[jj1];
02670                   w2 = w2+rhs[irhs]*a[jj2];
02671                   jj1 = jj1+1;
02672                   jj2 = jj2+1;
02673                }
02674             }
02675             rhs[i1rhs] = w1;
02676             rhs[i2rhs] = w2;
02677             jpos = jpos-2;
02678          }
02679       }
02680    }
02681 }
02682 
02683 //______________________________________________________________________________
02684 void TDecompSparse::Print(Option_t *opt) const
02685 {
02686 // Print class members
02687 
02688    TDecompBase::Print(opt);
02689 
02690    printf("fPrecision  = %.3f\n",fPrecision);
02691    printf("fIPessimism = %.3f\n",fIPessimism);
02692    printf("fRPessimism = %.3f\n",fRPessimism);
02693 
02694    TMatrixDSparse fact(0,fNrows-1,0,fNrows-1,fNnonZeros,
02695                        (Int_t*)fRowFact.GetArray(),(Int_t*)fColFact.GetArray(),(Double_t*)fFact.GetArray());
02696    fact.Print("fFact");
02697 }
02698 
02699 //______________________________________________________________________________
02700 TDecompSparse &TDecompSparse::operator=(const TDecompSparse &source)
02701 {
02702 // Assignment operator
02703 
02704    if (this != &source) {
02705       TDecompBase::operator=(source);
02706       memcpy(fIcntl,source.fIcntl,31*sizeof(Int_t));
02707       memcpy(fCntl,source.fCntl,6*sizeof(Double_t));
02708       memcpy(fInfo,source.fInfo,21*sizeof(Int_t));
02709       fVerbose    = source.fVerbose;
02710       fPrecision  = source.fPrecision;
02711       fIkeep      = source.fIkeep;
02712       fIw         = source.fIw;
02713       fIw1        = source.fIw1;
02714       fIw2        = source.fIw2;
02715       fNsteps     = source.fNsteps;
02716       fMaxfrt     = source.fMaxfrt;
02717       fW          = source.fW;
02718       fIPessimism = source.fIPessimism;
02719       fRPessimism = source.fRPessimism;
02720       if (fA.IsValid())
02721          fA.Use(*const_cast<TMatrixDSparse *>(&(source.fA)));
02722       fNrows      = source.fNrows;
02723       fNnonZeros  = source.fNnonZeros;
02724       fFact       = source.fFact;
02725       fRowFact    = source.fRowFact;
02726       fColFact    = source.fColFact;
02727    }
02728    return *this;
02729 }

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