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

BCTemplateFitter.cxx

Go to the documentation of this file.
00001 #include <iostream>
00002 
00003 #include <TROOT.h>
00004 #include <TH2D.h>
00005 #include <THStack.h>
00006 #include <TCanvas.h>
00007 #include <TLegend.h>
00008 #include <TMath.h>
00009 #include <TRandom3.h>
00010 #include <TGraphAsymmErrors.h>
00011 #include <TPostScript.h>
00012 #include <TPad.h>
00013 #include <TLine.h>
00014 #include <TDirectory.h>
00015 
00016 #include <BAT/BCMath.h>
00017 #include <BAT/BCLog.h>
00018 #include <BAT/BCH1D.h>
00019 
00020 #include "BCTemplateFitter.h"
00021 
00022 // ---------------------------------------------------------
00023 BCTemplateFitter::BCTemplateFitter()
00024    : BCModel()
00025    , fFlagFixNorm(false)
00026    , fFlagPhysicalLimits(true)
00027    , fUncertaintyHistogramExp(0)
00028    , fUncertaintyHistogramObsPosterior(0)
00029    , fNorm(-1)
00030    , fNBins(-1)
00031    , fXmin(1.)
00032    , fXmax(0.)
00033 {
00034 }
00035 
00036 // ---------------------------------------------------------
00037 BCTemplateFitter::BCTemplateFitter(const char * name)
00038    : BCModel(name)
00039    , fFlagFixNorm(false)
00040    , fFlagPhysicalLimits(true)
00041    , fUncertaintyHistogramExp(0)
00042    , fUncertaintyHistogramObsPosterior(0)
00043    , fNorm(-1)
00044    , fNBins(-1)
00045    , fXmin(1.)
00046    , fXmax(0.)
00047 {
00048 }
00049 
00050 // ---------------------------------------------------------
00051 BCTemplateFitter::~BCTemplateFitter()
00052 {
00053    if (fUncertaintyHistogramExp)
00054       delete fUncertaintyHistogramExp;
00055 
00056    if (fUncertaintyHistogramObsPosterior)
00057       delete fUncertaintyHistogramObsPosterior;
00058 }
00059 
00060 // ---------------------------------------------------------
00061 double BCTemplateFitter::LogLikelihood(std::vector <double> parameters)
00062 {
00063    double logprob = 0.;
00064 
00065    // fix the parameters of all functions
00066    FixTemplateFunctions(parameters);
00067 
00068    // loop over bins
00069    for (int ibin = 1; ibin <= fNBins; ++ibin) {
00070       // get expectation
00071       double nexp = Expectation(ibin, parameters);
00072 
00073       // get data
00074       double ndata = fHistData.GetBinContent(ibin);
00075 
00076       // add Poisson term
00077       logprob += BCMath::LogPoisson(ndata, nexp);
00078    }
00079 
00080    // return log likelihood
00081    return logprob;
00082 }
00083 
00084 // ---------------------------------------------------------
00085 double BCTemplateFitter::Expectation(int binnumber, std::vector<double>& parameters)
00086 {
00087    // initialize number of expected events with 0
00088    double nexp = 0;
00089 
00090    // get number of templates
00091    int ntemplates = GetNTemplates();
00092 
00093    // loop over all templates
00094    for (int itemp = 0; itemp < ntemplates; ++itemp) {
00095 
00096       // get template index
00097       int templateindex = fTemplateParIndexContainer.at(itemp);
00098 
00099       // calculate efficiency for this template and bin
00100       double efficiency = TemplateEfficiency(itemp, binnumber, parameters);
00101 
00102       // calculate probability for this template and bin
00103       double probability = TemplateProbability(itemp, binnumber, parameters);
00104 
00105       // calculate expectation value for this template
00106       nexp += parameters.at(templateindex)
00107          * efficiency
00108          * probability;
00109    }
00110    
00111    // check that expectation is larger or equal to zero
00112    if (nexp < 0) {
00113       BCLog::OutWarning("BCTemplateFitter::Expectation : Expectation value smaller than 0. Force it to be 0.");
00114       nexp = 0;
00115    }
00116    
00117    return nexp;
00118 }
00119 
00120 // ---------------------------------------------------------
00121 double BCTemplateFitter::TemplateEfficiency(int templatenumber, int binnumber, std::vector<double>& parameters)
00122 {
00123    // get number of sources of systematic uncertainties
00124    int nsyst = GetNSystErrors();
00125 
00126    // get efficiency index
00127    int effindex = fEffParIndexContainer.at(templatenumber);
00128    
00129    // get efficiency for the bin
00130    double efficiency = fEffHistogramContainer.at(templatenumber).GetBinContent(binnumber);
00131    
00132    // modify efficiency by uncertainty
00133    double efferr = fEffErrHistogramContainer.at(templatenumber).GetBinContent(binnumber);
00134    
00135    // check that efficiency error is larger than 0
00136    // add efficiency contribution from template
00137    if (efferr > 0)
00138       efficiency += parameters.at(effindex) * efferr;
00139    
00140    // loop over sources of systematic uncertainties
00141    for (int isyst = 0; isyst < nsyst; ++isyst) {
00142       // get parameter index
00143       int systindex = fSystErrorParIndexContainer.at(isyst);
00144       
00145       // add efficiency
00146       double deff = fSystErrorHistogramContainer.at(isyst).at(templatenumber).GetBinContent(binnumber);
00147       
00148       // check that efficiency error is larger than 0
00149       // add efficiency contribution from systematic
00150       if (deff > 0)
00151          efficiency += parameters.at(systindex) * deff;
00152    }
00153    
00154    // make sure efficiency is between 0 and 1
00155    if (efficiency < 0.)
00156       efficiency = 0.;
00157    if (efficiency > 1.)
00158       efficiency = 1.;
00159    
00160    // return efficiency
00161    return efficiency;
00162 }
00163 
00164 // ---------------------------------------------------------
00165 double BCTemplateFitter::TemplateProbability(int templatenumber, int binnumber, std::vector<double>& parameters)
00166 {
00167    // initialize probability
00168    double probability = 0;
00169 
00170    // get template type
00171    int templatetype = fTemplateTypeContainer[templatenumber];
00172 
00173    // calculate probability from template histogram
00174    if (templatetype == 0) {
00175       int index = GetHistogramIndex(templatenumber);
00176       
00177       probability = fTemplateHistogramContainer.at(index).GetBinContent(binnumber);
00178    }
00179    else if (templatetype == 1) {
00180       double bincenter = fHistData.GetBinCenter(binnumber);
00181       int index = GetFunctionIndex(templatenumber);
00182       probability = fTemplateFunctionContainer.at(index).Eval(bincenter);
00183    }
00184 
00185    // return probability
00186    return probability;
00187 }
00188 
00189 // ---------------------------------------------------------
00190 int BCTemplateFitter::FixTemplateFunctions(std::vector<double>& parameters)
00191 {
00192    // get number of templates
00193    int ntemplates = GetNTemplates();
00194 
00195    // fix all function parameters
00196    for (int i = 0; i < ntemplates; ++i) {
00197       if (fTemplateTypeContainer.at(i) == 1) {
00198          int index = GetFunctionIndex(i);
00199          TF1* func = &(fTemplateFunctionContainer.at(index));
00200          int npar = func->GetNpar();
00201          int tempindex = GetParIndexTemplate(i);
00202          for (int i = 0; i < npar; ++i) {
00203             func->FixParameter(i, parameters.at(tempindex+2+i));
00204          }
00205       }
00206    }
00207 
00208    // no errors
00209    return 1;
00210 }
00211 // ---------------------------------------------------------
00212 int BCTemplateFitter::SetData(const TH1D& hist)
00213 {
00214    // compare previous data set with this one
00215    if (fHistData.Integral() > 0.) {
00216       if (CompareHistogramProperties(fHistData, hist) != 1) {
00217          BCLog::OutError("BCTemplateFitter::SetData : Data and this histogram properties are incompatible.");
00218          return 0;
00219       }
00220    }
00221    
00222    // create histogram
00223    fNBins = hist.GetNbinsX();
00224    fXmin = (hist.GetXaxis())->GetXmin();
00225    fXmax = (hist.GetXaxis())->GetXmax();
00226    fHistData = TH1D(hist);
00227 
00228    // copy histogram content
00229    for (int i = 1; i <= fNBins; ++i)
00230       fHistData.SetBinContent(i, hist.GetBinContent(i));
00231 
00232    // set histogram style
00233    fHistData.SetXTitle((hist.GetXaxis())->GetTitle());
00234    fHistData.SetYTitle((hist.GetYaxis())->GetTitle());
00235    fHistData.SetMarkerStyle(20);
00236    fHistData.SetMarkerSize(1.1);
00237    fHistData.SetStats(kFALSE);
00238 
00239    // calculate norm
00240    fNorm = hist.Integral();
00241 
00242    // no errors
00243    return 1;
00244 }
00245 
00246 // ---------------------------------------------------------
00247 int BCTemplateFitter::AddTemplate(TH1D hist, const char * name, double Nmin, double Nmax)
00248 {
00249    // check if histogram if filled
00250    if (hist.Integral() <= 0.) {
00251       BCLog::OutError("BCTemplateFitter::AddTemplate : Normalization is zero or less than that.");
00252       return 0;
00253    }
00254 
00255    // compare template properties with data
00256    if (fHistData.Integral() > 0.) {
00257       if (CompareHistogramProperties(fHistData, hist) != 1) {
00258          BCLog::OutError("BCTemplateFitter::AddTemplate : Data and template histogram properties are incompatible.");
00259          return 0;
00260       }
00261    }
00262    else { 
00263       SetData( TH1D("", "", hist.GetNbinsX(), hist.GetXaxis()->GetXmin(), hist.GetXaxis()->GetXmax()) );
00264    } 
00265 
00266    // check if prior makes sense
00267    if (fFlagPhysicalLimits && Nmin < 0)
00268       Nmin = 0;
00269 
00270    if (Nmin > Nmax) {
00271       BCLog::OutError("BCTemplateFitter::AddTemplate : Lower limit exceeds upper limit.");
00272       return 0;
00273    }
00274 
00275    // get number of templates
00276    int ntemplates = GetNTemplates();
00277 
00278    // set histogram color and style
00279    hist.SetFillColor(2 + ntemplates);
00280    hist.SetFillStyle(1001);
00281    hist.SetStats(kFALSE);
00282 
00283    // scale histogram
00284    hist.Scale(1.0 / hist.Integral());
00285 
00286    // check if template is consistent with other templates
00287    if (ntemplates > 0)
00288       if (!CompareHistogramProperties(fTemplateHistogramContainer.at(0), hist)) {
00289          BCLog::OutError("BCTemplateFitter::AddTemplate : Properties of template histogram is not compatible with older template histograms.");
00290          return 0;
00291       }
00292 
00293    // histograms
00294    TH1D histsysterror = TH1D(fHistData);
00295 
00296    // fill histograms
00297    for (int i = 1; i <= fNBins; ++i)
00298       histsysterror.SetBinContent(i, 0.);
00299 
00300    // get parameter index
00301    int parindex = GetNParameters();
00302 
00303    // add a parameter for the expectation value
00304    AddParameter(name, Nmin, Nmax);
00305 
00306    // get efficiency parameter index
00307    int effindex = GetNParameters();
00308 
00309    // add a parameter for the efficiency
00310    AddParameter(Form("Efficiency_%s", name), -5.0, 5.0);
00311 
00312    // add histogram, name and index to containers
00313    fTemplateHistogramContainer.push_back(hist);
00314    fHistogramIndexContainer.push_back(GetNTemplatesType(0));
00315    fFunctionIndexContainer.push_back(-1);
00316    fTemplateTypeContainer.push_back(0);
00317    fTemplateNameContainer.push_back(name);
00318    fTemplateParIndexContainer.push_back(parindex);
00319    fEffParIndexContainer.push_back(effindex);
00320 
00321    // set efficiency histograms to one without uncertainty
00322    fEffHistogramContainer.push_back(TH1D());
00323    fEffErrHistogramContainer.push_back(TH1D());
00324    SetTemplateEfficiency(name, 1., 0.);
00325 
00326    // set prior for efficiencies
00327    SetPriorGauss(Form("Efficiency_%s", name), 0., 1.);
00328 
00329    // add systematic uncertainties
00330    for (int i = 0; i < GetNSystErrors(); ++i) {
00331       std::vector <TH1D> histvector = fSystErrorHistogramContainer.at(i);
00332       histvector.push_back(histsysterror);
00333    }
00334 
00335    // successfully added histogram to container
00336    return 1;
00337 }
00338 
00339 // ---------------------------------------------------------
00340 int BCTemplateFitter::AddTemplate(TF1 func, const char * name, double Nmin, double Nmax, const char* options)
00341 {
00342    // check if prior makes sense
00343    if (fFlagPhysicalLimits && Nmin < 0)
00344       Nmin = 0;
00345 
00346    if (Nmin > Nmax) {
00347       BCLog::OutError("BCTemplateFitter::AddTemplate : Lower limit exceeds upper limit.");
00348       return 0;
00349    }
00350 
00351    // get parameter index
00352    int parindex = GetNParameters();
00353 
00354    // add a parameter for the expectation value
00355    AddParameter(name, Nmin, Nmax);
00356 
00357    // get efficiency parameter index
00358    int effindex = GetNParameters();
00359 
00360    // add a parameter for the efficiency
00361    AddParameter(Form("Efficiency_%s", name), -5.0, 5.0);
00362 
00363    // add all additional parameters
00364    int npar = func.GetNpar(); 
00365 
00366    for (int i = 0; i < npar; ++i) {
00367       double parmin;
00368       double parmax;
00369       func.GetParLimits(i, parmin, parmax);
00370       AddParameter(func.GetParName(i), parmin, parmax);
00371    }
00372 
00373 
00374    // add histogram, name and index to containers
00375    fTemplateFunctionContainer.push_back(func);
00376    fHistogramIndexContainer.push_back(-1);
00377    fFunctionIndexContainer.push_back(GetNTemplatesType(1));
00378    fTemplateTypeContainer.push_back(1);
00379    fTemplateNameContainer.push_back(name);
00380    fTemplateParIndexContainer.push_back(parindex);
00381    fEffParIndexContainer.push_back(effindex);
00382 
00383    // set efficiency histograms to one without uncertainty
00384    fEffHistogramContainer.push_back(TH1D());
00385    fEffErrHistogramContainer.push_back(TH1D());
00386    SetTemplateEfficiency(name, 1., 0.);
00387 
00388    // set prior for efficiencies
00389    SetPriorGauss(Form("Efficiency_%s", name), 0., 1.);
00390 
00391    // successfully added histogram to container
00392    return 1;
00393 }
00394 
00395 // ---------------------------------------------------------
00396 int BCTemplateFitter::AddSystError(const char * errorname, const char * errtype)
00397 {
00398    // define parameter range
00399    double dx = 1.0;
00400 
00401    // check error type
00402    if (std::string(errtype) == std::string("gauss"))
00403       dx = 5.0;
00404    else if (std::string(errtype) == std::string("flat"))
00405       dx = 1.0;
00406    else {
00407       BCLog::OutError("BCTemplateFitter::AddSystError : Unknown error type.");
00408       return 0;
00409    }
00410 
00411    // add a parameter for the expectation value
00412    AddParameter(Form("%s", errorname), -dx, dx);
00413 
00414    // set prior for systematics
00415    SetPriorGauss(Form("%s", errorname), 0., 1.);
00416 
00417    // add name and index to containers
00418    fSystErrorNameContainer.push_back(errorname);
00419    fSystErrorParIndexContainer.push_back(GetNParameters()-1);
00420 
00421    // create histogram
00422    TH1D hist = TH1D(fHistData);
00423 
00424    // fill histograms
00425    for (int i = 1; i <= fNBins; ++i)
00426       hist.SetBinContent(i, 0.);
00427 
00428    // get number of templates
00429    int n = GetNTemplates();
00430 
00431    // define vector of histograms
00432    std::vector<TH1D> histvector;
00433 
00434    // add histograms
00435    for (int i = 0; i < n; ++i)
00436       histvector.push_back(hist);
00437 
00438    // add histogram vector
00439    fSystErrorHistogramContainer.push_back(histvector);
00440 
00441    // add error type to container
00442    fSystErrorTypeContainer.push_back(std::string(errtype));
00443 
00444    // no error
00445    return 1;
00446 }
00447 
00448 // ---------------------------------------------------------
00449 int BCTemplateFitter::SetTemplateSystError(const char * errorname, const char * templatename, TH1D parerror)
00450 {
00451    // get error index
00452    int errindex = GetIndexSystError(errorname);
00453 
00454    // check parameter index
00455    if (errindex < 0) {
00456       BCLog::OutError("BCTemplateFitter::SetTemplateSystError : Did not find parameter.");
00457       return 0;
00458    }
00459 
00460    // get template index
00461    int tempindex = GetIndexTemplate(templatename);
00462 
00463    // check index
00464    if (tempindex < 0) {
00465       BCLog::OutError("BCTemplateFitter::SetTemplateSystError : Could not find template.");
00466       return 0;
00467    }
00468 
00469    // set style
00470    parerror.SetStats(kFALSE);
00471 
00472    // set histogram
00473    (fSystErrorHistogramContainer.at(errindex))[tempindex] = parerror;
00474 
00475    // no error
00476    return 1;
00477 }
00478 
00479 // ---------------------------------------------------------
00480 int BCTemplateFitter::Initialize()
00481 {
00482    // check data integral
00483    if (fHistData.Integral() <= 0) {
00484       BCLog::OutError("BCTemplateFitter::Initialize : Normalization of data histogram is zero or less than that.");
00485       return 0;
00486    }
00487 
00488    // create histograms for uncertainty determination
00489    // double maximum = 1.5 * fHistData.GetMaximum();
00490    double minimum = TMath::Max(0., fHistData.GetMinimum() - 5.*sqrt(fHistData.GetMinimum()));
00491    double maximum = fHistData.GetMaximum() + 5.*sqrt(fHistData.GetMaximum());
00492    
00493    Double_t a[fHistData.GetNbinsX()+1];
00494    for (int i = 0; i < fHistData.GetNbinsX()+1; ++i) {
00495       a[i] = fHistData.GetXaxis()->GetBinLowEdge(i+1);
00496    }
00497 
00498    fUncertaintyHistogramExp = new TH2D(
00499          TString::Format("UncertaintyExp_%i", BCLog::GetHIndex()), "",
00500          //       fHistData.GetNbinsX(), fHistData.GetXaxis()->GetXbins()->GetArray(),
00501          fHistData.GetNbinsX(), a,
00502          //       fHistData.GetNbinsX(), fHistData.GetXaxis()->GetXmin(), fHistData.GetXaxis()->GetXmax(),
00503          100, minimum, maximum);
00504 
00505    fUncertaintyHistogramObsPosterior = new TH2D(
00506          TString::Format("UncertaintyObsPosterior_%i", BCLog::GetHIndex()), "",
00507          fHistData.GetNbinsX(), a,
00508          //       fHistData.GetNbinsX(), fHistData.GetXaxis()->GetXbins()->GetArray(),       //       fHistData.GetNbinsX(), fHistData.GetXaxis()->GetXmin(), fHistData.GetXaxis()->GetXmax(),
00509          int(maximum) + 1, -0.5, double(int(maximum))+0.5);
00510 
00511    // create histogram containing the normalization
00512    double xmin = 0;
00513    double xmax = 0;
00514    int ntemplates = GetNTemplates();
00515 
00516    // calculate the limits on the norm from the sum of all parameter
00517    // limits
00518    for (int i = 0; i < ntemplates; ++i) {
00519       // get parameter index
00520       int parindex = fTemplateParIndexContainer.at(i);
00521       int effindex = fEffParIndexContainer.at(i);
00522 
00523       // get parameter
00524       BCParameter * par = this->GetParameter(parindex);
00525       BCParameter * eff = this->GetParameter(effindex);
00526 
00527       // increate limits
00528       xmin += par->GetLowerLimit() * eff->GetLowerLimit();
00529       xmax += par->GetUpperLimit() * eff->GetUpperLimit();
00530    }
00531 
00532    // create new histogram for norm
00533    fHistNorm = TH1D("", ";N_{norm};dN/dN_{norm}", 100, xmin, xmax);
00534 
00535    // no error
00536    return 1;
00537 }
00538 
00539 // ---------------------------------------------------------
00540 int BCTemplateFitter::CalculateRatio(int index, std::vector<int> indices, double rmin, double rmax)
00541 {
00542    // get number of templates
00543    int ntemplates = GetNTemplates();
00544 
00545    // check index
00546    if (index < 0 || index >= ntemplates) {
00547       return 0;
00548    }
00549 
00550    // check indices
00551    for (int i = 0; i < int(indices.size()); ++i) {
00552       if (indices.at(i) < 0 || indices.at(i) >= ntemplates) {
00553          return 0;
00554       }
00555    }
00556 
00557    // create temporary vector
00558    std::vector<int> tempvec;
00559    tempvec.push_back(index);
00560    for (int i = 0; i < int(indices.size()); ++i)
00561       tempvec.push_back(indices.at(i));
00562 
00563    // add ratio
00564    fIndicesRatios1D.push_back(tempvec);
00565 
00566    // get number of ratios
00567    int nratios = int(fHistRatios1D.size());
00568 
00569    // create histogram
00570    double fmin = rmin;
00571    double fmax = rmax;
00572    if (fFlagPhysicalLimits) {
00573       fmin = TMath::Max(rmin, 0.0);
00574       fmax = TMath::Min(1.0, rmax);
00575    }
00576 
00577    TH1D hist_ratio1d(TString::Format("ratio %i", nratios), ";r;p(r|data)", 100, fmin, fmax);
00578    fHistRatios1D.push_back(hist_ratio1d);
00579 
00580    // no error
00581    return 1;
00582 }
00583 
00584 // ---------------------------------------------------------
00585 int BCTemplateFitter::CompareHistogramProperties(TH1D hist1, TH1D hist2)
00586 {
00587    // compare number of bins
00588    if (hist1.GetNbinsX() != hist2.GetNbinsX())
00589       return 0;
00590 
00591    // compare minimum x-values
00592    if (hist1.GetXaxis()->GetXmin() != hist2.GetXaxis()->GetXmin())
00593       return 0;
00594 
00595    // compare maximum x-values
00596    if (hist1.GetXaxis()->GetXmax() != hist2.GetXaxis()->GetXmax())
00597       return 0;
00598 
00599    // conclusion: they have the same properties
00600    return 1;
00601 }
00602 
00603 // ---------------------------------------------------------
00604 void BCTemplateFitter::PrintStack(const char * filename, const char * options)
00605 {
00606    int nbins = fHistData.GetNbinsX();
00607    int ntemplates = GetNTemplates();
00608 
00609    // check options
00610    bool flag_legend = false;
00611    bool flag_error0 = false; // symm. poisson error for data
00612    bool flag_error1 = false; // symm. poisson error for exp.
00613    bool flag_error2 = false; // asymm. poisson error of expectation value
00614    bool flag_error3 = false; // asymm. poisson error of expected no. of events
00615    bool flag_diff   = false; // plot difference between data and expectation below stack plot
00616    bool flag_logx   = false; // plot x-axis in log-scale
00617    bool flag_logy   = false; // plot y-axis in log-scale
00618 
00619    if (std::string(options).find("L") < std::string(options).size())
00620       flag_legend = true;
00621 
00622    if (std::string(options).find("E0") < std::string(options).size())
00623       flag_error0 = true;
00624 
00625    if (std::string(options).find("E1") < std::string(options).size())
00626       flag_error1 = true;
00627 
00628    if (std::string(options).find("E2") < std::string(options).size())
00629       flag_error2 = true;
00630 
00631    if (std::string(options).find("E3") < std::string(options).size())
00632       flag_error3 = true;
00633 
00634    if (std::string(options).find("D") < std::string(options).size())
00635       flag_diff = true;
00636 
00637    if (std::string(options).find("logx") < std::string(options).size())
00638       flag_logx = true;
00639 
00640    if (std::string(options).find("logy") < std::string(options).size())
00641       flag_logy = true;
00642 
00643    // create canvas
00644    TCanvas* c1 = new TCanvas("", "", 700, 700);
00645    c1->cd();
00646    TPad * pad1;
00647    TPad * pad2;
00648 
00649    double fraction_pads = 0.3;
00650 
00651    if(!flag_diff)
00652       fraction_pads=0.0;
00653 
00654    if (flag_diff) {
00655       pad1 = new TPad("pad1", "", 0.0, fraction_pads, 1.0, 1.0);
00656       pad1->SetTopMargin   (0.13/(1.0-fraction_pads));
00657       pad1->SetBottomMargin(0.0);
00658       pad1->SetLeftMargin  (0.15);
00659       pad1->SetRightMargin (0.13);
00660       pad1->SetFillColor(kWhite);
00661       pad2 = new TPad("pad2", "", 0.0, 0.0, 1.0, fraction_pads);
00662       pad2->SetTopMargin   (0.0);
00663       pad2->SetBottomMargin(0.15 / fraction_pads);
00664       pad2->SetLeftMargin  (0.15);
00665       pad2->SetRightMargin (0.13);
00666       pad2->SetFillColor(kWhite);
00667       if (flag_logx) {
00668          pad1->SetLogx();
00669          pad2->SetLogx();
00670       }
00671       if (flag_logy)
00672          pad1->SetLogy();
00673       pad1->Draw();
00674       pad2->Draw();
00675    }
00676    else {
00677       pad1 = new TPad("pad1", "",0.0, 0.0, 1.0, 1.0);
00678       pad1->SetFillColor(kWhite);
00679       pad2 = new TPad();
00680       if (flag_logx)
00681          pad1->SetLogx();
00682       if (flag_logy)
00683          pad1->SetLogy();
00684       pad1->Draw();
00685    }
00686 
00687    pad1->cd();
00688 
00689    // set style and draw data
00690    double ymin = 0.01;
00691    double ymax = 1.1 * (fHistData.GetMaximum() + sqrt(fHistData.GetMaximum()));
00692    fHistData.GetYaxis()->SetRangeUser(ymin, ymax);
00693    fHistData.GetXaxis()->SetNdivisions(505);
00694    if (flag_diff) {
00695       fHistData.GetXaxis()->SetLabelSize(fHistData.GetXaxis()->GetLabelSize()/(1.0-fraction_pads));
00696       fHistData.GetXaxis()->SetLabelOffset(fHistData.GetXaxis()->GetLabelOffset()*(1.0-fraction_pads));
00697       fHistData.GetXaxis()->SetTitleSize(fHistData.GetXaxis()->GetTitleSize()/(1.0-fraction_pads));
00698       fHistData.GetXaxis()->SetTitleOffset(fHistData.GetXaxis()->GetTitleOffset()*(1.0-fraction_pads));
00699       fHistData.GetYaxis()->SetLabelSize(fHistData.GetYaxis()->GetLabelSize()/(1.0-fraction_pads));
00700       fHistData.GetYaxis()->SetLabelOffset(fHistData.GetYaxis()->GetLabelOffset()/(fraction_pads));
00701       fHistData.GetYaxis()->SetTitleSize(fHistData.GetYaxis()->GetTitleSize()/(1.0-fraction_pads));
00702       fHistData.GetYaxis()->SetTitleOffset(fHistData.GetYaxis()->GetTitleOffset()*(1.0-fraction_pads));
00703    }
00704    fHistData.Draw("P");
00705 
00706    // create a histogram with the sum of all contributions
00707    //    TH1D * histsum = (TH1D*) fHistData.Clone("temp");
00708    TH1D * histsum = new TH1D(fHistData);
00709 
00710    // create stack
00711    THStack stack("histostack","");
00712 
00713    // create legends
00714    TLegend* legend1;
00715    TLegend* legend2;
00716 
00717    if (flag_diff)
00718       legend1 = new TLegend(0.15, (0.88-fraction_pads)/(1-fraction_pads), 0.50, 0.99);
00719    else
00720       legend1 = new TLegend(0.15, 0.88, 0.50, 0.99);
00721    legend1->SetBorderSize(0);
00722    legend1->SetFillColor(kWhite);
00723    legend1->AddEntry(&fHistData, "Data", "LEP");
00724    legend1->AddEntry(&fHistData, "Total expected uncertainty", "LE");
00725 
00726    double y = 0.99;
00727    if (ntemplates > 2 && ntemplates <7)
00728       y -= 0.11 / 4. * double(ntemplates - 2);
00729    legend2 = new TLegend(0.50,(y-fraction_pads)/(1-fraction_pads) , 0.85, 0.99);
00730    legend2->SetBorderSize(0);
00731    legend2->SetFillColor(kWhite);
00732 
00733 
00734    // create new histograms for plotting
00735    std::vector<TH1D*> histcontainer;
00736 
00737    // get best fit parameters
00738    std::vector<double> parameters = GetBestFitParameters(); 
00739    
00740    // fix the parameters of all functions
00741    FixTemplateFunctions(parameters);
00742 
00743    // scale histograms and add to stack and legend
00744    for (int itemp = 0; itemp < ntemplates; ++itemp) {
00745 
00746       // get template index
00747       int templateindex = fTemplateParIndexContainer.at(itemp);
00748 
00749       // create new histogram
00750       TH1D* hist;
00751       if (fTemplateTypeContainer[itemp] == 0) {
00752          int histogramindex = GetHistogramIndex(itemp);
00753          hist = new TH1D( fTemplateHistogramContainer.at(histogramindex) );
00754       }
00755       else if (fTemplateTypeContainer[itemp] == 1) {
00756          hist = new TH1D( fHistData );
00757          for (int i = 1; i <= fNBins; ++i) {
00758             double bincenter = fHistData.GetBinCenter(i);
00759             int index = GetFunctionIndex(itemp);
00760             double bincontent = fTemplateFunctionContainer.at(index).Eval(bincenter);
00761             hist->SetBinContent(i, bincontent);
00762          }
00763          // scale histogram
00764          hist->Scale(1.0/hist->Integral());
00765       }
00766       
00767       // set histogram color and style
00768       hist->SetFillColor(2 + itemp);
00769       hist->SetFillStyle(1001);
00770       hist->SetStats(kFALSE);
00771       
00772       // add histogram to list of histograms
00773       histcontainer.push_back(hist);
00774       
00775       // loop over bins and set them
00776       for (int ibin = 1; ibin <= fNBins; ++ibin) {
00777 
00778          // calculate efficiency for this template and bin
00779          double efficiency = TemplateEfficiency(itemp, ibin, parameters);
00780 
00781          // calculate probability for this template and bin
00782          double probability = TemplateProbability(itemp, ibin, parameters);
00783 
00784          // calculate expectation value for this template
00785          double nexp = parameters.at(templateindex)
00786             * efficiency
00787             * probability;
00788 
00789          // set bin content
00790          hist->SetBinContent(ibin, nexp);
00791       }
00792 
00793 
00794       // add histogram to stack
00795       stack.Add(hist);
00796       if (itemp < 2)
00797          legend1->AddEntry(hist, fTemplateNameContainer.at(itemp).data(), "F");
00798       else if (itemp < 6)
00799          legend2->AddEntry(hist, fTemplateNameContainer.at(itemp).data(), "F");
00800    }
00801 
00802    // loop over all bins
00803    for (int ibin = 1; ibin <= nbins; ++ibin) {
00804       // set bin content
00805       histsum->SetBinContent(ibin, Expectation(ibin, parameters));
00806    }
00807 
00808 
00809 
00810 
00811 
00812 
00813 
00814    // define error graph
00815    TGraphAsymmErrors * graph_error_exp = new TGraphAsymmErrors(nbins);
00816    graph_error_exp->SetLineWidth(2);
00817    graph_error_exp->SetMarkerStyle(0);
00818 
00819    TGraphAsymmErrors * graph_error_obs = new TGraphAsymmErrors(nbins);
00820    graph_error_obs->SetMarkerStyle(0);
00821 
00822    // calculate uncertainty
00823    if (flag_error1)
00824       for (int i = 1; i <= nbins; ++i) {
00825          double nexp = histsum->GetBinContent(i);
00826          histsum->SetBinError(i, sqrt(nexp));
00827          histsum->SetMarkerStyle(0);
00828       }
00829 
00830    if (flag_error2)
00831       for (int i = 1; i <= nbins; ++i) {
00832          TH1D * proj = fUncertaintyHistogramExp->ProjectionY("_py", i, i);
00833          if (proj->Integral() > 0)
00834             proj->Scale(1.0 / proj->Integral());
00835          double quantiles[3];
00836          double sums[3] = {0.16, 0.5, 0.84};
00837          proj->GetQuantiles(3, quantiles, sums);
00838          graph_error_exp->SetPoint(i-1, histsum->GetBinCenter(i), quantiles[1]);
00839          graph_error_exp->SetPointError(i-1, 0.0, 0.0, quantiles[1] - quantiles[0], quantiles[2]-quantiles[1]);
00840          delete proj;
00841       }
00842 
00843    if (flag_error3)
00844       for (int i = 1; i <= nbins; ++i) {
00845          TH1D * proj = fUncertaintyHistogramObsPosterior->ProjectionY("_py", i, i);
00846          if (proj->Integral() > 0)
00847             proj->Scale(1.0 / proj->Integral());
00848          double quantiles[3];
00849          double sums[3] = {0.16, 0.5, 0.84};
00850          proj->GetQuantiles(3, quantiles, sums);
00851          graph_error_obs->SetPoint(i-1, histsum->GetBinCenter(i), quantiles[1]);
00852          graph_error_obs->SetPointError(i-1, 0.0, 0.0, quantiles[1] - TMath::Floor(quantiles[0]), TMath::Ceil(quantiles[2])-quantiles[1]);
00853          delete proj;
00854       }
00855 
00856    // create difference histogram
00857    TH1D * hist_diff = 0;
00858 
00859    TGraphAsymmErrors * graph_diff_exp = 0;
00860 
00861    if (flag_diff) {
00862       ymin = 0;
00863       ymax = 0;
00864       //    hist_diff = new TH1D("hist_diff", "", nbins, histsum->GetXaxis()->GetXmin(), histsum->GetXaxis()->GetXmax() );
00865       //    hist_diff = new TH1D(*histsum);
00866       Double_t a[fHistData.GetNbinsX()+1];
00867       for (int i = 0; i < fHistData.GetNbinsX()+1; ++i) {
00868          a[i] = fHistData.GetXaxis()->GetBinLowEdge(i+1);
00869       }
00870 
00871       hist_diff = new TH1D("hist_diff", "", nbins, a);
00872       hist_diff->SetName("hist_diff");
00873       hist_diff->GetXaxis()->SetTitle(fHistData.GetXaxis()->GetTitle());
00874       hist_diff->GetYaxis()->SetTitle("#Delta N");
00875       hist_diff->GetXaxis()->SetNdivisions(505);
00876       hist_diff->GetXaxis()->SetLabelSize(hist_diff->GetXaxis()->GetLabelSize()/(fraction_pads));
00877       hist_diff->GetXaxis()->SetLabelOffset(hist_diff->GetXaxis()->GetLabelOffset()/fraction_pads*2.);
00878       hist_diff->GetXaxis()->SetTitleSize(hist_diff->GetXaxis()->GetTitleSize()/(fraction_pads));
00879       hist_diff->GetXaxis()->SetTitleOffset((hist_diff->GetXaxis()->GetTitleOffset()-(1.0-fraction_pads))/(fraction_pads));
00880       hist_diff->GetYaxis()->SetNdivisions(503);
00881       hist_diff->GetYaxis()->SetLabelSize(hist_diff->GetYaxis()->GetLabelSize()/(fraction_pads));
00882       hist_diff->GetYaxis()->SetLabelOffset(hist_diff->GetYaxis()->GetLabelOffset()/(fraction_pads));
00883       hist_diff->GetYaxis()->SetTitleSize(hist_diff->GetYaxis()->GetTitleSize()/(fraction_pads));
00884       hist_diff->GetYaxis()->SetTitleOffset(hist_diff->GetYaxis()->GetTitleOffset()*(fraction_pads));
00885       hist_diff->SetStats(kFALSE);
00886 
00887       graph_diff_exp = new TGraphAsymmErrors(nbins);
00888       graph_diff_exp->SetLineWidth(2);
00889       graph_diff_exp->SetMarkerStyle(0);
00890       graph_diff_exp->SetFillColor(kYellow);
00891       for (int i = 0; i < nbins; ++i) {
00892          hist_diff->SetBinContent(i+1, fHistData.GetBinContent(i+1)-histsum->GetBinContent(i+1));
00893          hist_diff->SetBinError(i+1, fHistData.GetBinError(i+1));
00894          graph_diff_exp->SetPoint(i, (graph_error_exp->GetX())[i], 0.0);
00895          graph_diff_exp->SetPointEXlow(i, 0.0);
00896          graph_diff_exp->SetPointEXhigh(i, 0.0);
00897          graph_diff_exp->SetPointEYlow(i, (graph_error_exp->GetEYlow())[i]);
00898          graph_diff_exp->SetPointEYhigh(i, (graph_error_exp->GetEYhigh())[i]);
00899 
00900          if (-(graph_error_exp->GetEYlow())[i] < ymin)
00901             ymin = -(graph_error_exp->GetEYlow())[i];
00902          if ((graph_error_exp->GetEYhigh())[i] > ymax)
00903             ymax = (graph_error_exp->GetEYhigh())[i];
00904       }
00905       if (ymax < (hist_diff->GetMaximum() + hist_diff->GetBinError(hist_diff->GetMaximumBin())))
00906          ymax = 1.1 * (hist_diff->GetMaximum() + hist_diff->GetBinError(hist_diff->GetMaximumBin()));
00907       if (ymin>(hist_diff->GetMinimum() - hist_diff->GetBinError(hist_diff->GetMaximumBin())))
00908          ymin = 1.1 * (hist_diff->GetMinimum() - hist_diff->GetBinError(hist_diff->GetMaximumBin()));
00909       (hist_diff->GetYaxis())->SetRangeUser(TMath::Floor(-1.1*TMath::Max(-ymin, ymax)), TMath::Ceil(1.1*TMath::Max(-ymin, ymax)));
00910 
00911    }
00912 
00913    // draw histograms
00914    stack.Draw("SAMEAF");
00915    stack.GetHistogram() -> SetXTitle("");
00916    stack.GetHistogram() -> SetYTitle("");
00917    stack.GetHistogram() -> GetXaxis() -> SetLabelSize(0);
00918    stack.GetHistogram() -> GetYaxis() -> SetLabelSize(0);
00919    stack.Draw("SAME");
00920    fHistData.Draw("SAMEP");
00921 
00922    if (flag_error0)
00923       fHistData.Draw("SAMEPE");
00924 
00925    if (flag_error1)
00926       histsum->Draw("SAMEE");
00927 
00928    if (flag_error3)
00929       graph_error_obs->Draw("SAMEZ");
00930 
00931    if (flag_error2)
00932       graph_error_exp->Draw("SAME||");
00933 
00934    if (flag_legend) {
00935       legend1->Draw();
00936       if (ntemplates>2)
00937       legend2->Draw();
00938    }
00939 
00940    TLine * line = 0;
00941    if (flag_diff) {
00942       pad2->cd();
00943       hist_diff->Draw("P");
00944       graph_diff_exp->Draw("SAME4");
00945       line = new TLine((hist_diff->GetXaxis())->GetXmin(), 0.0, (hist_diff->GetXaxis())->GetXmax(), 0.0);
00946       line->SetLineWidth(2);
00947       line->SetLineColor(kBlack);
00948       line->Draw("SAME");
00949       hist_diff->Draw("SAMEP");
00950    }
00951 
00952    c1->Print(filename);
00953 
00954    // delete temporary histograms
00955    delete pad1;
00956    delete pad2;
00957    delete c1;
00958    delete legend1;
00959    delete legend2;
00960    delete graph_error_exp;
00961    delete graph_error_obs;
00962    delete histsum;
00963    if (flag_diff) {
00964       delete hist_diff;
00965       delete graph_diff_exp;
00966       delete line;
00967    }
00968 }
00969 
00970 // ---------------------------------------------------------
00971 double BCTemplateFitter::CalculateChi2(std::vector<double> parameters)
00972 {
00973    double chi2 = 0;
00974 
00975    // loop over bins
00976    for (int ibin = 1; ibin <= fNBins; ++ibin) {
00977       // get expectation
00978       double nexp = Expectation(ibin, parameters);
00979 
00980       // get data
00981       double ndata = fHistData.GetBinContent(ibin);
00982 
00983       // add to chi2
00984       chi2 += (nexp - ndata) * (nexp - ndata) / nexp;
00985    }
00986 
00987    // return chi2
00988    return chi2;
00989 }
00990 
00991 // ---------------------------------------------------------
00992 double BCTemplateFitter::CalculateChi2Prob(std::vector<double> parameters)
00993 {
00994    double chi2 = CalculateChi2(parameters);
00995    int ndf = GetNDF();
00996 
00997    // return chi2 probability
00998    return TMath::Prob(chi2, ndf);
00999 }
01000 
01001 // ---------------------------------------------------------
01002 double BCTemplateFitter::CalculateMaxLike()
01003 {
01004    // return maximum likelihood
01005    return Eval( GetBestFitParameters() );
01006 }
01007 
01008 // ---------------------------------------------------------
01009 double BCTemplateFitter::CalculateKSProb()
01010 {
01011    // create a temporary histogram with the sum of all contributions
01012    TH1D * histsum = new TH1D(fHistData);
01013 
01014    // get best fit parameters
01015    std::vector<double> parameters = GetBestFitParameters();
01016 
01017    // loop over bins
01018    for (int ibin = 1; ibin <= fNBins; ++ibin) {
01019       // get expectation
01020       double nexp = Expectation(ibin, parameters);
01021 
01022       // set bin content
01023       histsum->SetBinContent(ibin, nexp);
01024    }
01025 
01026    // perform KS test
01027    double ksprob = histsum->KolmogorovTest(&fHistData);
01028 
01029    // delete histogram
01030    delete histsum;
01031 
01032    return ksprob;
01033 }
01034 
01035 // ---------------------------------------------------------
01036 double BCTemplateFitter::CalculatePValue()
01037 {
01038    // get best fit parameters
01039    std::vector<double> par = GetBestFitParameters();
01040 
01041    // check size of parameter vector
01042    if (par.size() != GetNParameters()) {
01043       BCLog::OutWarning("BCBCTemplateFitter::CalculatePValueFast : Number of parameters is inconsistent.");
01044       return -1;
01045    }
01046 
01047    std::vector <int> histogram;
01048    std::vector <double> expectation;
01049    histogram.assign(fNBins, 0);
01050    expectation.assign(fNBins, 0);
01051 
01052    double logp = 0;
01053    double logp_start = 0;
01054    int counter_pvalue = 0;
01055 
01056    // define starting distribution
01057    for (int ibin = 1; ibin <= fNBins; ++ibin) {
01058 
01059       // get expectation
01060       double nexp = Expectation(ibin, par);
01061       
01062       /*
01063       double nexp = 0;
01064 
01065       // loop over all templates
01066       for (int itemp = 0; itemp < ntemplates; ++itemp) {
01067          int tempindex = fTemplateParIndexContainer.at(itemp);
01068          int effindex = fEffParIndexContainer.at(itemp);
01069 
01070          // get efficiency for the bin
01071          double efficiency = fEffHistogramContainer.at(itemp).GetBinContent(ibin);
01072 
01073          // modify efficiency by uncertainty
01074          double efferr = fEffErrHistogramContainer.at(itemp).GetBinContent(ibin);
01075 
01076          // check efficiency error
01077          if (efferr > 0)
01078             efficiency = TMath::Max(0., efficiency + par.at(effindex) * efferr);
01079 
01080          // add expectation from bin
01081          nexp += par.at(tempindex) * efficiency * fTemplateHistogramContainer.at(itemp).GetBinContent(ibin);
01082       }
01083       */
01084 
01085 
01086       histogram[ibin-1]   = int(nexp);
01087       expectation[ibin-1] = nexp;
01088 
01089       // calculate p;
01090       logp += BCMath::LogPoisson(double(int(nexp)), nexp);
01091       logp_start += BCMath::LogPoisson(fHistData.GetBinContent(ibin), nexp);
01092    }
01093 
01094    int niter = 100000;
01095 
01096    // loop over iterations
01097    for (int iiter = 0; iiter < niter; ++iiter) {
01098       // loop over bins
01099       for (int ibin = 0; ibin < fNBins; ++ibin) {
01100          // random step up or down in statistics for this bin
01101          double ptest = fRandom->Rndm() - 0.5;
01102 
01103          // increase statistics by 1
01104          if (ptest > 0) {
01105             // calculate factor of probability
01106             double r = expectation.at(ibin) / double(histogram.at(ibin) + 1);
01107 
01108             // walk, or don't (this is the Metropolis part)
01109             if (fRandom->Rndm() < r) {
01110                histogram[ibin] = histogram.at(ibin) + 1;
01111                logp += log(r);
01112             }
01113          }
01114 
01115          // decrease statistics by 1
01116          else if (ptest <= 0 && histogram[ibin] > 0) {
01117             // calculate factor of probability
01118             double r = double(histogram.at(ibin)) / expectation.at(ibin);
01119 
01120             // walk, or don't (this is the Metropolis part)
01121             if (fRandom->Rndm() < r) {
01122                histogram[ibin] = histogram.at(ibin) - 1;
01123                logp += log(r);
01124             }
01125          }
01126       } // end of looping over bins
01127 
01128       // increase counter
01129       if (logp < logp_start)
01130          counter_pvalue++;
01131 
01132    } // end of looping over iterations
01133 
01134    // calculate p-value
01135    return double(counter_pvalue) / double(niter);
01136 }
01137 
01138 // ---------------------------------------------------------
01139 void BCTemplateFitter::MCMCUserIterationInterface()
01140 {
01141    // loop over all bins
01142    for (int ibin = 1; ibin <= fNBins; ++ibin) {
01143       // get bin content
01144       double bincontent = Expectation(ibin, fMCMCx);
01145 
01146       // set bin content
01147       fUncertaintyHistogramExp->Fill(fHistData.GetBinCenter(ibin), bincontent);
01148 
01149       // loop over bins in the other direction
01150       int nbinsy = fUncertaintyHistogramObsPosterior->GetNbinsY();
01151       for (int jbin = 1; jbin <= nbinsy; ++jbin) {
01152          int n = jbin - 1;
01153          if (fabs(n - bincontent) < 2*sqrt(bincontent))
01154             fUncertaintyHistogramObsPosterior->Fill(fHistData.GetBinCenter(ibin), n, TMath::Poisson(bincontent, n));
01155       }
01156    }
01157 
01158    // fill normalization
01159    fHistNorm.Fill(fNorm);
01160 
01161    // fill ratios
01162    int nratios = int( fIndicesRatios1D.size() );
01163 
01164    // loop over fractions to fill
01165    for (int i = 0; i < nratios; ++i) {
01166       int nsum = int( (fIndicesRatios1D.at(i)).size() ) - 1;
01167       double sum = 0;
01168       for (int j = 1; j <= nsum; ++j) {
01169          int indexsum = fIndicesRatios1D.at(i).at(j);
01170          sum += fMCMCx.at(fTemplateParIndexContainer.at(indexsum));
01171       }
01172 
01173       fHistRatios1D.at(i).Fill(fMCMCx.at(fTemplateParIndexContainer.at(fIndicesRatios1D.at(i).at(0)))/sum);
01174    }
01175 
01176 }
01177 
01178 // ---------------------------------------------------------
01179 void BCTemplateFitter::PrintRatios(const char * filename, int options, double ovalue)
01180 {
01181    int nratios = int(fHistRatios1D.size());
01182 
01183    TCanvas * c1 = new TCanvas("");
01184 
01185    TPostScript * ps = new TPostScript(filename, 112);
01186    ps->NewPage();
01187 
01188    c1->cd();
01189    BCH1D * h1temp = new BCH1D(&fHistNorm);
01190    h1temp->Draw();
01191    c1->Update();
01192    ps->NewPage();
01193    for (int i = 0; i < nratios; ++i) {
01194       c1->Update();
01195       ps->NewPage();
01196       c1->cd();
01197       BCH1D* h1temp = new BCH1D(&fHistRatios1D.at(i));
01198       h1temp->Draw(options, ovalue);
01199    }
01200    c1->Update();
01201    ps->Close();
01202 
01203    // free memory
01204    delete c1;
01205    delete ps;
01206 }
01207 
01208 // ---------------------------------------------------------
01209 int BCTemplateFitter::SetTemplateEfficiency(const char * name, TH1D eff, TH1D efferr)
01210 {
01211    // get index
01212    int index = GetIndexTemplate(name);
01213 
01214    // check index
01215    if (index < 0) {
01216       BCLog::OutError("BCTemplateFitter::SetTemplateEfficiency : Could not find template.");
01217       return 0;
01218    }
01219 
01220    // check efficiency histogram
01221    if (fTemplateTypeContainer[index]==0)
01222       if (CompareHistogramProperties(fTemplateHistogramContainer.at(index), eff) != 1) {
01223          BCLog::OutError("BCTemplateFitter::SetTemplate efficiency : Template and efficiency histogram properties are incompatible.");
01224          return 0;
01225       }
01226    
01227    // set histogram style
01228    eff.SetXTitle((fHistData.GetXaxis())->GetTitle());
01229    eff.SetYTitle("Efficiency");
01230 
01231    efferr.SetXTitle((fHistData.GetXaxis())->GetTitle());
01232    efferr.SetYTitle("Efficiency uncertainty");
01233 
01234    // set efficiency histogram
01235    fEffHistogramContainer[index] = eff;
01236    fEffErrHistogramContainer[index] = efferr;
01237 
01238    // no error
01239    return 1;
01240 }
01241 
01242 // ---------------------------------------------------------
01243 int BCTemplateFitter::SetTemplateEfficiency(const char * name, double effmean, double effsigma)
01244 {
01245    // get index
01246    int index = GetIndexTemplate(name);
01247 
01248    // check index
01249    if (index < 0) {
01250       BCLog::OutError("BCTemplateFitter::SetTemplateEfficiency : Could not find template.");
01251       return 0;
01252    }
01253 
01254    // create histograms
01255    TH1D histeff = TH1D(fHistData);
01256    TH1D histefferr = TH1D(fHistData);
01257 
01258    // fill histograms
01259    for (int i = 1; i <= fNBins; ++i) {
01260       histeff.SetBinContent(i, effmean);
01261       histefferr.SetBinContent(i, effsigma);
01262    }
01263 
01264    // set histograms
01265    int err = SetTemplateEfficiency(name, histeff, histefferr);
01266 
01267    // return error code
01268    return err;
01269 }
01270 
01271 // ---------------------------------------------------------
01272 int BCTemplateFitter::ConstrainSum(std::vector <int> indices, double mean, double rms)
01273 {
01274    // add contraint to container(s)
01275    fConstraintSumIndices.push_back(indices);
01276    fConstraintSumMean.push_back(mean);
01277    fConstraintSumRMS.push_back(rms);
01278 
01279    // no error
01280    return 1;
01281 }
01282 
01283 // ---------------------------------------------------------
01284 void BCTemplateFitter::PrintTemp()
01285 {
01286    TCanvas * c1 = new TCanvas("");
01287 
01288    c1->cd();
01289    fUncertaintyHistogramExp->Draw("COL");
01290    c1->Print("uncertainty_exp.eps");
01291 
01292    c1->cd();
01293    fUncertaintyHistogramObsPosterior->Draw("COL");
01294    c1->Print("uncertainty_obs.eps");
01295 
01296    delete c1;
01297 }
01298 
01299 // ---------------------------------------------------------
01300 int BCTemplateFitter::PrintTemplate(const char * name, const char * filename)
01301 {
01302    // get number of sources of systematic uncertainty
01303    int nsyst = GetNSystErrors();
01304 
01305    // get index
01306    int index = GetIndexTemplate(name);
01307 
01308    // check index
01309    if (index < 0) {
01310       BCLog::OutError("BCTemplateFitter::PrintTemplate : Could not find template.");
01311       return 0;
01312    }
01313 
01314    // remember old directory
01315    TDirectory* dir = gDirectory;
01316 
01317    // create postscript
01318    TPostScript * ps = new TPostScript(filename);
01319 
01320    // create new canvas
01321    TCanvas * c1 = new TCanvas("", "", 700, 700);
01322 
01323    c1->Update();
01324    ps->NewPage();
01325    c1->cd();
01326 
01327    // create legend
01328    TLegend l1(0.18, 0.75, 0.85, 0.85);
01329    l1.SetBorderSize(0);
01330    l1.SetFillColor(kWhite);
01331 
01332    // draw histogram and uncertainties
01333    TH1D hist_template = fTemplateHistogramContainer.at(index);
01334    hist_template.SetFillColor(kWhite);
01335    hist_template.SetFillStyle(0);
01336    hist_template.SetMarkerSize(0);
01337    hist_template.SetLineWidth(0);
01338    l1.AddEntry(&hist_template, name, "L");
01339    TH1D hist_totalerr = CombineUncertainties(name);
01340    TH1D hist_template_totalerr = CreateErrorHist(hist_template, hist_totalerr);
01341    hist_template_totalerr.SetFillColor(kYellow);
01342    hist_template_totalerr.SetFillStyle(1001);
01343    hist_template_totalerr.SetMarkerSize(0);
01344    l1.AddEntry(&hist_template_totalerr, "Systematic uncertainties", "F");
01345    TH1D hist_efferr = fEffErrHistogramContainer.at(index);
01346    TH1D hist_template_efferr = CreateErrorHist(hist_template, hist_efferr);
01347    hist_template_totalerr.SetFillColor(kRed);
01348    hist_template_totalerr.SetFillStyle(1001);
01349    hist_template_efferr.SetMarkerSize(0);
01350    l1.AddEntry(&hist_template_efferr, "Efficiency uncertainties", "F");
01351    int binmax = hist_template.GetMaximumBin();
01352    double ymax = hist_template.GetBinContent(binmax) + 2.0 * hist_template_totalerr.GetBinError(binmax);
01353    hist_template_totalerr.GetYaxis()->SetRangeUser(0.0, 1.25 * ymax);
01354    hist_template_totalerr.Draw("E2");
01355    hist_template_efferr.Draw("SAMEE2");
01356    hist_template.Draw("SAME");
01357 
01358    // draw legend
01359    l1.Draw();
01360 
01361    // update ps
01362    c1->Update();
01363    ps->NewPage();
01364    c1->cd();
01365 
01366    // create legend
01367    TLegend l2(0.18, 0.75, 0.85, 0.85);
01368    l2.SetBorderSize(0);
01369    l2.SetFillColor(kWhite);
01370 
01371    // print uncertainties
01372    c1->cd(2);
01373    hist_efferr = fEffErrHistogramContainer.at(index);
01374    double ymin = hist_efferr.GetMinimum();
01375    ymax = hist_efferr.GetMaximum();
01376    l2.AddEntry(&hist_efferr, "Efficiency", "L");
01377    hist_efferr.SetStats(kFALSE);
01378    hist_efferr.Draw();
01379 
01380    // loop over all uncertainties
01381    for (int i = 0; i < nsyst; ++i) {
01382       TH1D * hist = new TH1D(fSystErrorHistogramContainer.at(i).at(index));
01383       hist->SetLineColor(2 + i);
01384       if (hist->GetMaximum()>ymax)
01385          ymax = hist->GetMaximum();
01386       if (hist->GetMinimum()<ymin)
01387          ymin = hist->GetMinimum();
01388       l2.AddEntry(hist, fSystErrorNameContainer.at(i).c_str(), "L");
01389       hist->Draw("SAME");
01390    }
01391    if (ymin < 0)
01392       ymin = 1.25*ymin;
01393    else
01394       ymin = 0.8*ymin;
01395 
01396    if (ymax > 0)
01397       ymax = 1.25*ymax;
01398    else
01399       ymax = 0.8*ymax;
01400 
01401    hist_efferr.GetYaxis()->SetRangeUser(ymin, ymax);
01402 
01403    // draw legend
01404    l2.Draw();
01405 
01406    // close ps
01407    c1->Update();
01408    ps->Close();
01409 
01410    // change to old directory
01411    dir->cd();
01412 
01413    // free memory
01414    delete c1;
01415    delete ps;
01416 
01417    // no error
01418    return 1;
01419 }
01420 
01421 // ---------------------------------------------------------
01422 TH1D BCTemplateFitter::CreateErrorHist(TH1D hist, TH1D histerr)
01423 {
01424    // check histogram properties
01425    if (CompareHistogramProperties(fHistData, hist) != 1) {
01426       BCLog::OutError("BCTemplateFitter::CreateErrorHist : Histograms are incompatible.");
01427       return hist;
01428    }
01429 
01430    // copy histogram
01431    TH1D h = hist;
01432 
01433    // set style
01434    h.SetStats(kFALSE);
01435    h.SetFillColor(kYellow);
01436    h.SetFillStyle(1001);
01437 
01438    // get number of bins
01439    int nbins = hist.GetNbinsX();
01440 
01441    // loop over bins
01442    for (int i = 1; i <= nbins; ++i)
01443       h.SetBinError(i, histerr.GetBinContent(i) * hist.GetBinContent(i));
01444 
01445    // return histogram
01446    return h;
01447 }
01448 
01449 // ---------------------------------------------------------
01450 TH1D BCTemplateFitter::CombineUncertainties(const char * name)
01451 {
01452    // get number of sources of systematic uncertainty
01453    int nsyst = GetNSystErrors();
01454 
01455    // get template index
01456    int tempindex = GetIndexTemplate(name);
01457 
01458    // create new histogram
01459    //    TH1D hist = TH1D("", "", fNBins, fXmin, fXmax);
01460    TH1D hist = TH1D(fHistData);
01461 
01462    // fill histogram
01463    for (int ibin = 1; ibin <= fNBins; ++ibin) {
01464 
01465       // define total uncertainty
01466       double err = 0;
01467 
01468       // add efficiency uncertainty squared
01469       double erreff = fEffErrHistogramContainer.at(tempindex).GetBinContent(ibin);
01470       err += erreff * erreff;
01471 
01472       // add systematic uncertainty squared
01473       for (int isyst = 0; isyst < nsyst; ++isyst) {
01474          double errsyst = fSystErrorHistogramContainer.at(isyst).at(tempindex).GetBinContent(ibin);
01475          err += errsyst*errsyst;
01476       }
01477 
01478       // take square root
01479       err = sqrt(err);
01480 
01481       // set bin content
01482       hist.SetBinContent(ibin, err);
01483    }
01484 
01485    // return histogram
01486    return hist;
01487 }
01488 
01489 // ---------------------------------------------------------
01490 int BCTemplateFitter::GetIndexTemplate(const char * name)
01491 {
01492    int index = -1;
01493    int n = GetNTemplates();
01494 
01495    for (int i = 0; i < n; i++)
01496       if (name == fTemplateNameContainer.at(i))
01497          index = i;
01498 
01499    if (index < 0) {
01500       BCLog::OutWarning("BCTemplateFitter::GetIndexTemplate : Template does not exist.");
01501       return 0;
01502    }
01503 
01504    // return index
01505    return index;
01506 }
01507 
01508 // ---------------------------------------------------------
01509 int BCTemplateFitter::GetIndexSystError(const char * name)
01510 {
01511    int index = -1;
01512    int n = GetNSystErrors();
01513 
01514    for (int i = 0; i < n; i++)
01515       if (name == fSystErrorNameContainer.at(i))
01516          index = i;
01517 
01518    if (index < 0) {
01519       BCLog::OutWarning("BCTemplateFitter::GetIndexSystError : Template does not exist.");
01520       return 0;
01521    }
01522 
01523    // return index
01524    return index;
01525 }
01526 
01527 // ---------------------------------------------------------
01528 int BCTemplateFitter::GetParIndexTemplate(const char * name)
01529 {
01530    int index = -1;
01531    int n = GetNTemplates();
01532 
01533    for (int i = 0; i < n; i++)
01534       if (name == fTemplateNameContainer.at(i))
01535          index = fTemplateParIndexContainer.at(i);
01536 
01537    if (index < 0) {
01538       BCLog::OutWarning("BCTemplateFitter::GetParIndexTemplate : Template does not exist.");
01539       return 0;
01540    }
01541 
01542    // return index
01543    return index;
01544 }
01545 
01546 // ---------------------------------------------------------
01547 int BCTemplateFitter::GetParIndexTemplate(int index)
01548 {
01549    // get number of templates
01550    int n = GetNTemplates();
01551 
01552    if (index < 0 || index > n) {
01553       BCLog::OutError("BCTemplateFitter::GetParIndexTemplate : Index out of range.");
01554       return -1;
01555    }
01556 
01557    // return index
01558    return fTemplateParIndexContainer.at(index);
01559 }
01560 
01561 // ---------------------------------------------------------
01562 int BCTemplateFitter::GetParIndexEff(const char * name)
01563 {
01564    int index = -1;
01565    int n = GetNTemplates();
01566 
01567    for (int i = 0; i < n; i++)
01568       if (name == fTemplateNameContainer.at(i))
01569          index = fTemplateParIndexContainer.at(i);
01570 
01571    if (index < 0) {
01572       BCLog::OutWarning("BCTemplateFitter::GetParIndexEff : Template does not exist.");
01573       return 0;
01574    }
01575 
01576    // return index
01577    return index;
01578 }
01579 
01580 // ---------------------------------------------------------
01581 int BCTemplateFitter::GetParIndexSystError(const char * name)
01582 {
01583    int index = -1;
01584    int n = GetNTemplates();
01585 
01586    for (int i = 0; i < n; i++)
01587       if (name == fSystErrorNameContainer.at(i))
01588          index = fSystErrorParIndexContainer.at(i);
01589 
01590    if (index < 0) {
01591       BCLog::OutWarning("BCTemplateFitter::GetParIndexStatError : Systematic error does not exist.");
01592       return 0;
01593    }
01594 
01595    // return index
01596    return index;
01597 }
01598 
01599 // ---------------------------------------------------------
01600 int  BCTemplateFitter::PerformFit()
01601 {
01602    // initialize
01603    if (!Initialize()) {
01604       BCLog::OutError("BCTemplateFitter::PerformFit : Could not initialize template fitter.");
01605       return 0;
01606    }
01607 
01608    // run Markov Chains
01609    MarginalizeAll();
01610 
01611    // find global mode
01612    FindMode();
01613 
01614    // no error
01615    return 1;
01616 }
01617 // ---------------------------------------------------------
01618 
01619 int BCTemplateFitter::GetNTemplatesType(int templatetype)
01620 {
01621    // get number of templates
01622    int ntemplates = GetNTemplates();
01623 
01624    // initialize counter to zero
01625    int counter = 0;
01626 
01627    // loop over all templates
01628    for (int i = 0; i < ntemplates; ++i) {
01629       // increase counter if template of same type
01630       if (templatetype == fTemplateTypeContainer.at(i))
01631          counter++;
01632    }
01633 
01634    // return counter
01635    return counter;
01636 }
01637 
01638 // ---------------------------------------------------------

Generated on Fri Nov 19 2010 17:02:53 for Bayesian Analysis Toolkit by  doxygen 1.7.1