TRobustEstimator.cxx

Go to the documentation of this file.
00001 // @(#)root/physics:$Id: TRobustEstimator.cxx 30749 2009-10-15 16:33:04Z brun $
00002 // Author: Anna Kreshuk  08/10/2004
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //////////////////////////////////////////////////////////////////////////////
00013 //
00014 //  TRobustEstimator
00015 //
00016 // Minimum Covariance Determinant Estimator - a Fast Algorithm
00017 // invented by Peter J.Rousseeuw and Katrien Van Dreissen
00018 // "A Fast Algorithm for the Minimum covariance Determinant Estimator"
00019 // Technometrics, August 1999, Vol.41, NO.3
00020 //
00021 // What are robust estimators?
00022 // "An important property of an estimator is its robustness. An estimator
00023 // is called robust if it is insensitive to measurements that deviate
00024 // from the expected behaviour. There are 2 ways to treat such deviating
00025 // measurements: one may either try to recongize them and then remove
00026 // them from the data sample; or one may leave them in the sample, taking
00027 // care that they do not influence the estimate unduly. In both cases robust
00028 // estimators are needed...Robust procedures compensate for systematic errors
00029 // as much as possible, and indicate any situation in which a danger of not being
00030 // able to operate reliably is detected."
00031 // R.Fruhwirth, M.Regler, R.K.Bock, H.Grote, D.Notz
00032 // "Data Analysis Techniques for High-Energy Physics", 2nd edition
00033 //
00034 // What does this algorithm do?
00035 // It computes a highly robust estimator of multivariate location and scatter.
00036 // Then, it takes those estimates to compute robust distances of all the
00037 // data vectors. Those with large robust distances are considered outliers.
00038 // Robust distances can then be plotted for better visualization of the data.
00039 //
00040 // How does this algorithm do it?
00041 // The MCD objective is to find h observations(out of n) whose classical
00042 // covariance matrix has the lowest determinant. The MCD estimator of location
00043 // is then the average of those h points and the MCD estimate of scatter
00044 // is their covariance matrix. The minimum(and default) h = (n+nvariables+1)/2
00045 // so the algorithm is effective when less than (n+nvar+1)/2 variables are outliers.
00046 // The algorithm also allows for exact fit situations - that is, when h or more
00047 // observations lie on a hyperplane. Then the algorithm still yields the MCD location T
00048 // and scatter matrix S, the latter being singular as it should be. From (T,S) the
00049 // program then computes the equation of the hyperplane.
00050 //
00051 // How can this algorithm be used?
00052 // In any case, when contamination of data is suspected, that might influence
00053 // the classical estimates.
00054 // Also, robust estimation of location and scatter is a tool to robustify
00055 // other multivariate techniques such as, for example, principal-component analysis
00056 // and discriminant analysis.
00057 //
00058 //
00059 //
00060 //
00061 // Technical details of the algorithm:
00062 // 0.The default h = (n+nvariables+1)/2, but the user may choose any interger h with
00063 //   (n+nvariables+1)/2<=h<=n. The program then reports the MCD's breakdown value
00064 //   (n-h+1)/n. If you are sure that the dataset contains less than 25% contamination
00065 //   which is usually the case, a good compromise between breakdown value and
00066 //  efficiency is obtained by putting h=[.75*n].
00067 // 1.If h=n,the MCD location estimate is the average of the whole dataset, and
00068 //   the MCD scatter estimate is its covariance matrix. Report this and stop
00069 // 2.If nvariables=1 (univariate data), compute the MCD estimate by the exact
00070 //   algorithm of Rousseeuw and Leroy (1987, pp.171-172) in O(nlogn)time and stop
00071 // 3.From here on, h<n and nvariables>=2.
00072 //   3a.If n is small:
00073 //    - repeat (say) 500 times:
00074 //    -- construct an initial h-subset, starting from a random (nvar+1)-subset
00075 //    -- carry out 2 C-steps (described in the comments of CStep function)
00076 //    - for the 10 results with lowest det(S):
00077 //    -- carry out C-steps until convergence
00078 //    - report the solution (T, S) with the lowest det(S)
00079 //   3b.If n is larger (say, n>600), then
00080 //    - construct up to 5 disjoint random subsets of size nsub (say, nsub=300)
00081 //    - inside each subset repeat 500/5 times:
00082 //    -- construct an initial subset of size hsub=[nsub*h/n]
00083 //    -- carry out 2 C-steps
00084 //    -- keep the best 10 results (Tsub, Ssub)
00085 //    - pool the subsets, yielding the merged set (say, of size nmerged=1500)
00086 //    - in the merged set, repeat for each of the 50 solutions (Tsub, Ssub)
00087 //    -- carry out 2 C-steps
00088 //    -- keep the 10 best results
00089 //    - in the full dataset, repeat for those best results:
00090 //    -- take several C-steps, using n and h
00091 //    -- report the best final result (T, S)
00092 // 4.To obtain consistency when the data comes from a multivariate normal
00093 //   distribution, covariance matrix is multiplied by a correction factor
00094 // 5.Robust distances for all elements, using the final (T, S) are calculated
00095 //   Then the very final mean and covariance estimates are calculated only for
00096 //   values, whose robust distances are less than a cutoff value (0.975 quantile
00097 //   of chi2 distribution with nvariables degrees of freedom)
00098 //
00099 //////////////////////////////////////////////////////////////////////////////
00100 
00101 #include "TRobustEstimator.h"
00102 #include "TRandom.h"
00103 #include "TMath.h"
00104 #include "TH1D.h"
00105 #include "TPaveLabel.h"
00106 #include "TDecompChol.h"
00107 
00108 ClassImp(TRobustEstimator)
00109 
00110 const Double_t kChiMedian[50]= {
00111          0.454937, 1.38629, 2.36597, 3.35670, 4.35146, 5.34812, 6.34581, 7.34412, 8.34283,
00112          9.34182, 10.34, 11.34, 12.34, 13.34, 14.34, 15.34, 16.34, 17.34, 18.34, 19.34,
00113         20.34, 21.34, 22.34, 23.34, 24.34, 25.34, 26.34, 27.34, 28.34, 29.34, 30.34,
00114         31.34, 32.34, 33.34, 34.34, 35.34, 36.34, 37.34, 38.34, 39.34, 40.34,
00115         41.34, 42.34, 43.34, 44.34, 45.34, 46.34, 47.34, 48.34, 49.33};
00116 
00117 const Double_t kChiQuant[50]={
00118          5.02389, 7.3776,9.34840,11.1433,12.8325,
00119         14.4494,16.0128,17.5346,19.0228,20.4831,21.920,23.337,
00120         24.736,26.119,27.488,28.845,30.191,31.526,32.852,34.170,
00121         35.479,36.781,38.076,39.364,40.646,41.923,43.194,44.461,
00122         45.722,46.979,48.232,49.481,50.725,51.966,53.203,54.437,
00123         55.668,56.896,58.120,59.342,60.561,61.777,62.990,64.201,
00124         65.410,66.617,67.821,69.022,70.222,71.420};
00125 
00126 //_____________________________________________________________________________
00127 TRobustEstimator::TRobustEstimator(){
00128   //this constructor should be used in a univariate case:
00129   //first call this constructor, then - the EvaluateUni(..) function
00130 
00131 }
00132 
00133 //______________________________________________________________________________
00134 TRobustEstimator::TRobustEstimator(Int_t nvectors, Int_t nvariables, Int_t hh)
00135    :fMean(nvariables),
00136     fCovariance(nvariables),
00137     fInvcovariance(nvariables),
00138     fCorrelation(nvariables),
00139     fRd(nvectors),
00140     fSd(nvectors),
00141     fOut(1),
00142     fHyperplane(nvariables),
00143     fData(nvectors, nvariables)
00144 {
00145    //constructor
00146 
00147    if ((nvectors<=1)||(nvariables<=0)){
00148       Error("TRobustEstimator","Not enough vectors or variables");
00149       return;
00150    }
00151    if (nvariables==1){
00152       Error("TRobustEstimator","For the univariate case, use the default constructor and EvaluateUni() function");
00153       return;
00154    }
00155    
00156    fN=nvectors;
00157    fNvar=nvariables;
00158    if (hh<(fN+fNvar+1)/2){
00159       if (hh>0)
00160          Warning("TRobustEstimator","chosen h is too small, default h is taken instead");
00161       fH=(fN+fNvar+1)/2;
00162    } else
00163       fH=hh;
00164    
00165    fVarTemp=0;
00166    fVecTemp=0;
00167    fExact=0;
00168 }
00169 
00170 //_____________________________________________________________________________
00171 void TRobustEstimator::AddColumn(Double_t *col)
00172 {
00173    //adds a column to the data matrix
00174    //it is assumed that the column has size fN
00175    //variable fVarTemp keeps the number of columns l
00176    //already added
00177    
00178    if (fVarTemp==fNvar) {
00179       fNvar++;
00180       fCovariance.ResizeTo(fNvar, fNvar);
00181       fInvcovariance.ResizeTo(fNvar, fNvar);
00182       fCorrelation.ResizeTo(fNvar, fNvar);
00183       fMean.ResizeTo(fNvar);
00184       fHyperplane.ResizeTo(fNvar);
00185       fData.ResizeTo(fN, fNvar);
00186    }
00187    for (Int_t i=0; i<fN; i++) {
00188       fData(i, fVarTemp)=col[i];
00189    }
00190    fVarTemp++;
00191 }
00192 
00193 //_______________________________________________________________________________
00194 void TRobustEstimator::AddRow(Double_t *row)
00195 {
00196    //adds a vector to the data matrix
00197    //it is supposed that the vector is of size fNvar
00198 
00199    if(fVecTemp==fN) {
00200       fN++;
00201       fRd.ResizeTo(fN);
00202       fSd.ResizeTo(fN);
00203       fData.ResizeTo(fN, fNvar);
00204    }
00205    for (Int_t i=0; i<fNvar; i++)
00206       fData(fVecTemp, i)=row[i];
00207    
00208    fVecTemp++;
00209 }
00210 
00211 //_______________________________________________________________________________
00212 void TRobustEstimator::Evaluate()
00213 {
00214 //Finds the estimate of multivariate mean and variance
00215 
00216    Double_t kEps=1e-14;
00217    
00218    if (fH==fN){
00219       Warning("Evaluate","Chosen h = #observations, so classic estimates of location and scatter will be calculated");
00220       Classic();
00221       return;
00222    }
00223    
00224    Int_t i, j, k;
00225    Int_t ii, jj;
00226    Int_t nmini = 300;
00227    Int_t k1=500;
00228    Int_t nbest=10;
00229    TMatrixD sscp(fNvar+1, fNvar+1);
00230    TVectorD vec(fNvar);
00231    
00232    Int_t *index = new Int_t[fN];
00233    Double_t *ndist = new Double_t[fN];
00234    Double_t det;
00235    Double_t *deti=new Double_t[nbest];
00236    for (i=0; i<nbest; i++)
00237       deti[i]=1e16;
00238    
00239    for (i=0; i<fN; i++)
00240       fRd(i)=0;
00241    ////////////////////////////
00242    //for small n
00243    ////////////////////////////
00244    if (fN<nmini*2) {
00245       //for storing the best fMeans and covariances
00246       
00247       TMatrixD mstock(nbest, fNvar);
00248       TMatrixD cstock(fNvar, fNvar*nbest);
00249       
00250       for (k=0; k<k1; k++) {
00251          CreateSubset(fN, fH, fNvar, index, fData, sscp, ndist);
00252          //calculate the mean and covariance of the created subset
00253          ClearSscp(sscp);
00254          for (i=0; i<fH; i++) {
00255             for(j=0; j<fNvar; j++)
00256                vec(j)=fData[index[i]][j];
00257             AddToSscp(sscp, vec);
00258          }
00259          Covar(sscp, fMean, fCovariance, fSd, fH);
00260          det = fCovariance.Determinant();
00261          if (det < kEps) {
00262             fExact = Exact(ndist);
00263             delete [] index;
00264             delete [] ndist;
00265             delete [] deti;
00266             return;
00267          }
00268          //make 2 CSteps
00269          det = CStep(fN, fH, index, fData, sscp, ndist);
00270          if (det < kEps) {
00271             fExact = Exact(ndist);
00272             delete [] index;
00273             delete [] ndist;
00274             delete [] deti;
00275             return;
00276          }
00277          det = CStep(fN, fH, index, fData, sscp, ndist);
00278          if (det < kEps) {
00279             fExact = Exact(ndist);
00280             delete [] index;
00281             delete [] ndist;
00282             delete [] deti;
00283             return;
00284          } else {
00285             Int_t maxind=TMath::LocMax(nbest, deti);
00286             if(det<deti[maxind]) {
00287                deti[maxind]=det;
00288                for(ii=0; ii<fNvar; ii++) {
00289                   mstock(maxind, ii)=fMean(ii);
00290                   for(jj=0; jj<fNvar; jj++)
00291                      cstock(ii, jj+maxind*fNvar)=fCovariance(ii, jj);
00292                }
00293             }
00294          }
00295       }
00296       
00297       //now for nbest best results perform CSteps until convergence
00298       
00299       for (i=0; i<nbest; i++) {
00300          for(ii=0; ii<fNvar; ii++) {
00301             fMean(ii)=mstock(i, ii);
00302             for (jj=0; jj<fNvar; jj++)
00303                fCovariance(ii, jj)=cstock(ii, jj+i*fNvar);
00304          }
00305          
00306          det=1;
00307          while (det>kEps) {
00308             det=CStep(fN, fH, index, fData, sscp, ndist);
00309             if(TMath::Abs(det-deti[i])<kEps)
00310                break;
00311             else
00312                deti[i]=det;
00313          }
00314          for(ii=0; ii<fNvar; ii++) {
00315             mstock(i,ii)=fMean(ii);
00316             for (jj=0; jj<fNvar; jj++)
00317                cstock(ii,jj+i*fNvar)=fCovariance(ii, jj);
00318          }
00319       }
00320       
00321       Int_t detind=TMath::LocMin(nbest, deti);
00322       for(ii=0; ii<fNvar; ii++) {
00323          fMean(ii)=mstock(detind,ii);
00324          
00325          for(jj=0; jj<fNvar; jj++)
00326             fCovariance(ii, jj)=cstock(ii,jj+detind*fNvar);
00327       }
00328       
00329       if (deti[detind]!=0) {
00330          //calculate robust distances and throw out the bad points
00331          Int_t nout = RDist(sscp);
00332          Double_t cutoff=kChiQuant[fNvar-1];
00333          
00334          fOut.Set(nout);
00335          
00336          j=0;
00337          for (i=0; i<fN; i++) {
00338             if(fRd(i)>cutoff) {
00339                fOut[j]=i;
00340                j++;
00341             }
00342          }
00343          
00344       } else {
00345          fExact=Exact(ndist);
00346       }
00347       delete [] index;
00348       delete [] ndist;
00349       delete [] deti;
00350       return;
00351       
00352    }
00353    /////////////////////////////////////////////////
00354   //if n>nmini, the dataset should be partitioned
00355   //partitioning
00356   ////////////////////////////////////////////////
00357    Int_t indsubdat[5];
00358    Int_t nsub;
00359    for (ii=0; ii<5; ii++)
00360       indsubdat[ii]=0;
00361    
00362    nsub = Partition(nmini, indsubdat);
00363    
00364    Int_t sum=0;
00365    for (ii=0; ii<5; ii++)
00366       sum+=indsubdat[ii];
00367    Int_t *subdat=new Int_t[sum];
00368    RDraw(subdat, nsub, indsubdat);
00369    //now the indexes of selected cases are in the array subdat
00370    //matrices to store best means and covariances
00371    Int_t nbestsub=nbest*nsub;
00372    TMatrixD mstockbig(nbestsub, fNvar);
00373    TMatrixD cstockbig(fNvar, fNvar*nbestsub);
00374    TMatrixD hyperplane(nbestsub, fNvar);
00375    for (i=0; i<nbestsub; i++) {
00376       for(j=0; j<fNvar; j++)
00377          hyperplane(i,j)=0;
00378    }
00379    Double_t *detibig = new Double_t[nbestsub];
00380    Int_t maxind;
00381    maxind=TMath::LocMax(5, indsubdat);
00382    TMatrixD dattemp(indsubdat[maxind], fNvar);
00383    
00384    Int_t k2=Int_t(k1/nsub);
00385    //construct h-subsets and perform 2 CSteps in subgroups
00386    
00387    for (Int_t kgroup=0; kgroup<nsub; kgroup++) {
00388       //printf("group #%d\n", kgroup);
00389       Int_t ntemp=indsubdat[kgroup];
00390       Int_t temp=0;
00391       for (i=0; i<kgroup; i++)
00392          temp+=indsubdat[i];
00393       Int_t par;
00394       
00395       for(i=0; i<ntemp; i++) {
00396          for (j=0; j<fNvar; j++) {
00397             dattemp(i,j)=fData[subdat[temp+i]][j];
00398          }
00399       }
00400       Int_t htemp=Int_t(fH*ntemp/fN);
00401       
00402       for (i=0; i<nbest; i++)
00403          deti[i]=1e16;
00404       
00405       for(k=0; k<k2; k++) {
00406          CreateSubset(ntemp, htemp, fNvar, index, dattemp, sscp, ndist);
00407          ClearSscp(sscp);
00408          for (i=0; i<htemp; i++) {
00409             for(j=0; j<fNvar; j++) {
00410                vec(j)=dattemp(index[i],j);
00411             }
00412             AddToSscp(sscp, vec);
00413          }
00414          Covar(sscp, fMean, fCovariance, fSd, htemp);
00415          det = fCovariance.Determinant();
00416          if (det<kEps) {
00417             par =Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
00418             if(par==nbest+1) {
00419                
00420                delete [] detibig;
00421                delete [] deti;
00422                delete [] subdat;
00423                delete [] ndist;
00424                delete [] index;
00425                return;
00426             } else
00427                deti[par]=det;
00428          } else {
00429             det = CStep(ntemp, htemp, index, dattemp, sscp, ndist);
00430             if (det<kEps) {
00431                par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp, ndist);
00432                if(par==nbest+1) {
00433                   
00434                   delete [] detibig;
00435                   delete [] deti;
00436                   delete [] subdat;
00437                   delete [] ndist;
00438                   delete [] index;
00439                   return;
00440                } else
00441                   deti[par]=det;
00442             } else {
00443                det=CStep(ntemp,htemp, index, dattemp, sscp, ndist);
00444                if(det<kEps){
00445                   par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
00446                   if(par==nbest+1) {
00447                      
00448                      delete [] detibig;
00449                      delete [] deti;
00450                      delete [] subdat;
00451                      delete [] ndist;
00452                      delete [] index;
00453                      return;
00454                   } else {
00455                      deti[par]=det;
00456                   }
00457                } else {
00458                   maxind=TMath::LocMax(nbest, deti);
00459                   if(det<deti[maxind]) {
00460                      deti[maxind]=det;
00461                      for(i=0; i<fNvar; i++) {
00462                         mstockbig(nbest*kgroup+maxind,i)=fMean(i);
00463                         for(j=0; j<fNvar; j++) {
00464                            cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
00465                            
00466                         }
00467                      }
00468                   }
00469                   
00470                }
00471             }
00472          }
00473          
00474          maxind=TMath::LocMax(nbest, deti);
00475          if (deti[maxind]<kEps)
00476             break;
00477       }
00478       
00479       
00480       for(i=0; i<nbest; i++) {
00481          detibig[kgroup*nbest + i]=deti[i];
00482          
00483       }
00484       
00485    }
00486    
00487    //now the arrays mstockbig and cstockbig store nbest*nsub best means and covariances
00488    //detibig stores nbest*nsub their determinants
00489    //merge the subsets and carry out 2 CSteps on the merged set for all 50 best solutions
00490    
00491    TMatrixD datmerged(sum, fNvar);
00492    for(i=0; i<sum; i++) {
00493       for (j=0; j<fNvar; j++)
00494          datmerged(i,j)=fData[subdat[i]][j];
00495    }
00496    //  printf("performing calculations for merged set\n");
00497    Int_t hmerged=Int_t(sum*fH/fN);
00498    
00499    Int_t nh;
00500    for(k=0; k<nbestsub; k++) {
00501       //for all best solutions perform 2 CSteps and then choose the very best
00502       for(ii=0; ii<fNvar; ii++) {
00503          fMean(ii)=mstockbig(k,ii);
00504          for(jj=0; jj<fNvar; jj++)
00505             fCovariance(ii, jj)=cstockbig(ii,k*fNvar+jj);
00506       }
00507       if(detibig[k]==0) {
00508          for(i=0; i<fNvar; i++)
00509             fHyperplane(i)=hyperplane(k,i);
00510          CreateOrtSubset(datmerged,index, hmerged, sum, sscp, ndist);
00511          
00512       }
00513       det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
00514       if (det<kEps) {
00515          nh= Exact(ndist);
00516          if (nh>=fH) {
00517             fExact = nh;
00518             
00519             delete [] detibig;
00520             delete [] deti;
00521             delete [] subdat;
00522             delete [] ndist;
00523             delete [] index;
00524             return;
00525          } else {
00526             CreateOrtSubset(datmerged, index, hmerged, sum, sscp, ndist);
00527          }
00528       }
00529       
00530       det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
00531       if (det<kEps) {
00532          nh=Exact(ndist);
00533          if (nh>=fH) {
00534             fExact = nh;
00535             delete [] detibig;
00536             delete [] deti;
00537             delete [] subdat;
00538             delete [] ndist;
00539             delete [] index;
00540             return;
00541          }
00542       }
00543       detibig[k]=det;
00544       for(i=0; i<fNvar; i++) {
00545          mstockbig(k,i)=fMean(i);
00546          for(j=0; j<fNvar; j++) {
00547             cstockbig(i,k*fNvar+j)=fCovariance(i, j);
00548          }
00549       }
00550    }
00551    //now for the subset with the smallest determinant
00552    //repeat CSteps until convergence
00553    Int_t minind=TMath::LocMin(nbestsub, detibig);
00554    det=detibig[minind];
00555    for(i=0; i<fNvar; i++) {
00556       fMean(i)=mstockbig(minind,i);
00557       fHyperplane(i)=hyperplane(minind,i);
00558       for(j=0; j<fNvar; j++)
00559          fCovariance(i, j)=cstockbig(i,minind*fNvar + j);
00560    }
00561    if(det<kEps)
00562       CreateOrtSubset(fData, index, fH, fN, sscp, ndist);
00563    det=1;
00564    while (det>kEps) {
00565       det=CStep(fN, fH, index, fData, sscp, ndist);
00566       if(TMath::Abs(det-detibig[minind])<kEps) {
00567          break;
00568       } else {
00569          detibig[minind]=det;
00570       }
00571    }
00572    if(det<kEps) {
00573       Exact(ndist);
00574       fExact=kTRUE;
00575    }
00576    Int_t nout = RDist(sscp);
00577    Double_t cutoff=kChiQuant[fNvar-1];
00578    
00579    fOut.Set(nout);
00580    
00581    j=0;
00582    for (i=0; i<fN; i++) {
00583       if(fRd(i)>cutoff) {
00584          fOut[j]=i;
00585          j++;
00586       }
00587    }
00588    
00589    delete [] detibig;
00590    delete [] deti;
00591    delete [] subdat;
00592    delete [] ndist;
00593    delete [] index;
00594    return;
00595 }
00596 
00597 //___________________________________________________________________________________________
00598 void TRobustEstimator::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh)
00599 {
00600    //for the univariate case
00601    //estimates of location and scatter are returned in mean and sigma parameters
00602    //the algorithm works on the same principle as in multivariate case -
00603    //it finds a subset of size hh with smallest sigma, and then returns mean and
00604    //sigma of this subset
00605    
00606    if (hh==0)
00607       hh=(nvectors+2)/2;
00608    Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
00609    Int_t *index=new Int_t[nvectors];
00610    TMath::Sort(nvectors, data, index, kFALSE);
00611    
00612    Int_t nquant;
00613    nquant=TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
00614    Double_t factor=faclts[nquant-1];
00615    
00616    Double_t *aw=new Double_t[nvectors];
00617    Double_t *aw2=new Double_t[nvectors];
00618    Double_t sq=0;
00619    Double_t sqmin=0;
00620    Int_t ndup=0;
00621    Int_t len=nvectors-hh;
00622    Double_t *slutn=new Double_t[len];
00623    for(Int_t i=0; i<len; i++)
00624       slutn[i]=0;
00625    for(Int_t jint=0; jint<len; jint++) {
00626       aw[jint]=0;
00627       for (Int_t j=0; j<hh; j++) {
00628          aw[jint]+=data[index[j+jint]];
00629          if(jint==0)
00630             sq+=data[index[j]]*data[index[j]];
00631       }
00632       aw2[jint]=aw[jint]*aw[jint]/hh;
00633       
00634       if(jint==0) {
00635          sq=sq-aw2[jint];
00636          sqmin=sq;
00637          slutn[ndup]=aw[jint];
00638          
00639       } else {
00640          sq=sq - data[index[jint-1]]*data[index[jint-1]]+
00641             data[index[jint+hh]]*data[index[jint+hh]]-
00642             aw2[jint]+aw2[jint-1];
00643          if(sq<sqmin) {
00644             ndup=0;
00645             sqmin=sq;
00646             slutn[ndup]=aw[jint];
00647             
00648          } else {
00649             if(sq==sqmin) {
00650                ndup++;
00651                slutn[ndup]=aw[jint];
00652             }
00653          }
00654       }
00655    }
00656    
00657    slutn[0]=slutn[Int_t((ndup)/2)]/hh;
00658    Double_t bstd=factor*TMath::Sqrt(sqmin/hh);
00659    mean=slutn[0];
00660    sigma=bstd;
00661    delete [] aw;
00662    delete [] aw2;
00663    delete [] slutn;
00664    delete [] index;
00665 }
00666 
00667 //_______________________________________________________________________
00668 Int_t TRobustEstimator::GetBDPoint()
00669 {
00670   //returns the breakdown point of the algorithm
00671 
00672    Int_t n;
00673    n=(fN-fH+1)/fN;
00674    return n;
00675 }
00676 
00677 //_______________________________________________________________________
00678 Double_t TRobustEstimator::GetChiQuant(Int_t i) const
00679 {
00680    //returns the chi2 quantiles
00681 
00682    if (i < 0 || i >= 50) return 0;
00683    return kChiQuant[i];
00684 }
00685 
00686 //_______________________________________________________________________
00687 void TRobustEstimator::GetCovariance(TMatrixDSym &matr)
00688 {
00689    //returns the covariance matrix
00690 
00691    if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar){
00692       Warning("GetCovariance","provided matrix is of the wrong size, it will be resized");
00693       matr.ResizeTo(fNvar, fNvar);
00694    }
00695    matr=fCovariance;
00696 }
00697 
00698 //_______________________________________________________________________
00699 void TRobustEstimator::GetCorrelation(TMatrixDSym &matr)
00700 {
00701    //returns the correlation matrix
00702    
00703    if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar) {
00704       Warning("GetCorrelation","provided matrix is of the wrong size, it will be resized");
00705       matr.ResizeTo(fNvar, fNvar);
00706    }
00707    matr=fCorrelation;
00708 }
00709 
00710 //____________________________________________________________________
00711 const TVectorD* TRobustEstimator::GetHyperplane() const
00712 {
00713    //if the points are on a hyperplane, returns this hyperplane
00714 
00715    if (fExact==0) {
00716       Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
00717       return 0;
00718    } else {
00719       return &fHyperplane;
00720    }
00721 }
00722 
00723 //______________________________________________________________________
00724 void TRobustEstimator::GetHyperplane(TVectorD &vec)
00725 {
00726    //if the points are on a hyperplane, returns this hyperplane
00727 
00728    if (fExact==0){
00729       Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
00730       return;
00731    }
00732    if (vec.GetNoElements()!=fNvar) {
00733       Warning("GetHyperPlane","provided vector is of the wrong size, it will be resized");
00734       vec.ResizeTo(fNvar);
00735    }
00736    vec=fHyperplane;
00737 }
00738 
00739 //________________________________________________________________________
00740 void TRobustEstimator::GetMean(TVectorD &means)
00741 {
00742    //return the estimate of the mean
00743 
00744    if (means.GetNoElements()!=fNvar) {
00745       Warning("GetMean","provided vector is of the wrong size, it will be resized");
00746       means.ResizeTo(fNvar);
00747    }
00748    means=fMean;
00749 }
00750 
00751 //_________________________________________________________________________
00752 void TRobustEstimator::GetRDistances(TVectorD &rdist)
00753 {
00754    //returns the robust distances (helps to find outliers)
00755 
00756    if (rdist.GetNoElements()!=fN) {
00757       Warning("GetRDistances","provided vector is of the wrong size, it will be resized");
00758       rdist.ResizeTo(fN);
00759    }
00760    rdist=fRd;
00761 }
00762 
00763 //__________________________________________________________________________
00764 Int_t TRobustEstimator::GetNOut()
00765 {
00766    //returns the number of outliers
00767 
00768    return fOut.GetSize();
00769 }
00770 
00771 //_________________________________________________________________________
00772 void TRobustEstimator::AddToSscp(TMatrixD &sscp, TVectorD &vec)
00773 {
00774   //update the sscp matrix with vector vec
00775 
00776    Int_t i, j;
00777    for (j=1; j<fNvar+1; j++) {
00778       sscp(0, j) +=vec(j-1);
00779       sscp(j, 0) = sscp(0, j);
00780    }
00781    for (i=1; i<fNvar+1; i++) {
00782       for (j=1; j<fNvar+1; j++) {
00783          sscp(i, j) += vec(i-1)*vec(j-1);
00784       }
00785    }
00786 }
00787 
00788 //__________________________________________________________________________
00789 void TRobustEstimator::ClearSscp(TMatrixD &sscp)
00790 {
00791    //clear the sscp matrix, used for covariance and mean calculation
00792 
00793    for (Int_t i=0; i<fNvar+1; i++) {
00794       for (Int_t j=0; j<fNvar+1; j++) {
00795          sscp(i, j)=0;
00796       }
00797    }
00798 }
00799 
00800 //_______________________________________________________________
00801 void TRobustEstimator::Classic()
00802 {
00803   //called when h=n. Returns classic covariance matrix
00804   //and mean
00805    TMatrixD sscp(fNvar+1, fNvar+1);
00806    TVectorD temp(fNvar);
00807    ClearSscp(sscp);
00808    for (Int_t i=0; i<fN; i++) {
00809       for (Int_t j=0; j<fNvar; j++)
00810          temp(j)=fData(i, j);
00811       AddToSscp(sscp, temp);
00812    }
00813    Covar(sscp, fMean, fCovariance, fSd, fN);
00814    Correl();
00815    
00816 }
00817 
00818 //____________________________________________________________________
00819 void TRobustEstimator::Covar(TMatrixD &sscp, TVectorD &m, TMatrixDSym &cov, TVectorD &sd, Int_t nvec)
00820 {
00821    //calculates mean and covariance
00822 
00823    Int_t i, j;
00824    Double_t f;
00825    for (i=0; i<fNvar; i++) {
00826       m(i)=sscp(0, i+1);
00827       sd[i]=sscp(i+1, i+1);
00828       f=(sd[i]-m(i)*m(i)/nvec)/(nvec-1);
00829       if (f>1e-14) sd[i]=TMath::Sqrt(f);
00830       else sd[i]=0;
00831       m(i)/=nvec;
00832    }
00833    for (i=0; i<fNvar; i++) {
00834       for (j=0; j<fNvar; j++) {
00835          cov(i, j)=sscp(i+1, j+1)-nvec*m(i)*m(j);
00836       cov(i, j)/=nvec-1;
00837       }
00838    }
00839 }
00840 
00841 //____________________________________________________________________
00842 void TRobustEstimator::Correl()
00843 {
00844   //transforms covariance matrix into correlation matrix
00845 
00846    Int_t i, j;
00847    Double_t *sd=new Double_t[fNvar];
00848    for(j=0; j<fNvar; j++)
00849       sd[j]=1./TMath::Sqrt(fCovariance(j, j));
00850    for(i=0; i<fNvar; i++) {
00851       for (j=0; j<fNvar; j++) {
00852          if (i==j)
00853             fCorrelation(i, j)=1.;
00854          else
00855             fCorrelation(i, j)=fCovariance(i, j)*sd[i]*sd[j];
00856       }
00857    }
00858    delete [] sd;
00859 }
00860 
00861 //____________________________________________________________________
00862 void TRobustEstimator::CreateSubset(Int_t ntotal, Int_t htotal, Int_t p, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
00863 {
00864   //creates a subset of htotal elements from ntotal elements
00865   //first, p+1 elements are drawn randomly(without repetitions)
00866   //if their covariance matrix is singular, more elements are
00867   //added one by one, until their covariance matrix becomes regular
00868   //or it becomes clear that htotal observations lie on a hyperplane
00869   //If covariance matrix determinant!=0, distances of all ntotal elements
00870   //are calculated, using formula d_i=Sqrt((x_i-M)*S_inv*(x_i-M)), where
00871   //M is mean and S_inv is the inverse of the covariance matrix
00872   //htotal points with smallest distances are included in the returned subset.
00873 
00874    Double_t kEps = 1e-14;
00875    Int_t i, j;
00876    Bool_t repeat=kFALSE;
00877    Int_t nindex=0;
00878    Int_t num;
00879    for(i=0; i<ntotal; i++)
00880       index[i]=ntotal+1;
00881    
00882    for (i=0; i<p+1; i++) {
00883       num=Int_t(gRandom->Uniform(0, 1)*(ntotal-1));
00884       if (i>0){
00885          for(j=0; j<=i-1; j++) {
00886             if(index[j]==num)
00887             repeat=kTRUE;
00888          }
00889       }
00890       if(repeat==kTRUE) {
00891          i--;
00892          repeat=kFALSE;
00893       } else {
00894          index[i]=num;
00895          nindex++;
00896       }
00897    }
00898    
00899    ClearSscp(sscp);
00900    
00901    TVectorD vec(fNvar);
00902    Double_t det;
00903    for (i=0; i<p+1; i++) {
00904       for (j=0; j<fNvar; j++) {
00905          vec[j]=data[index[i]][j];
00906          
00907       }
00908       AddToSscp(sscp, vec);
00909    }
00910    
00911    Covar(sscp, fMean, fCovariance, fSd, p+1);
00912    det=fCovariance.Determinant();
00913    while((det<kEps)&&(nindex < htotal)) {
00914     //if covariance matrix is singular,another vector is added until
00915     //the matrix becomes regular or it becomes clear that all
00916     //vectors of the group lie on a hyperplane
00917       repeat=kFALSE;
00918       do{
00919          num=Int_t(gRandom->Uniform(0,1)*(ntotal-1));
00920          repeat=kFALSE;
00921          for(i=0; i<nindex; i++) {
00922             if(index[i]==num) {
00923                repeat=kTRUE;
00924                break;
00925             }
00926          }
00927       }while(repeat==kTRUE);
00928       
00929       index[nindex]=num;
00930       nindex++;
00931       //check if covariance matrix is singular
00932       for(j=0; j<fNvar; j++)
00933          vec[j]=data[index[nindex-1]][j];
00934       AddToSscp(sscp, vec);
00935       Covar(sscp, fMean, fCovariance, fSd, nindex);
00936       det=fCovariance.Determinant();
00937    }
00938    
00939    if(nindex!=htotal) {
00940       TDecompChol chol(fCovariance);
00941       fInvcovariance = chol.Invert();
00942       
00943       TVectorD temp(fNvar);
00944       for(j=0; j<ntotal; j++) {
00945          ndist[j]=0;
00946          for(i=0; i<fNvar; i++)
00947             temp[i]=data[j][i] - fMean(i);
00948          temp*=fInvcovariance;
00949          for(i=0; i<fNvar; i++)
00950             ndist[j]+=(data[j][i]-fMean(i))*temp[i];
00951       }
00952       KOrdStat(ntotal, ndist, htotal-1,index);
00953    }
00954    
00955 }
00956 
00957 //_____________________________________________________________________________
00958 void TRobustEstimator::CreateOrtSubset(TMatrixD &dat,Int_t *index, Int_t hmerged, Int_t nmerged, TMatrixD &sscp, Double_t *ndist)
00959 {
00960   //creates a subset of hmerged vectors with smallest orthogonal distances to the hyperplane
00961   //hyp[1]*(x1-mean[1])+...+hyp[nvar]*(xnvar-mean[nvar])=0
00962   //This function is called in case when less than fH samples lie on a hyperplane.
00963 
00964    Int_t i, j;
00965       TVectorD vec(fNvar);
00966    for (i=0; i<nmerged; i++) {
00967       ndist[i]=0;
00968       for(j=0; j<fNvar; j++) {
00969          ndist[i]+=fHyperplane[j]*(dat[i][j]-fMean[j]);
00970          ndist[i]=TMath::Abs(ndist[i]);
00971       }
00972    }
00973    KOrdStat(nmerged, ndist, hmerged-1, index);
00974    ClearSscp(sscp);
00975    for (i=0; i<hmerged; i++) {
00976       for(j=0; j<fNvar; j++)
00977          vec[j]=dat[index[i]][j];
00978       AddToSscp(sscp, vec);
00979    }
00980    Covar(sscp, fMean, fCovariance, fSd, hmerged);
00981 }
00982 
00983 //__________________________________________________________________________
00984 Double_t TRobustEstimator::CStep(Int_t ntotal, Int_t htotal, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
00985 {
00986   //from the input htotal-subset constructs another htotal subset with lower determinant
00987   //
00988   //As proven by Peter J.Rousseeuw and Katrien Van Driessen, if distances for all elements
00989   //are calculated, using the formula:d_i=Sqrt((x_i-M)*S_inv*(x_i-M)), where M is the mean
00990   //of the input htotal-subset, and S_inv - the inverse of its covariance matrix, then
00991   //htotal elements with smallest distances will have covariance matrix with determinant
00992   //less or equal to the determinant of the input subset covariance matrix.
00993   //
00994   //determinant for this htotal-subset with smallest distances is returned
00995 
00996    Int_t i, j;
00997    TVectorD vec(fNvar);
00998    Double_t det;
00999    
01000    TDecompChol chol(fCovariance);
01001    fInvcovariance = chol.Invert();
01002 
01003    TVectorD temp(fNvar);
01004    for(j=0; j<ntotal; j++) {
01005       ndist[j]=0;
01006       for(i=0; i<fNvar; i++)
01007          temp[i]=data[j][i]-fMean[i];
01008       temp*=fInvcovariance;
01009       for(i=0; i<fNvar; i++)
01010          ndist[j]+=(data[j][i]-fMean[i])*temp[i];
01011    }
01012    
01013    //taking h smallest
01014    KOrdStat(ntotal, ndist, htotal-1, index);
01015    //writing their mean and covariance
01016    ClearSscp(sscp);
01017    for (i=0; i<htotal; i++) {
01018       for (j=0; j<fNvar; j++)
01019          temp[j]=data[index[i]][j];
01020       AddToSscp(sscp, temp);
01021    }
01022    Covar(sscp, fMean, fCovariance, fSd, htotal);
01023    det = fCovariance.Determinant();
01024    return det;
01025 }
01026 
01027 //_______________________________________________________________
01028 Int_t TRobustEstimator::Exact(Double_t *ndist)
01029 {
01030   //for the exact fit situaions
01031   //returns number of observations on the hyperplane
01032 
01033    Int_t i, j;
01034    
01035    TMatrixDSymEigen eigen(fCovariance);
01036    TVectorD eigenValues=eigen.GetEigenValues();
01037    TMatrixD eigenMatrix=eigen.GetEigenVectors();
01038    
01039    for (j=0; j<fNvar; j++) {
01040       fHyperplane[j]=eigenMatrix(j,fNvar-1);
01041    }
01042    //calculate and return how many observations lie on the hyperplane
01043    for (i=0; i<fN; i++) {
01044       ndist[i]=0;
01045       for(j=0; j<fNvar; j++) {
01046          ndist[i]+=fHyperplane[j]*(fData[i][j]-fMean[j]);
01047          ndist[i]=TMath::Abs(ndist[i]);
01048       }
01049    }
01050    Int_t nhyp=0;
01051    
01052    for (i=0; i<fN; i++) {
01053       if(ndist[i] < 1e-14) nhyp++;
01054    }
01055    return nhyp;
01056    
01057 }
01058 
01059 //____________________________________________________________________________
01060 Int_t TRobustEstimator::Exact2(TMatrixD &mstockbig, TMatrixD &cstockbig, TMatrixD &hyperplane,
01061                              Double_t *deti, Int_t nbest, Int_t kgroup,
01062                              TMatrixD &sscp, Double_t *ndist)
01063 {
01064   //This function is called if determinant of the covariance matrix of a subset=0.
01065   //
01066   //If there are more then fH vectors on a hyperplane,
01067   //returns this hyperplane and stops
01068   //else stores the hyperplane coordinates in hyperplane matrix
01069 
01070    Int_t i, j;
01071    
01072    TVectorD vec(fNvar);
01073    Int_t maxind = TMath::LocMax(nbest, deti);
01074    Int_t nh=Exact(ndist);
01075    //now nh is the number of observation on the hyperplane
01076    //ndist stores distances of observation from this hyperplane
01077    if(nh>=fH) {
01078       ClearSscp(sscp);
01079       for (i=0; i<fN; i++) {
01080          if(ndist[i]<1e-14) {
01081             for (j=0; j<fNvar; j++)
01082                vec[j]=fData[i][j];
01083             AddToSscp(sscp, vec);
01084          }
01085       }
01086       Covar(sscp, fMean, fCovariance, fSd, nh);
01087       
01088       fExact=nh;
01089       return nbest+1;
01090       
01091    } else {
01092       //if less than fH observations lie on a hyperplane,
01093       //mean and covariance matrix are stored in mstockbig
01094       //and cstockbig in place of the previous maximum determinant
01095       //mean and covariance
01096       for(i=0; i<fNvar; i++) {
01097          mstockbig(nbest*kgroup+maxind,i)=fMean(i);
01098          hyperplane(nbest*kgroup+maxind,i)=fHyperplane(i);
01099          for(j=0; j<fNvar; j++) {
01100             cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
01101          }   
01102       }      
01103       return maxind;
01104    }
01105 }
01106 
01107 
01108 //_____________________________________________________________________________
01109 Int_t TRobustEstimator::Partition(Int_t nmini, Int_t *indsubdat)
01110 {
01111   //divides the elements into approximately equal subgroups
01112   //number of elements in each subgroup is stored in indsubdat
01113   //number of subgroups is returned
01114 
01115    Int_t nsub;
01116    if ((fN>=2*nmini) && (fN<=(3*nmini-1))) {
01117       if (fN%2==1){
01118          indsubdat[0]=Int_t(fN*0.5);
01119       indsubdat[1]=Int_t(fN*0.5)+1;
01120       } else
01121          indsubdat[0]=indsubdat[1]=Int_t(fN/2);
01122       nsub=2;
01123    }
01124    else{
01125       if((fN>=3*nmini) && (fN<(4*nmini -1))) {
01126          if(fN%3==0){
01127             indsubdat[0]=indsubdat[1]=indsubdat[2]=Int_t(fN/3);
01128          } else {
01129             indsubdat[0]=Int_t(fN/3);
01130             indsubdat[1]=Int_t(fN/3)+1;
01131             if (fN%3==1) indsubdat[2]=Int_t(fN/3);
01132             else indsubdat[2]=Int_t(fN/3)+1;
01133          }
01134          nsub=3;
01135       }
01136       else{
01137          if((fN>=4*nmini)&&(fN<=(5*nmini-1))){
01138             if (fN%4==0) indsubdat[0]=indsubdat[1]=indsubdat[2]=indsubdat[3]=Int_t(fN/4);
01139             else {
01140                indsubdat[0]=Int_t(fN/4);
01141                indsubdat[1]=Int_t(fN/4)+1;
01142                if(fN%4==1) indsubdat[2]=indsubdat[3]=Int_t(fN/4);
01143                if(fN%4==2) {
01144                   indsubdat[2]=Int_t(fN/4)+1;
01145                   indsubdat[3]=Int_t(fN/4);
01146                }
01147                if(fN%4==3) indsubdat[2]=indsubdat[3]=Int_t(fN/4)+1;
01148             }
01149             nsub=4;
01150          } else {
01151             for(Int_t i=0; i<5; i++)
01152                indsubdat[i]=nmini;
01153             nsub=5;
01154          }
01155       }
01156    }
01157    return nsub;
01158 }
01159 
01160 //___________________________________________________________________________
01161 Int_t TRobustEstimator::RDist(TMatrixD &sscp)
01162 {
01163   //Calculates robust distances.Then the samples with robust distances
01164   //greater than a cutoff value (0.975 quantile of chi2 distribution with
01165   //fNvar degrees of freedom, multiplied by a correction factor), are given
01166   //weiht=0, and new, reweighted estimates of location and scatter are calculated
01167   //The function returns the number of outliers.
01168 
01169    Int_t i, j;
01170    Int_t nout=0;
01171    
01172    TVectorD temp(fNvar);
01173    TDecompChol chol(fCovariance);
01174    fInvcovariance = chol.Invert();
01175    
01176 
01177    for (i=0; i<fN; i++) {
01178       fRd[i]=0;
01179       for(j=0; j<fNvar; j++) {
01180          temp[j]=fData[i][j]-fMean[j];
01181       }
01182       temp*=fInvcovariance;
01183       for(j=0; j<fNvar; j++) {
01184          fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
01185       }
01186    }
01187    
01188    Double_t med;
01189    Double_t chi = kChiMedian[fNvar-1];
01190    
01191    med=TMath::Median(fN, fRd.GetMatrixArray());
01192    med/=chi;
01193    fCovariance*=med;
01194    TDecompChol chol2(fCovariance);
01195    fInvcovariance = chol2.Invert();
01196    
01197    for (i=0; i<fN; i++) {
01198       fRd[i]=0;
01199       for(j=0; j<fNvar; j++) {
01200          temp[j]=fData[i][j]-fMean[j];
01201    }
01202       
01203       temp*=fInvcovariance;
01204       for(j=0; j<fNvar; j++) {
01205          fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
01206       }
01207    }
01208    
01209    Double_t cutoff = kChiQuant[fNvar-1];
01210    
01211    ClearSscp(sscp);
01212    for(i=0; i<fN; i++) {
01213       if (fRd[i]<=cutoff) {
01214          for(j=0; j<fNvar; j++)
01215             temp[j]=fData[i][j];
01216          AddToSscp(sscp,temp);
01217       } else {
01218          nout++;
01219       }
01220    }
01221    
01222    Covar(sscp, fMean, fCovariance, fSd, fN-nout);
01223    return nout;
01224 }
01225 
01226 //____________________________________________________________________________
01227 void TRobustEstimator::RDraw(Int_t *subdat, Int_t ngroup, Int_t *indsubdat)
01228 {
01229   //Draws ngroup nonoverlapping subdatasets out of a dataset of size n
01230   //such that the selected case numbers are uniformly distributed from 1 to n
01231 
01232    Int_t jndex = 0;
01233    Int_t nrand;
01234    Int_t i, k, m, j;
01235    for (k=1; k<=ngroup; k++) {
01236       for (m=1; m<=indsubdat[k-1]; m++) {
01237 
01238          nrand = Int_t(gRandom->Uniform(0, 1) * (fN-jndex))+1;
01239          
01240          jndex++;
01241          if (jndex==1) {
01242             subdat[0]=nrand;
01243          } else {
01244             subdat[jndex-1]=nrand+jndex-2;
01245             for (i=1; i<=jndex-1; i++) {
01246                if(subdat[i-1] > nrand+i-2) {
01247                   for(j=jndex; j>=i+1; j--) {
01248                      subdat[j-1]=subdat[j-2];
01249                   }
01250                   subdat[i-1]=nrand+i-2;
01251                   break;  //breaking the loop for(i=1...
01252                }
01253             }
01254          }
01255       }
01256    }
01257 }
01258 
01259 //_____________________________________________________________________________
01260 Double_t TRobustEstimator::KOrdStat(Int_t ntotal, Double_t *a, Int_t k, Int_t *work){
01261   //because I need an Int_t work array
01262    Bool_t isAllocated = kFALSE;
01263    const Int_t kWorkMax=100;
01264    Int_t i, ir, j, l, mid;
01265    Int_t arr;
01266    Int_t *ind;
01267    Int_t workLocal[kWorkMax];
01268    Int_t temp;
01269 
01270 
01271    if (work) {
01272       ind = work;
01273    } else {
01274       ind = workLocal;
01275       if (ntotal > kWorkMax) {
01276          isAllocated = kTRUE;
01277          ind = new Int_t[ntotal];
01278       }
01279    }
01280 
01281    for (Int_t ii=0; ii<ntotal; ii++) {
01282       ind[ii]=ii;
01283    }
01284    Int_t rk = k;
01285    l=0;
01286    ir = ntotal-1;
01287    for(;;) {
01288       if (ir<=l+1) { //active partition contains 1 or 2 elements
01289          if (ir == l+1 && a[ind[ir]]<a[ind[l]])
01290             {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
01291          Double_t tmp = a[ind[rk]];
01292          if (isAllocated)
01293             delete [] ind;
01294          return tmp;
01295       } else {
01296          mid = (l+ir) >> 1; //choose median of left, center and right
01297          {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}//elements as partitioning element arr.
01298          if (a[ind[l]]>a[ind[ir]])  //also rearrange so that a[l]<=a[l+1]
01299             {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
01300 
01301          if (a[ind[l+1]]>a[ind[ir]])
01302             {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}
01303 
01304          if (a[ind[l]]>a[ind[l+1]])
01305                 {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}
01306 
01307          i=l+1;        //initialize pointers for partitioning
01308          j=ir;
01309          arr = ind[l+1];
01310          for (;;) {
01311             do i++; while (a[ind[i]]<a[arr]);
01312             do j--; while (a[ind[j]]>a[arr]);
01313             if (j<i) break;  //pointers crossed, partitioning complete
01314                {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
01315          }
01316          ind[l+1]=ind[j];
01317          ind[j]=arr;
01318          if (j>=rk) ir = j-1; //keep active the partition that
01319          if (j<=rk) l=i;      //contains the k_th element
01320       }
01321    }
01322 }

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