00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 #include "TRobustEstimator.h"
00102 #include "TRandom.h"
00103 #include "TMath.h"
00104 #include "TH1D.h"
00105 #include "TPaveLabel.h"
00106 #include "TDecompChol.h"
00107
00108 ClassImp(TRobustEstimator)
00109
00110 const Double_t kChiMedian[50]= {
00111 0.454937, 1.38629, 2.36597, 3.35670, 4.35146, 5.34812, 6.34581, 7.34412, 8.34283,
00112 9.34182, 10.34, 11.34, 12.34, 13.34, 14.34, 15.34, 16.34, 17.34, 18.34, 19.34,
00113 20.34, 21.34, 22.34, 23.34, 24.34, 25.34, 26.34, 27.34, 28.34, 29.34, 30.34,
00114 31.34, 32.34, 33.34, 34.34, 35.34, 36.34, 37.34, 38.34, 39.34, 40.34,
00115 41.34, 42.34, 43.34, 44.34, 45.34, 46.34, 47.34, 48.34, 49.33};
00116
00117 const Double_t kChiQuant[50]={
00118 5.02389, 7.3776,9.34840,11.1433,12.8325,
00119 14.4494,16.0128,17.5346,19.0228,20.4831,21.920,23.337,
00120 24.736,26.119,27.488,28.845,30.191,31.526,32.852,34.170,
00121 35.479,36.781,38.076,39.364,40.646,41.923,43.194,44.461,
00122 45.722,46.979,48.232,49.481,50.725,51.966,53.203,54.437,
00123 55.668,56.896,58.120,59.342,60.561,61.777,62.990,64.201,
00124 65.410,66.617,67.821,69.022,70.222,71.420};
00125
00126
00127 TRobustEstimator::TRobustEstimator(){
00128
00129
00130
00131 }
00132
00133
00134 TRobustEstimator::TRobustEstimator(Int_t nvectors, Int_t nvariables, Int_t hh)
00135 :fMean(nvariables),
00136 fCovariance(nvariables),
00137 fInvcovariance(nvariables),
00138 fCorrelation(nvariables),
00139 fRd(nvectors),
00140 fSd(nvectors),
00141 fOut(1),
00142 fHyperplane(nvariables),
00143 fData(nvectors, nvariables)
00144 {
00145
00146
00147 if ((nvectors<=1)||(nvariables<=0)){
00148 Error("TRobustEstimator","Not enough vectors or variables");
00149 return;
00150 }
00151 if (nvariables==1){
00152 Error("TRobustEstimator","For the univariate case, use the default constructor and EvaluateUni() function");
00153 return;
00154 }
00155
00156 fN=nvectors;
00157 fNvar=nvariables;
00158 if (hh<(fN+fNvar+1)/2){
00159 if (hh>0)
00160 Warning("TRobustEstimator","chosen h is too small, default h is taken instead");
00161 fH=(fN+fNvar+1)/2;
00162 } else
00163 fH=hh;
00164
00165 fVarTemp=0;
00166 fVecTemp=0;
00167 fExact=0;
00168 }
00169
00170
00171 void TRobustEstimator::AddColumn(Double_t *col)
00172 {
00173
00174
00175
00176
00177
00178 if (fVarTemp==fNvar) {
00179 fNvar++;
00180 fCovariance.ResizeTo(fNvar, fNvar);
00181 fInvcovariance.ResizeTo(fNvar, fNvar);
00182 fCorrelation.ResizeTo(fNvar, fNvar);
00183 fMean.ResizeTo(fNvar);
00184 fHyperplane.ResizeTo(fNvar);
00185 fData.ResizeTo(fN, fNvar);
00186 }
00187 for (Int_t i=0; i<fN; i++) {
00188 fData(i, fVarTemp)=col[i];
00189 }
00190 fVarTemp++;
00191 }
00192
00193
00194 void TRobustEstimator::AddRow(Double_t *row)
00195 {
00196
00197
00198
00199 if(fVecTemp==fN) {
00200 fN++;
00201 fRd.ResizeTo(fN);
00202 fSd.ResizeTo(fN);
00203 fData.ResizeTo(fN, fNvar);
00204 }
00205 for (Int_t i=0; i<fNvar; i++)
00206 fData(fVecTemp, i)=row[i];
00207
00208 fVecTemp++;
00209 }
00210
00211
00212 void TRobustEstimator::Evaluate()
00213 {
00214
00215
00216 Double_t kEps=1e-14;
00217
00218 if (fH==fN){
00219 Warning("Evaluate","Chosen h = #observations, so classic estimates of location and scatter will be calculated");
00220 Classic();
00221 return;
00222 }
00223
00224 Int_t i, j, k;
00225 Int_t ii, jj;
00226 Int_t nmini = 300;
00227 Int_t k1=500;
00228 Int_t nbest=10;
00229 TMatrixD sscp(fNvar+1, fNvar+1);
00230 TVectorD vec(fNvar);
00231
00232 Int_t *index = new Int_t[fN];
00233 Double_t *ndist = new Double_t[fN];
00234 Double_t det;
00235 Double_t *deti=new Double_t[nbest];
00236 for (i=0; i<nbest; i++)
00237 deti[i]=1e16;
00238
00239 for (i=0; i<fN; i++)
00240 fRd(i)=0;
00241
00242
00243
00244 if (fN<nmini*2) {
00245
00246
00247 TMatrixD mstock(nbest, fNvar);
00248 TMatrixD cstock(fNvar, fNvar*nbest);
00249
00250 for (k=0; k<k1; k++) {
00251 CreateSubset(fN, fH, fNvar, index, fData, sscp, ndist);
00252
00253 ClearSscp(sscp);
00254 for (i=0; i<fH; i++) {
00255 for(j=0; j<fNvar; j++)
00256 vec(j)=fData[index[i]][j];
00257 AddToSscp(sscp, vec);
00258 }
00259 Covar(sscp, fMean, fCovariance, fSd, fH);
00260 det = fCovariance.Determinant();
00261 if (det < kEps) {
00262 fExact = Exact(ndist);
00263 delete [] index;
00264 delete [] ndist;
00265 delete [] deti;
00266 return;
00267 }
00268
00269 det = CStep(fN, fH, index, fData, sscp, ndist);
00270 if (det < kEps) {
00271 fExact = Exact(ndist);
00272 delete [] index;
00273 delete [] ndist;
00274 delete [] deti;
00275 return;
00276 }
00277 det = CStep(fN, fH, index, fData, sscp, ndist);
00278 if (det < kEps) {
00279 fExact = Exact(ndist);
00280 delete [] index;
00281 delete [] ndist;
00282 delete [] deti;
00283 return;
00284 } else {
00285 Int_t maxind=TMath::LocMax(nbest, deti);
00286 if(det<deti[maxind]) {
00287 deti[maxind]=det;
00288 for(ii=0; ii<fNvar; ii++) {
00289 mstock(maxind, ii)=fMean(ii);
00290 for(jj=0; jj<fNvar; jj++)
00291 cstock(ii, jj+maxind*fNvar)=fCovariance(ii, jj);
00292 }
00293 }
00294 }
00295 }
00296
00297
00298
00299 for (i=0; i<nbest; i++) {
00300 for(ii=0; ii<fNvar; ii++) {
00301 fMean(ii)=mstock(i, ii);
00302 for (jj=0; jj<fNvar; jj++)
00303 fCovariance(ii, jj)=cstock(ii, jj+i*fNvar);
00304 }
00305
00306 det=1;
00307 while (det>kEps) {
00308 det=CStep(fN, fH, index, fData, sscp, ndist);
00309 if(TMath::Abs(det-deti[i])<kEps)
00310 break;
00311 else
00312 deti[i]=det;
00313 }
00314 for(ii=0; ii<fNvar; ii++) {
00315 mstock(i,ii)=fMean(ii);
00316 for (jj=0; jj<fNvar; jj++)
00317 cstock(ii,jj+i*fNvar)=fCovariance(ii, jj);
00318 }
00319 }
00320
00321 Int_t detind=TMath::LocMin(nbest, deti);
00322 for(ii=0; ii<fNvar; ii++) {
00323 fMean(ii)=mstock(detind,ii);
00324
00325 for(jj=0; jj<fNvar; jj++)
00326 fCovariance(ii, jj)=cstock(ii,jj+detind*fNvar);
00327 }
00328
00329 if (deti[detind]!=0) {
00330
00331 Int_t nout = RDist(sscp);
00332 Double_t cutoff=kChiQuant[fNvar-1];
00333
00334 fOut.Set(nout);
00335
00336 j=0;
00337 for (i=0; i<fN; i++) {
00338 if(fRd(i)>cutoff) {
00339 fOut[j]=i;
00340 j++;
00341 }
00342 }
00343
00344 } else {
00345 fExact=Exact(ndist);
00346 }
00347 delete [] index;
00348 delete [] ndist;
00349 delete [] deti;
00350 return;
00351
00352 }
00353
00354
00355
00356
00357 Int_t indsubdat[5];
00358 Int_t nsub;
00359 for (ii=0; ii<5; ii++)
00360 indsubdat[ii]=0;
00361
00362 nsub = Partition(nmini, indsubdat);
00363
00364 Int_t sum=0;
00365 for (ii=0; ii<5; ii++)
00366 sum+=indsubdat[ii];
00367 Int_t *subdat=new Int_t[sum];
00368 RDraw(subdat, nsub, indsubdat);
00369
00370
00371 Int_t nbestsub=nbest*nsub;
00372 TMatrixD mstockbig(nbestsub, fNvar);
00373 TMatrixD cstockbig(fNvar, fNvar*nbestsub);
00374 TMatrixD hyperplane(nbestsub, fNvar);
00375 for (i=0; i<nbestsub; i++) {
00376 for(j=0; j<fNvar; j++)
00377 hyperplane(i,j)=0;
00378 }
00379 Double_t *detibig = new Double_t[nbestsub];
00380 Int_t maxind;
00381 maxind=TMath::LocMax(5, indsubdat);
00382 TMatrixD dattemp(indsubdat[maxind], fNvar);
00383
00384 Int_t k2=Int_t(k1/nsub);
00385
00386
00387 for (Int_t kgroup=0; kgroup<nsub; kgroup++) {
00388
00389 Int_t ntemp=indsubdat[kgroup];
00390 Int_t temp=0;
00391 for (i=0; i<kgroup; i++)
00392 temp+=indsubdat[i];
00393 Int_t par;
00394
00395 for(i=0; i<ntemp; i++) {
00396 for (j=0; j<fNvar; j++) {
00397 dattemp(i,j)=fData[subdat[temp+i]][j];
00398 }
00399 }
00400 Int_t htemp=Int_t(fH*ntemp/fN);
00401
00402 for (i=0; i<nbest; i++)
00403 deti[i]=1e16;
00404
00405 for(k=0; k<k2; k++) {
00406 CreateSubset(ntemp, htemp, fNvar, index, dattemp, sscp, ndist);
00407 ClearSscp(sscp);
00408 for (i=0; i<htemp; i++) {
00409 for(j=0; j<fNvar; j++) {
00410 vec(j)=dattemp(index[i],j);
00411 }
00412 AddToSscp(sscp, vec);
00413 }
00414 Covar(sscp, fMean, fCovariance, fSd, htemp);
00415 det = fCovariance.Determinant();
00416 if (det<kEps) {
00417 par =Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
00418 if(par==nbest+1) {
00419
00420 delete [] detibig;
00421 delete [] deti;
00422 delete [] subdat;
00423 delete [] ndist;
00424 delete [] index;
00425 return;
00426 } else
00427 deti[par]=det;
00428 } else {
00429 det = CStep(ntemp, htemp, index, dattemp, sscp, ndist);
00430 if (det<kEps) {
00431 par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp, ndist);
00432 if(par==nbest+1) {
00433
00434 delete [] detibig;
00435 delete [] deti;
00436 delete [] subdat;
00437 delete [] ndist;
00438 delete [] index;
00439 return;
00440 } else
00441 deti[par]=det;
00442 } else {
00443 det=CStep(ntemp,htemp, index, dattemp, sscp, ndist);
00444 if(det<kEps){
00445 par=Exact2(mstockbig, cstockbig, hyperplane, deti, nbest, kgroup, sscp,ndist);
00446 if(par==nbest+1) {
00447
00448 delete [] detibig;
00449 delete [] deti;
00450 delete [] subdat;
00451 delete [] ndist;
00452 delete [] index;
00453 return;
00454 } else {
00455 deti[par]=det;
00456 }
00457 } else {
00458 maxind=TMath::LocMax(nbest, deti);
00459 if(det<deti[maxind]) {
00460 deti[maxind]=det;
00461 for(i=0; i<fNvar; i++) {
00462 mstockbig(nbest*kgroup+maxind,i)=fMean(i);
00463 for(j=0; j<fNvar; j++) {
00464 cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
00465
00466 }
00467 }
00468 }
00469
00470 }
00471 }
00472 }
00473
00474 maxind=TMath::LocMax(nbest, deti);
00475 if (deti[maxind]<kEps)
00476 break;
00477 }
00478
00479
00480 for(i=0; i<nbest; i++) {
00481 detibig[kgroup*nbest + i]=deti[i];
00482
00483 }
00484
00485 }
00486
00487
00488
00489
00490
00491 TMatrixD datmerged(sum, fNvar);
00492 for(i=0; i<sum; i++) {
00493 for (j=0; j<fNvar; j++)
00494 datmerged(i,j)=fData[subdat[i]][j];
00495 }
00496
00497 Int_t hmerged=Int_t(sum*fH/fN);
00498
00499 Int_t nh;
00500 for(k=0; k<nbestsub; k++) {
00501
00502 for(ii=0; ii<fNvar; ii++) {
00503 fMean(ii)=mstockbig(k,ii);
00504 for(jj=0; jj<fNvar; jj++)
00505 fCovariance(ii, jj)=cstockbig(ii,k*fNvar+jj);
00506 }
00507 if(detibig[k]==0) {
00508 for(i=0; i<fNvar; i++)
00509 fHyperplane(i)=hyperplane(k,i);
00510 CreateOrtSubset(datmerged,index, hmerged, sum, sscp, ndist);
00511
00512 }
00513 det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
00514 if (det<kEps) {
00515 nh= Exact(ndist);
00516 if (nh>=fH) {
00517 fExact = nh;
00518
00519 delete [] detibig;
00520 delete [] deti;
00521 delete [] subdat;
00522 delete [] ndist;
00523 delete [] index;
00524 return;
00525 } else {
00526 CreateOrtSubset(datmerged, index, hmerged, sum, sscp, ndist);
00527 }
00528 }
00529
00530 det=CStep(sum, hmerged, index, datmerged, sscp, ndist);
00531 if (det<kEps) {
00532 nh=Exact(ndist);
00533 if (nh>=fH) {
00534 fExact = nh;
00535 delete [] detibig;
00536 delete [] deti;
00537 delete [] subdat;
00538 delete [] ndist;
00539 delete [] index;
00540 return;
00541 }
00542 }
00543 detibig[k]=det;
00544 for(i=0; i<fNvar; i++) {
00545 mstockbig(k,i)=fMean(i);
00546 for(j=0; j<fNvar; j++) {
00547 cstockbig(i,k*fNvar+j)=fCovariance(i, j);
00548 }
00549 }
00550 }
00551
00552
00553 Int_t minind=TMath::LocMin(nbestsub, detibig);
00554 det=detibig[minind];
00555 for(i=0; i<fNvar; i++) {
00556 fMean(i)=mstockbig(minind,i);
00557 fHyperplane(i)=hyperplane(minind,i);
00558 for(j=0; j<fNvar; j++)
00559 fCovariance(i, j)=cstockbig(i,minind*fNvar + j);
00560 }
00561 if(det<kEps)
00562 CreateOrtSubset(fData, index, fH, fN, sscp, ndist);
00563 det=1;
00564 while (det>kEps) {
00565 det=CStep(fN, fH, index, fData, sscp, ndist);
00566 if(TMath::Abs(det-detibig[minind])<kEps) {
00567 break;
00568 } else {
00569 detibig[minind]=det;
00570 }
00571 }
00572 if(det<kEps) {
00573 Exact(ndist);
00574 fExact=kTRUE;
00575 }
00576 Int_t nout = RDist(sscp);
00577 Double_t cutoff=kChiQuant[fNvar-1];
00578
00579 fOut.Set(nout);
00580
00581 j=0;
00582 for (i=0; i<fN; i++) {
00583 if(fRd(i)>cutoff) {
00584 fOut[j]=i;
00585 j++;
00586 }
00587 }
00588
00589 delete [] detibig;
00590 delete [] deti;
00591 delete [] subdat;
00592 delete [] ndist;
00593 delete [] index;
00594 return;
00595 }
00596
00597
00598 void TRobustEstimator::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh)
00599 {
00600
00601
00602
00603
00604
00605
00606 if (hh==0)
00607 hh=(nvectors+2)/2;
00608 Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473};
00609 Int_t *index=new Int_t[nvectors];
00610 TMath::Sort(nvectors, data, index, kFALSE);
00611
00612 Int_t nquant;
00613 nquant=TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11);
00614 Double_t factor=faclts[nquant-1];
00615
00616 Double_t *aw=new Double_t[nvectors];
00617 Double_t *aw2=new Double_t[nvectors];
00618 Double_t sq=0;
00619 Double_t sqmin=0;
00620 Int_t ndup=0;
00621 Int_t len=nvectors-hh;
00622 Double_t *slutn=new Double_t[len];
00623 for(Int_t i=0; i<len; i++)
00624 slutn[i]=0;
00625 for(Int_t jint=0; jint<len; jint++) {
00626 aw[jint]=0;
00627 for (Int_t j=0; j<hh; j++) {
00628 aw[jint]+=data[index[j+jint]];
00629 if(jint==0)
00630 sq+=data[index[j]]*data[index[j]];
00631 }
00632 aw2[jint]=aw[jint]*aw[jint]/hh;
00633
00634 if(jint==0) {
00635 sq=sq-aw2[jint];
00636 sqmin=sq;
00637 slutn[ndup]=aw[jint];
00638
00639 } else {
00640 sq=sq - data[index[jint-1]]*data[index[jint-1]]+
00641 data[index[jint+hh]]*data[index[jint+hh]]-
00642 aw2[jint]+aw2[jint-1];
00643 if(sq<sqmin) {
00644 ndup=0;
00645 sqmin=sq;
00646 slutn[ndup]=aw[jint];
00647
00648 } else {
00649 if(sq==sqmin) {
00650 ndup++;
00651 slutn[ndup]=aw[jint];
00652 }
00653 }
00654 }
00655 }
00656
00657 slutn[0]=slutn[Int_t((ndup)/2)]/hh;
00658 Double_t bstd=factor*TMath::Sqrt(sqmin/hh);
00659 mean=slutn[0];
00660 sigma=bstd;
00661 delete [] aw;
00662 delete [] aw2;
00663 delete [] slutn;
00664 delete [] index;
00665 }
00666
00667
00668 Int_t TRobustEstimator::GetBDPoint()
00669 {
00670
00671
00672 Int_t n;
00673 n=(fN-fH+1)/fN;
00674 return n;
00675 }
00676
00677
00678 Double_t TRobustEstimator::GetChiQuant(Int_t i) const
00679 {
00680
00681
00682 if (i < 0 || i >= 50) return 0;
00683 return kChiQuant[i];
00684 }
00685
00686
00687 void TRobustEstimator::GetCovariance(TMatrixDSym &matr)
00688 {
00689
00690
00691 if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar){
00692 Warning("GetCovariance","provided matrix is of the wrong size, it will be resized");
00693 matr.ResizeTo(fNvar, fNvar);
00694 }
00695 matr=fCovariance;
00696 }
00697
00698
00699 void TRobustEstimator::GetCorrelation(TMatrixDSym &matr)
00700 {
00701
00702
00703 if (matr.GetNrows()!=fNvar || matr.GetNcols()!=fNvar) {
00704 Warning("GetCorrelation","provided matrix is of the wrong size, it will be resized");
00705 matr.ResizeTo(fNvar, fNvar);
00706 }
00707 matr=fCorrelation;
00708 }
00709
00710
00711 const TVectorD* TRobustEstimator::GetHyperplane() const
00712 {
00713
00714
00715 if (fExact==0) {
00716 Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
00717 return 0;
00718 } else {
00719 return &fHyperplane;
00720 }
00721 }
00722
00723
00724 void TRobustEstimator::GetHyperplane(TVectorD &vec)
00725 {
00726
00727
00728 if (fExact==0){
00729 Error("GetHyperplane","the data doesn't lie on a hyperplane!\n");
00730 return;
00731 }
00732 if (vec.GetNoElements()!=fNvar) {
00733 Warning("GetHyperPlane","provided vector is of the wrong size, it will be resized");
00734 vec.ResizeTo(fNvar);
00735 }
00736 vec=fHyperplane;
00737 }
00738
00739
00740 void TRobustEstimator::GetMean(TVectorD &means)
00741 {
00742
00743
00744 if (means.GetNoElements()!=fNvar) {
00745 Warning("GetMean","provided vector is of the wrong size, it will be resized");
00746 means.ResizeTo(fNvar);
00747 }
00748 means=fMean;
00749 }
00750
00751
00752 void TRobustEstimator::GetRDistances(TVectorD &rdist)
00753 {
00754
00755
00756 if (rdist.GetNoElements()!=fN) {
00757 Warning("GetRDistances","provided vector is of the wrong size, it will be resized");
00758 rdist.ResizeTo(fN);
00759 }
00760 rdist=fRd;
00761 }
00762
00763
00764 Int_t TRobustEstimator::GetNOut()
00765 {
00766
00767
00768 return fOut.GetSize();
00769 }
00770
00771
00772 void TRobustEstimator::AddToSscp(TMatrixD &sscp, TVectorD &vec)
00773 {
00774
00775
00776 Int_t i, j;
00777 for (j=1; j<fNvar+1; j++) {
00778 sscp(0, j) +=vec(j-1);
00779 sscp(j, 0) = sscp(0, j);
00780 }
00781 for (i=1; i<fNvar+1; i++) {
00782 for (j=1; j<fNvar+1; j++) {
00783 sscp(i, j) += vec(i-1)*vec(j-1);
00784 }
00785 }
00786 }
00787
00788
00789 void TRobustEstimator::ClearSscp(TMatrixD &sscp)
00790 {
00791
00792
00793 for (Int_t i=0; i<fNvar+1; i++) {
00794 for (Int_t j=0; j<fNvar+1; j++) {
00795 sscp(i, j)=0;
00796 }
00797 }
00798 }
00799
00800
00801 void TRobustEstimator::Classic()
00802 {
00803
00804
00805 TMatrixD sscp(fNvar+1, fNvar+1);
00806 TVectorD temp(fNvar);
00807 ClearSscp(sscp);
00808 for (Int_t i=0; i<fN; i++) {
00809 for (Int_t j=0; j<fNvar; j++)
00810 temp(j)=fData(i, j);
00811 AddToSscp(sscp, temp);
00812 }
00813 Covar(sscp, fMean, fCovariance, fSd, fN);
00814 Correl();
00815
00816 }
00817
00818
00819 void TRobustEstimator::Covar(TMatrixD &sscp, TVectorD &m, TMatrixDSym &cov, TVectorD &sd, Int_t nvec)
00820 {
00821
00822
00823 Int_t i, j;
00824 Double_t f;
00825 for (i=0; i<fNvar; i++) {
00826 m(i)=sscp(0, i+1);
00827 sd[i]=sscp(i+1, i+1);
00828 f=(sd[i]-m(i)*m(i)/nvec)/(nvec-1);
00829 if (f>1e-14) sd[i]=TMath::Sqrt(f);
00830 else sd[i]=0;
00831 m(i)/=nvec;
00832 }
00833 for (i=0; i<fNvar; i++) {
00834 for (j=0; j<fNvar; j++) {
00835 cov(i, j)=sscp(i+1, j+1)-nvec*m(i)*m(j);
00836 cov(i, j)/=nvec-1;
00837 }
00838 }
00839 }
00840
00841
00842 void TRobustEstimator::Correl()
00843 {
00844
00845
00846 Int_t i, j;
00847 Double_t *sd=new Double_t[fNvar];
00848 for(j=0; j<fNvar; j++)
00849 sd[j]=1./TMath::Sqrt(fCovariance(j, j));
00850 for(i=0; i<fNvar; i++) {
00851 for (j=0; j<fNvar; j++) {
00852 if (i==j)
00853 fCorrelation(i, j)=1.;
00854 else
00855 fCorrelation(i, j)=fCovariance(i, j)*sd[i]*sd[j];
00856 }
00857 }
00858 delete [] sd;
00859 }
00860
00861
00862 void TRobustEstimator::CreateSubset(Int_t ntotal, Int_t htotal, Int_t p, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
00863 {
00864
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874 Double_t kEps = 1e-14;
00875 Int_t i, j;
00876 Bool_t repeat=kFALSE;
00877 Int_t nindex=0;
00878 Int_t num;
00879 for(i=0; i<ntotal; i++)
00880 index[i]=ntotal+1;
00881
00882 for (i=0; i<p+1; i++) {
00883 num=Int_t(gRandom->Uniform(0, 1)*(ntotal-1));
00884 if (i>0){
00885 for(j=0; j<=i-1; j++) {
00886 if(index[j]==num)
00887 repeat=kTRUE;
00888 }
00889 }
00890 if(repeat==kTRUE) {
00891 i--;
00892 repeat=kFALSE;
00893 } else {
00894 index[i]=num;
00895 nindex++;
00896 }
00897 }
00898
00899 ClearSscp(sscp);
00900
00901 TVectorD vec(fNvar);
00902 Double_t det;
00903 for (i=0; i<p+1; i++) {
00904 for (j=0; j<fNvar; j++) {
00905 vec[j]=data[index[i]][j];
00906
00907 }
00908 AddToSscp(sscp, vec);
00909 }
00910
00911 Covar(sscp, fMean, fCovariance, fSd, p+1);
00912 det=fCovariance.Determinant();
00913 while((det<kEps)&&(nindex < htotal)) {
00914
00915
00916
00917 repeat=kFALSE;
00918 do{
00919 num=Int_t(gRandom->Uniform(0,1)*(ntotal-1));
00920 repeat=kFALSE;
00921 for(i=0; i<nindex; i++) {
00922 if(index[i]==num) {
00923 repeat=kTRUE;
00924 break;
00925 }
00926 }
00927 }while(repeat==kTRUE);
00928
00929 index[nindex]=num;
00930 nindex++;
00931
00932 for(j=0; j<fNvar; j++)
00933 vec[j]=data[index[nindex-1]][j];
00934 AddToSscp(sscp, vec);
00935 Covar(sscp, fMean, fCovariance, fSd, nindex);
00936 det=fCovariance.Determinant();
00937 }
00938
00939 if(nindex!=htotal) {
00940 TDecompChol chol(fCovariance);
00941 fInvcovariance = chol.Invert();
00942
00943 TVectorD temp(fNvar);
00944 for(j=0; j<ntotal; j++) {
00945 ndist[j]=0;
00946 for(i=0; i<fNvar; i++)
00947 temp[i]=data[j][i] - fMean(i);
00948 temp*=fInvcovariance;
00949 for(i=0; i<fNvar; i++)
00950 ndist[j]+=(data[j][i]-fMean(i))*temp[i];
00951 }
00952 KOrdStat(ntotal, ndist, htotal-1,index);
00953 }
00954
00955 }
00956
00957
00958 void TRobustEstimator::CreateOrtSubset(TMatrixD &dat,Int_t *index, Int_t hmerged, Int_t nmerged, TMatrixD &sscp, Double_t *ndist)
00959 {
00960
00961
00962
00963
00964 Int_t i, j;
00965 TVectorD vec(fNvar);
00966 for (i=0; i<nmerged; i++) {
00967 ndist[i]=0;
00968 for(j=0; j<fNvar; j++) {
00969 ndist[i]+=fHyperplane[j]*(dat[i][j]-fMean[j]);
00970 ndist[i]=TMath::Abs(ndist[i]);
00971 }
00972 }
00973 KOrdStat(nmerged, ndist, hmerged-1, index);
00974 ClearSscp(sscp);
00975 for (i=0; i<hmerged; i++) {
00976 for(j=0; j<fNvar; j++)
00977 vec[j]=dat[index[i]][j];
00978 AddToSscp(sscp, vec);
00979 }
00980 Covar(sscp, fMean, fCovariance, fSd, hmerged);
00981 }
00982
00983
00984 Double_t TRobustEstimator::CStep(Int_t ntotal, Int_t htotal, Int_t *index, TMatrixD &data, TMatrixD &sscp, Double_t *ndist)
00985 {
00986
00987
00988
00989
00990
00991
00992
00993
00994
00995
00996 Int_t i, j;
00997 TVectorD vec(fNvar);
00998 Double_t det;
00999
01000 TDecompChol chol(fCovariance);
01001 fInvcovariance = chol.Invert();
01002
01003 TVectorD temp(fNvar);
01004 for(j=0; j<ntotal; j++) {
01005 ndist[j]=0;
01006 for(i=0; i<fNvar; i++)
01007 temp[i]=data[j][i]-fMean[i];
01008 temp*=fInvcovariance;
01009 for(i=0; i<fNvar; i++)
01010 ndist[j]+=(data[j][i]-fMean[i])*temp[i];
01011 }
01012
01013
01014 KOrdStat(ntotal, ndist, htotal-1, index);
01015
01016 ClearSscp(sscp);
01017 for (i=0; i<htotal; i++) {
01018 for (j=0; j<fNvar; j++)
01019 temp[j]=data[index[i]][j];
01020 AddToSscp(sscp, temp);
01021 }
01022 Covar(sscp, fMean, fCovariance, fSd, htotal);
01023 det = fCovariance.Determinant();
01024 return det;
01025 }
01026
01027
01028 Int_t TRobustEstimator::Exact(Double_t *ndist)
01029 {
01030
01031
01032
01033 Int_t i, j;
01034
01035 TMatrixDSymEigen eigen(fCovariance);
01036 TVectorD eigenValues=eigen.GetEigenValues();
01037 TMatrixD eigenMatrix=eigen.GetEigenVectors();
01038
01039 for (j=0; j<fNvar; j++) {
01040 fHyperplane[j]=eigenMatrix(j,fNvar-1);
01041 }
01042
01043 for (i=0; i<fN; i++) {
01044 ndist[i]=0;
01045 for(j=0; j<fNvar; j++) {
01046 ndist[i]+=fHyperplane[j]*(fData[i][j]-fMean[j]);
01047 ndist[i]=TMath::Abs(ndist[i]);
01048 }
01049 }
01050 Int_t nhyp=0;
01051
01052 for (i=0; i<fN; i++) {
01053 if(ndist[i] < 1e-14) nhyp++;
01054 }
01055 return nhyp;
01056
01057 }
01058
01059
01060 Int_t TRobustEstimator::Exact2(TMatrixD &mstockbig, TMatrixD &cstockbig, TMatrixD &hyperplane,
01061 Double_t *deti, Int_t nbest, Int_t kgroup,
01062 TMatrixD &sscp, Double_t *ndist)
01063 {
01064
01065
01066
01067
01068
01069
01070 Int_t i, j;
01071
01072 TVectorD vec(fNvar);
01073 Int_t maxind = TMath::LocMax(nbest, deti);
01074 Int_t nh=Exact(ndist);
01075
01076
01077 if(nh>=fH) {
01078 ClearSscp(sscp);
01079 for (i=0; i<fN; i++) {
01080 if(ndist[i]<1e-14) {
01081 for (j=0; j<fNvar; j++)
01082 vec[j]=fData[i][j];
01083 AddToSscp(sscp, vec);
01084 }
01085 }
01086 Covar(sscp, fMean, fCovariance, fSd, nh);
01087
01088 fExact=nh;
01089 return nbest+1;
01090
01091 } else {
01092
01093
01094
01095
01096 for(i=0; i<fNvar; i++) {
01097 mstockbig(nbest*kgroup+maxind,i)=fMean(i);
01098 hyperplane(nbest*kgroup+maxind,i)=fHyperplane(i);
01099 for(j=0; j<fNvar; j++) {
01100 cstockbig(i,nbest*kgroup*fNvar+maxind*fNvar+j)=fCovariance(i,j);
01101 }
01102 }
01103 return maxind;
01104 }
01105 }
01106
01107
01108
01109 Int_t TRobustEstimator::Partition(Int_t nmini, Int_t *indsubdat)
01110 {
01111
01112
01113
01114
01115 Int_t nsub;
01116 if ((fN>=2*nmini) && (fN<=(3*nmini-1))) {
01117 if (fN%2==1){
01118 indsubdat[0]=Int_t(fN*0.5);
01119 indsubdat[1]=Int_t(fN*0.5)+1;
01120 } else
01121 indsubdat[0]=indsubdat[1]=Int_t(fN/2);
01122 nsub=2;
01123 }
01124 else{
01125 if((fN>=3*nmini) && (fN<(4*nmini -1))) {
01126 if(fN%3==0){
01127 indsubdat[0]=indsubdat[1]=indsubdat[2]=Int_t(fN/3);
01128 } else {
01129 indsubdat[0]=Int_t(fN/3);
01130 indsubdat[1]=Int_t(fN/3)+1;
01131 if (fN%3==1) indsubdat[2]=Int_t(fN/3);
01132 else indsubdat[2]=Int_t(fN/3)+1;
01133 }
01134 nsub=3;
01135 }
01136 else{
01137 if((fN>=4*nmini)&&(fN<=(5*nmini-1))){
01138 if (fN%4==0) indsubdat[0]=indsubdat[1]=indsubdat[2]=indsubdat[3]=Int_t(fN/4);
01139 else {
01140 indsubdat[0]=Int_t(fN/4);
01141 indsubdat[1]=Int_t(fN/4)+1;
01142 if(fN%4==1) indsubdat[2]=indsubdat[3]=Int_t(fN/4);
01143 if(fN%4==2) {
01144 indsubdat[2]=Int_t(fN/4)+1;
01145 indsubdat[3]=Int_t(fN/4);
01146 }
01147 if(fN%4==3) indsubdat[2]=indsubdat[3]=Int_t(fN/4)+1;
01148 }
01149 nsub=4;
01150 } else {
01151 for(Int_t i=0; i<5; i++)
01152 indsubdat[i]=nmini;
01153 nsub=5;
01154 }
01155 }
01156 }
01157 return nsub;
01158 }
01159
01160
01161 Int_t TRobustEstimator::RDist(TMatrixD &sscp)
01162 {
01163
01164
01165
01166
01167
01168
01169 Int_t i, j;
01170 Int_t nout=0;
01171
01172 TVectorD temp(fNvar);
01173 TDecompChol chol(fCovariance);
01174 fInvcovariance = chol.Invert();
01175
01176
01177 for (i=0; i<fN; i++) {
01178 fRd[i]=0;
01179 for(j=0; j<fNvar; j++) {
01180 temp[j]=fData[i][j]-fMean[j];
01181 }
01182 temp*=fInvcovariance;
01183 for(j=0; j<fNvar; j++) {
01184 fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
01185 }
01186 }
01187
01188 Double_t med;
01189 Double_t chi = kChiMedian[fNvar-1];
01190
01191 med=TMath::Median(fN, fRd.GetMatrixArray());
01192 med/=chi;
01193 fCovariance*=med;
01194 TDecompChol chol2(fCovariance);
01195 fInvcovariance = chol2.Invert();
01196
01197 for (i=0; i<fN; i++) {
01198 fRd[i]=0;
01199 for(j=0; j<fNvar; j++) {
01200 temp[j]=fData[i][j]-fMean[j];
01201 }
01202
01203 temp*=fInvcovariance;
01204 for(j=0; j<fNvar; j++) {
01205 fRd[i]+=(fData[i][j]-fMean[j])*temp[j];
01206 }
01207 }
01208
01209 Double_t cutoff = kChiQuant[fNvar-1];
01210
01211 ClearSscp(sscp);
01212 for(i=0; i<fN; i++) {
01213 if (fRd[i]<=cutoff) {
01214 for(j=0; j<fNvar; j++)
01215 temp[j]=fData[i][j];
01216 AddToSscp(sscp,temp);
01217 } else {
01218 nout++;
01219 }
01220 }
01221
01222 Covar(sscp, fMean, fCovariance, fSd, fN-nout);
01223 return nout;
01224 }
01225
01226
01227 void TRobustEstimator::RDraw(Int_t *subdat, Int_t ngroup, Int_t *indsubdat)
01228 {
01229
01230
01231
01232 Int_t jndex = 0;
01233 Int_t nrand;
01234 Int_t i, k, m, j;
01235 for (k=1; k<=ngroup; k++) {
01236 for (m=1; m<=indsubdat[k-1]; m++) {
01237
01238 nrand = Int_t(gRandom->Uniform(0, 1) * (fN-jndex))+1;
01239
01240 jndex++;
01241 if (jndex==1) {
01242 subdat[0]=nrand;
01243 } else {
01244 subdat[jndex-1]=nrand+jndex-2;
01245 for (i=1; i<=jndex-1; i++) {
01246 if(subdat[i-1] > nrand+i-2) {
01247 for(j=jndex; j>=i+1; j--) {
01248 subdat[j-1]=subdat[j-2];
01249 }
01250 subdat[i-1]=nrand+i-2;
01251 break;
01252 }
01253 }
01254 }
01255 }
01256 }
01257 }
01258
01259
01260 Double_t TRobustEstimator::KOrdStat(Int_t ntotal, Double_t *a, Int_t k, Int_t *work){
01261
01262 Bool_t isAllocated = kFALSE;
01263 const Int_t kWorkMax=100;
01264 Int_t i, ir, j, l, mid;
01265 Int_t arr;
01266 Int_t *ind;
01267 Int_t workLocal[kWorkMax];
01268 Int_t temp;
01269
01270
01271 if (work) {
01272 ind = work;
01273 } else {
01274 ind = workLocal;
01275 if (ntotal > kWorkMax) {
01276 isAllocated = kTRUE;
01277 ind = new Int_t[ntotal];
01278 }
01279 }
01280
01281 for (Int_t ii=0; ii<ntotal; ii++) {
01282 ind[ii]=ii;
01283 }
01284 Int_t rk = k;
01285 l=0;
01286 ir = ntotal-1;
01287 for(;;) {
01288 if (ir<=l+1) {
01289 if (ir == l+1 && a[ind[ir]]<a[ind[l]])
01290 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
01291 Double_t tmp = a[ind[rk]];
01292 if (isAllocated)
01293 delete [] ind;
01294 return tmp;
01295 } else {
01296 mid = (l+ir) >> 1;
01297 {temp = ind[mid]; ind[mid]=ind[l+1]; ind[l+1]=temp;}
01298 if (a[ind[l]]>a[ind[ir]])
01299 {temp = ind[l]; ind[l]=ind[ir]; ind[ir]=temp;}
01300
01301 if (a[ind[l+1]]>a[ind[ir]])
01302 {temp=ind[l+1]; ind[l+1]=ind[ir]; ind[ir]=temp;}
01303
01304 if (a[ind[l]]>a[ind[l+1]])
01305 {temp = ind[l]; ind[l]=ind[l+1]; ind[l+1]=temp;}
01306
01307 i=l+1;
01308 j=ir;
01309 arr = ind[l+1];
01310 for (;;) {
01311 do i++; while (a[ind[i]]<a[arr]);
01312 do j--; while (a[ind[j]]>a[arr]);
01313 if (j<i) break;
01314 {temp=ind[i]; ind[i]=ind[j]; ind[j]=temp;}
01315 }
01316 ind[l+1]=ind[j];
01317 ind[j]=arr;
01318 if (j>=rk) ir = j-1;
01319 if (j<=rk) l=i;
01320 }
01321 }
01322 }