102 ClassImp(TRobustEstimator);
104 const Double_t kChiMedian[50]= {
105 0.454937, 1.38629, 2.36597, 3.35670, 4.35146, 5.34812, 6.34581, 7.34412, 8.34283,
106 9.34182, 10.34, 11.34, 12.34, 13.34, 14.34, 15.34, 16.34, 17.34, 18.34, 19.34,
107 20.34, 21.34, 22.34, 23.34, 24.34, 25.34, 26.34, 27.34, 28.34, 29.34, 30.34,
108 31.34, 32.34, 33.34, 34.34, 35.34, 36.34, 37.34, 38.34, 39.34, 40.34,
109 41.34, 42.34, 43.34, 44.34, 45.34, 46.34, 47.34, 48.34, 49.33};
111 const Double_t kChiQuant[50]={
112 5.02389, 7.3776,9.34840,11.1433,12.8325,
113 14.4494,16.0128,17.5346,19.0228,20.4831,21.920,23.337,
114 24.736,26.119,27.488,28.845,30.191,31.526,32.852,34.170,
115 35.479,36.781,38.076,39.364,40.646,41.923,43.194,44.461,
116 45.722,46.979,48.232,49.481,50.725,51.966,53.203,54.437,
117 55.668,56.896,58.120,59.342,60.561,61.777,62.990,64.201,
118 65.410,66.617,67.821,69.022,70.222,71.420};
124 TRobustEstimator::TRobustEstimator(){
130 TRobustEstimator::TRobustEstimator(Int_t nvectors, Int_t nvariables, Int_t hh)
132 fCovariance(nvariables),
133 fInvcovariance(nvariables),
134 fCorrelation(nvariables),
138 fHyperplane(nvariables),
139 fData(nvectors, nvariables)
141 if ((nvectors<=1)||(nvariables<=0)){
142 Error(
"TRobustEstimator",
"Not enough vectors or variables");
146 Error(
"TRobustEstimator",
"For the univariate case, use the default constructor and EvaluateUni() function");
152 if (hh<(fN+fNvar+1)/2){
154 Warning(
"TRobustEstimator",
"chosen h is too small, default h is taken instead");
170 void TRobustEstimator::AddColumn(Double_t *col)
172 if (fVarTemp==fNvar) {
174 fCovariance.ResizeTo(fNvar, fNvar);
175 fInvcovariance.ResizeTo(fNvar, fNvar);
176 fCorrelation.ResizeTo(fNvar, fNvar);
177 fMean.ResizeTo(fNvar);
178 fHyperplane.ResizeTo(fNvar);
179 fData.ResizeTo(fN, fNvar);
181 for (Int_t i=0; i<fN; i++) {
182 fData(i, fVarTemp)=col[i];
191 void TRobustEstimator::AddRow(Double_t *row)
197 fData.ResizeTo(fN, fNvar);
199 for (Int_t i=0; i<fNvar; i++)
200 fData(fVecTemp, i)=row[i];
208 void TRobustEstimator::Evaluate()
213 Warning(
"Evaluate",
"Chosen h = #observations, so classic estimates of location and scatter will be calculated");
223 TMatrixD sscp(fNvar+1, fNvar+1);
226 Int_t *index =
new Int_t[fN];
227 Double_t *ndist =
new Double_t[fN];
229 Double_t *deti=
new Double_t[nbest];
230 for (i=0; i<nbest; i++)
241 TMatrixD mstock(nbest, fNvar);
242 TMatrixD cstock(fNvar, fNvar*nbest);
244 for (k=0; k<k1; k++) {
245 CreateSubset(fN, fH, fNvar, index, fData, sscp, ndist);
248 for (i=0; i<fH; i++) {
249 for(j=0; j<fNvar; j++)
250 vec(j)=fData[index[i]][j];
251 AddToSscp(sscp, vec);
253 Covar(sscp, fMean, fCovariance, fSd, fH);
254 det = fCovariance.Determinant();
256 fExact = Exact(ndist);
263 det = CStep(fN, fH, index, fData, sscp, ndist);
265 fExact = Exact(ndist);
271 det = CStep(fN, fH, index, fData, sscp, ndist);
273 fExact = Exact(ndist);
279 Int_t maxind=TMath::LocMax(nbest, deti);
280 if(det<deti[maxind]) {
282 for(ii=0; ii<fNvar; ii++) {
283 mstock(maxind, ii)=fMean(ii);
284 for(jj=0; jj<fNvar; jj++)
285 cstock(ii, jj+maxind*fNvar)=fCovariance(ii, jj);
293 for (i=0; i<nbest; i++) {
294 for(ii=0; ii<fNvar; ii++) {
295 fMean(ii)=mstock(i, ii);
296 for (jj=0; jj<fNvar; jj++)
297 fCovariance(ii, jj)=cstock(ii, jj+i*fNvar);
302 det=CStep(fN, fH, index, fData, sscp, ndist);
303 if(TMath::Abs(det-deti[i])<kEps)
308 for(ii=0; ii<fNvar; ii++) {
309 mstock(i,ii)=fMean(ii);
310 for (jj=0; jj<fNvar; jj++)
311 cstock(ii,jj+i*fNvar)=fCovariance(ii, jj);
315 Int_t detind=TMath::LocMin(nbest, deti);
316 for(ii=0; ii<fNvar; ii++) {
317 fMean(ii)=mstock(detind,ii);
319 for(jj=0; jj<fNvar; jj++)
320 fCovariance(ii, jj)=cstock(ii,jj+detind*fNvar);
323 if (deti[detind]!=0) {
325 Int_t nout = RDist(sscp);
326 Double_t cutoff=kChiQuant[fNvar-1];
331 for (i=0; i<fN; i++) {
353 for (ii=0; ii<5; ii++)
356 nsub = Partition(nmini, indsubdat);
359 for (ii=0; ii<5; ii++)
361 Int_t *subdat=
new Int_t[sum];
364 for (
int iii = 0; iii < sum; ++iii) subdat[iii] = -999;
365 RDraw(subdat, nsub, indsubdat);
366 for (
int iii = 0; iii < sum; ++iii) {
367 if (subdat[iii] < 0 || subdat[iii] >= fN ) {
368 Error(
"Evaluate",
"subdat index is invalid subdat[%d] = %d",iii, subdat[iii] );
374 Int_t nbestsub=nbest*nsub;
375 TMatrixD mstockbig(nbestsub, fNvar);
376 TMatrixD cstockbig(fNvar, fNvar*nbestsub);
377 TMatrixD hyperplane(nbestsub, fNvar);
378 for (i=0; i<nbestsub; i++) {
379 for(j=0; j<fNvar; j++)
382 Double_t *detibig =
new Double_t[nbestsub];
384 maxind=TMath::LocMax(5, indsubdat);
385 TMatrixD dattemp(indsubdat[maxind], fNvar);
387 Int_t k2=Int_t(k1/nsub);
390 for (Int_t kgroup=0; kgroup<nsub; kgroup++) {
392 Int_t ntemp=indsubdat[kgroup];
394 for (i=0; i<kgroup; i++)
399 for(i=0; i<ntemp; i++) {
400 for (j=0; j<fNvar; j++) {
401 dattemp(i,j)=fData[subdat[temp+i]][j];
404 Int_t htemp=Int_t(fH*ntemp/fN);
406 for (i=0; i<nbest; i++)
409 for(k=0; k<k2; k++) {
410 CreateSubset(ntemp, htemp, fNvar, index, dattemp, sscp, ndist);
412 for (i=0; i<htemp; i++) {
413 for(j=0; j<fNvar; j++) {
414 vec(j)=dattemp(index[i],j);
416 AddToSscp(sscp, vec);
418 Covar(sscp, fMean, fCovariance, fSd, htemp);
419 det = fCovariance.Determinant();
421 par =Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
433 det = CStep(ntemp, htemp, index, dattemp, sscp, ndist);
435 par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp, ndist);
447 det=CStep(ntemp,htemp, index, dattemp, sscp, ndist);
449 par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
462 maxind=TMath::LocMax(nbest, deti);
463 if(det<deti[maxind]) {
465 for(i=0; i<fNvar; i++) {
466 mstockbig(nbest*kgroup+maxind,i)=fMean(i);
467 for(j=0; j<fNvar; j++) {
468 cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
478 maxind=TMath::LocMax(nbest, deti);
479 if (deti[maxind]<kEps)
484 for(i=0; i<nbest; i++) {
485 detibig[kgroup*nbest + i]=deti[i];
495 TMatrixD datmerged(sum, fNvar);
496 for(i=0; i<sum; i++) {
497 for (j=0; j<fNvar; j++)
498 datmerged(i,j)=fData[subdat[i]][j];
501 Int_t hmerged=Int_t(sum*fH/fN);
504 for(k=0; k<nbestsub; k++) {
506 for(ii=0; ii<fNvar; ii++) {
507 fMean(ii)=mstockbig(k,ii);
508 for(jj=0; jj<fNvar; jj++)
509 fCovariance(ii, jj)=cstockbig(ii,k*fNvar+jj);
512 for(i=0; i<fNvar; i++)
513 fHyperplane(i)=hyperplane(k,i);
514 CreateOrtSubset(datmerged,index, hmerged, sum, sscp, ndist);
517 det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
530 CreateOrtSubset(datmerged, index, hmerged, sum, sscp, ndist);
534 det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
548 for(i=0; i<fNvar; i++) {
549 mstockbig(k,i)=fMean(i);
550 for(j=0; j<fNvar; j++) {
551 cstockbig(i,k*fNvar+j)=fCovariance(i, j);
557 Int_t minind=TMath::LocMin(nbestsub, detibig);
559 for(i=0; i<fNvar; i++) {
560 fMean(i)=mstockbig(minind,i);
561 fHyperplane(i)=hyperplane(minind,i);
562 for(j=0; j<fNvar; j++)
563 fCovariance(i, j)=cstockbig(i,minind*fNvar + j);
566 CreateOrtSubset(fData, index, fH, fN, sscp, ndist);
569 det=CStep(fN, fH, index, fData, sscp, ndist);
570 if(TMath::Abs(det-detibig[minind])<kEps) {
580 Int_t nout = RDist(sscp);
581 Double_t cutoff=kChiQuant[fNvar-1];
586 for (i=0; i<fN; i++) {
608 void TRobustEstimator::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh)
612 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};
613 Int_t *index=
new Int_t[nvectors];
614 TMath::Sort(nvectors, data, index, kFALSE);
617 nquant=TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
618 Double_t factor=faclts[nquant-1];
620 Double_t *aw=
new Double_t[nvectors];
621 Double_t *aw2=
new Double_t[nvectors];
625 Int_t len=nvectors-hh;
626 Double_t *slutn=
new Double_t[len];
627 for(Int_t i=0; i<len; i++)
629 for(Int_t jint=0; jint<len; jint++) {
631 for (Int_t j=0; j<hh; j++) {
632 aw[jint]+=data[index[j+jint]];
634 sq+=data[index[j]]*data[index[j]];
636 aw2[jint]=aw[jint]*aw[jint]/hh;
641 slutn[ndup]=aw[jint];
644 sq=sq - data[index[jint-1]]*data[index[jint-1]]+
645 data[index[jint+hh]]*data[index[jint+hh]]-
646 aw2[jint]+aw2[jint-1];
650 slutn[ndup]=aw[jint];
655 slutn[ndup]=aw[jint];
661 slutn[0]=slutn[Int_t((ndup)/2)]/hh;
662 Double_t bstd=factor*TMath::Sqrt(sqmin/hh);
674 Int_t TRobustEstimator::GetBDPoint()
684 Double_t TRobustEstimator::GetChiQuant(Int_t i)
const
686 if (i < 0 || i >= 50)
return 0;
693 void TRobustEstimator::GetCovariance(TMatrixDSym &matr)
695 if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar){
696 Warning(
"GetCovariance",
"provided matrix is of the wrong size, it will be resized");
697 matr.ResizeTo(fNvar, fNvar);
705 void TRobustEstimator::GetCorrelation(TMatrixDSym &matr)
707 if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar) {
708 Warning(
"GetCorrelation",
"provided matrix is of the wrong size, it will be resized");
709 matr.ResizeTo(fNvar, fNvar);
717 const TVectorD* TRobustEstimator::GetHyperplane()
const
720 Error(
"GetHyperplane",
"the data doesn't lie on a hyperplane!\n");
730 void TRobustEstimator::GetHyperplane(TVectorD &vec)
733 Error(
"GetHyperplane",
"the data doesn't lie on a hyperplane!\n");
736 if (vec.GetNoElements()!=fNvar) {
737 Warning(
"GetHyperPlane",
"provided vector is of the wrong size, it will be resized");
746 void TRobustEstimator::GetMean(TVectorD &means)
748 if (means.GetNoElements()!=fNvar) {
749 Warning(
"GetMean",
"provided vector is of the wrong size, it will be resized");
750 means.ResizeTo(fNvar);
758 void TRobustEstimator::GetRDistances(TVectorD &rdist)
760 if (rdist.GetNoElements()!=fN) {
761 Warning(
"GetRDistances",
"provided vector is of the wrong size, it will be resized");
770 Int_t TRobustEstimator::GetNOut()
772 return fOut.GetSize();
778 void TRobustEstimator::AddToSscp(TMatrixD &sscp, TVectorD &vec)
781 for (j=1; j<fNvar+1; j++) {
782 sscp(0, j) +=vec(j-1);
783 sscp(j, 0) = sscp(0, j);
785 for (i=1; i<fNvar+1; i++) {
786 for (j=1; j<fNvar+1; j++) {
787 sscp(i, j) += vec(i-1)*vec(j-1);
795 void TRobustEstimator::ClearSscp(TMatrixD &sscp)
797 for (Int_t i=0; i<fNvar+1; i++) {
798 for (Int_t j=0; j<fNvar+1; j++) {
808 void TRobustEstimator::Classic()
810 TMatrixD sscp(fNvar+1, fNvar+1);
811 TVectorD temp(fNvar);
813 for (Int_t i=0; i<fN; i++) {
814 for (Int_t j=0; j<fNvar; j++)
816 AddToSscp(sscp, temp);
818 Covar(sscp, fMean, fCovariance, fSd, fN);
826 void TRobustEstimator::Covar(TMatrixD &sscp, TVectorD &m, TMatrixDSym &cov, TVectorD &sd, Int_t nvec)
830 for (i=0; i<fNvar; i++) {
832 sd[i]=sscp(i+1, i+1);
833 f=(sd[i]-m(i)*m(i)/nvec)/(nvec-1);
834 if (f>1e-14) sd[i]=TMath::Sqrt(f);
838 for (i=0; i<fNvar; i++) {
839 for (j=0; j<fNvar; j++) {
840 cov(i, j)=sscp(i+1, j+1)-nvec*m(i)*m(j);
849 void TRobustEstimator::Correl()
852 Double_t *sd=
new Double_t[fNvar];
853 for(j=0; j<fNvar; j++)
854 sd[j]=1./TMath::Sqrt(fCovariance(j, j));
855 for(i=0; i<fNvar; i++) {
856 for (j=0; j<fNvar; j++) {
858 fCorrelation(i, j)=1.;
860 fCorrelation(i, j)=fCovariance(i, j)*sd[i]*sd[j];
877 void TRobustEstimator::CreateSubset(Int_t ntotal, Int_t htotal, Int_t p, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
879 Double_t kEps = 1e-14;
881 Bool_t repeat=kFALSE;
884 for(i=0; i<ntotal; i++)
887 for (i=0; i<p+1; i++) {
888 num=Int_t(gRandom->Uniform(0, 1)*(ntotal-1));
890 for(j=0; j<=i-1; j++) {
908 for (i=0; i<p+1; i++) {
909 for (j=0; j<fNvar; j++) {
910 vec[j]=data[index[i]][j];
913 AddToSscp(sscp, vec);
916 Covar(sscp, fMean, fCovariance, fSd, p+1);
917 det=fCovariance.Determinant();
918 while((det<kEps)&&(nindex < htotal)) {
924 num=Int_t(gRandom->Uniform(0,1)*(ntotal-1));
926 for(i=0; i<nindex; i++) {
932 }
while(repeat==kTRUE);
937 for(j=0; j<fNvar; j++)
938 vec[j]=data[index[nindex-1]][j];
939 AddToSscp(sscp, vec);
940 Covar(sscp, fMean, fCovariance, fSd, nindex);
941 det=fCovariance.Determinant();
945 TDecompChol chol(fCovariance);
946 fInvcovariance = chol.Invert();
948 TVectorD temp(fNvar);
949 for(j=0; j<ntotal; j++) {
951 for(i=0; i<fNvar; i++)
952 temp[i]=data[j][i] - fMean(i);
953 temp*=fInvcovariance;
954 for(i=0; i<fNvar; i++)
955 ndist[j]+=(data[j][i]-fMean(i))*temp[i];
957 KOrdStat(ntotal, ndist, htotal-1,index);
967 void TRobustEstimator::CreateOrtSubset(TMatrixD &dat,Int_t *index, Int_t hmerged, Int_t nmerged, TMatrixD &sscp, Double_t *ndist)
971 for (i=0; i<nmerged; i++) {
973 for(j=0; j<fNvar; j++) {
974 ndist[i]+=fHyperplane[j]*(dat[i][j]-fMean[j]);
975 ndist[i]=TMath::Abs(ndist[i]);
978 KOrdStat(nmerged, ndist, hmerged-1, index);
980 for (i=0; i<hmerged; i++) {
981 for(j=0; j<fNvar; j++)
982 vec[j]=dat[index[i]][j];
983 AddToSscp(sscp, vec);
985 Covar(sscp, fMean, fCovariance, fSd, hmerged);
999 Double_t TRobustEstimator::CStep(Int_t ntotal, Int_t htotal, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
1002 TVectorD vec(fNvar);
1005 TDecompChol chol(fCovariance);
1006 fInvcovariance = chol.Invert();
1008 TVectorD temp(fNvar);
1009 for(j=0; j<ntotal; j++) {
1011 for(i=0; i<fNvar; i++)
1012 temp[i]=data[j][i]-fMean[i];
1013 temp*=fInvcovariance;
1014 for(i=0; i<fNvar; i++)
1015 ndist[j]+=(data[j][i]-fMean[i])*temp[i];
1019 KOrdStat(ntotal, ndist, htotal-1, index);
1022 for (i=0; i<htotal; i++) {
1023 for (j=0; j<fNvar; j++)
1024 temp[j]=data[index[i]][j];
1025 AddToSscp(sscp, temp);
1027 Covar(sscp, fMean, fCovariance, fSd, htotal);
1028 det = fCovariance.Determinant();
1036 Int_t TRobustEstimator::Exact(Double_t *ndist)
1040 TMatrixDSymEigen eigen(fCovariance);
1041 TVectorD eigenValues=eigen.GetEigenValues();
1042 TMatrixD eigenMatrix=eigen.GetEigenVectors();
1044 for (j=0; j<fNvar; j++) {
1045 fHyperplane[j]=eigenMatrix(j,fNvar-1);
1048 for (i=0; i<fN; i++) {
1050 for(j=0; j<fNvar; j++) {
1051 ndist[i]+=fHyperplane[j]*(fData[i][j]-fMean[j]);
1052 ndist[i]=TMath::Abs(ndist[i]);
1057 for (i=0; i<fN; i++) {
1058 if(ndist[i] < 1e-14) nhyp++;
1071 Int_t TRobustEstimator::Exact2(TMatrixD &mstockbig, TMatrixD &cstockbig, TMatrixD &hyperplane,
1072 Double_t *deti, Int_t nbest, Int_t kgroup,
1073 TMatrixD &sscp, Double_t *ndist)
1077 TVectorD vec(fNvar);
1078 Int_t maxind = TMath::LocMax(nbest, deti);
1079 Int_t nh=Exact(ndist);
1084 for (i=0; i<fN; i++) {
1085 if(ndist[i]<1e-14) {
1086 for (j=0; j<fNvar; j++)
1088 AddToSscp(sscp, vec);
1091 Covar(sscp, fMean, fCovariance, fSd, nh);
1101 for(i=0; i<fNvar; i++) {
1102 mstockbig(nbest*kgroup+maxind,i)=fMean(i);
1103 hyperplane(nbest*kgroup+maxind,i)=fHyperplane(i);
1104 for(j=0; j<fNvar; j++) {
1105 cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
1118 Int_t TRobustEstimator::Partition(Int_t nmini, Int_t *indsubdat)
1121 if ((fN>=2*nmini) && (fN<=(3*nmini-1))) {
1123 indsubdat[0]=Int_t(fN*0.5);
1124 indsubdat[1]=Int_t(fN*0.5)+1;
1126 indsubdat[0]=indsubdat[1]=Int_t(fN/2);
1130 if((fN>=3*nmini) && (fN<(4*nmini -1))) {
1132 indsubdat[0]=indsubdat[1]=indsubdat[2]=Int_t(fN/3);
1134 indsubdat[0]=Int_t(fN/3);
1135 indsubdat[1]=Int_t(fN/3)+1;
1136 if (fN%3==1) indsubdat[2]=Int_t(fN/3);
1137 else indsubdat[2]=Int_t(fN/3)+1;
1142 if((fN>=4*nmini)&&(fN<=(5*nmini-1))){
1143 if (fN%4==0) indsubdat[0]=indsubdat[1]=indsubdat[2]=indsubdat[3]=Int_t(fN/4);
1145 indsubdat[0]=Int_t(fN/4);
1146 indsubdat[1]=Int_t(fN/4)+1;
1147 if(fN%4==1) indsubdat[2]=indsubdat[3]=Int_t(fN/4);
1149 indsubdat[2]=Int_t(fN/4)+1;
1150 indsubdat[3]=Int_t(fN/4);
1152 if(fN%4==3) indsubdat[2]=indsubdat[3]=Int_t(fN/4)+1;
1156 for(Int_t i=0; i<5; i++)
1172 Int_t TRobustEstimator::RDist(TMatrixD &sscp)
1177 TVectorD temp(fNvar);
1178 TDecompChol chol(fCovariance);
1179 fInvcovariance = chol.Invert();
1182 for (i=0; i<fN; i++) {
1184 for(j=0; j<fNvar; j++) {
1185 temp[j]=fData[i][j]-fMean[j];
1187 temp*=fInvcovariance;
1188 for(j=0; j<fNvar; j++) {
1189 fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
1194 Double_t chi = kChiMedian[fNvar-1];
1196 med=TMath::Median(fN, fRd.GetMatrixArray());
1199 TDecompChol chol2(fCovariance);
1200 fInvcovariance = chol2.Invert();
1202 for (i=0; i<fN; i++) {
1204 for(j=0; j<fNvar; j++) {
1205 temp[j]=fData[i][j]-fMean[j];
1208 temp*=fInvcovariance;
1209 for(j=0; j<fNvar; j++) {
1210 fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
1214 Double_t cutoff = kChiQuant[fNvar-1];
1217 for(i=0; i<fN; i++) {
1218 if (fRd[i]<=cutoff) {
1219 for(j=0; j<fNvar; j++)
1220 temp[j]=fData[i][j];
1221 AddToSscp(sscp,temp);
1227 Covar(sscp, fMean, fCovariance, fSd, fN-nout);
1235 void TRobustEstimator::RDraw(Int_t *subdat, Int_t ngroup, Int_t *indsubdat)
1240 for (k=1; k<=ngroup; k++) {
1241 for (m=1; m<=indsubdat[k-1]; m++) {
1242 nrand = Int_t(gRandom->Uniform(0, 1) * double(fN-jndex))+1;
1248 subdat[jndex-1]=nrand+jndex-2;
1249 for (i=1; i<=jndex-1; i++) {
1250 if(subdat[i-1] > nrand+i-2) {
1251 for(j=jndex; j>=i+1; j--) {
1252 subdat[j-1]=subdat[j-2];
1255 subdat[i-1]=nrand+i-2;
1267 Double_t TRobustEstimator::KOrdStat(Int_t ntotal, Double_t *a, Int_t k, Int_t *work){
1268 Bool_t isAllocated = kFALSE;
1269 const Int_t kWorkMax=100;
1270 Int_t i, ir, j, l, mid;
1273 Int_t workLocal[kWorkMax];
1281 if (ntotal > kWorkMax) {
1282 isAllocated = kTRUE;
1283 ind =
new Int_t[ntotal];
1287 for (Int_t ii=0; ii<ntotal; ii++) {
1295 if (ir == l+1 && a[ind[ir]]<a[ind[l]])
1296 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1297 Double_t tmp = a[ind[rk]];
1303 {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}
1304 if (a[ind[l]]>a[ind[ir]])
1305 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
1307 if (a[ind[l+1]]>a[ind[ir]])
1308 {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}
1310 if (a[ind[l]]>a[ind[l+1]])
1311 {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}
1317 do i++;
while (a[ind[i]]<a[arr]);
1318 do j--;
while (a[ind[j]]>a[arr]);
1320 {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
1324 if (j>=rk) ir = j-1;