• Main Page
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

BCH1D.cxx

Go to the documentation of this file.
00001 /*
00002  * Copyright (C) 2008-2010, Daniel Kollar and Kevin Kroeninger.
00003  * All rights reserved.
00004  *
00005  * For the licensing terms see doc/COPYING.
00006  */
00007 
00008 // ---------------------------------------------------------
00009 
00010 #include "BAT/BCH1D.h"
00011 
00012 #include "BAT/BCLog.h"
00013 #include "BAT/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.; // in percent
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    // use ROOT function to calculat quantile.
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    // use ROOT function to calculate integral.
00079    integral = fHistogram -> Integral(binmin, binmax);
00080 
00081    return integral;
00082 }
00083 
00084 // ---------------------------------------------------------
00085 
00086 double BCH1D::GetPValue(double probability)
00087 {
00088    // use ROOT function to calculate the integral from 0 to "probability".
00089    return fHistogram -> Integral(1, fHistogram -> FindBin(probability));
00090 }
00091 
00092 // ---------------------------------------------------------
00093 
00094 void BCH1D::SetDefaultCLLimit(double limit)
00095 {
00096    // check if limit is between 68% and 100%. Give a warning if not ...
00097    if(limit>=100. || limit<68.)
00098       BCLog::Out(BCLog::warning,BCLog::warning,
00099          Form("BCH1D::SetDefaultCLLimit : Value %.1f out of range. Keeping default %.1f%% CL limit.",limit,fDefaultCLLimit));
00100 
00101    // ... or set value if true.
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 // temporary
00117    fHistogram -> SetLineWidth(1);
00118 
00119    // create temporary canvas
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    this->Draw(options, ovalue);
00128 
00129    gPad->RedrawAxis();
00130 
00131    // print to file.
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 = this->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    // define temporary line.
00153    TLine * line;
00154 
00155    //control if legend is required
00156    bool flagLegend=false;
00157    char confidenceLabel[25];
00158 
00159    // check drawing option.
00160    switch(options)
00161    {
00162       // Draw a band between 16% and 84% probability.
00163       // If the mode is outside the band only draw a limit.
00164       case 0:
00165          if (fabs(ovalue) >= 100 || ovalue==0.)
00166          {//default case if no args to Draw() supplied
00167 
00168             min = this -> GetQuantile(.16);
00169             max = this -> GetQuantile(.84);
00170 
00171             //draw a legend later
00172             flagLegend = true;
00173             sprintf(confidenceLabel, "Central 68%%");
00174 
00175 //          if ( thismode > max || fHistogram -> FindBin(thismode) == fHistogram -> GetNbinsX() )
00176             if ( fHistogram -> FindBin(thismode) == fHistogram -> GetNbinsX() )
00177             {
00178                min = this -> GetQuantile(1.-(double)fDefaultCLLimit/100.);
00179                max = fHistogram->GetXaxis()->GetXmax();
00180                ovalue = fDefaultCLLimit;
00181                sprintf(confidenceLabel, "%g%% region", fDefaultCLLimit);
00182             }
00183 //          else if (thismode < min || fHistogram->FindBin(thismode)==1)
00184             else if ( fHistogram->FindBin(thismode)==1)
00185             {
00186                min = fHistogram->GetXaxis()->GetXmin();
00187                max = this -> GetQuantile((double)fDefaultCLLimit/100.);
00188                ovalue = -fDefaultCLLimit;
00189                sprintf(confidenceLabel, "%g%% region", fDefaultCLLimit);
00190             }
00191          }
00192 
00193          else if(ovalue < 0)
00194          {
00195             min = fHistogram->GetXaxis()->GetXmin();
00196             max = this -> GetQuantile(-ovalue/100.);
00197          }
00198          else
00199          {
00200             min = this -> GetQuantile(1. - ovalue/100.);
00201             max = fHistogram->GetXaxis()->GetXmax();
00202          }
00203 
00204          // do the drawing
00205          this->DrawShadedLimits(mode, min, max, ovalue);
00206 
00207          // add legend for the symbols mean, mode, median, confidence band
00208          if(flagLegend)
00209             this ->DrawLegend(confidenceLabel);
00210 
00211 
00212          break;
00213 
00214       // Draw a line at "ovalue".
00215       case 1:
00216 
00217          fHistogram -> Draw();
00218          min = fHistogram->GetBinLowEdge(1);
00219          max = fHistogram->GetBinLowEdge(nbins+1);
00220          if(min<=ovalue && ovalue<=max)
00221          {
00222             line = new TLine();
00223             line -> SetLineColor(kRed);
00224             line -> DrawLine(ovalue, 0., ovalue, 1.05 * fHistogram -> GetMaximum());
00225          }
00226 
00227          break;
00228 
00229       // Draw a shaded band at the smallest interval.
00230       case 2:
00231 
00232          if(ovalue<50) // just to ensure there's some sense in the number
00233             ovalue = 68.; // default is 68%
00234 
00235          this->GetSmallestInterval(min, max, ovalue/100.);
00236          this->DrawShadedLimits(mode, min, max, 0.);
00237 
00238 //       this -> DrawLegend(TString::Format("Smallest %g%",ovalue));
00239 
00240          break;
00241 
00242       // Draw a shaded band at the smallest intervals
00243       case 3:
00244 
00245          if(ovalue<50) // just to ensure there's some sense in the number
00246             ovalue = 68.; // default is 68%
00247 
00248          this -> DrawSmallest(mode,ovalue);
00249 
00250          break;
00251 
00252       // Sort out bad options and warn.
00253       default:
00254 
00255          BCLog::Out(BCLog::warning, BCLog::warning, Form("BCH1D::Draw : Invalid option %d",options));
00256          break;
00257    }
00258 }
00259 // ---------------------------------------------------------
00260 
00261 void BCH1D::DrawLegend(const char* text)
00262 {
00263    //draw on top right corner
00264 
00265    TLegend* legend = new TLegend(0.73, 0.72, 0.86, 0.85);
00266    legend -> SetFillColor(kWhite);
00267    legend -> SetBorderSize(1);
00268 
00269    TMarker* triangle = new TMarker(0, 0, 23);
00270    triangle -> SetMarkerColor(kRed);
00271    legend -> AddEntry(triangle, "Global mode", "P");
00272 
00273    TMarker* diamond = new TMarker(0, 0, 27);
00274    diamond -> SetMarkerColor(kBlue);
00275    legend -> AddEntry(diamond, "Mean", "P");
00276 
00277    TLine * line;
00278    line = new TLine();
00279    line -> SetLineStyle(2);
00280    line -> SetLineColor(kGreen + 2);
00281    legend -> AddEntry(line, "Median", "l");
00282 
00283    TLegend* band = new TLegend(0, 0, 1, 1);
00284    band -> SetFillColor(kYellow);
00285    legend -> AddEntry(band, text, "F");
00286 
00287    legend -> SetTextAlign(12);
00288 
00289    //fine tuned by hand so text legible even with several plots in a row
00290    legend -> SetTextSize(0.016);
00291 
00292    legend -> Draw();
00293 }
00294 
00295 // ---------------------------------------------------------
00296 
00297 //TODO Are graphics objects ever deleted from the heap?
00298 
00299 void BCH1D::DrawShadedLimits(double mode, double min, double max, double limit)
00300 {
00301    double maximum = fHistogram -> GetMaximum();
00302 
00303    double x0 = mode;
00304    double y0 = fHistogram->GetBinContent( fHistogram->FindBin(mode) );
00305 
00306    double x1 = fHistogram->GetMean();
00307    double y1 = fHistogram->GetBinContent( fHistogram->FindBin(x1) );
00308 
00309    double x2 = this -> GetQuantile(.5); // median
00310    double y2 = fHistogram->GetBinContent( fHistogram->FindBin(x2) );
00311 
00312    double ysize = maximum*1.3;
00313 
00314    double xmin = fHistogram->GetXaxis()->GetXmin();
00315    double xmax = fHistogram->GetXaxis()->GetXmax();
00316 
00317    // draw histogram with axes first
00318    TH2D * hsc = new TH2D(
00319          TString::Format("h2scale_%s_%d",fHistogram->GetName(),BCLog::GetHIndex()),"",
00320          10, xmin, xmax, 10, 0., ysize);
00321    hsc -> SetStats(0);
00322    hsc -> SetXTitle(fHistogram->GetXaxis()->GetTitle());
00323    hsc -> SetYTitle(fHistogram->GetYaxis()->GetTitle());
00324    hsc -> Draw();
00325 
00326    // draw histogram
00327    fHistogram -> SetLineWidth(1);
00328    fHistogram -> Draw("same");
00329 
00330    // draw yellow shaded region between min and max
00331    TH1D * hist_shaded = this->GetSubHisto(min,max,TString::Format("%s_sub_%d",fHistogram->GetName(),BCLog::GetHIndex()));
00332    hist_shaded -> SetFillStyle(1001);
00333    hist_shaded -> SetFillColor(kYellow);
00334 
00335    // draw shaded histogram
00336    hist_shaded -> Draw("same");
00337 
00338    gPad->RedrawAxis();
00339 
00340    // draw triangle for mode
00341    TPolyLine * tmax;
00342 
00343    double dx = 0.01*(xmax-xmin);
00344    double dy = 0.04*(ysize);
00345    y0+=dy/5.;
00346    double tmax_x[] = {x0, x0 + dx, x0 - dx, x0};
00347    double tmax_y[] = {y0, y0 + dy, y0 + dy, y0};
00348    tmax = new TPolyLine(4,tmax_x,tmax_y);
00349    tmax->SetLineColor(kRed);
00350    tmax->SetLineWidth(1);
00351    tmax->SetFillColor(kRed);
00352    tmax->Draw();
00353    tmax->Draw("f");
00354 
00355    // draw diamond for mean
00356    TPolyLine * tmean;
00357 
00358    y1+=dy/5.;
00359 // double tmean_x[] = {x1, x1 + dx, x1 - dx, x1};
00360 // double tmean_y[] = {y1, y1 + dy, y1 + dy, y1};
00361    double tmean_x[] = {x1, x1 + dx, x1 , x1 - dx, x1};
00362    double tmean_y[] = {y1, y1 + dy/2., y1 + dy, y1 + dy/2., y1};
00363    tmean = new TPolyLine(5,tmean_x,tmean_y);
00364    tmean->SetLineColor(kBlue);
00365 // tmean->SetFillColor(kWhite);
00366    tmean->SetLineWidth(1);
00367 // tmean->SetLineStyle(1);
00368    tmean->Draw();
00369 
00370 /*
00371    // draw arrow for median
00372    TPolyLine * tmed;
00373    TPolyLine * tmed2;
00374 
00375 // y2+=dy/5.;
00376    y2=0.+dy/5.;
00377    double tmed_x[] = {x2 + dx, x2, x2 - dx};
00378    double tmed_y[] = {y2 + dy, y2, y2 + dy};
00379    double tmed2_x[] = {x2, x2};
00380    double tmed2_y[] = {y2, y2 + dy*2.};
00381    tmed = new TPolyLine(3,tmed_x,tmed_y);
00382    tmed2 = new TPolyLine(2,tmed2_x,tmed2_y);
00383    tmed->SetLineColor(kGreen+2);
00384    tmed->SetLineWidth(1);
00385    tmed2->SetLineColor(kGreen+2);
00386    tmed2->SetLineWidth(1);
00387    tmed->Draw();
00388    tmed2->Draw();
00389 */
00390    // draw dashed line for median
00391    TLine * line;
00392    line = new TLine();
00393    line -> SetLineStyle(2);
00394    line -> SetLineColor(kGreen+2);
00395    line -> DrawLine(x2, 0., x2, y2);
00396 
00397 
00398    // write mode location and shaded band
00399 
00400    // format of the number
00401    double delta_max = fmax(fabs(max-x1),fabs(x1-min));
00402 
00403    int sd = 2 + (int)log10(fabs(x1/delta_max));
00404 
00405    if( (int)log10(x1) > (int)log10(delta_max) )
00406       sd++;
00407 
00408    TLatex * tmax_text = new TLatex();
00409    tmax_text->SetTextSize(0.045);
00410    tmax_text->SetTextFont(62);
00411    tmax_text->SetTextAlign(22); // center
00412 // tmax_text->SetTextAlign(13); // top-left
00413 
00414    double xprint=(xmax+xmin)/2.;
00415 // double xprint=xmin+(xmax-xmin)/20.;
00416    double yprint=ysize*(1-1.4*tmax_text->GetTextSize());
00417 
00418    if(fabs(limit)<50) // just to ensure there's some sense in the number
00419       tmax_text->DrawLatex(xprint,yprint,
00420          TString::Format( TString::Format("%%s^{med} = %%.%dg ^{+%%.2g}_{ -%%.2g}",sd),
00421             fHistogram->GetXaxis()->GetTitle(), x2, max-x2, x2-min));
00422 
00423    else if (limit<0)
00424       tmax_text->DrawLatex(xprint,yprint,
00425          TString::Format( TString::Format("%%s (%d%%%% prob.) < %%.4g",-(int)limit),
00426             fHistogram->GetXaxis()->GetTitle(), max));
00427 
00428    else if (limit>0)
00429       tmax_text->DrawLatex(xprint,yprint,
00430          TString::Format( TString::Format("%%s (%d%%%% prob.) > %%.4g",(int)limit),
00431             fHistogram->GetXaxis()->GetTitle(), min));
00432 
00433 /*
00434    TLegend * leg = new TLegend(.61,.7,.9,.88);
00435    leg -> SetFillColor(kWhite);
00436    leg -> SetFillStyle(0);
00437    leg -> SetBorderSize(0);
00438    leg -> AddEntry(line,"Median", "l");
00439    TH1D * hh0 = new TH1D("hh0","",10,0,10);
00440    hh0->SetMarkerStyle(23);
00441    hh0->SetMarkerColor(kBlue);
00442    hh0->SetMarkerSize(1);
00443    leg -> AddEntry(hh0,"Mean","p");
00444    TH1D * hh1 = new TH1D("hh1","",10,0,10);
00445    hh1->SetMarkerStyle(23);
00446    hh1->SetMarkerColor(kRed);
00447    hh1->SetMarkerSize(1);
00448    leg -> AddEntry(hh1,"Mode","p");
00449    leg -> AddEntry(hist_shaded, "Central 68%", "f");
00450    leg ->Draw();
00451 */
00452 
00453 }
00454 
00455 // ---------------------------------------------------------
00456 
00457 void BCH1D::DrawSmallest(double mode, double prob, bool drawmean)
00458 {
00459    // prepare triangle for mode
00460    TPolyLine * tmax;
00461 
00462    double x0 = mode;
00463    double y0 = fHistogram->GetBinContent( fHistogram->FindBin(mode) );
00464    double xmin = fHistogram->GetXaxis()->GetXmin();
00465    double xmax = fHistogram->GetXaxis()->GetXmax();
00466    double ysize = 1.3 * fHistogram -> GetMaximum();
00467 
00468    double x1 = fHistogram->GetMean();
00469    double y1 = fHistogram->GetBinContent( fHistogram->FindBin(x1) );
00470 
00471    double x2 = this -> GetQuantile(.5); // median
00472    double y2 = fHistogram->GetBinContent( fHistogram->FindBin(x2) );
00473 
00474    double dx = 0.01*(xmax-xmin);
00475    double dy = 0.04*(ysize);
00476    double tmax_x[] = {x0, x0 + dx, x0 - dx, x0};
00477    double tmax_y[] = {y0, y0 + dy, y0 + dy, y0};
00478    tmax = new TPolyLine(4,tmax_x,tmax_y);
00479    tmax->SetLineColor(kRed);
00480    tmax->SetFillColor(kRed);
00481 
00482    // draw histogram with axes first
00483    TH2D * hsc = new TH2D(
00484          TString::Format("h2scale_%s_%d",fHistogram->GetName(),BCLog::GetHIndex()),"",
00485          10, xmin, xmax, 10, 0., ysize);
00486    hsc -> SetStats(0);
00487    hsc -> SetXTitle(fHistogram->GetXaxis()->GetTitle());
00488    hsc -> SetYTitle(fHistogram->GetYaxis()->GetTitle());
00489    hsc -> Draw();
00490 
00491    // histogram to be filled with band
00492    TH1D * hist_temp1 = (TH1D*) fHistogram -> Clone();
00493    hist_temp1 -> Scale(1.0/fHistogram -> Integral("width"));
00494    hist_temp1 -> SetFillColor(kYellow);
00495    hist_temp1 -> SetFillStyle(1001);
00496 
00497    // temporary histogram
00498    TH1D * hist_temp2 = (TH1D*) fHistogram -> Clone();
00499    hist_temp2 -> Scale(1.0/fHistogram -> Integral("width"));
00500 
00501    // clear content
00502    hist_temp1 -> Reset();
00503 
00504    // loop over original histogram and copy bin untils integral equals
00505    // "prob"
00506    prob /= 100.;
00507 
00508    // integral
00509    double sum = 0.0;
00510 
00511    int n=0;
00512    int nbins=fHistogram->GetNbinsX();
00513 
00514    // loop
00515    while (sum < prob && n < nbins)
00516    {
00517       n++;
00518 
00519       // find bin with maximum
00520       int bin = hist_temp2 -> GetMaximumBin();
00521 
00522       // copy bin to new histogram
00523       double val = hist_temp2 -> GetBinContent(bin);
00524       hist_temp1 -> SetBinContent(bin, val);
00525 
00526       // remove maximum from temporary histogram
00527       hist_temp2 -> SetBinContent(bin, 0.0);
00528 
00529       // integrate by summing
00530       sum += val * hist_temp2->GetBinWidth(bin);
00531    }
00532 
00533    // scale histogram
00534    hist_temp1 -> Scale(fHistogram -> Integral("width"));
00535 
00536    // draw histograms
00537    fHistogram -> Draw("same");
00538    hist_temp1 -> Draw("same");
00539 
00540    // draw triangle for mode
00541    tmax->Draw("f");
00542 
00543    if(drawmean)
00544    {
00545       // draw triangle for mean
00546       // draw diamond for mean
00547       TPolyLine * tmean;
00548 
00549       y1+=dy/5.;
00550 //    double tmean_x[] = {x1, x1 + dx, x1 - dx, x1};
00551 //    double tmean_y[] = {y1, y1 + dy, y1 + dy, y1};
00552       double tmean_x[] = {x1, x1 + dx, x1 , x1 - dx, x1};
00553       double tmean_y[] = {y1, y1 + dy/2., y1 + dy, y1 + dy/2., y1};
00554       tmean = new TPolyLine(5,tmean_x,tmean_y);
00555       tmean->SetLineColor(kBlue);
00556 //    tmean->SetFillColor(kWhite);
00557       tmean->SetLineWidth(1);
00558 //    tmean->SetLineStyle(1);
00559       tmean->Draw();
00560 
00561       // draw dashed line for median
00562       TLine * line;
00563       line = new TLine();
00564       line -> SetLineStyle(2);
00565       line -> SetLineColor(kGreen+2);
00566       line -> DrawLine(x2, 0., x2, y2);
00567    }
00568 
00569    // free memory
00570    delete hist_temp2;
00571 }
00572 
00573 // ---------------------------------------------------------
00574 
00575 double BCH1D::GetSmallestInterval(double & min, double & max, double content)
00576 {
00577    if(content<=0. || content>= 1.)
00578       return 0.;
00579 
00580    int nbins = fHistogram -> GetNbinsX();
00581 
00582    double factor = fHistogram->Integral("width");
00583    if(factor==0)
00584       return 0.;
00585 
00586    // normalize if not yet done
00587    fHistogram->Scale(1./factor);
00588 
00589    double xmin = fHistogram->GetXaxis()->GetXmin();
00590    double xmax = fHistogram->GetXaxis()->GetXmax();
00591    double width = xmax - xmin;
00592 
00593    double xdown=xmin;
00594    double xup=xmax;
00595 
00596    int ndiv = 100;
00597    if(nbins<100)
00598       ndiv = 1000;
00599    if(nbins>1000)
00600       ndiv = 10;
00601 
00602    int warn=0;
00603 
00604    double integral_out=0.;
00605 
00606    // loop through the bins
00607    for(int i=1;i<nbins+1;i++)
00608    {
00609       if(fHistogram->Integral(i,nbins,"width") < content)
00610          break;
00611 
00612       // get width of the starting bin
00613       double firstbinwidth = fHistogram->GetBinWidth(i);
00614 
00615       // divide i-th bin in ndiv divisions
00616       // integrate starting at each of these divisions
00617       for(int j=0;j<ndiv;j++)
00618       {
00619          double dxdown = (double)(ndiv-j)/(double)ndiv * firstbinwidth;
00620          xdown = fHistogram->GetBinLowEdge(i) + firstbinwidth - dxdown;
00621          double integral = dxdown * fHistogram->GetBinContent(i);
00622 
00623          if(integral>content)
00624          {
00625             // content fits within one bin
00626             xup = xdown + content / fHistogram->GetBinContent(i);
00627             warn = 1;
00628          }
00629          else
00630             for(int k=i+1;k<nbins+1;k++)
00631             {
00632 
00633                double thisbin = fHistogram->GetBinContent(k) * fHistogram->GetBinWidth(k);
00634 
00635                if(integral+thisbin==content)
00636                {
00637                   xup = fHistogram->GetBinLowEdge(k+1);
00638                   break;
00639                }
00640 
00641                if(integral+thisbin>content)
00642                {
00643                   xup = fHistogram->GetBinLowEdge(k) + (content-integral)/thisbin * fHistogram->GetBinWidth(k);
00644                   integral += thisbin * (xup-fHistogram->GetBinLowEdge(k))/fHistogram->GetBinWidth(k);
00645                   break;
00646                }
00647 
00648                integral += thisbin;
00649             }
00650 
00651          if(integral < content)
00652             continue;
00653 
00654          if(xup - xdown < width)
00655          {
00656             // store necessary information
00657             width = xup - xdown;
00658             xmin  = xdown;
00659             xmax  = xup;
00660             integral_out = integral;
00661          }
00662       }
00663    }
00664 
00665    if(warn)
00666    {
00667       BCLog::Out(BCLog::warning,BCLog::warning,
00668             Form("BCH1D::GetSmallestInterval : The requested content of %d%% fits within one bin.",(int)(content*100)));
00669       BCLog::Out(BCLog::warning,BCLog::warning,
00670             "BCH1D::GetSmallestInterval : MAKE FINER BINNING (OR CHANGE RANGE) !!!");
00671    }
00672 
00673    // restore normalization to state before calling this method
00674    fHistogram->Scale(factor);
00675 
00676    min=xmin;
00677    max=xmax;
00678 
00679    return integral_out;
00680 }
00681 
00682 // ---------------------------------------------------------
00683 
00684 double BCH1D::IntegralWidth(double min, double max)
00685 {
00686    int imin = fHistogram->FindBin(min);
00687    int imax = fHistogram->FindBin(max);
00688 
00689    int nbins = fHistogram->GetNbinsX();
00690 
00691    // if outside of histogram range, return -1.
00692    if ( imin<1 || imin>nbins || imax<1 || imax>nbins )
00693       return -1.;
00694 
00695    if ( imin==imax )
00696       return -1.;
00697 
00698    // swap values if necessary
00699    if (imin>imax)
00700    {
00701       int i=imin;
00702       double x=min;
00703       imin=imax, min=max;
00704       imax=i, max=x;
00705    }
00706 
00707    // calculate first bin
00708    double first = ( fHistogram->GetBinLowEdge(imin+1) - min ) * fHistogram->GetBinContent(imin);
00709 
00710    // calculate last bin
00711    double last = ( max - fHistogram->GetBinLowEdge(imax) ) * fHistogram->GetBinContent(imax);
00712 
00713    // calculate inbetween
00714    double inbetween=0.;
00715    if(imax-imin>1)
00716       inbetween = fHistogram->Integral(imin+1, imax-1, "width");
00717 
00718    return first + inbetween + last;
00719 }
00720 
00721 // ---------------------------------------------------------
00722 
00723 TH1D * BCH1D::GetSubHisto(double min, double max, const char * name)
00724 {
00725    if(min>max)
00726    {
00727       double t=min;
00728       min=max;
00729       max=t;
00730    }
00731 
00732    int imin = fHistogram->FindBin(min);
00733    int imax = fHistogram->FindBin(max);
00734 
00735    int nbins = fHistogram->GetNbinsX();
00736    double xmin = fHistogram->GetXaxis()->GetXmin();
00737    double xmax = fHistogram->GetXaxis()->GetXmax();
00738 
00739    if( min==max || (min<=xmin && max>=xmax) )
00740    {
00741       TH1D * h0 = (TH1D*) fHistogram->Clone(name);
00742       return h0;
00743    }
00744 
00745    if (imin<1)
00746    {
00747       imin=1;
00748       min=xmin;
00749    }
00750    if (imax>nbins)
00751    {
00752       imax=nbins;
00753       max=xmax;
00754    }
00755 
00756    double * xb = new double[nbins+3]; // nbins+1 original bin edges + 2 new bins
00757    int n=0; // counter
00758 
00759    int domin=1;
00760    int domax=1;
00761 
00762    for (int i=1;i<nbins+2;i++)
00763    {
00764       double x0 = fHistogram->GetBinLowEdge(i);
00765 
00766       if (min<x0 && domin)
00767       {
00768          xb[n++]=min;
00769          domin=0;
00770       }
00771       else if (min==x0)
00772          domin=0;
00773 
00774       if (max<x0 && domax)
00775       {
00776          xb[n++]=max;
00777          domax=0;
00778       }
00779       else if (max==x0)
00780          domax=0;
00781 
00782       xb[n++]=x0;
00783    }
00784 
00785    // now define the new histogram
00786    TH1D * h0 = new TH1D(name,"",n-1,xb);
00787    for(int i=1;i<n;i++)
00788    {
00789       double x0 = h0->GetBinCenter(i);
00790       if(x0<min || x0>max)
00791          continue;
00792 
00793       int bin=fHistogram->FindBin(x0);
00794       h0->SetBinContent(i, fHistogram->GetBinContent(bin));
00795    }
00796 
00797    return h0;
00798 }
00799 
00800 // ---------------------------------------------------------
00801 
00802 TH1D * BCH1D::GetSmallestIntervalHistogram(double level)
00803 {
00804    // create a new histogram which will be returned and all yellow
00805    TH1D * hist_yellow = (TH1D*) fHistogram -> Clone();
00806    hist_yellow -> Reset();
00807    hist_yellow -> SetFillStyle(1001);
00808    hist_yellow -> SetFillColor(kYellow);
00809 
00810    // copy a temporary histogram
00811    TH1D * hist_temp = (TH1D*) fHistogram -> Clone(TString::Format("%s_%i",fHistogram->GetName(),BCLog::GetHIndex()));
00812    double factor = hist_temp -> Integral("");
00813 
00814    if(factor == 0)
00815       return 0;
00816 
00817    hist_temp -> Scale(1. / factor);
00818 
00819    // here's the algorithm:
00820    // 1. find the maximum bin in the temporary histogram and copy it to
00821    // the yellow histogram.
00822    // 2. remove this bin from the temporary histogram.
00823    // 3. repeat this until a total of "level" probability is copied to
00824    // the yellow histogram.
00825 
00826    double sumprob = 0.0;
00827 
00828    while (sumprob < level)
00829    {
00830       // find maximum bin and value
00831       int bin = hist_temp -> GetMaximumBin();
00832       double value = hist_temp -> GetMaximum();
00833 
00834       // copy "1" into the corresponding bin in the yellow histogram
00835       hist_yellow -> SetBinContent(bin, 1.0);
00836 
00837       // set the bin value in the temporary histogram to zero
00838       hist_temp -> SetBinContent(bin, 0.0);
00839 
00840       // increase probability sum
00841       sumprob += value;
00842    }
00843 
00844    delete hist_temp;
00845 
00846    return hist_yellow;
00847 }
00848 
00849 // ---------------------------------------------------------
00850 
00851 std::vector <double> BCH1D::GetSmallestIntervals(double content)
00852 {
00853    std::vector <double> v;
00854 
00855    TH1D * hist = this -> GetSmallestIntervalHistogram(content);
00856 
00857    int nbins = hist -> GetNbinsX();
00858    int ninter = 0;
00859    int lastbin = -1;
00860 
00861    double max = -1;
00862    double localmax = -1;
00863    double localmaxpos = -1;
00864    double localint = 0;
00865    bool flag = false;
00866 
00867    for (int i = 1; i <= nbins; ++i)
00868    {
00869       // interval starts here
00870       if (!flag && hist -> GetBinContent(i) > 0.)
00871       {
00872          flag = true;
00873          v.push_back(hist -> GetBinLowEdge(i));
00874 
00875          // remember start position of the interval
00876          lastbin = i;
00877 
00878          // increase number of intervals
00879          ninter++;
00880 
00881          // reset local maximum
00882          localmax = fHistogram -> GetBinContent(i);
00883          localmaxpos = hist -> GetBinLowEdge(i);
00884 
00885          // reset local integral
00886          localint = 0;
00887       }
00888 
00889       // interval stops here
00890       if ((flag && !(hist -> GetBinContent(i) > 0.)) || (flag && i == nbins))
00891       {
00892          flag = false;
00893          v.push_back(hist -> GetBinLowEdge(i) + hist -> GetBinWidth(i));
00894 
00895          // set right bin to maximum if on edge
00896          if (i == nbins && localmax < fHistogram -> GetBinContent(i))
00897             localmaxpos = hist -> GetBinCenter(i) + 0.5 * hist -> GetBinWidth(i);
00898 
00899          // find the absolute maximum
00900          if (localmax > max)
00901             max = localmax;
00902 
00903          // save local maximum
00904          v.push_back(localmax);
00905          v.push_back(localmaxpos);
00906 
00907          // save local integral
00908          v.push_back(localint);
00909       }
00910 
00911       // find local maximum
00912       if (i < nbins && localmax < fHistogram -> GetBinContent(i))
00913       {
00914          localmax = fHistogram -> GetBinContent(i);
00915          localmaxpos = hist -> GetBinCenter(i);
00916       }
00917 
00918       // increase area
00919       localint += fHistogram -> GetBinContent(i) / fHistogram -> Integral();
00920    }
00921 
00922    // rescale absolute heights to relative heights
00923    for (int i = 0; i < ninter; ++i)
00924       v[i*5+2] = v.at(i*5+2) / max;
00925 
00926    return v;
00927 }
00928 
00929 // ---------------------------------------------------------

Generated on Mon Aug 30 2010 22:14:54 for Bayesian Analysis Toolkit by  doxygen 1.7.1