ROOT logo
// @(#)root/physics:$Id: TRobustEstimator.cxx 26746 2008-12-09 08:58:27Z brun $
// Author: Anna Kreshuk  08/10/2004

/*************************************************************************
 * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//////////////////////////////////////////////////////////////////////////////
//
//  TRobustEstimator
//
// Minimum Covariance Determinant Estimator - a Fast Algorithm
// invented by Peter J.Rousseeuw and Katrien Van Dreissen
// "A Fast Algorithm for the Minimum covariance Determinant Estimator"
// Technometrics, August 1999, Vol.41, NO.3
//
// What are robust estimators?
// "An important property of an estimator is its robustness. An estimator
// is called robust if it is insensitive to measurements that deviate
// from the expected behaviour. There are 2 ways to treat such deviating
// measurements: one may either try to recongize them and then remove
// them from the data sample; or one may leave them in the sample, taking
// care that they do not influence the estimate unduly. In both cases robust
// estimators are needed...Robust procedures compensate for systematic errors
// as much as possible, and indicate any situation in which a danger of not being
// able to operate reliably is detected."
// R.Fruhwirth, M.Regler, R.K.Bock, H.Grote, D.Notz
// "Data Analysis Techniques for High-Energy Physics", 2nd edition
//
// What does this algorithm do?
// It computes a highly robust estimator of multivariate location and scatter.
// Then, it takes those estimates to compute robust distances of all the
// data vectors. Those with large robust distances are considered outliers.
// Robust distances can then be plotted for better visualization of the data.
//
// How does this algorithm do it?
// The MCD objective is to find h observations(out of n) whose classical
// covariance matrix has the lowest determinant. The MCD estimator of location
// is then the average of those h points and the MCD estimate of scatter
// is their covariance matrix. The minimum(and default) h = (n+nvariables+1)/2
// so the algorithm is effective when less than (n+nvar+1)/2 variables are outliers.
// The algorithm also allows for exact fit situations - that is, when h or more
// observations lie on a hyperplane. Then the algorithm still yields the MCD location T
// and scatter matrix S, the latter being singular as it should be. From (T,S) the
// program then computes the equation of the hyperplane.
//
// How can this algorithm be used?
// In any case, when contamination of data is suspected, that might influence
// the classical estimates.
// Also, robust estimation of location and scatter is a tool to robustify
// other multivariate techniques such as, for example, principal-component analysis
// and discriminant analysis.
//
//
//
//
// Technical details of the algorithm:
// 0.The default h = (n+nvariables+1)/2, but the user may choose any interger h with
//   (n+nvariables+1)/2<=h<=n. The program then reports the MCD's breakdown value
//   (n-h+1)/n. If you are sure that the dataset contains less than 25% contamination
//   which is usually the case, a good compromise between breakdown value and
//  efficiency is obtained by putting h=[.75*n].
// 1.If h=n,the MCD location estimate is the average of the whole dataset, and
//   the MCD scatter estimate is its covariance matrix. Report this and stop
// 2.If nvariables=1 (univariate data), compute the MCD estimate by the exact
//   algorithm of Rousseeuw and Leroy (1987, pp.171-172) in O(nlogn)time and stop
// 3.From here on, h<n and nvariables>=2.
//   3a.If n is small:
//    - repeat (say) 500 times:
//    -- construct an initial h-subset, starting from a random (nvar+1)-subset
//    -- carry out 2 C-steps (described in the comments of CStep funtion)
//    - for the 10 results with lowest det(S):
//    -- carry out C-steps until convergence
//    - report the solution (T, S) with the lowest det(S)
//   3b.If n is larger (say, n>600), then
//    - construct up to 5 disjoint random subsets of size nsub (say, nsub=300)
//    - inside each subset repeat 500/5 times:
//    -- construct an initial subset of size hsub=[nsub*h/n]
//    -- carry out 2 C-steps
//    -- keep the best 10 results (Tsub, Ssub)
//    - pool the subsets, yielding the merged set (say, of size nmerged=1500)
//    - in the merged set, repeat for each of the 50 solutions (Tsub, Ssub)
//    -- carry out 2 C-steps
//    -- keep the 10 best results
//    - in the full dataset, repeat for those best results:
//    -- take several C-steps, using n and h
//    -- report the best final result (T, S)
// 4.To obtain consistency when the data comes from a multivariate normal
//   distribution, covariance matrix is multiplied by a correction factor
// 5.Robust distances for all elements, using the final (T, S) are calculated
//   Then the very final mean and covariance estimates are calculated only for
//   values, whose robust distances are less than a cutoff value (0.975 quantile
//   of chi2 distribution with nvariables degrees of freedom)
//
//////////////////////////////////////////////////////////////////////////////

#include "TRobustEstimator.h"
#include "TRandom.h"
#include "TMath.h"
#include "TH1D.h"
#include "TPaveLabel.h"
#include "TDecompChol.h"

ClassImp(TRobustEstimator)

const Double_t kChiMedian[50]= {
         0.454937, 1.38629, 2.36597, 3.35670, 4.35146, 5.34812, 6.34581, 7.34412, 8.34283,
         9.34182, 10.34, 11.34, 12.34, 13.34, 14.34, 15.34, 16.34, 17.34, 18.34, 19.34,
        20.34, 21.34, 22.34, 23.34, 24.34, 25.34, 26.34, 27.34, 28.34, 29.34, 30.34,
        31.34, 32.34, 33.34, 34.34, 35.34, 36.34, 37.34, 38.34, 39.34, 40.34,
        41.34, 42.34, 43.34, 44.34, 45.34, 46.34, 47.34, 48.34, 49.33};

const Double_t kChiQuant[50]={
         5.02389, 7.3776,9.34840,11.1433,12.8325,
        14.4494,16.0128,17.5346,19.0228,20.4831,21.920,23.337,
        24.736,26.119,27.488,28.845,30.191,31.526,32.852,34.170,
        35.479,36.781,38.076,39.364,40.646,41.923,43.194,44.461,
        45.722,46.979,48.232,49.481,50.725,51.966,53.203,54.437,
        55.668,56.896,58.120,59.342,60.561,61.777,62.990,64.201,
        65.410,66.617,67.821,69.022,70.222,71.420};

//_____________________________________________________________________________
TRobustEstimator::TRobustEstimator(){
  //this constructor should be used in a univariate case:
  //first call this constructor, then - the EvaluateUni(..) fucntion

}

//______________________________________________________________________________
TRobustEstimator::TRobustEstimator(Int_t nvectors, Int_t nvariables, Int_t hh)
   :fMean(nvariables),
    fCovariance(nvariables),
    fInvcovariance(nvariables),
    fCorrelation(nvariables),
    fRd(nvectors),
    fSd(nvectors),
    fOut(1),
    fHyperplane(nvariables),
    fData(nvectors, nvariables)
{
   //constructor

   if ((nvectors<=1)||(nvariables<=0)){
      Error("TRobustEstimator","Not enough vectors or variables");
      return;
   }
   if (nvariables==1){
      Error("TRobustEstimator","For the univariate case, use the default constructor and EvaluateUni() function");
      return;
   }
   
   fN=nvectors;
   fNvar=nvariables;
   if (hh<(fN+fNvar+1)/2){
      if (hh>0)
         Warning("TRobustEstimator","chosen h is too small, default h is taken instead");
      fH=(fN+fNvar+1)/2;
   } else
      fH=hh;
   
   fVarTemp=0;
   fVecTemp=0;
   fExact=0;
}

//_____________________________________________________________________________
void TRobustEstimator::AddColumn(Double_t *col)
{
   //adds a column to the data matrix
   //it is assumed that the column has size fN
   //variable fVarTemp keeps the number of columns l
   //already added
   
   if (fVarTemp==fNvar) {
      fNvar++;
      fCovariance.ResizeTo(fNvar, fNvar);
      fInvcovariance.ResizeTo(fNvar, fNvar);
      fCorrelation.ResizeTo(fNvar, fNvar);
      fMean.ResizeTo(fNvar);
      fHyperplane.ResizeTo(fNvar);
      fData.ResizeTo(fN, fNvar);
   }
   for (Int_t i=0; i<fN; i++) {
      fData(i, fVarTemp)=col[i];
   }
   fVarTemp++;
}

//_______________________________________________________________________________
void TRobustEstimator::AddRow(Double_t *row)
{
   //adds a vector to the data matrix
   //it is supposed that the vector is of size fNvar

   if(fVecTemp==fN) {
      fN++;
      fRd.ResizeTo(fN);
      fSd.ResizeTo(fN);
      fData.ResizeTo(fN, fNvar);
   }
   for (Int_t i=0; i<fNvar; i++)
      fData(fVecTemp, i)=row[i];
   
   fVecTemp++;
}

//_______________________________________________________________________________
void TRobustEstimator::Evaluate()
{
//Finds the estimate of multivariate mean and variance

   Double_t kEps=1e-14;
   
   if (fH==fN){
      Warning("Evaluate","Chosen h = #observations, so classic estimates of location and scatter will be calculated");
      Classic();
      return;
   }
   
   Int_t i, j, k;
   Int_t ii, jj;
   Int_t nmini = 300;
   Int_t k1=500;
   Int_t nbest=10;
   TMatrixD sscp(fNvar+1, fNvar+1);
   TVectorD vec(fNvar);
   
   Int_t *index = new Int_t[fN];
   Double_t *ndist = new Double_t[fN];
   Double_t det;
   Double_t *deti=new Double_t[nbest];
   for (i=0; i<nbest; i++)
      deti[i]=1e16;
   
   for (i=0; i<fN; i++)
      fRd(i)=0;
   ////////////////////////////
   //for small n
   ////////////////////////////
   if (fN<nmini*2) {
      //for storing the best fMeans and covariances
      
      TMatrixD mstock(nbest, fNvar);
      TMatrixD cstock(fNvar, fNvar*nbest);
      
      for (k=0; k<k1; k++) {
         CreateSubset(fN, fH, fNvar, index, fData, sscp, ndist);
         //calculate the mean and covariance of the created subset
         ClearSscp(sscp);
         for (i=0; i<fH; i++) {
            for(j=0; j<fNvar; j++)
               vec(j)=fData[index[i]][j];
            AddToSscp(sscp, vec);
         }
         Covar(sscp, fMean, fCovariance, fSd, fH);
         det = fCovariance.Determinant();
         if (det < kEps) {
            fExact = Exact(ndist);
            delete [] index;
            delete [] ndist;
            delete [] deti;
            return;
         }
         //make 2 CSteps
         det = CStep(fN, fH, index, fData, sscp, ndist);
         if (det < kEps) {
            fExact = Exact(ndist);
            delete [] index;
            delete [] ndist;
            delete [] deti;
            return;
         }
         det = CStep(fN, fH, index, fData, sscp, ndist);
         if (det < kEps) {
            fExact = Exact(ndist);
            delete [] index;
            delete [] ndist;
            delete [] deti;
            return;
         } else {
            Int_t maxind=TMath::LocMax(nbest, deti);
            if(det<deti[maxind]) {
               deti[maxind]=det;
               for(ii=0; ii<fNvar; ii++) {
                  mstock(maxind, ii)=fMean(ii);
                  for(jj=0; jj<fNvar; jj++)
                     cstock(ii, jj+maxind*fNvar)=fCovariance(ii, jj);
               }
            }
         }
      }
      
      //now for nbest best results perform CSteps until convergence
      
      for (i=0; i<nbest; i++) {
         for(ii=0; ii<fNvar; ii++) {
            fMean(ii)=mstock(i, ii);
            for (jj=0; jj<fNvar; jj++)
               fCovariance(ii, jj)=cstock(ii, jj+i*fNvar);
         }
         
         det=1;
         while (det>kEps) {
            det=CStep(fN, fH, index, fData, sscp, ndist);
            if(TMath::Abs(det-deti[i])<kEps)
               break;
            else
               deti[i]=det;
         }
         for(ii=0; ii<fNvar; ii++) {
            mstock(i,ii)=fMean(ii);
            for (jj=0; jj<fNvar; jj++)
               cstock(ii,jj+i*fNvar)=fCovariance(ii, jj);
         }
      }
      
      Int_t detind=TMath::LocMin(nbest, deti);
      for(ii=0; ii<fNvar; ii++) {
         fMean(ii)=mstock(detind,ii);
         
         for(jj=0; jj<fNvar; jj++)
            fCovariance(ii, jj)=cstock(ii,jj+detind*fNvar);
      }
      
      if (deti[detind]!=0) {
         //calculate robust distances and throw out the bad points
         Int_t nout = RDist(sscp);
         Double_t cutoff=kChiQuant[fNvar-1];
         
         fOut.Set(nout);
         
         j=0;
         for (i=0; i<fN; i++) {
            if(fRd(i)>cutoff) {
               fOut[j]=i;
               j++;
            }
         }
         
      } else {
         fExact=Exact(ndist);
      }
      delete [] index;
      delete [] ndist;
      delete [] deti;
      return;
      
   }
   /////////////////////////////////////////////////
  //if n>nmini, the dataset should be partitioned
  //partitioning
  ////////////////////////////////////////////////
   Int_t indsubdat[5];
   Int_t nsub;
   for (ii=0; ii<5; ii++)
      indsubdat[ii]=0;
   
   nsub = Partition(nmini, indsubdat);
   
   Int_t sum=0;
   for (ii=0; ii<5; ii++)
      sum+=indsubdat[ii];
   Int_t *subdat=new Int_t[sum];
   RDraw(subdat, nsub, indsubdat);
   //now the indexes of selected cases are in the array subdat
   //matrices to store best means and covariances
   Int_t nbestsub=nbest*nsub;
   TMatrixD mstockbig(nbestsub, fNvar);
   TMatrixD cstockbig(fNvar, fNvar*nbestsub);
   TMatrixD hyperplane(nbestsub, fNvar);
   for (i=0; i<nbestsub; i++) {
      for(j=0; j<fNvar; j++)
         hyperplane(i,j)=0;
   }
   Double_t *detibig = new Double_t[nbestsub];
   Int_t maxind;
   maxind=TMath::LocMax(5, indsubdat);
   TMatrixD dattemp(indsubdat[maxind], fNvar);
   
   Int_t k2=Int_t(k1/nsub);
   //construct h-subsets and perform 2 CSteps in subgroups
   
   for (Int_t kgroup=0; kgroup<nsub; kgroup++) {
      //printf("group #%d\n", kgroup);
      Int_t ntemp=indsubdat[kgroup];
      Int_t temp=0;
      for (i=0; i<kgroup; i++)
         temp+=indsubdat[i];
      Int_t par;
      
      for(i=0; i<ntemp; i++) {
         for (j=0; j<fNvar; j++) {
            dattemp(i,j)=fData[subdat[temp+i]][j];
         }
      }
      Int_t htemp=Int_t(fH*ntemp/fN);
      
      for (i=0; i<nbest; i++)
         deti[i]=1e16;
      
      for(k=0; k<k2; k++) {
         CreateSubset(ntemp, htemp, fNvar, index, dattemp, sscp, ndist);
         ClearSscp(sscp);
         for (i=0; i<htemp; i++) {
            for(j=0; j<fNvar; j++) {
               vec(j)=dattemp(index[i],j);
            }
            AddToSscp(sscp, vec);
         }
         Covar(sscp, fMean, fCovariance, fSd, htemp);
         det = fCovariance.Determinant();
         if (det<kEps) {
            par =Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
            if(par==nbest+1) {
               
               delete [] detibig;
               delete [] deti;
               delete [] subdat;
               delete [] ndist;
               delete [] index;
               return;
            } else
               deti[par]=det;
         } else {
            det = CStep(ntemp, htemp, index, dattemp, sscp, ndist);
            if (det<kEps) {
               par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp, ndist);
               if(par==nbest+1) {
                  
                  delete [] detibig;
                  delete [] deti;
                  delete [] subdat;
                  delete [] ndist;
                  delete [] index;
                  return;
               } else
                  deti[par]=det;
            } else {
               det=CStep(ntemp,htemp, index, dattemp, sscp, ndist);
               if(det<kEps){
                  par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
                  if(par==nbest+1) {
                     
                     delete [] detibig;
                     delete [] deti;
                     delete [] subdat;
                     delete [] ndist;
                     delete [] index;
                     return;
                  } else {
                     deti[par]=det;
                  }
               } else {
                  maxind=TMath::LocMax(nbest, deti);
                  if(det<deti[maxind]) {
                     deti[maxind]=det;
                     for(i=0; i<fNvar; i++) {
                        mstockbig(nbest*kgroup+maxind,i)=fMean(i);
                        for(j=0; j<fNvar; j++) {
                           cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
                           
                        }
                     }
                  }
                  
               }
            }
         }
         
         maxind=TMath::LocMax(nbest, deti);
         if (deti[maxind]<kEps)
            break;
      }
      
      
      for(i=0; i<nbest; i++) {
         detibig[kgroup*nbest + i]=deti[i];
         
      }
      
   }
   
   //now the arrays mstockbig and cstockbig store nbest*nsub best means and covariances
   //detibig stores nbest*nsub their determinants
   //merge the subsets and carry out 2 CSteps on the merged set for all 50 best solutions
   
   TMatrixD datmerged(sum, fNvar);
   for(i=0; i<sum; i++) {
      for (j=0; j<fNvar; j++)
         datmerged(i,j)=fData[subdat[i]][j];
   }
   //  printf("performing calculations for merged set\n");
   Int_t hmerged=Int_t(sum*fH/fN);
   
   Int_t nh;
   for(k=0; k<nbestsub; k++) {
      //for all best solutions perform 2 CSteps and then choose the very best
      for(ii=0; ii<fNvar; ii++) {
         fMean(ii)=mstockbig(k,ii);
         for(jj=0; jj<fNvar; jj++)
            fCovariance(ii, jj)=cstockbig(ii,k*fNvar+jj);
      }
      if(detibig[k]==0) {
         for(i=0; i<fNvar; i++)
            fHyperplane(i)=hyperplane(k,i);
         CreateOrtSubset(datmerged,index, hmerged, sum, sscp, ndist);
         
      }
      det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
      if (det<kEps) {
         nh= Exact(ndist);
         if (nh>=fH) {
            fExact = nh;
            
            delete [] detibig;
            delete [] deti;
            delete [] subdat;
            delete [] ndist;
            delete [] index;
            return;
         } else {
            CreateOrtSubset(datmerged, index, hmerged, sum, sscp, ndist);
         }
      }
      
      det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
      if (det<kEps) {
         nh=Exact(ndist);
         if (nh>=fH) {
            fExact = nh;
            delete [] detibig;
            delete [] deti;
            delete [] subdat;
            delete [] ndist;
            delete [] index;
            return;
         }
      }
      detibig[k]=det;
      for(i=0; i<fNvar; i++) {
         mstockbig(k,i)=fMean(i);
         for(j=0; j<fNvar; j++) {
            cstockbig(i,k*fNvar+j)=fCovariance(i, j);
         }
      }
   }
   //now for the subset with the smallest determinant
   //repeat CSteps until convergence
   Int_t minind=TMath::LocMin(nbestsub, detibig);
   det=detibig[minind];
   for(i=0; i<fNvar; i++) {
      fMean(i)=mstockbig(minind,i);
      fHyperplane(i)=hyperplane(minind,i);
      for(j=0; j<fNvar; j++)
         fCovariance(i, j)=cstockbig(i,minind*fNvar + j);
   }
   if(det<kEps)
      CreateOrtSubset(fData, index, fH, fN, sscp, ndist);
   det=1;
   while (det>kEps) {
      det=CStep(fN, fH, index, fData, sscp, ndist);
      if(TMath::Abs(det-detibig[minind])<kEps) {
         break;
      } else {
         detibig[minind]=det;
      }
   }
   if(det<kEps) {
      Exact(ndist);
      fExact=kTRUE;
   }
   Int_t nout = RDist(sscp);
   Double_t cutoff=kChiQuant[fNvar-1];
   
   fOut.Set(nout);
   
   j=0;
   for (i=0; i<fN; i++) {
      if(fRd(i)>cutoff) {
         fOut[j]=i;
         j++;
      }
   }
   
   delete [] detibig;
   delete [] deti;
   delete [] subdat;
   delete [] ndist;
   delete [] index;
   return;
}

//___________________________________________________________________________________________
void TRobustEstimator::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh)
{
   //for the univariate case
   //estimates of location and scatter are returned in mean and sigma parameters
   //the algorithm works on the same principle as in multivariate case -
   //it finds a subset of size hh with smallest sigma, and then returns mean and
   //sigma of this subset
   
   if (hh==0)
      hh=(nvectors+2)/2;
   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};
   Int_t *index=new Int_t[nvectors];
   TMath::Sort(nvectors, data, index, kFALSE);
   
   Int_t nquant;
   nquant=TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
   Double_t factor=faclts[nquant-1];
   
   Double_t *aw=new Double_t[nvectors];
   Double_t *aw2=new Double_t[nvectors];
   Double_t sq=0;
   Double_t sqmin=0;
   Int_t ndup=0;
   Int_t len=nvectors-hh;
   Double_t *slutn=new Double_t[len];
   for(Int_t i=0; i<len; i++)
      slutn[i]=0;
   for(Int_t jint=0; jint<len; jint++) {
      aw[jint]=0;
      for (Int_t j=0; j<hh; j++) {
         aw[jint]+=data[index[j+jint]];
         if(jint==0)
            sq+=data[index[j]]*data[index[j]];
      }
      aw2[jint]=aw[jint]*aw[jint]/hh;
      
      if(jint==0) {
         sq=sq-aw2[jint];
         sqmin=sq;
         slutn[ndup]=aw[jint];
         
      } else {
         sq=sq - data[index[jint-1]]*data[index[jint-1]]+
            data[index[jint+hh]]*data[index[jint+hh]]-
            aw2[jint]+aw2[jint-1];
         if(sq<sqmin) {
            ndup=0;
            sqmin=sq;
            slutn[ndup]=aw[jint];
            
         } else {
            if(sq==sqmin) {
               ndup++;
               slutn[ndup]=aw[jint];
            }
         }
      }
   }
   
   slutn[0]=slutn[Int_t((ndup)/2)]/hh;
   Double_t bstd=factor*TMath::Sqrt(sqmin/hh);
   mean=slutn[0];
   sigma=bstd;
   delete [] aw;
   delete [] aw2;
   delete [] slutn;
   delete [] index;
}

//_______________________________________________________________________
Int_t TRobustEstimator::GetBDPoint()
{
  //returns the breakdown point of the algorithm

   Int_t n;
   n=(fN-fH+1)/fN;
   return n;
}

//_______________________________________________________________________
Double_t TRobustEstimator::GetChiQuant(Int_t i) const
{
   //returns the chi2 quantiles

   if (i < 0 || i >= 50) return 0;
   return kChiQuant[i];
}

//_______________________________________________________________________
void TRobustEstimator::GetCovariance(TMatrixDSym &matr)
{
   //returns the covariance matrix

   if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar){
      Warning("GetCovariance","provided matrix is of the wrong size, it will be resized");
      matr.ResizeTo(fNvar, fNvar);
   }
   matr=fCovariance;
}

//_______________________________________________________________________
void TRobustEstimator::GetCorrelation(TMatrixDSym &matr)
{
   //returns the correlation matrix
   
   if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar) {
      Warning("GetCorrelation","provided matrix is of the wrong size, it will be resized");
      matr.ResizeTo(fNvar, fNvar);
   }
   matr=fCorrelation;
}

//____________________________________________________________________
const TVectorD* TRobustEstimator::GetHyperplane() const
{
   //if the points are on a hyperplane, returns this hyperplane

   if (fExact==0) {
      Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
      return 0;
   } else {
      return &fHyperplane;
   }
}

//______________________________________________________________________
void TRobustEstimator::GetHyperplane(TVectorD &vec)
{
   //if the points are on a hyperplane, returns this hyperplane

   if (fExact==0){
      Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
      return;
   }
   if (vec.GetNoElements()!=fNvar) {
      Warning("GetHyperPlane","provided vector is of the wrong size, it will be resized");
      vec.ResizeTo(fNvar);
   }
   vec=fHyperplane;
}

//________________________________________________________________________
void TRobustEstimator::GetMean(TVectorD &means)
{
   //return the estimate of the mean

   if (means.GetNoElements()!=fNvar) {
      Warning("GetMean","provided vector is of the wrong size, it will be resized");
      means.ResizeTo(fNvar);
   }
   means=fMean;
}

//_________________________________________________________________________
void TRobustEstimator::GetRDistances(TVectorD &rdist)
{
   //returns the robust distances (helps to find outliers)

   if (rdist.GetNoElements()!=fN) {
      Warning("GetRDistances","provided vector is of the wrong size, it will be resized");
      rdist.ResizeTo(fN);
   }
   rdist=fRd;
}

//__________________________________________________________________________
Int_t TRobustEstimator::GetNOut()
{
   //returns the number of outliers

   return fOut.GetSize();
}

//_________________________________________________________________________
void TRobustEstimator::AddToSscp(TMatrixD &sscp, TVectorD &vec)
{
  //update the sscp matrix with vector vec

   Int_t i, j;
   for (j=1; j<fNvar+1; j++) {
      sscp(0, j) +=vec(j-1);
      sscp(j, 0) = sscp(0, j);
   }
   for (i=1; i<fNvar+1; i++) {
      for (j=1; j<fNvar+1; j++) {
         sscp(i, j) += vec(i-1)*vec(j-1);
      }
   }
}

//__________________________________________________________________________
void TRobustEstimator::ClearSscp(TMatrixD &sscp)
{
   //clear the sscp matrix, used for covariance and mean calculation

   for (Int_t i=0; i<fNvar+1; i++) {
      for (Int_t j=0; j<fNvar+1; j++) {
         sscp(i, j)=0;
      }
   }
}

//_______________________________________________________________
void TRobustEstimator::Classic()
{
  //called when h=n. Returns classic covariance matrix
  //and mean
   TMatrixD sscp(fNvar+1, fNvar+1);
   TVectorD temp(fNvar);
   ClearSscp(sscp);
   for (Int_t i=0; i<fN; i++) {
      for (Int_t j=0; j<fNvar; j++)
         temp(j)=fData(i, j);
      AddToSscp(sscp, temp);
   }
   Covar(sscp, fMean, fCovariance, fSd, fN);
   Correl();
   
}

//____________________________________________________________________
void TRobustEstimator::Covar(TMatrixD &sscp, TVectorD &m, TMatrixDSym &cov, TVectorD &sd, Int_t nvec)
{
   //calculates mean and covariance

   Int_t i, j;
   Double_t f;
   for (i=0; i<fNvar; i++) {
      m(i)=sscp(0, i+1);
      sd[i]=sscp(i+1, i+1);
      f=(sd[i]-m(i)*m(i)/nvec)/(nvec-1);
      if (f>1e-14) sd[i]=TMath::Sqrt(f);
      else sd[i]=0;
      m(i)/=nvec;
   }
   for (i=0; i<fNvar; i++) {
      for (j=0; j<fNvar; j++) {
         cov(i, j)=sscp(i+1, j+1)-nvec*m(i)*m(j);
      cov(i, j)/=nvec-1;
      }
   }
}

//____________________________________________________________________
void TRobustEstimator::Correl()
{
  //transforms covariance matrix into correlation matrix

   Int_t i, j;
   Double_t *sd=new Double_t[fNvar];
   for(j=0; j<fNvar; j++)
      sd[j]=1./TMath::Sqrt(fCovariance(j, j));
   for(i=0; i<fNvar; i++) {
      for (j=0; j<fNvar; j++) {
         if (i==j)
            fCorrelation(i, j)=1.;
         else
            fCorrelation(i, j)=fCovariance(i, j)*sd[i]*sd[j];
      }
   }
   delete [] sd;
}

//____________________________________________________________________
void TRobustEstimator::CreateSubset(Int_t ntotal, Int_t htotal, Int_t p, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
{
  //creates a subset of htotal elements from ntotal elements
  //first, p+1 elements are drawn randomly(without repetitions)
  //if their covariance matrix is singular, more elements are
  //added one by one, until their covariance matrix becomes regular
  //or it becomes clear that htotal observations lie on a hyperplane
  //If covariance matrix determinant!=0, distances of all ntotal elements
  //are calculated, using formula d_i=Sqrt((x_i-M)*S_inv*(x_i-M)), where
  //M is mean and S_inv is the inverse of the covariance matrix
  //htotal points with smallest distances are included in the returned subset.

   Double_t kEps = 1e-14;
   Int_t i, j;
   Bool_t repeat=kFALSE;
   Int_t nindex=0;
   Int_t num;
   for(i=0; i<ntotal; i++)
      index[i]=ntotal+1;
   
   for (i=0; i<p+1; i++) {
      num=Int_t(gRandom->Uniform(0, 1)*(ntotal-1));
      if (i>0){
         for(j=0; j<=i-1; j++) {
            if(index[j]==num)
            repeat=kTRUE;
         }
      }
      if(repeat==kTRUE) {
         i--;
         repeat=kFALSE;
      } else {
         index[i]=num;
         nindex++;
      }
   }
   
   ClearSscp(sscp);
   
   TVectorD vec(fNvar);
   Double_t det;
   for (i=0; i<p+1; i++) {
      for (j=0; j<fNvar; j++) {
         vec[j]=data[index[i]][j];
         
      }
      AddToSscp(sscp, vec);
   }
   
   Covar(sscp, fMean, fCovariance, fSd, p+1);
   det=fCovariance.Determinant();
   while((det<kEps)&&(nindex < htotal)) {
    //if covariance matrix is singular,another vector is added until
    //the matrix becomes regular or it becomes clear that all
    //vectors of the group lie on a hyperplane
      repeat=kFALSE;
      do{
         num=Int_t(gRandom->Uniform(0,1)*(ntotal-1));
         repeat=kFALSE;
         for(i=0; i<nindex; i++) {
            if(index[i]==num) {
               repeat=kTRUE;
               break;
            }
         }
      }while(repeat==kTRUE);
      
      index[nindex]=num;
      nindex++;
      //check if covariance matrix is singular
      for(j=0; j<fNvar; j++)
         vec[j]=data[index[nindex-1]][j];
      AddToSscp(sscp, vec);
      Covar(sscp, fMean, fCovariance, fSd, nindex);
      det=fCovariance.Determinant();
   }
   
   if(nindex!=htotal) {
      TDecompChol chol(fCovariance);
      fInvcovariance = chol.Invert();
      
      TVectorD temp(fNvar);
      for(j=0; j<ntotal; j++) {
         ndist[j]=0;
         for(i=0; i<fNvar; i++)
            temp[i]=data[j][i] - fMean(i);
         temp*=fInvcovariance;
         for(i=0; i<fNvar; i++)
            ndist[j]+=(data[j][i]-fMean(i))*temp[i];
      }
      KOrdStat(ntotal, ndist, htotal-1,index);
   }
   
}

//_____________________________________________________________________________
void TRobustEstimator::CreateOrtSubset(TMatrixD &dat,Int_t *index, Int_t hmerged, Int_t nmerged, TMatrixD &sscp, Double_t *ndist)
{
  //creates a subset of hmerged vectors with smallest orthogonal distances to the hyperplane
  //hyp[1]*(x1-mean[1])+...+hyp[nvar]*(xnvar-mean[nvar])=0
  //This function is called in case when less than fH samples lie on a hyperplane.

   Int_t i, j;
      TVectorD vec(fNvar);
   for (i=0; i<nmerged; i++) {
      ndist[i]=0;
      for(j=0; j<fNvar; j++) {
         ndist[i]+=fHyperplane[j]*(dat[i][j]-fMean[j]);
         ndist[i]=TMath::Abs(ndist[i]);
      }
   }
   KOrdStat(nmerged, ndist, hmerged-1, index);
   ClearSscp(sscp);
   for (i=0; i<hmerged; i++) {
      for(j=0; j<fNvar; j++)
         vec[j]=dat[index[i]][j];
      AddToSscp(sscp, vec);
   }
   Covar(sscp, fMean, fCovariance, fSd, hmerged);
}

//__________________________________________________________________________
Double_t TRobustEstimator::CStep(Int_t ntotal, Int_t htotal, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
{
  //from the input htotal-subset constructs another htotal subset with lower determinant
  //
  //As proven by Peter J.Rousseeuw and Katrien Van Driessen, if distances for all elements
  //are calculated, using the formula:d_i=Sqrt((x_i-M)*S_inv*(x_i-M)), where M is the mean
  //of the input htotal-subset, and S_inv - the inverse of its covariance matrix, then
  //htotal elements with smallest distances will have covariance matrix with determinant
  //less or equal to the determinant of the input subset covariance matrix.
  //
  //determinant for this htotal-subset with smallest distances is returned

   Int_t i, j;
   TVectorD vec(fNvar);
   Double_t det;
   
   TDecompChol chol(fCovariance);
   fInvcovariance = chol.Invert();

   TVectorD temp(fNvar);
   for(j=0; j<ntotal; j++) {
      ndist[j]=0;
      for(i=0; i<fNvar; i++)
         temp[i]=data[j][i]-fMean[i];
      temp*=fInvcovariance;
      for(i=0; i<fNvar; i++)
         ndist[j]+=(data[j][i]-fMean[i])*temp[i];
   }
   
   //taking h smallest
   KOrdStat(ntotal, ndist, htotal-1, index);
   //writing their mean and covariance
   ClearSscp(sscp);
   for (i=0; i<htotal; i++) {
      for (j=0; j<fNvar; j++)
         temp[j]=data[index[i]][j];
      AddToSscp(sscp, temp);
   }
   Covar(sscp, fMean, fCovariance, fSd, htotal);
   det = fCovariance.Determinant();
   return det;
}

//_______________________________________________________________
Int_t TRobustEstimator::Exact(Double_t *ndist)
{
  //for the exact fit situaions
  //returns number of observations on the hyperplane

   Int_t i, j;
   
   TMatrixDSymEigen eigen(fCovariance);
   TVectorD eigenValues=eigen.GetEigenValues();
   TMatrixD eigenMatrix=eigen.GetEigenVectors();
   
   for (j=0; j<fNvar; j++) {
      fHyperplane[j]=eigenMatrix(j,fNvar-1);
   }
   //calculate and return how many observations lie on the hyperplane
   for (i=0; i<fN; i++) {
      ndist[i]=0;
      for(j=0; j<fNvar; j++) {
         ndist[i]+=fHyperplane[j]*(fData[i][j]-fMean[j]);
         ndist[i]=TMath::Abs(ndist[i]);
      }
   }
   Int_t nhyp=0;
   
   for (i=0; i<fN; i++) {
      if(ndist[i] < 1e-14) nhyp++;
   }
   return nhyp;
   
}

//____________________________________________________________________________
Int_t TRobustEstimator::Exact2(TMatrixD &mstockbig, TMatrixD &cstockbig, TMatrixD &hyperplane,
                             Double_t *deti, Int_t nbest, Int_t kgroup,
                             TMatrixD &sscp, Double_t *ndist)
{
  //This function is called if determinant of the covariance matrix of a subset=0.
  //
  //If there are more then fH vectors on a hyperplane,
  //returns this hyperplane and stops
  //else stores the hyperplane coordinates in hyperplane matrix

   Int_t i, j;
   
   TVectorD vec(fNvar);
   Int_t maxind = TMath::LocMax(nbest, deti);
   Int_t nh=Exact(ndist);
   //now nh is the number of observation on the hyperplane
   //ndist stores distances of observation from this hyperplane
   if(nh>=fH) {
      ClearSscp(sscp);
      for (i=0; i<fN; i++) {
         if(ndist[i]<1e-14) {
            for (j=0; j<fNvar; j++)
               vec[j]=fData[i][j];
            AddToSscp(sscp, vec);
         }
      }
      Covar(sscp, fMean, fCovariance, fSd, nh);
      
      fExact=nh;
      return nbest+1;
      
   } else {
      //if less than fH observations lie on a hyperplane,
      //mean and covariance matrix are stored in mstockbig
      //and cstockbig in place of the previous maximum determinant
      //mean and covariance
      for(i=0; i<fNvar; i++) {
         mstockbig(nbest*kgroup+maxind,i)=fMean(i);
         hyperplane(nbest*kgroup+maxind,i)=fHyperplane(i);
         for(j=0; j<fNvar; j++) {
            cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
         }   
      }      
      return maxind;
   }
}


//_____________________________________________________________________________
Int_t TRobustEstimator::Partition(Int_t nmini, Int_t *indsubdat)
{
  //divides the elements into approximately equal subgroups
  //number of elements in each subgroup is stored in indsubdat
  //number of subgroups is returned

   Int_t nsub;
   if ((fN>=2*nmini) && (fN<=(3*nmini-1))) {
      if (fN%2==1){
         indsubdat[0]=Int_t(fN*0.5);
      indsubdat[1]=Int_t(fN*0.5)+1;
      } else
         indsubdat[0]=indsubdat[1]=Int_t(fN/2);
      nsub=2;
   }
   else{
      if((fN>=3*nmini) && (fN<(4*nmini -1))) {
         if(fN%3==0){
            indsubdat[0]=indsubdat[1]=indsubdat[2]=Int_t(fN/3);
         } else {
            indsubdat[0]=Int_t(fN/3);
            indsubdat[1]=Int_t(fN/3)+1;
            if (fN%3==1) indsubdat[2]=Int_t(fN/3);
            else indsubdat[2]=Int_t(fN/3)+1;
         }
         nsub=3;
      }
      else{
         if((fN>=4*nmini)&&(fN<=(5*nmini-1))){
            if (fN%4==0) indsubdat[0]=indsubdat[1]=indsubdat[2]=indsubdat[3]=Int_t(fN/4);
            else {
               indsubdat[0]=Int_t(fN/4);
               indsubdat[1]=Int_t(fN/4)+1;
               if(fN%4==1) indsubdat[2]=indsubdat[3]=Int_t(fN/4);
               if(fN%4==2) {
                  indsubdat[2]=Int_t(fN/4)+1;
                  indsubdat[3]=Int_t(fN/4);
               }
               if(fN%4==3) indsubdat[2]=indsubdat[3]=Int_t(fN/4)+1;
            }
            nsub=4;
         } else {
            for(Int_t i=0; i<5; i++)
               indsubdat[i]=nmini;
            nsub=5;
         }
      }
   }
   return nsub;
}

//___________________________________________________________________________
Int_t TRobustEstimator::RDist(TMatrixD &sscp)
{
  //Calculates robust distances.Then the samples with robust distances
  //greater than a cutoff value (0.975 quantile of chi2 distribution with
  //fNvar degrees of freedom, multiplied by a correction factor), are given
  //weiht=0, and new, reweighted estimates of location and scatter are calculated
  //The function returns the number of outliers.

   Int_t i, j;
   Int_t nout=0;
   
   TVectorD temp(fNvar);
   TDecompChol chol(fCovariance);
   fInvcovariance = chol.Invert();
   

   for (i=0; i<fN; i++) {
      fRd[i]=0;
      for(j=0; j<fNvar; j++) {
         temp[j]=fData[i][j]-fMean[j];
      }
      temp*=fInvcovariance;
      for(j=0; j<fNvar; j++) {
         fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
      }
   }
   
   Double_t med;
   Double_t chi = kChiMedian[fNvar-1];
   
   med=TMath::Median(fN, fRd.GetMatrixArray());
   med/=chi;
   fCovariance*=med;
   TDecompChol chol2(fCovariance);
   fInvcovariance = chol2.Invert();
   
   for (i=0; i<fN; i++) {
      fRd[i]=0;
      for(j=0; j<fNvar; j++) {
         temp[j]=fData[i][j]-fMean[j];
   }
      
      temp*=fInvcovariance;
      for(j=0; j<fNvar; j++) {
         fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
      }
   }
   
   Double_t cutoff = kChiQuant[fNvar-1];
   
   ClearSscp(sscp);
   for(i=0; i<fN; i++) {
      if (fRd[i]<=cutoff) {
         for(j=0; j<fNvar; j++)
            temp[j]=fData[i][j];
         AddToSscp(sscp,temp);
      } else {
         nout++;
      }
   }
   
   Covar(sscp, fMean, fCovariance, fSd, fN-nout);
   return nout;
}

//____________________________________________________________________________
void TRobustEstimator::RDraw(Int_t *subdat, Int_t ngroup, Int_t *indsubdat)
{
  //Draws ngroup nonoverlapping subdatasets out of a dataset of size n
  //such that the selected case numbers are uniformly distributed from 1 to n

   Int_t jndex = 0;
   Int_t nrand;
   Int_t i, k, m, j;
   for (k=1; k<=ngroup; k++) {
      for (m=1; m<=indsubdat[k-1]; m++) {

         nrand = Int_t(gRandom->Uniform(0, 1) * (fN-jndex))+1;
         
         jndex++;
         if (jndex==1) {
            subdat[0]=nrand;
         } else {
            subdat[jndex-1]=nrand+jndex-2;
            for (i=1; i<=jndex-1; i++) {
               if(subdat[i-1] > nrand+i-2) {
                  for(j=jndex; j>=i+1; j--) {
                     subdat[j-1]=subdat[j-2];
                  }
                  subdat[i-1]=nrand+i-2;
                  break;  //breaking the loop for(i=1...
               }
            }
         }
      }
   }
}

//_____________________________________________________________________________
Double_t TRobustEstimator::KOrdStat(Int_t ntotal, Double_t *a, Int_t k, Int_t *work){
  //because I need an Int_t work array
   Bool_t isAllocated = kFALSE;
   const Int_t kWorkMax=100;
   Int_t i, ir, j, l, mid;
   Int_t arr;
   Int_t *ind;
   Int_t workLocal[kWorkMax];
   Int_t temp;


   if (work) {
      ind = work;
   } else {
      ind = workLocal;
      if (ntotal > kWorkMax) {
         isAllocated = kTRUE;
         ind = new Int_t[ntotal];
      }
   }

   for (Int_t ii=0; ii<ntotal; ii++) {
      ind[ii]=ii;
   }
   Int_t rk = k;
   l=0;
   ir = ntotal-1;
   for(;;) {
      if (ir<=l+1) { //active partition contains 1 or 2 elements
         if (ir == l+1 && a[ind[ir]]<a[ind[l]])
            {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
         Double_t tmp = a[ind[rk]];
         if (isAllocated)
            delete [] ind;
         return tmp;
      } else {
         mid = (l+ir) >> 1; //choose median of left, center and right
         {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}//elements as partitioning element arr.
         if (a[ind[l]]>a[ind[ir]])  //also rearrange so that a[l]<=a[l+1]
            {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}

         if (a[ind[l+1]]>a[ind[ir]])
            {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}

         if (a[ind[l]]>a[ind[l+1]])
                {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}

         i=l+1;        //initialize pointers for partitioning
         j=ir;
         arr = ind[l+1];
         for (;;) {
            do i++; while (a[ind[i]]<a[arr]);
            do j--; while (a[ind[j]]>a[arr]);
            if (j<i) break;  //pointers crossed, partitioning complete
               {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
         }
         ind[l+1]=ind[j];
         ind[j]=arr;
         if (j>=rk) ir = j-1; //keep active the partition that
         if (j<=rk) l=i;      //contains the k_th element
      }
   }
}
 TRobustEstimator.cxx:1
 TRobustEstimator.cxx:2
 TRobustEstimator.cxx:3
 TRobustEstimator.cxx:4
 TRobustEstimator.cxx:5
 TRobustEstimator.cxx:6
 TRobustEstimator.cxx:7
 TRobustEstimator.cxx:8
 TRobustEstimator.cxx:9
 TRobustEstimator.cxx:10
 TRobustEstimator.cxx:11
 TRobustEstimator.cxx:12
 TRobustEstimator.cxx:13
 TRobustEstimator.cxx:14
 TRobustEstimator.cxx:15
 TRobustEstimator.cxx:16
 TRobustEstimator.cxx:17
 TRobustEstimator.cxx:18
 TRobustEstimator.cxx:19
 TRobustEstimator.cxx:20
 TRobustEstimator.cxx:21
 TRobustEstimator.cxx:22
 TRobustEstimator.cxx:23
 TRobustEstimator.cxx:24
 TRobustEstimator.cxx:25
 TRobustEstimator.cxx:26
 TRobustEstimator.cxx:27
 TRobustEstimator.cxx:28
 TRobustEstimator.cxx:29
 TRobustEstimator.cxx:30
 TRobustEstimator.cxx:31
 TRobustEstimator.cxx:32
 TRobustEstimator.cxx:33
 TRobustEstimator.cxx:34
 TRobustEstimator.cxx:35
 TRobustEstimator.cxx:36
 TRobustEstimator.cxx:37
 TRobustEstimator.cxx:38
 TRobustEstimator.cxx:39
 TRobustEstimator.cxx:40
 TRobustEstimator.cxx:41
 TRobustEstimator.cxx:42
 TRobustEstimator.cxx:43
 TRobustEstimator.cxx:44
 TRobustEstimator.cxx:45
 TRobustEstimator.cxx:46
 TRobustEstimator.cxx:47
 TRobustEstimator.cxx:48
 TRobustEstimator.cxx:49
 TRobustEstimator.cxx:50
 TRobustEstimator.cxx:51
 TRobustEstimator.cxx:52
 TRobustEstimator.cxx:53
 TRobustEstimator.cxx:54
 TRobustEstimator.cxx:55
 TRobustEstimator.cxx:56
 TRobustEstimator.cxx:57
 TRobustEstimator.cxx:58
 TRobustEstimator.cxx:59
 TRobustEstimator.cxx:60
 TRobustEstimator.cxx:61
 TRobustEstimator.cxx:62
 TRobustEstimator.cxx:63
 TRobustEstimator.cxx:64
 TRobustEstimator.cxx:65
 TRobustEstimator.cxx:66
 TRobustEstimator.cxx:67
 TRobustEstimator.cxx:68
 TRobustEstimator.cxx:69
 TRobustEstimator.cxx:70
 TRobustEstimator.cxx:71
 TRobustEstimator.cxx:72
 TRobustEstimator.cxx:73
 TRobustEstimator.cxx:74
 TRobustEstimator.cxx:75
 TRobustEstimator.cxx:76
 TRobustEstimator.cxx:77
 TRobustEstimator.cxx:78
 TRobustEstimator.cxx:79
 TRobustEstimator.cxx:80
 TRobustEstimator.cxx:81
 TRobustEstimator.cxx:82
 TRobustEstimator.cxx:83
 TRobustEstimator.cxx:84
 TRobustEstimator.cxx:85
 TRobustEstimator.cxx:86
 TRobustEstimator.cxx:87
 TRobustEstimator.cxx:88
 TRobustEstimator.cxx:89
 TRobustEstimator.cxx:90
 TRobustEstimator.cxx:91
 TRobustEstimator.cxx:92
 TRobustEstimator.cxx:93
 TRobustEstimator.cxx:94
 TRobustEstimator.cxx:95
 TRobustEstimator.cxx:96
 TRobustEstimator.cxx:97
 TRobustEstimator.cxx:98
 TRobustEstimator.cxx:99
 TRobustEstimator.cxx:100
 TRobustEstimator.cxx:101
 TRobustEstimator.cxx:102
 TRobustEstimator.cxx:103
 TRobustEstimator.cxx:104
 TRobustEstimator.cxx:105
 TRobustEstimator.cxx:106
 TRobustEstimator.cxx:107
 TRobustEstimator.cxx:108
 TRobustEstimator.cxx:109
 TRobustEstimator.cxx:110
 TRobustEstimator.cxx:111
 TRobustEstimator.cxx:112
 TRobustEstimator.cxx:113
 TRobustEstimator.cxx:114
 TRobustEstimator.cxx:115
 TRobustEstimator.cxx:116
 TRobustEstimator.cxx:117
 TRobustEstimator.cxx:118
 TRobustEstimator.cxx:119
 TRobustEstimator.cxx:120
 TRobustEstimator.cxx:121
 TRobustEstimator.cxx:122
 TRobustEstimator.cxx:123
 TRobustEstimator.cxx:124
 TRobustEstimator.cxx:125
 TRobustEstimator.cxx:126
 TRobustEstimator.cxx:127
 TRobustEstimator.cxx:128
 TRobustEstimator.cxx:129
 TRobustEstimator.cxx:130
 TRobustEstimator.cxx:131
 TRobustEstimator.cxx:132
 TRobustEstimator.cxx:133
 TRobustEstimator.cxx:134
 TRobustEstimator.cxx:135
 TRobustEstimator.cxx:136
 TRobustEstimator.cxx:137
 TRobustEstimator.cxx:138
 TRobustEstimator.cxx:139
 TRobustEstimator.cxx:140
 TRobustEstimator.cxx:141
 TRobustEstimator.cxx:142
 TRobustEstimator.cxx:143
 TRobustEstimator.cxx:144
 TRobustEstimator.cxx:145
 TRobustEstimator.cxx:146
 TRobustEstimator.cxx:147
 TRobustEstimator.cxx:148
 TRobustEstimator.cxx:149
 TRobustEstimator.cxx:150
 TRobustEstimator.cxx:151
 TRobustEstimator.cxx:152
 TRobustEstimator.cxx:153
 TRobustEstimator.cxx:154
 TRobustEstimator.cxx:155
 TRobustEstimator.cxx:156
 TRobustEstimator.cxx:157
 TRobustEstimator.cxx:158
 TRobustEstimator.cxx:159
 TRobustEstimator.cxx:160
 TRobustEstimator.cxx:161
 TRobustEstimator.cxx:162
 TRobustEstimator.cxx:163
 TRobustEstimator.cxx:164
 TRobustEstimator.cxx:165
 TRobustEstimator.cxx:166
 TRobustEstimator.cxx:167
 TRobustEstimator.cxx:168
 TRobustEstimator.cxx:169
 TRobustEstimator.cxx:170
 TRobustEstimator.cxx:171
 TRobustEstimator.cxx:172
 TRobustEstimator.cxx:173
 TRobustEstimator.cxx:174
 TRobustEstimator.cxx:175
 TRobustEstimator.cxx:176
 TRobustEstimator.cxx:177
 TRobustEstimator.cxx:178
 TRobustEstimator.cxx:179
 TRobustEstimator.cxx:180
 TRobustEstimator.cxx:181
 TRobustEstimator.cxx:182
 TRobustEstimator.cxx:183
 TRobustEstimator.cxx:184
 TRobustEstimator.cxx:185
 TRobustEstimator.cxx:186
 TRobustEstimator.cxx:187
 TRobustEstimator.cxx:188
 TRobustEstimator.cxx:189
 TRobustEstimator.cxx:190
 TRobustEstimator.cxx:191
 TRobustEstimator.cxx:192
 TRobustEstimator.cxx:193
 TRobustEstimator.cxx:194
 TRobustEstimator.cxx:195
 TRobustEstimator.cxx:196
 TRobustEstimator.cxx:197
 TRobustEstimator.cxx:198
 TRobustEstimator.cxx:199
 TRobustEstimator.cxx:200
 TRobustEstimator.cxx:201
 TRobustEstimator.cxx:202
 TRobustEstimator.cxx:203
 TRobustEstimator.cxx:204
 TRobustEstimator.cxx:205
 TRobustEstimator.cxx:206
 TRobustEstimator.cxx:207
 TRobustEstimator.cxx:208
 TRobustEstimator.cxx:209
 TRobustEstimator.cxx:210
 TRobustEstimator.cxx:211
 TRobustEstimator.cxx:212
 TRobustEstimator.cxx:213
 TRobustEstimator.cxx:214
 TRobustEstimator.cxx:215
 TRobustEstimator.cxx:216
 TRobustEstimator.cxx:217
 TRobustEstimator.cxx:218
 TRobustEstimator.cxx:219
 TRobustEstimator.cxx:220
 TRobustEstimator.cxx:221
 TRobustEstimator.cxx:222
 TRobustEstimator.cxx:223
 TRobustEstimator.cxx:224
 TRobustEstimator.cxx:225
 TRobustEstimator.cxx:226
 TRobustEstimator.cxx:227
 TRobustEstimator.cxx:228
 TRobustEstimator.cxx:229
 TRobustEstimator.cxx:230
 TRobustEstimator.cxx:231
 TRobustEstimator.cxx:232
 TRobustEstimator.cxx:233
 TRobustEstimator.cxx:234
 TRobustEstimator.cxx:235
 TRobustEstimator.cxx:236
 TRobustEstimator.cxx:237
 TRobustEstimator.cxx:238
 TRobustEstimator.cxx:239
 TRobustEstimator.cxx:240
 TRobustEstimator.cxx:241
 TRobustEstimator.cxx:242
 TRobustEstimator.cxx:243
 TRobustEstimator.cxx:244
 TRobustEstimator.cxx:245
 TRobustEstimator.cxx:246
 TRobustEstimator.cxx:247
 TRobustEstimator.cxx:248
 TRobustEstimator.cxx:249
 TRobustEstimator.cxx:250
 TRobustEstimator.cxx:251
 TRobustEstimator.cxx:252
 TRobustEstimator.cxx:253
 TRobustEstimator.cxx:254
 TRobustEstimator.cxx:255
 TRobustEstimator.cxx:256
 TRobustEstimator.cxx:257
 TRobustEstimator.cxx:258
 TRobustEstimator.cxx:259
 TRobustEstimator.cxx:260
 TRobustEstimator.cxx:261
 TRobustEstimator.cxx:262
 TRobustEstimator.cxx:263
 TRobustEstimator.cxx:264
 TRobustEstimator.cxx:265
 TRobustEstimator.cxx:266
 TRobustEstimator.cxx:267
 TRobustEstimator.cxx:268
 TRobustEstimator.cxx:269
 TRobustEstimator.cxx:270
 TRobustEstimator.cxx:271
 TRobustEstimator.cxx:272
 TRobustEstimator.cxx:273
 TRobustEstimator.cxx:274
 TRobustEstimator.cxx:275
 TRobustEstimator.cxx:276
 TRobustEstimator.cxx:277
 TRobustEstimator.cxx:278
 TRobustEstimator.cxx:279
 TRobustEstimator.cxx:280
 TRobustEstimator.cxx:281
 TRobustEstimator.cxx:282
 TRobustEstimator.cxx:283
 TRobustEstimator.cxx:284
 TRobustEstimator.cxx:285
 TRobustEstimator.cxx:286
 TRobustEstimator.cxx:287
 TRobustEstimator.cxx:288
 TRobustEstimator.cxx:289
 TRobustEstimator.cxx:290
 TRobustEstimator.cxx:291
 TRobustEstimator.cxx:292
 TRobustEstimator.cxx:293
 TRobustEstimator.cxx:294
 TRobustEstimator.cxx:295
 TRobustEstimator.cxx:296
 TRobustEstimator.cxx:297
 TRobustEstimator.cxx:298
 TRobustEstimator.cxx:299
 TRobustEstimator.cxx:300
 TRobustEstimator.cxx:301
 TRobustEstimator.cxx:302
 TRobustEstimator.cxx:303
 TRobustEstimator.cxx:304
 TRobustEstimator.cxx:305
 TRobustEstimator.cxx:306
 TRobustEstimator.cxx:307
 TRobustEstimator.cxx:308
 TRobustEstimator.cxx:309
 TRobustEstimator.cxx:310
 TRobustEstimator.cxx:311
 TRobustEstimator.cxx:312
 TRobustEstimator.cxx:313
 TRobustEstimator.cxx:314
 TRobustEstimator.cxx:315
 TRobustEstimator.cxx:316
 TRobustEstimator.cxx:317
 TRobustEstimator.cxx:318
 TRobustEstimator.cxx:319
 TRobustEstimator.cxx:320
 TRobustEstimator.cxx:321
 TRobustEstimator.cxx:322
 TRobustEstimator.cxx:323
 TRobustEstimator.cxx:324
 TRobustEstimator.cxx:325
 TRobustEstimator.cxx:326
 TRobustEstimator.cxx:327
 TRobustEstimator.cxx:328
 TRobustEstimator.cxx:329
 TRobustEstimator.cxx:330
 TRobustEstimator.cxx:331
 TRobustEstimator.cxx:332
 TRobustEstimator.cxx:333
 TRobustEstimator.cxx:334
 TRobustEstimator.cxx:335
 TRobustEstimator.cxx:336
 TRobustEstimator.cxx:337
 TRobustEstimator.cxx:338
 TRobustEstimator.cxx:339
 TRobustEstimator.cxx:340
 TRobustEstimator.cxx:341
 TRobustEstimator.cxx:342
 TRobustEstimator.cxx:343
 TRobustEstimator.cxx:344
 TRobustEstimator.cxx:345
 TRobustEstimator.cxx:346
 TRobustEstimator.cxx:347
 TRobustEstimator.cxx:348
 TRobustEstimator.cxx:349
 TRobustEstimator.cxx:350
 TRobustEstimator.cxx:351
 TRobustEstimator.cxx:352
 TRobustEstimator.cxx:353
 TRobustEstimator.cxx:354
 TRobustEstimator.cxx:355
 TRobustEstimator.cxx:356
 TRobustEstimator.cxx:357
 TRobustEstimator.cxx:358
 TRobustEstimator.cxx:359
 TRobustEstimator.cxx:360
 TRobustEstimator.cxx:361
 TRobustEstimator.cxx:362
 TRobustEstimator.cxx:363
 TRobustEstimator.cxx:364
 TRobustEstimator.cxx:365
 TRobustEstimator.cxx:366
 TRobustEstimator.cxx:367
 TRobustEstimator.cxx:368
 TRobustEstimator.cxx:369
 TRobustEstimator.cxx:370
 TRobustEstimator.cxx:371
 TRobustEstimator.cxx:372
 TRobustEstimator.cxx:373
 TRobustEstimator.cxx:374
 TRobustEstimator.cxx:375
 TRobustEstimator.cxx:376
 TRobustEstimator.cxx:377
 TRobustEstimator.cxx:378
 TRobustEstimator.cxx:379
 TRobustEstimator.cxx:380
 TRobustEstimator.cxx:381
 TRobustEstimator.cxx:382
 TRobustEstimator.cxx:383
 TRobustEstimator.cxx:384
 TRobustEstimator.cxx:385
 TRobustEstimator.cxx:386
 TRobustEstimator.cxx:387
 TRobustEstimator.cxx:388
 TRobustEstimator.cxx:389
 TRobustEstimator.cxx:390
 TRobustEstimator.cxx:391
 TRobustEstimator.cxx:392
 TRobustEstimator.cxx:393
 TRobustEstimator.cxx:394
 TRobustEstimator.cxx:395
 TRobustEstimator.cxx:396
 TRobustEstimator.cxx:397
 TRobustEstimator.cxx:398
 TRobustEstimator.cxx:399
 TRobustEstimator.cxx:400
 TRobustEstimator.cxx:401
 TRobustEstimator.cxx:402
 TRobustEstimator.cxx:403
 TRobustEstimator.cxx:404
 TRobustEstimator.cxx:405
 TRobustEstimator.cxx:406
 TRobustEstimator.cxx:407
 TRobustEstimator.cxx:408
 TRobustEstimator.cxx:409
 TRobustEstimator.cxx:410
 TRobustEstimator.cxx:411
 TRobustEstimator.cxx:412
 TRobustEstimator.cxx:413
 TRobustEstimator.cxx:414
 TRobustEstimator.cxx:415
 TRobustEstimator.cxx:416
 TRobustEstimator.cxx:417
 TRobustEstimator.cxx:418
 TRobustEstimator.cxx:419
 TRobustEstimator.cxx:420
 TRobustEstimator.cxx:421
 TRobustEstimator.cxx:422
 TRobustEstimator.cxx:423
 TRobustEstimator.cxx:424
 TRobustEstimator.cxx:425
 TRobustEstimator.cxx:426
 TRobustEstimator.cxx:427
 TRobustEstimator.cxx:428
 TRobustEstimator.cxx:429
 TRobustEstimator.cxx:430
 TRobustEstimator.cxx:431
 TRobustEstimator.cxx:432
 TRobustEstimator.cxx:433
 TRobustEstimator.cxx:434
 TRobustEstimator.cxx:435
 TRobustEstimator.cxx:436
 TRobustEstimator.cxx:437
 TRobustEstimator.cxx:438
 TRobustEstimator.cxx:439
 TRobustEstimator.cxx:440
 TRobustEstimator.cxx:441
 TRobustEstimator.cxx:442
 TRobustEstimator.cxx:443
 TRobustEstimator.cxx:444
 TRobustEstimator.cxx:445
 TRobustEstimator.cxx:446
 TRobustEstimator.cxx:447
 TRobustEstimator.cxx:448
 TRobustEstimator.cxx:449
 TRobustEstimator.cxx:450
 TRobustEstimator.cxx:451
 TRobustEstimator.cxx:452
 TRobustEstimator.cxx:453
 TRobustEstimator.cxx:454
 TRobustEstimator.cxx:455
 TRobustEstimator.cxx:456
 TRobustEstimator.cxx:457
 TRobustEstimator.cxx:458
 TRobustEstimator.cxx:459
 TRobustEstimator.cxx:460
 TRobustEstimator.cxx:461
 TRobustEstimator.cxx:462
 TRobustEstimator.cxx:463
 TRobustEstimator.cxx:464
 TRobustEstimator.cxx:465
 TRobustEstimator.cxx:466
 TRobustEstimator.cxx:467
 TRobustEstimator.cxx:468
 TRobustEstimator.cxx:469
 TRobustEstimator.cxx:470
 TRobustEstimator.cxx:471
 TRobustEstimator.cxx:472
 TRobustEstimator.cxx:473
 TRobustEstimator.cxx:474
 TRobustEstimator.cxx:475
 TRobustEstimator.cxx:476
 TRobustEstimator.cxx:477
 TRobustEstimator.cxx:478
 TRobustEstimator.cxx:479
 TRobustEstimator.cxx:480
 TRobustEstimator.cxx:481
 TRobustEstimator.cxx:482
 TRobustEstimator.cxx:483
 TRobustEstimator.cxx:484
 TRobustEstimator.cxx:485
 TRobustEstimator.cxx:486
 TRobustEstimator.cxx:487
 TRobustEstimator.cxx:488
 TRobustEstimator.cxx:489
 TRobustEstimator.cxx:490
 TRobustEstimator.cxx:491
 TRobustEstimator.cxx:492
 TRobustEstimator.cxx:493
 TRobustEstimator.cxx:494
 TRobustEstimator.cxx:495
 TRobustEstimator.cxx:496
 TRobustEstimator.cxx:497
 TRobustEstimator.cxx:498
 TRobustEstimator.cxx:499
 TRobustEstimator.cxx:500
 TRobustEstimator.cxx:501
 TRobustEstimator.cxx:502
 TRobustEstimator.cxx:503
 TRobustEstimator.cxx:504
 TRobustEstimator.cxx:505
 TRobustEstimator.cxx:506
 TRobustEstimator.cxx:507
 TRobustEstimator.cxx:508
 TRobustEstimator.cxx:509
 TRobustEstimator.cxx:510
 TRobustEstimator.cxx:511
 TRobustEstimator.cxx:512
 TRobustEstimator.cxx:513
 TRobustEstimator.cxx:514
 TRobustEstimator.cxx:515
 TRobustEstimator.cxx:516
 TRobustEstimator.cxx:517
 TRobustEstimator.cxx:518
 TRobustEstimator.cxx:519
 TRobustEstimator.cxx:520
 TRobustEstimator.cxx:521
 TRobustEstimator.cxx:522
 TRobustEstimator.cxx:523
 TRobustEstimator.cxx:524
 TRobustEstimator.cxx:525
 TRobustEstimator.cxx:526
 TRobustEstimator.cxx:527
 TRobustEstimator.cxx:528
 TRobustEstimator.cxx:529
 TRobustEstimator.cxx:530
 TRobustEstimator.cxx:531
 TRobustEstimator.cxx:532
 TRobustEstimator.cxx:533
 TRobustEstimator.cxx:534
 TRobustEstimator.cxx:535
 TRobustEstimator.cxx:536
 TRobustEstimator.cxx:537
 TRobustEstimator.cxx:538
 TRobustEstimator.cxx:539
 TRobustEstimator.cxx:540
 TRobustEstimator.cxx:541
 TRobustEstimator.cxx:542
 TRobustEstimator.cxx:543
 TRobustEstimator.cxx:544
 TRobustEstimator.cxx:545
 TRobustEstimator.cxx:546
 TRobustEstimator.cxx:547
 TRobustEstimator.cxx:548
 TRobustEstimator.cxx:549
 TRobustEstimator.cxx:550
 TRobustEstimator.cxx:551
 TRobustEstimator.cxx:552
 TRobustEstimator.cxx:553
 TRobustEstimator.cxx:554
 TRobustEstimator.cxx:555
 TRobustEstimator.cxx:556
 TRobustEstimator.cxx:557
 TRobustEstimator.cxx:558
 TRobustEstimator.cxx:559
 TRobustEstimator.cxx:560
 TRobustEstimator.cxx:561
 TRobustEstimator.cxx:562
 TRobustEstimator.cxx:563
 TRobustEstimator.cxx:564
 TRobustEstimator.cxx:565
 TRobustEstimator.cxx:566
 TRobustEstimator.cxx:567
 TRobustEstimator.cxx:568
 TRobustEstimator.cxx:569
 TRobustEstimator.cxx:570
 TRobustEstimator.cxx:571
 TRobustEstimator.cxx:572
 TRobustEstimator.cxx:573
 TRobustEstimator.cxx:574
 TRobustEstimator.cxx:575
 TRobustEstimator.cxx:576
 TRobustEstimator.cxx:577
 TRobustEstimator.cxx:578
 TRobustEstimator.cxx:579
 TRobustEstimator.cxx:580
 TRobustEstimator.cxx:581
 TRobustEstimator.cxx:582
 TRobustEstimator.cxx:583
 TRobustEstimator.cxx:584
 TRobustEstimator.cxx:585
 TRobustEstimator.cxx:586
 TRobustEstimator.cxx:587
 TRobustEstimator.cxx:588
 TRobustEstimator.cxx:589
 TRobustEstimator.cxx:590
 TRobustEstimator.cxx:591
 TRobustEstimator.cxx:592
 TRobustEstimator.cxx:593
 TRobustEstimator.cxx:594
 TRobustEstimator.cxx:595
 TRobustEstimator.cxx:596
 TRobustEstimator.cxx:597
 TRobustEstimator.cxx:598
 TRobustEstimator.cxx:599
 TRobustEstimator.cxx:600
 TRobustEstimator.cxx:601
 TRobustEstimator.cxx:602
 TRobustEstimator.cxx:603
 TRobustEstimator.cxx:604
 TRobustEstimator.cxx:605
 TRobustEstimator.cxx:606
 TRobustEstimator.cxx:607
 TRobustEstimator.cxx:608
 TRobustEstimator.cxx:609
 TRobustEstimator.cxx:610
 TRobustEstimator.cxx:611
 TRobustEstimator.cxx:612
 TRobustEstimator.cxx:613
 TRobustEstimator.cxx:614
 TRobustEstimator.cxx:615
 TRobustEstimator.cxx:616
 TRobustEstimator.cxx:617
 TRobustEstimator.cxx:618
 TRobustEstimator.cxx:619
 TRobustEstimator.cxx:620
 TRobustEstimator.cxx:621
 TRobustEstimator.cxx:622
 TRobustEstimator.cxx:623
 TRobustEstimator.cxx:624
 TRobustEstimator.cxx:625
 TRobustEstimator.cxx:626
 TRobustEstimator.cxx:627
 TRobustEstimator.cxx:628
 TRobustEstimator.cxx:629
 TRobustEstimator.cxx:630
 TRobustEstimator.cxx:631
 TRobustEstimator.cxx:632
 TRobustEstimator.cxx:633
 TRobustEstimator.cxx:634
 TRobustEstimator.cxx:635
 TRobustEstimator.cxx:636
 TRobustEstimator.cxx:637
 TRobustEstimator.cxx:638
 TRobustEstimator.cxx:639
 TRobustEstimator.cxx:640
 TRobustEstimator.cxx:641
 TRobustEstimator.cxx:642
 TRobustEstimator.cxx:643
 TRobustEstimator.cxx:644
 TRobustEstimator.cxx:645
 TRobustEstimator.cxx:646
 TRobustEstimator.cxx:647
 TRobustEstimator.cxx:648
 TRobustEstimator.cxx:649
 TRobustEstimator.cxx:650
 TRobustEstimator.cxx:651
 TRobustEstimator.cxx:652
 TRobustEstimator.cxx:653
 TRobustEstimator.cxx:654
 TRobustEstimator.cxx:655
 TRobustEstimator.cxx:656
 TRobustEstimator.cxx:657
 TRobustEstimator.cxx:658
 TRobustEstimator.cxx:659
 TRobustEstimator.cxx:660
 TRobustEstimator.cxx:661
 TRobustEstimator.cxx:662
 TRobustEstimator.cxx:663
 TRobustEstimator.cxx:664
 TRobustEstimator.cxx:665
 TRobustEstimator.cxx:666
 TRobustEstimator.cxx:667
 TRobustEstimator.cxx:668
 TRobustEstimator.cxx:669
 TRobustEstimator.cxx:670
 TRobustEstimator.cxx:671
 TRobustEstimator.cxx:672
 TRobustEstimator.cxx:673
 TRobustEstimator.cxx:674
 TRobustEstimator.cxx:675
 TRobustEstimator.cxx:676
 TRobustEstimator.cxx:677
 TRobustEstimator.cxx:678
 TRobustEstimator.cxx:679
 TRobustEstimator.cxx:680
 TRobustEstimator.cxx:681
 TRobustEstimator.cxx:682
 TRobustEstimator.cxx:683
 TRobustEstimator.cxx:684
 TRobustEstimator.cxx:685
 TRobustEstimator.cxx:686
 TRobustEstimator.cxx:687
 TRobustEstimator.cxx:688
 TRobustEstimator.cxx:689
 TRobustEstimator.cxx:690
 TRobustEstimator.cxx:691
 TRobustEstimator.cxx:692
 TRobustEstimator.cxx:693
 TRobustEstimator.cxx:694
 TRobustEstimator.cxx:695
 TRobustEstimator.cxx:696
 TRobustEstimator.cxx:697
 TRobustEstimator.cxx:698
 TRobustEstimator.cxx:699
 TRobustEstimator.cxx:700
 TRobustEstimator.cxx:701
 TRobustEstimator.cxx:702
 TRobustEstimator.cxx:703
 TRobustEstimator.cxx:704
 TRobustEstimator.cxx:705
 TRobustEstimator.cxx:706
 TRobustEstimator.cxx:707
 TRobustEstimator.cxx:708
 TRobustEstimator.cxx:709
 TRobustEstimator.cxx:710
 TRobustEstimator.cxx:711
 TRobustEstimator.cxx:712
 TRobustEstimator.cxx:713
 TRobustEstimator.cxx:714
 TRobustEstimator.cxx:715
 TRobustEstimator.cxx:716
 TRobustEstimator.cxx:717
 TRobustEstimator.cxx:718
 TRobustEstimator.cxx:719
 TRobustEstimator.cxx:720
 TRobustEstimator.cxx:721
 TRobustEstimator.cxx:722
 TRobustEstimator.cxx:723
 TRobustEstimator.cxx:724
 TRobustEstimator.cxx:725
 TRobustEstimator.cxx:726
 TRobustEstimator.cxx:727
 TRobustEstimator.cxx:728
 TRobustEstimator.cxx:729
 TRobustEstimator.cxx:730
 TRobustEstimator.cxx:731
 TRobustEstimator.cxx:732
 TRobustEstimator.cxx:733
 TRobustEstimator.cxx:734
 TRobustEstimator.cxx:735
 TRobustEstimator.cxx:736
 TRobustEstimator.cxx:737
 TRobustEstimator.cxx:738
 TRobustEstimator.cxx:739
 TRobustEstimator.cxx:740
 TRobustEstimator.cxx:741
 TRobustEstimator.cxx:742
 TRobustEstimator.cxx:743
 TRobustEstimator.cxx:744
 TRobustEstimator.cxx:745
 TRobustEstimator.cxx:746
 TRobustEstimator.cxx:747
 TRobustEstimator.cxx:748
 TRobustEstimator.cxx:749
 TRobustEstimator.cxx:750
 TRobustEstimator.cxx:751
 TRobustEstimator.cxx:752
 TRobustEstimator.cxx:753
 TRobustEstimator.cxx:754
 TRobustEstimator.cxx:755
 TRobustEstimator.cxx:756
 TRobustEstimator.cxx:757
 TRobustEstimator.cxx:758
 TRobustEstimator.cxx:759
 TRobustEstimator.cxx:760
 TRobustEstimator.cxx:761
 TRobustEstimator.cxx:762
 TRobustEstimator.cxx:763
 TRobustEstimator.cxx:764
 TRobustEstimator.cxx:765
 TRobustEstimator.cxx:766
 TRobustEstimator.cxx:767
 TRobustEstimator.cxx:768
 TRobustEstimator.cxx:769
 TRobustEstimator.cxx:770
 TRobustEstimator.cxx:771
 TRobustEstimator.cxx:772
 TRobustEstimator.cxx:773
 TRobustEstimator.cxx:774
 TRobustEstimator.cxx:775
 TRobustEstimator.cxx:776
 TRobustEstimator.cxx:777
 TRobustEstimator.cxx:778
 TRobustEstimator.cxx:779
 TRobustEstimator.cxx:780
 TRobustEstimator.cxx:781
 TRobustEstimator.cxx:782
 TRobustEstimator.cxx:783
 TRobustEstimator.cxx:784
 TRobustEstimator.cxx:785
 TRobustEstimator.cxx:786
 TRobustEstimator.cxx:787
 TRobustEstimator.cxx:788
 TRobustEstimator.cxx:789
 TRobustEstimator.cxx:790
 TRobustEstimator.cxx:791
 TRobustEstimator.cxx:792
 TRobustEstimator.cxx:793
 TRobustEstimator.cxx:794
 TRobustEstimator.cxx:795
 TRobustEstimator.cxx:796
 TRobustEstimator.cxx:797
 TRobustEstimator.cxx:798
 TRobustEstimator.cxx:799
 TRobustEstimator.cxx:800
 TRobustEstimator.cxx:801
 TRobustEstimator.cxx:802
 TRobustEstimator.cxx:803
 TRobustEstimator.cxx:804
 TRobustEstimator.cxx:805
 TRobustEstimator.cxx:806
 TRobustEstimator.cxx:807
 TRobustEstimator.cxx:808
 TRobustEstimator.cxx:809
 TRobustEstimator.cxx:810
 TRobustEstimator.cxx:811
 TRobustEstimator.cxx:812
 TRobustEstimator.cxx:813
 TRobustEstimator.cxx:814
 TRobustEstimator.cxx:815
 TRobustEstimator.cxx:816
 TRobustEstimator.cxx:817
 TRobustEstimator.cxx:818
 TRobustEstimator.cxx:819
 TRobustEstimator.cxx:820
 TRobustEstimator.cxx:821
 TRobustEstimator.cxx:822
 TRobustEstimator.cxx:823
 TRobustEstimator.cxx:824
 TRobustEstimator.cxx:825
 TRobustEstimator.cxx:826
 TRobustEstimator.cxx:827
 TRobustEstimator.cxx:828
 TRobustEstimator.cxx:829
 TRobustEstimator.cxx:830
 TRobustEstimator.cxx:831
 TRobustEstimator.cxx:832
 TRobustEstimator.cxx:833
 TRobustEstimator.cxx:834
 TRobustEstimator.cxx:835
 TRobustEstimator.cxx:836
 TRobustEstimator.cxx:837
 TRobustEstimator.cxx:838
 TRobustEstimator.cxx:839
 TRobustEstimator.cxx:840
 TRobustEstimator.cxx:841
 TRobustEstimator.cxx:842
 TRobustEstimator.cxx:843
 TRobustEstimator.cxx:844
 TRobustEstimator.cxx:845
 TRobustEstimator.cxx:846
 TRobustEstimator.cxx:847
 TRobustEstimator.cxx:848
 TRobustEstimator.cxx:849
 TRobustEstimator.cxx:850
 TRobustEstimator.cxx:851
 TRobustEstimator.cxx:852
 TRobustEstimator.cxx:853
 TRobustEstimator.cxx:854
 TRobustEstimator.cxx:855
 TRobustEstimator.cxx:856
 TRobustEstimator.cxx:857
 TRobustEstimator.cxx:858
 TRobustEstimator.cxx:859
 TRobustEstimator.cxx:860
 TRobustEstimator.cxx:861
 TRobustEstimator.cxx:862
 TRobustEstimator.cxx:863
 TRobustEstimator.cxx:864
 TRobustEstimator.cxx:865
 TRobustEstimator.cxx:866
 TRobustEstimator.cxx:867
 TRobustEstimator.cxx:868
 TRobustEstimator.cxx:869
 TRobustEstimator.cxx:870
 TRobustEstimator.cxx:871
 TRobustEstimator.cxx:872
 TRobustEstimator.cxx:873
 TRobustEstimator.cxx:874
 TRobustEstimator.cxx:875
 TRobustEstimator.cxx:876
 TRobustEstimator.cxx:877
 TRobustEstimator.cxx:878
 TRobustEstimator.cxx:879
 TRobustEstimator.cxx:880
 TRobustEstimator.cxx:881
 TRobustEstimator.cxx:882
 TRobustEstimator.cxx:883
 TRobustEstimator.cxx:884
 TRobustEstimator.cxx:885
 TRobustEstimator.cxx:886
 TRobustEstimator.cxx:887
 TRobustEstimator.cxx:888
 TRobustEstimator.cxx:889
 TRobustEstimator.cxx:890
 TRobustEstimator.cxx:891
 TRobustEstimator.cxx:892
 TRobustEstimator.cxx:893
 TRobustEstimator.cxx:894
 TRobustEstimator.cxx:895
 TRobustEstimator.cxx:896
 TRobustEstimator.cxx:897
 TRobustEstimator.cxx:898
 TRobustEstimator.cxx:899
 TRobustEstimator.cxx:900
 TRobustEstimator.cxx:901
 TRobustEstimator.cxx:902
 TRobustEstimator.cxx:903
 TRobustEstimator.cxx:904
 TRobustEstimator.cxx:905
 TRobustEstimator.cxx:906
 TRobustEstimator.cxx:907
 TRobustEstimator.cxx:908
 TRobustEstimator.cxx:909
 TRobustEstimator.cxx:910
 TRobustEstimator.cxx:911
 TRobustEstimator.cxx:912
 TRobustEstimator.cxx:913
 TRobustEstimator.cxx:914
 TRobustEstimator.cxx:915
 TRobustEstimator.cxx:916
 TRobustEstimator.cxx:917
 TRobustEstimator.cxx:918
 TRobustEstimator.cxx:919
 TRobustEstimator.cxx:920
 TRobustEstimator.cxx:921
 TRobustEstimator.cxx:922
 TRobustEstimator.cxx:923
 TRobustEstimator.cxx:924
 TRobustEstimator.cxx:925
 TRobustEstimator.cxx:926
 TRobustEstimator.cxx:927
 TRobustEstimator.cxx:928
 TRobustEstimator.cxx:929
 TRobustEstimator.cxx:930
 TRobustEstimator.cxx:931
 TRobustEstimator.cxx:932
 TRobustEstimator.cxx:933
 TRobustEstimator.cxx:934
 TRobustEstimator.cxx:935
 TRobustEstimator.cxx:936
 TRobustEstimator.cxx:937
 TRobustEstimator.cxx:938
 TRobustEstimator.cxx:939
 TRobustEstimator.cxx:940
 TRobustEstimator.cxx:941
 TRobustEstimator.cxx:942
 TRobustEstimator.cxx:943
 TRobustEstimator.cxx:944
 TRobustEstimator.cxx:945
 TRobustEstimator.cxx:946
 TRobustEstimator.cxx:947
 TRobustEstimator.cxx:948
 TRobustEstimator.cxx:949
 TRobustEstimator.cxx:950
 TRobustEstimator.cxx:951
 TRobustEstimator.cxx:952
 TRobustEstimator.cxx:953
 TRobustEstimator.cxx:954
 TRobustEstimator.cxx:955
 TRobustEstimator.cxx:956
 TRobustEstimator.cxx:957
 TRobustEstimator.cxx:958
 TRobustEstimator.cxx:959
 TRobustEstimator.cxx:960
 TRobustEstimator.cxx:961
 TRobustEstimator.cxx:962
 TRobustEstimator.cxx:963
 TRobustEstimator.cxx:964
 TRobustEstimator.cxx:965
 TRobustEstimator.cxx:966
 TRobustEstimator.cxx:967
 TRobustEstimator.cxx:968
 TRobustEstimator.cxx:969
 TRobustEstimator.cxx:970
 TRobustEstimator.cxx:971
 TRobustEstimator.cxx:972
 TRobustEstimator.cxx:973
 TRobustEstimator.cxx:974
 TRobustEstimator.cxx:975
 TRobustEstimator.cxx:976
 TRobustEstimator.cxx:977
 TRobustEstimator.cxx:978
 TRobustEstimator.cxx:979
 TRobustEstimator.cxx:980
 TRobustEstimator.cxx:981
 TRobustEstimator.cxx:982
 TRobustEstimator.cxx:983
 TRobustEstimator.cxx:984
 TRobustEstimator.cxx:985
 TRobustEstimator.cxx:986
 TRobustEstimator.cxx:987
 TRobustEstimator.cxx:988
 TRobustEstimator.cxx:989
 TRobustEstimator.cxx:990
 TRobustEstimator.cxx:991
 TRobustEstimator.cxx:992
 TRobustEstimator.cxx:993
 TRobustEstimator.cxx:994
 TRobustEstimator.cxx:995
 TRobustEstimator.cxx:996
 TRobustEstimator.cxx:997
 TRobustEstimator.cxx:998
 TRobustEstimator.cxx:999
 TRobustEstimator.cxx:1000
 TRobustEstimator.cxx:1001
 TRobustEstimator.cxx:1002
 TRobustEstimator.cxx:1003
 TRobustEstimator.cxx:1004
 TRobustEstimator.cxx:1005
 TRobustEstimator.cxx:1006
 TRobustEstimator.cxx:1007
 TRobustEstimator.cxx:1008
 TRobustEstimator.cxx:1009
 TRobustEstimator.cxx:1010
 TRobustEstimator.cxx:1011
 TRobustEstimator.cxx:1012
 TRobustEstimator.cxx:1013
 TRobustEstimator.cxx:1014
 TRobustEstimator.cxx:1015
 TRobustEstimator.cxx:1016
 TRobustEstimator.cxx:1017
 TRobustEstimator.cxx:1018
 TRobustEstimator.cxx:1019
 TRobustEstimator.cxx:1020
 TRobustEstimator.cxx:1021
 TRobustEstimator.cxx:1022
 TRobustEstimator.cxx:1023
 TRobustEstimator.cxx:1024
 TRobustEstimator.cxx:1025
 TRobustEstimator.cxx:1026
 TRobustEstimator.cxx:1027
 TRobustEstimator.cxx:1028
 TRobustEstimator.cxx:1029
 TRobustEstimator.cxx:1030
 TRobustEstimator.cxx:1031
 TRobustEstimator.cxx:1032
 TRobustEstimator.cxx:1033
 TRobustEstimator.cxx:1034
 TRobustEstimator.cxx:1035
 TRobustEstimator.cxx:1036
 TRobustEstimator.cxx:1037
 TRobustEstimator.cxx:1038
 TRobustEstimator.cxx:1039
 TRobustEstimator.cxx:1040
 TRobustEstimator.cxx:1041
 TRobustEstimator.cxx:1042
 TRobustEstimator.cxx:1043
 TRobustEstimator.cxx:1044
 TRobustEstimator.cxx:1045
 TRobustEstimator.cxx:1046
 TRobustEstimator.cxx:1047
 TRobustEstimator.cxx:1048
 TRobustEstimator.cxx:1049
 TRobustEstimator.cxx:1050
 TRobustEstimator.cxx:1051
 TRobustEstimator.cxx:1052
 TRobustEstimator.cxx:1053
 TRobustEstimator.cxx:1054
 TRobustEstimator.cxx:1055
 TRobustEstimator.cxx:1056
 TRobustEstimator.cxx:1057
 TRobustEstimator.cxx:1058
 TRobustEstimator.cxx:1059
 TRobustEstimator.cxx:1060
 TRobustEstimator.cxx:1061
 TRobustEstimator.cxx:1062
 TRobustEstimator.cxx:1063
 TRobustEstimator.cxx:1064
 TRobustEstimator.cxx:1065
 TRobustEstimator.cxx:1066
 TRobustEstimator.cxx:1067
 TRobustEstimator.cxx:1068
 TRobustEstimator.cxx:1069
 TRobustEstimator.cxx:1070
 TRobustEstimator.cxx:1071
 TRobustEstimator.cxx:1072
 TRobustEstimator.cxx:1073
 TRobustEstimator.cxx:1074
 TRobustEstimator.cxx:1075
 TRobustEstimator.cxx:1076
 TRobustEstimator.cxx:1077
 TRobustEstimator.cxx:1078
 TRobustEstimator.cxx:1079
 TRobustEstimator.cxx:1080
 TRobustEstimator.cxx:1081
 TRobustEstimator.cxx:1082
 TRobustEstimator.cxx:1083
 TRobustEstimator.cxx:1084
 TRobustEstimator.cxx:1085
 TRobustEstimator.cxx:1086
 TRobustEstimator.cxx:1087
 TRobustEstimator.cxx:1088
 TRobustEstimator.cxx:1089
 TRobustEstimator.cxx:1090
 TRobustEstimator.cxx:1091
 TRobustEstimator.cxx:1092
 TRobustEstimator.cxx:1093
 TRobustEstimator.cxx:1094
 TRobustEstimator.cxx:1095
 TRobustEstimator.cxx:1096
 TRobustEstimator.cxx:1097
 TRobustEstimator.cxx:1098
 TRobustEstimator.cxx:1099
 TRobustEstimator.cxx:1100
 TRobustEstimator.cxx:1101
 TRobustEstimator.cxx:1102
 TRobustEstimator.cxx:1103
 TRobustEstimator.cxx:1104
 TRobustEstimator.cxx:1105
 TRobustEstimator.cxx:1106
 TRobustEstimator.cxx:1107
 TRobustEstimator.cxx:1108
 TRobustEstimator.cxx:1109
 TRobustEstimator.cxx:1110
 TRobustEstimator.cxx:1111
 TRobustEstimator.cxx:1112
 TRobustEstimator.cxx:1113
 TRobustEstimator.cxx:1114
 TRobustEstimator.cxx:1115
 TRobustEstimator.cxx:1116
 TRobustEstimator.cxx:1117
 TRobustEstimator.cxx:1118
 TRobustEstimator.cxx:1119
 TRobustEstimator.cxx:1120
 TRobustEstimator.cxx:1121
 TRobustEstimator.cxx:1122
 TRobustEstimator.cxx:1123
 TRobustEstimator.cxx:1124
 TRobustEstimator.cxx:1125
 TRobustEstimator.cxx:1126
 TRobustEstimator.cxx:1127
 TRobustEstimator.cxx:1128
 TRobustEstimator.cxx:1129
 TRobustEstimator.cxx:1130
 TRobustEstimator.cxx:1131
 TRobustEstimator.cxx:1132
 TRobustEstimator.cxx:1133
 TRobustEstimator.cxx:1134
 TRobustEstimator.cxx:1135
 TRobustEstimator.cxx:1136
 TRobustEstimator.cxx:1137
 TRobustEstimator.cxx:1138
 TRobustEstimator.cxx:1139
 TRobustEstimator.cxx:1140
 TRobustEstimator.cxx:1141
 TRobustEstimator.cxx:1142
 TRobustEstimator.cxx:1143
 TRobustEstimator.cxx:1144
 TRobustEstimator.cxx:1145
 TRobustEstimator.cxx:1146
 TRobustEstimator.cxx:1147
 TRobustEstimator.cxx:1148
 TRobustEstimator.cxx:1149
 TRobustEstimator.cxx:1150
 TRobustEstimator.cxx:1151
 TRobustEstimator.cxx:1152
 TRobustEstimator.cxx:1153
 TRobustEstimator.cxx:1154
 TRobustEstimator.cxx:1155
 TRobustEstimator.cxx:1156
 TRobustEstimator.cxx:1157
 TRobustEstimator.cxx:1158
 TRobustEstimator.cxx:1159
 TRobustEstimator.cxx:1160
 TRobustEstimator.cxx:1161
 TRobustEstimator.cxx:1162
 TRobustEstimator.cxx:1163
 TRobustEstimator.cxx:1164
 TRobustEstimator.cxx:1165
 TRobustEstimator.cxx:1166
 TRobustEstimator.cxx:1167
 TRobustEstimator.cxx:1168
 TRobustEstimator.cxx:1169
 TRobustEstimator.cxx:1170
 TRobustEstimator.cxx:1171
 TRobustEstimator.cxx:1172
 TRobustEstimator.cxx:1173
 TRobustEstimator.cxx:1174
 TRobustEstimator.cxx:1175
 TRobustEstimator.cxx:1176
 TRobustEstimator.cxx:1177
 TRobustEstimator.cxx:1178
 TRobustEstimator.cxx:1179
 TRobustEstimator.cxx:1180
 TRobustEstimator.cxx:1181
 TRobustEstimator.cxx:1182
 TRobustEstimator.cxx:1183
 TRobustEstimator.cxx:1184
 TRobustEstimator.cxx:1185
 TRobustEstimator.cxx:1186
 TRobustEstimator.cxx:1187
 TRobustEstimator.cxx:1188
 TRobustEstimator.cxx:1189
 TRobustEstimator.cxx:1190
 TRobustEstimator.cxx:1191
 TRobustEstimator.cxx:1192
 TRobustEstimator.cxx:1193
 TRobustEstimator.cxx:1194
 TRobustEstimator.cxx:1195
 TRobustEstimator.cxx:1196
 TRobustEstimator.cxx:1197
 TRobustEstimator.cxx:1198
 TRobustEstimator.cxx:1199
 TRobustEstimator.cxx:1200
 TRobustEstimator.cxx:1201
 TRobustEstimator.cxx:1202
 TRobustEstimator.cxx:1203
 TRobustEstimator.cxx:1204
 TRobustEstimator.cxx:1205
 TRobustEstimator.cxx:1206
 TRobustEstimator.cxx:1207
 TRobustEstimator.cxx:1208
 TRobustEstimator.cxx:1209
 TRobustEstimator.cxx:1210
 TRobustEstimator.cxx:1211
 TRobustEstimator.cxx:1212
 TRobustEstimator.cxx:1213
 TRobustEstimator.cxx:1214
 TRobustEstimator.cxx:1215
 TRobustEstimator.cxx:1216
 TRobustEstimator.cxx:1217
 TRobustEstimator.cxx:1218
 TRobustEstimator.cxx:1219
 TRobustEstimator.cxx:1220
 TRobustEstimator.cxx:1221
 TRobustEstimator.cxx:1222
 TRobustEstimator.cxx:1223
 TRobustEstimator.cxx:1224
 TRobustEstimator.cxx:1225
 TRobustEstimator.cxx:1226
 TRobustEstimator.cxx:1227
 TRobustEstimator.cxx:1228
 TRobustEstimator.cxx:1229
 TRobustEstimator.cxx:1230
 TRobustEstimator.cxx:1231
 TRobustEstimator.cxx:1232
 TRobustEstimator.cxx:1233
 TRobustEstimator.cxx:1234
 TRobustEstimator.cxx:1235
 TRobustEstimator.cxx:1236
 TRobustEstimator.cxx:1237
 TRobustEstimator.cxx:1238
 TRobustEstimator.cxx:1239
 TRobustEstimator.cxx:1240
 TRobustEstimator.cxx:1241
 TRobustEstimator.cxx:1242
 TRobustEstimator.cxx:1243
 TRobustEstimator.cxx:1244
 TRobustEstimator.cxx:1245
 TRobustEstimator.cxx:1246
 TRobustEstimator.cxx:1247
 TRobustEstimator.cxx:1248
 TRobustEstimator.cxx:1249
 TRobustEstimator.cxx:1250
 TRobustEstimator.cxx:1251
 TRobustEstimator.cxx:1252
 TRobustEstimator.cxx:1253
 TRobustEstimator.cxx:1254
 TRobustEstimator.cxx:1255
 TRobustEstimator.cxx:1256
 TRobustEstimator.cxx:1257
 TRobustEstimator.cxx:1258
 TRobustEstimator.cxx:1259
 TRobustEstimator.cxx:1260
 TRobustEstimator.cxx:1261
 TRobustEstimator.cxx:1262
 TRobustEstimator.cxx:1263
 TRobustEstimator.cxx:1264
 TRobustEstimator.cxx:1265
 TRobustEstimator.cxx:1266
 TRobustEstimator.cxx:1267
 TRobustEstimator.cxx:1268
 TRobustEstimator.cxx:1269
 TRobustEstimator.cxx:1270
 TRobustEstimator.cxx:1271
 TRobustEstimator.cxx:1272
 TRobustEstimator.cxx:1273
 TRobustEstimator.cxx:1274
 TRobustEstimator.cxx:1275
 TRobustEstimator.cxx:1276
 TRobustEstimator.cxx:1277
 TRobustEstimator.cxx:1278
 TRobustEstimator.cxx:1279
 TRobustEstimator.cxx:1280
 TRobustEstimator.cxx:1281
 TRobustEstimator.cxx:1282
 TRobustEstimator.cxx:1283
 TRobustEstimator.cxx:1284
 TRobustEstimator.cxx:1285
 TRobustEstimator.cxx:1286
 TRobustEstimator.cxx:1287
 TRobustEstimator.cxx:1288
 TRobustEstimator.cxx:1289
 TRobustEstimator.cxx:1290
 TRobustEstimator.cxx:1291
 TRobustEstimator.cxx:1292
 TRobustEstimator.cxx:1293
 TRobustEstimator.cxx:1294
 TRobustEstimator.cxx:1295
 TRobustEstimator.cxx:1296
 TRobustEstimator.cxx:1297
 TRobustEstimator.cxx:1298
 TRobustEstimator.cxx:1299
 TRobustEstimator.cxx:1300
 TRobustEstimator.cxx:1301
 TRobustEstimator.cxx:1302
 TRobustEstimator.cxx:1303
 TRobustEstimator.cxx:1304
 TRobustEstimator.cxx:1305
 TRobustEstimator.cxx:1306
 TRobustEstimator.cxx:1307
 TRobustEstimator.cxx:1308
 TRobustEstimator.cxx:1309
 TRobustEstimator.cxx:1310
 TRobustEstimator.cxx:1311
 TRobustEstimator.cxx:1312
 TRobustEstimator.cxx:1313
 TRobustEstimator.cxx:1314
 TRobustEstimator.cxx:1315
 TRobustEstimator.cxx:1316
 TRobustEstimator.cxx:1317
 TRobustEstimator.cxx:1318
 TRobustEstimator.cxx:1319
 TRobustEstimator.cxx:1320
 TRobustEstimator.cxx:1321
 TRobustEstimator.cxx:1322