00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "BCH1D.h"
00011
00012 #include "BCLog.h"
00013 #include "BCMath.h"
00014
00015 #include <TH1D.h>
00016 #include <TH2.h>
00017 #include <TLine.h>
00018 #include <TPolyLine.h>
00019 #include <TPaveLabel.h>
00020 #include <TLatex.h>
00021 #include <TError.h>
00022 #include <TCanvas.h>
00023 #include <TMarker.h>
00024 #include <TLegend.h>
00025 #include <TString.h>
00026
00027 #include <math.h>
00028
00029
00030
00031 BCH1D::BCH1D()
00032 {
00033 fHistogram = 0;
00034 fDefaultCLLimit = 95.;
00035
00036 fModeFlag = 0;
00037 }
00038
00039
00040
00041 BCH1D::~BCH1D()
00042 {
00043 if (fHistogram) delete fHistogram;
00044 }
00045
00046
00047
00048 double BCH1D::GetMode()
00049 {
00050 return fHistogram->GetBinCenter(fHistogram->GetMaximumBin());
00051 }
00052
00053
00054
00055 double BCH1D::GetQuantile(double probability)
00056 {
00057 int nquantiles = 1;
00058 double quantiles[1];
00059 double probsum[1];
00060
00061 probsum[0] = probability;
00062
00063
00064 fHistogram->GetQuantiles(nquantiles, quantiles, probsum);
00065
00066 return quantiles[0];
00067 }
00068
00069
00070
00071 double BCH1D::GetIntegral(double valuemin, double valuemax)
00072 {
00073 double integral = 0;
00074
00075 int binmin = fHistogram->FindBin(valuemin);
00076 int binmax = fHistogram->FindBin(valuemax);
00077
00078
00079 integral = fHistogram->Integral(binmin, binmax);
00080
00081 return integral;
00082 }
00083
00084
00085
00086 double BCH1D::GetPValue(double probability)
00087 {
00088
00089 return fHistogram->Integral(1, fHistogram->FindBin(probability));
00090 }
00091
00092
00093
00094 void BCH1D::SetDefaultCLLimit(double limit)
00095 {
00096
00097 if(limit>=100. || limit<68.)
00098 BCLog::OutWarning(
00099 Form("BCH1D::SetDefaultCLLimit : Value %.1f out of range. Keeping default %.1f%% CL limit.",limit,fDefaultCLLimit));
00100
00101
00102 else
00103 fDefaultCLLimit = limit;
00104 }
00105
00106
00107
00108 void BCH1D::Print(const char * filename, int options, double ovalue, int ww, int wh)
00109 {
00110 char file[256];
00111 int i=0;
00112 while(i<255 && *filename!='\0')
00113 file[i++]=*filename++;
00114 file[i]='\0';
00115
00116
00117 fHistogram->SetLineWidth(1);
00118
00119
00120 TCanvas * c;
00121 if(ww > 0 && wh > 0)
00122 c = new TCanvas("c","c",ww,wh);
00123 else
00124 c = new TCanvas("c","c");
00125
00126 c->cd();
00127 Draw(options, ovalue);
00128
00129 gPad->RedrawAxis();
00130
00131
00132 c->Print(file);
00133 }
00134
00135
00136
00137 void BCH1D::Draw(int options, double ovalue)
00138 {
00139 double min, max;
00140 double mode;
00141 double thismode = GetMode();
00142
00143 int nbins = fHistogram->GetNbinsX();
00144
00145 fHistogram->Scale(1./fHistogram->Integral("width"));
00146
00147 if(fModeFlag)
00148 mode=fMode;
00149 else
00150 mode=thismode;
00151
00152
00153 TLine * line;
00154
00155
00156 bool flagLegend=false;
00157 char confidenceLabel[25];
00158
00159
00160 fHistogram->SetLineColor(kBlack);
00161 fHistogram->SetLineWidth(1);
00162 fHistogram->SetFillStyle(0);
00163
00164
00165 switch(options)
00166 {
00167
00168
00169 case 0:
00170 if (fabs(ovalue) >= 100 || ovalue==0.)
00171 {
00172
00173 min = GetQuantile(.16);
00174 max = GetQuantile(.84);
00175
00176
00177 flagLegend = true;
00178 sprintf(confidenceLabel, "Central 68%%");
00179
00180 if ( fHistogram->FindBin(thismode) == fHistogram->GetNbinsX() )
00181 {
00182 min = GetQuantile(1.-(double)fDefaultCLLimit/100.);
00183 max = fHistogram->GetXaxis()->GetXmax();
00184 ovalue = fDefaultCLLimit;
00185 sprintf(confidenceLabel, "%g%% region", fDefaultCLLimit);
00186 }
00187 else if ( fHistogram->FindBin(thismode)==1)
00188 {
00189 min = fHistogram->GetXaxis()->GetXmin();
00190 max = GetQuantile((double)fDefaultCLLimit/100.);
00191 ovalue = -fDefaultCLLimit;
00192 sprintf(confidenceLabel, "%g%% region", fDefaultCLLimit);
00193 }
00194 }
00195
00196 else if(ovalue < 0)
00197 {
00198 min = fHistogram->GetXaxis()->GetXmin();
00199 max = GetQuantile(-ovalue/100.);
00200 }
00201 else
00202 {
00203 min = GetQuantile(1. - ovalue/100.);
00204 max = fHistogram->GetXaxis()->GetXmax();
00205 }
00206
00207
00208 DrawShadedLimits(mode, min, max, ovalue);
00209
00210
00211 if(flagLegend)
00212 this ->DrawLegend(confidenceLabel);
00213
00214 break;
00215
00216
00217 case 1:
00218
00219 fHistogram->Draw();
00220 min = fHistogram->GetBinLowEdge(1);
00221 max = fHistogram->GetBinLowEdge(nbins+1);
00222 if(min<=ovalue && ovalue<=max)
00223 {
00224 line = new TLine();
00225 line->SetLineColor(kRed);
00226 line->DrawLine(ovalue, 0., ovalue, 1.05 * fHistogram->GetMaximum());
00227 }
00228
00229 break;
00230
00231
00232 case 2:
00233
00234 if(ovalue<50)
00235 ovalue = 68.;
00236
00237 GetSmallestInterval(min, max, ovalue/100.);
00238 DrawShadedLimits(mode, min, max, 0.);
00239
00240 break;
00241
00242
00243 case 3:
00244
00245 if(ovalue<50)
00246 ovalue = 68.;
00247
00248 DrawSmallest(mode,ovalue);
00249
00250 break;
00251
00252
00253 case 4:
00254
00255 DrawDelta(ovalue);
00256
00257 break;
00258
00259
00260 default:
00261
00262 BCLog::OutError(Form("BCH1D::Draw : Invalid option %d",options));
00263 break;
00264 }
00265 }
00266
00267 void BCH1D::DrawDelta(double value) const
00268 {
00269
00270 double xmin = fHistogram->GetXaxis()->GetXmin();
00271 double xmax = fHistogram->GetXaxis()->GetXmax();
00272 double ysize = 1.3 * fHistogram->GetMaximum();
00273
00274 TH2D* hsc = new TH2D(
00275 TString::Format("h2scale_%s_%d",fHistogram->GetName(),BCLog::GetHIndex()),"",
00276 50, xmin, xmax, 10, 0., ysize);
00277 hsc->SetStats(0);
00278 hsc->SetXTitle(fHistogram->GetXaxis()->GetTitle());
00279 hsc->SetYTitle(fHistogram->GetYaxis()->GetTitle());
00280 hsc->Draw();
00281
00282
00283
00284 TLatex * tmax_text = new TLatex();
00285 tmax_text->SetTextSize(0.035);
00286 tmax_text->SetTextFont(62);
00287 tmax_text->SetTextAlign(22);
00288
00289 double xprint=(xmax+xmin)/2.;
00290 double yprint=ysize*(1-1.4*tmax_text->GetTextSize());
00291
00292 tmax_text->DrawLatex(xprint,yprint, TString::Format("%s = %g", fHistogram->GetXaxis()->GetTitle(), value));
00293 delete tmax_text;
00294
00295
00296 TH1D * hist_temp = new TH1D(TString::Format("h1scale_%s_%d", fHistogram->GetName(), BCLog::GetHIndex()), "", 100, xmin, xmax);
00297 hist_temp->SetBinContent(hist_temp->FindBin(value), fHistogram->GetMaximum());
00298 hist_temp->Draw("same");
00299 }
00300
00301
00302 void BCH1D::DrawLegend(const char* text)
00303 {
00304
00305
00306 TLegend* legend = new TLegend(0.73, 0.72, 0.86, 0.85);
00307 legend->SetFillColor(kWhite);
00308 legend->SetBorderSize(1);
00309
00310 TMarker* triangle = new TMarker(0, 0, 23);
00311 triangle->SetMarkerColor(kRed);
00312 legend->AddEntry(triangle, "Global mode", "P");
00313
00314 TMarker* diamond = new TMarker(0, 0, 27);
00315 diamond->SetMarkerColor(kBlue);
00316 legend->AddEntry(diamond, "Mean", "P");
00317
00318 TLine * line;
00319 line = new TLine();
00320 line->SetLineStyle(2);
00321 line->SetLineColor(kGreen + 2);
00322 legend->AddEntry(line, "Median", "l");
00323
00324 TLegend* band = new TLegend(0, 0, 1, 1);
00325 band->SetFillColor(kYellow);
00326 legend->AddEntry(band, text, "F");
00327
00328 legend->SetTextAlign(12);
00329
00330
00331 legend->SetTextSize(0.016);
00332
00333 legend->Draw();
00334 }
00335
00336
00337
00338
00339
00340
00341
00342 void BCH1D::DrawShadedLimits(double mode, double min, double max, double limit)
00343 {
00344 double maximum = fHistogram->GetMaximum();
00345
00346 double x0 = mode;
00347 double y0 = fHistogram->GetBinContent( fHistogram->FindBin(mode) );
00348
00349 double x1 = fHistogram->GetMean();
00350 double y1 = fHistogram->GetBinContent( fHistogram->FindBin(x1) );
00351
00352 double x2 = GetQuantile(.5);
00353 double y2 = fHistogram->GetBinContent( fHistogram->FindBin(x2) );
00354
00355 double ysize = maximum*1.3;
00356
00357 double xmin = fHistogram->GetXaxis()->GetXmin();
00358 double xmax = fHistogram->GetXaxis()->GetXmax();
00359
00360
00361 TH2D * hsc = new TH2D(
00362 TString::Format("h2scale_%s_%d",fHistogram->GetName(),BCLog::GetHIndex()),"",
00363 10, xmin, xmax, 10, 0., ysize);
00364 hsc->SetStats(0);
00365 hsc->SetXTitle(fHistogram->GetXaxis()->GetTitle());
00366 hsc->SetYTitle(fHistogram->GetYaxis()->GetTitle());
00367 hsc->Draw();
00368
00369
00370 fHistogram->SetLineWidth(1);
00371 fHistogram->Draw("same");
00372
00373
00374 TH1D * hist_shaded = GetSubHisto(min,max,TString::Format("%s_sub_%d",fHistogram->GetName(),BCLog::GetHIndex()));
00375 hist_shaded->SetFillStyle(1001);
00376 hist_shaded->SetFillColor(kYellow);
00377
00378
00379 hist_shaded->Draw("same");
00380
00381 gPad->RedrawAxis();
00382
00383
00384 TPolyLine * tmax;
00385
00386 double dx = 0.01*(xmax-xmin);
00387 double dy = 0.04*(ysize);
00388 y0+=dy/5.;
00389 double tmax_x[] = {x0, x0 + dx, x0 - dx, x0};
00390 double tmax_y[] = {y0, y0 + dy, y0 + dy, y0};
00391 tmax = new TPolyLine(4,tmax_x,tmax_y);
00392 tmax->SetLineColor(kRed);
00393 tmax->SetLineWidth(1);
00394 tmax->SetFillColor(kRed);
00395 tmax->Draw();
00396 tmax->Draw("f");
00397
00398
00399 TPolyLine * tmean;
00400
00401 y1+=dy/5.;
00402
00403
00404 double tmean_x[] = {x1, x1 + dx, x1 , x1 - dx, x1};
00405 double tmean_y[] = {y1, y1 + dy/2., y1 + dy, y1 + dy/2., y1};
00406 tmean = new TPolyLine(5,tmean_x,tmean_y);
00407 tmean->SetLineColor(kBlue);
00408
00409 tmean->SetLineWidth(1);
00410
00411 tmean->Draw();
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434 TLine * line;
00435 line = new TLine();
00436 line->SetLineStyle(2);
00437 line->SetLineColor(kGreen+2);
00438 line->DrawLine(x2, 0., x2, y2);
00439
00440
00441
00442
00443
00444 double delta_max = fmax(fabs(max-x1),fabs(x1-min));
00445
00446 int sd = 2 + (int)log10(fabs(x1/delta_max));
00447
00448 if( (int)log10(x1) > (int)log10(delta_max) )
00449 sd++;
00450
00451 TLatex * tmax_text = new TLatex();
00452 tmax_text->SetTextSize(0.035);
00453 tmax_text->SetTextFont(62);
00454 tmax_text->SetTextAlign(22);
00455
00456
00457 double xprint=(xmax+xmin)/2.;
00458
00459 double yprint=ysize*(1-1.4*tmax_text->GetTextSize());
00460
00461 if(fabs(limit)<50)
00462 tmax_text->DrawLatex(xprint,yprint,
00463 TString::Format( TString::Format("%%s^{med} = %%.%dg ^{+%%.2g}_{ -%%.2g}",sd),
00464 fHistogram->GetXaxis()->GetTitle(), x2, max-x2, x2-min));
00465
00466 else if (limit<0)
00467 tmax_text->DrawLatex(xprint,yprint,
00468 TString::Format( TString::Format("%%s (%d%%%% prob.) < %%.4g",-(int)limit),
00469 fHistogram->GetXaxis()->GetTitle(), max));
00470
00471 else if (limit>0)
00472 tmax_text->DrawLatex(xprint,yprint,
00473 TString::Format( TString::Format("%%s (%d%%%% prob.) > %%.4g",(int)limit),
00474 fHistogram->GetXaxis()->GetTitle(), min));
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485
00486
00487
00488
00489
00490
00491
00492
00493
00494
00495
00496 }
00497
00498
00499
00500 void BCH1D::DrawSmallest(double mode, double prob, bool drawmean)
00501 {
00502
00503 TPolyLine * tmax;
00504
00505 double x0 = mode;
00506 double y0 = fHistogram->GetBinContent( fHistogram->FindBin(mode) );
00507 double xmin = fHistogram->GetXaxis()->GetXmin();
00508 double xmax = fHistogram->GetXaxis()->GetXmax();
00509 double ysize = 1.3 * fHistogram->GetMaximum();
00510
00511 double x1 = fHistogram->GetMean();
00512 double y1 = fHistogram->GetBinContent( fHistogram->FindBin(x1) );
00513
00514 double x2 = GetQuantile(.5);
00515 double y2 = fHistogram->GetBinContent( fHistogram->FindBin(x2) );
00516
00517 double dx = 0.01*(xmax-xmin);
00518 double dy = 0.04*(ysize);
00519 double tmax_x[] = {x0, x0 + dx, x0 - dx, x0};
00520 double tmax_y[] = {y0, y0 + dy, y0 + dy, y0};
00521 tmax = new TPolyLine(4,tmax_x,tmax_y);
00522 tmax->SetLineColor(kRed);
00523 tmax->SetFillColor(kRed);
00524
00525
00526 TH2D * hsc = new TH2D(
00527 TString::Format("h2scale_%s_%d",fHistogram->GetName(),BCLog::GetHIndex()),"",
00528 10, xmin, xmax, 10, 0., ysize);
00529 hsc->SetStats(0);
00530 hsc->SetXTitle(fHistogram->GetXaxis()->GetTitle());
00531 hsc->SetYTitle(fHistogram->GetYaxis()->GetTitle());
00532 hsc->Draw();
00533
00534
00535 TH1D * hist_temp1 = (TH1D*) fHistogram->Clone();
00536 hist_temp1->Scale(1.0/fHistogram->Integral("width"));
00537 hist_temp1->SetFillColor(kYellow);
00538 hist_temp1->SetFillStyle(1001);
00539
00540
00541 TH1D * hist_temp2 = (TH1D*) fHistogram->Clone();
00542 hist_temp2->Scale(1.0/fHistogram->Integral("width"));
00543
00544
00545 hist_temp1->Reset();
00546
00547
00548
00549 prob /= 100.;
00550
00551
00552 double sum = 0.0;
00553
00554 int n=0;
00555 int nbins=fHistogram->GetNbinsX();
00556
00557
00558 while (sum < prob && n < nbins)
00559 {
00560 n++;
00561
00562
00563 int bin = hist_temp2->GetMaximumBin();
00564
00565
00566 double val = hist_temp2->GetBinContent(bin);
00567 hist_temp1->SetBinContent(bin, val);
00568
00569
00570 hist_temp2->SetBinContent(bin, 0.0);
00571
00572
00573 sum += val * hist_temp2->GetBinWidth(bin);
00574 }
00575
00576
00577 hist_temp1->Scale(fHistogram->Integral("width"));
00578
00579
00580 fHistogram->Draw("same");
00581 hist_temp1->Draw("same");
00582
00583
00584 tmax->Draw("f");
00585
00586 if(drawmean)
00587 {
00588
00589
00590 TPolyLine * tmean;
00591
00592 y1+=dy/5.;
00593
00594
00595 double tmean_x[] = {x1, x1 + dx, x1 , x1 - dx, x1};
00596 double tmean_y[] = {y1, y1 + dy/2., y1 + dy, y1 + dy/2., y1};
00597 tmean = new TPolyLine(5,tmean_x,tmean_y);
00598 tmean->SetLineColor(kBlue);
00599
00600 tmean->SetLineWidth(1);
00601
00602 tmean->Draw();
00603
00604
00605 TLine * line;
00606 line = new TLine();
00607 line->SetLineStyle(2);
00608 line->SetLineColor(kGreen+2);
00609 line->DrawLine(x2, 0., x2, y2);
00610 }
00611
00612
00613 delete hist_temp2;
00614 }
00615
00616
00617
00618 double BCH1D::GetSmallestInterval(double & min, double & max, double content)
00619 {
00620 if(content<=0. || content>= 1.)
00621 return 0.;
00622
00623 int nbins = fHistogram->GetNbinsX();
00624
00625 double factor = fHistogram->Integral("width");
00626 if(factor==0)
00627 return 0.;
00628
00629
00630 fHistogram->Scale(1./factor);
00631
00632 double xmin = fHistogram->GetXaxis()->GetXmin();
00633 double xmax = fHistogram->GetXaxis()->GetXmax();
00634 double width = xmax - xmin;
00635
00636 double xdown=xmin;
00637 double xup=xmax;
00638
00639 int ndiv = 10;
00640
00641
00642
00643
00644
00645 int warn=0;
00646
00647 double integral_out=0.;
00648
00649
00650 for(int i=1;i<nbins+1;i++)
00651 {
00652 if(fHistogram->Integral(i,nbins,"width") < content)
00653 break;
00654
00655
00656 double firstbinwidth = fHistogram->GetBinWidth(i);
00657
00658
00659
00660 for(int j=0;j<ndiv;j++)
00661 {
00662 double dxdown = (double)(ndiv-j)/(double)ndiv * firstbinwidth;
00663 xdown = fHistogram->GetBinLowEdge(i) + firstbinwidth - dxdown;
00664 double integral = dxdown * fHistogram->GetBinContent(i);
00665
00666 if(integral>content)
00667 {
00668
00669 xup = xdown + content / fHistogram->GetBinContent(i);
00670 warn = 1;
00671 }
00672 else
00673 for(int k=i+1;k<nbins+1;k++)
00674 {
00675
00676 double thisbin = fHistogram->GetBinContent(k) * fHistogram->GetBinWidth(k);
00677
00678 if(integral+thisbin==content)
00679 {
00680 xup = fHistogram->GetBinLowEdge(k+1);
00681 break;
00682 }
00683
00684 if(integral+thisbin>content)
00685 {
00686 xup = fHistogram->GetBinLowEdge(k) + (content-integral)/thisbin * fHistogram->GetBinWidth(k);
00687 integral += thisbin * (xup-fHistogram->GetBinLowEdge(k))/fHistogram->GetBinWidth(k);
00688 break;
00689 }
00690
00691 integral += thisbin;
00692 }
00693
00694 if(integral < content)
00695 continue;
00696
00697 if(xup - xdown < width)
00698 {
00699
00700 width = xup - xdown;
00701 xmin = xdown;
00702 xmax = xup;
00703 integral_out = integral;
00704 }
00705 }
00706 }
00707
00708 if(warn)
00709 {
00710 BCLog::OutWarning(
00711 Form("BCH1D::GetSmallestInterval : The requested content of %d%% fits within one bin.",(int)(content*100)));
00712 BCLog::OutWarning("BCH1D::GetSmallestInterval : MAKE FINER BINNING (OR CHANGE RANGE) !!!");
00713 }
00714
00715
00716 fHistogram->Scale(factor);
00717
00718 min=xmin;
00719 max=xmax;
00720
00721 return integral_out;
00722 }
00723
00724
00725
00726 double BCH1D::IntegralWidth(double min, double max)
00727 {
00728 int imin = fHistogram->FindBin(min);
00729 int imax = fHistogram->FindBin(max);
00730
00731 int nbins = fHistogram->GetNbinsX();
00732
00733
00734 if ( imin<1 || imin>nbins || imax<1 || imax>nbins )
00735 return -1.;
00736
00737 if ( imin==imax )
00738 return -1.;
00739
00740
00741 if (imin>imax)
00742 {
00743 int i=imin;
00744 double x=min;
00745 imin=imax, min=max;
00746 imax=i, max=x;
00747 }
00748
00749
00750 double first = ( fHistogram->GetBinLowEdge(imin+1) - min ) * fHistogram->GetBinContent(imin);
00751
00752
00753 double last = ( max - fHistogram->GetBinLowEdge(imax) ) * fHistogram->GetBinContent(imax);
00754
00755
00756 double inbetween=0.;
00757 if(imax-imin>1)
00758 inbetween = fHistogram->Integral(imin+1, imax-1, "width");
00759
00760 return first + inbetween + last;
00761 }
00762
00763
00764
00765 TH1D * BCH1D::GetSubHisto(double min, double max, const char * name)
00766 {
00767 if(min>max)
00768 {
00769 double t=min;
00770 min=max;
00771 max=t;
00772 }
00773
00774 int imin = fHistogram->FindBin(min);
00775 int imax = fHistogram->FindBin(max);
00776
00777 int nbins = fHistogram->GetNbinsX();
00778 double xmin = fHistogram->GetXaxis()->GetXmin();
00779 double xmax = fHistogram->GetXaxis()->GetXmax();
00780
00781 if( min==max || (min<=xmin && max>=xmax) )
00782 {
00783 TH1D * h0 = (TH1D*) fHistogram->Clone(name);
00784 return h0;
00785 }
00786
00787 if (imin<1)
00788 {
00789 imin=1;
00790 min=xmin;
00791 }
00792 if (imax>nbins)
00793 {
00794 imax=nbins;
00795 max=xmax;
00796 }
00797
00798 double * xb = new double[nbins+3];
00799 int n=0;
00800
00801 int domin=1;
00802 int domax=1;
00803
00804 for (int i=1;i<nbins+2;i++)
00805 {
00806 double x0 = fHistogram->GetBinLowEdge(i);
00807
00808 if (min<x0 && domin)
00809 {
00810 xb[n++]=min;
00811 domin=0;
00812 }
00813 else if (min==x0)
00814 domin=0;
00815
00816 if (max<x0 && domax)
00817 {
00818 xb[n++]=max;
00819 domax=0;
00820 }
00821 else if (max==x0)
00822 domax=0;
00823
00824 xb[n++]=x0;
00825 }
00826
00827
00828 TH1D * h0 = new TH1D(name,"",n-1,xb);
00829 for(int i=1;i<n;i++)
00830 {
00831 double x0 = h0->GetBinCenter(i);
00832 if(x0<min || x0>max)
00833 continue;
00834
00835 int bin=fHistogram->FindBin(x0);
00836 h0->SetBinContent(i, fHistogram->GetBinContent(bin));
00837 }
00838
00839 return h0;
00840 }
00841
00842
00843
00844 TH1D * BCH1D::GetSmallestIntervalHistogram(double level)
00845 {
00846
00847 TH1D * hist_yellow = (TH1D*) fHistogram->Clone();
00848 hist_yellow->Reset();
00849 hist_yellow->SetFillStyle(1001);
00850 hist_yellow->SetFillColor(kYellow);
00851
00852
00853 TH1D * hist_temp = (TH1D*) fHistogram->Clone(TString::Format("%s_%i",fHistogram->GetName(),BCLog::GetHIndex()));
00854 double factor = hist_temp->Integral("");
00855
00856 if(factor == 0)
00857 return 0;
00858
00859 hist_temp->Scale(1. / factor);
00860
00861
00862
00863
00864
00865
00866
00867
00868 double sumprob = 0.0;
00869
00870 while (sumprob < level)
00871 {
00872
00873 int bin = hist_temp->GetMaximumBin();
00874 double value = hist_temp->GetMaximum();
00875
00876
00877 hist_yellow->SetBinContent(bin, 1.0);
00878
00879
00880 hist_temp->SetBinContent(bin, 0.0);
00881
00882
00883 sumprob += value;
00884 }
00885
00886 delete hist_temp;
00887
00888 return hist_yellow;
00889 }
00890
00891
00892
00893 std::vector<double> BCH1D::GetSmallestIntervals(double content)
00894 {
00895 std::vector<double> v;
00896
00897 TH1D * hist = GetSmallestIntervalHistogram(content);
00898
00899 int nbins = hist->GetNbinsX();
00900 int ninter = 0;
00901 int lastbin = -1;
00902
00903 double max = -1;
00904 double localmax = -1;
00905 double localmaxpos = -1;
00906 double localint = 0;
00907 bool flag = false;
00908
00909 for (int i = 1; i <= nbins; ++i)
00910 {
00911
00912 if (!flag && hist->GetBinContent(i) > 0.)
00913 {
00914 flag = true;
00915 v.push_back(hist->GetBinLowEdge(i));
00916
00917
00918 lastbin = i;
00919
00920
00921 ninter++;
00922
00923
00924 localmax = fHistogram->GetBinContent(i);
00925 localmaxpos = hist->GetBinLowEdge(i);
00926
00927
00928 localint = 0;
00929 }
00930
00931
00932 if ((flag && !(hist->GetBinContent(i) > 0.)) || (flag && i == nbins))
00933 {
00934 flag = false;
00935 v.push_back(hist->GetBinLowEdge(i) + hist->GetBinWidth(i));
00936
00937
00938 if (i == nbins && localmax < fHistogram->GetBinContent(i))
00939 localmaxpos = hist->GetBinCenter(i) + 0.5 * hist->GetBinWidth(i);
00940
00941
00942 if (localmax > max)
00943 max = localmax;
00944
00945
00946 v.push_back(localmax);
00947 v.push_back(localmaxpos);
00948
00949
00950 v.push_back(localint);
00951 }
00952
00953
00954 if (i < nbins && localmax < fHistogram->GetBinContent(i))
00955 {
00956 localmax = fHistogram->GetBinContent(i);
00957 localmaxpos = hist->GetBinCenter(i);
00958 }
00959
00960
00961 localint += fHistogram->GetBinContent(i) / fHistogram->Integral();
00962 }
00963
00964
00965 for (int i = 0; i < ninter; ++i)
00966 v[i*5+2] = v.at(i*5+2) / max;
00967
00968 return v;
00969 }
00970
00971