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