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

BCMTF.cxx

Go to the documentation of this file.
00001 /*
00002  * Copyright (C) 2008-2012, Daniel Kollar and Kevin Kroeninger.
00003  * All rights reserved.
00004  *
00005  * For the licensing terms see doc/COPYING.
00006  */
00007 
00008 // ---------------------------------------------------------
00009 
00010 #include <iostream>
00011 #include <fstream>
00012 
00013 #include <TCanvas.h>
00014 #include <THStack.h>
00015 #include <TH1D.h>
00016 #include <TAxis.h>
00017 #include <TF1.h>
00018 
00019 #include "../../BAT/BCMath.h"
00020 #include "../../BAT/BCLog.h"
00021 
00022 #include "BCMTFChannel.h"
00023 #include "BCMTFProcess.h"
00024 #include "BCMTFTemplate.h"
00025 #include "BCMTFSystematic.h"
00026 #include "BCMTFSystematicVariation.h"
00027 
00028 #include "BCMTF.h"
00029 
00030 // ---------------------------------------------------------
00031 BCMTF::BCMTF()
00032  : BCModel("Multi-template Fitter")
00033  , fNChannels(0)
00034  , fNProcesses(0)
00035  , fNSystematics(0)
00036  , fFlagEfficiencyConstraint(false)
00037 {}
00038 
00039 // ---------------------------------------------------------
00040 BCMTF::BCMTF(const char * name)
00041  : BCModel(name)
00042  , fNChannels(0)
00043  , fNProcesses(0)
00044  , fNSystematics(0)
00045  , fFlagEfficiencyConstraint(false)
00046 {}
00047 
00048 // ---------------------------------------------------------
00049 BCMTF::~BCMTF()
00050    // default destructor
00051 {
00052    for (int i = 0; i < fNChannels; ++i)
00053       delete fChannelContainer.at(i);
00054 }
00055 
00056 // ---------------------------------------------------------
00057 int BCMTF::GetChannelIndex(const char * name)
00058 {
00059    // loop over all channels and compare names
00060    for (int i = 0; i < fNChannels; ++i) {
00061       // get channel
00062       BCMTFChannel * channel = GetChannel(i);
00063 
00064       // compare names
00065       if (!channel->GetName().compare(name))
00066          return i;
00067    }
00068 
00069    // if channel does not exist, return -1
00070    return -1;
00071 }
00072 
00073 // ---------------------------------------------------------
00074 int BCMTF::GetProcessIndex(const char * name)
00075 {
00076    // loop over all processs and compare names
00077    for (int i = 0; i < fNProcesses; ++i) {
00078       // get process
00079       BCMTFProcess * process = GetProcess(i);
00080 
00081       // compare names
00082       if (!process->GetName().compare(name))
00083          return i;
00084    }
00085 
00086    // if process does not exist, return -1
00087    return -1;
00088 }
00089 
00090 // ---------------------------------------------------------
00091 int BCMTF::GetSystematicIndex(const char * name)
00092 {
00093    // loop over all systematics and compare names
00094    for (int i = 0; i < fNSystematics; ++i) {
00095       // get systematic
00096       BCMTFSystematic * systematic = GetSystematic(i);
00097 
00098       // compare names
00099       if (!systematic->GetName().compare(name))
00100          return i;
00101    }
00102 
00103    // if process does not exist, return -1
00104    return -1;
00105 }
00106 
00107 // ---------------------------------------------------------
00108 int BCMTF::SetTemplate(const char * channelname, const char * processname, TH1D hist, double efficiency)
00109 {
00110    // get channel index
00111    int channelindex = GetChannelIndex(channelname);
00112 
00113    // check if channel exists
00114    if (channelindex < 0) {
00115       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Channel does not exist.");
00116          return -1;
00117    }
00118 
00119    // get process index
00120    int processindex = GetProcessIndex(processname);
00121 
00122    // check if process exists
00123    if (processindex < 0) {
00124       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Process does not exist.");
00125          return -1;
00126    }
00127 
00128    // get channel
00129    BCMTFChannel * channel = GetChannel(channelindex);
00130 
00131    // get template
00132    BCMTFTemplate * bctemplate = channel->GetTemplate(processindex);
00133 
00134    // normalize histogram
00135    if (hist.Integral())
00136       hist.Scale(1.0 / hist.Integral());
00137 
00138    // remove statistics box
00139    hist.SetStats(kFALSE);
00140 
00141    // set color and fill style
00142    hist.SetFillColor(2 + processindex);
00143    hist.SetFillStyle(1001);
00144 
00145    // create new histogram
00146    TH1D * temphist = new TH1D(hist);
00147 
00148    // set histogram
00149    bctemplate->SetHistogram(temphist);
00150 
00151    // set efficiency
00152    bctemplate->SetEfficiency(efficiency);
00153 
00154    // no error
00155    return 1;
00156 }
00157 
00158 // ---------------------------------------------------------
00159 int BCMTF::SetTemplate(const char * channelname, const char * processname, std::vector<TF1 *> * funccont, int nbins, double efficiency)
00160 {
00161    // get channel index
00162    int channelindex = GetChannelIndex(channelname);
00163 
00164    // check if channel exists
00165    if (channelindex < 0) {
00166       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Channel does not exist.");
00167          return -1;
00168    }
00169 
00170    // get process index
00171    int processindex = GetProcessIndex(processname);
00172 
00173    // check if process exists
00174    if (processindex < 0) {
00175       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Process does not exist.");
00176          return -1;
00177    }
00178 
00179    // get channel
00180    BCMTFChannel * channel = GetChannel(channelindex);
00181 
00182    // get template
00183    BCMTFTemplate * bctemplate = channel->GetTemplate(processindex);
00184 
00185    // set histogram
00186    bctemplate->SetFunctionContainer(funccont, nbins);
00187 
00188    // set efficiency
00189    bctemplate->SetEfficiency(efficiency);
00190 
00191    // no error
00192    return 1;
00193 }
00194 
00195 // ---------------------------------------------------------
00196 int BCMTF::SetData(const char * channelname, TH1D hist)
00197 {
00198    int channelindex = GetChannelIndex(channelname);
00199 
00200    // check if channel exists
00201    if (channelindex < 0) {
00202       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Channel does not exist.");
00203          return -1;
00204    }
00205 
00206    // get channel
00207    BCMTFChannel * channel = GetChannel(channelindex);
00208 
00209    // get template
00210    BCMTFTemplate * data = channel->GetData();
00211 
00212    // remove statistics box
00213    hist.SetStats(kFALSE);
00214 
00215    // set marker
00216    hist.SetMarkerStyle(20);
00217    hist.SetMarkerSize(1.1);
00218 
00219    // remove old data set if it exists
00220    if (data->GetHistogram()) {
00221       delete data->GetHistogram();
00222       data->SetHistogram(0);
00223    }
00224 
00225    // set histogram
00226    data->SetHistogram(new TH1D(hist));
00227 
00228    // no error
00229    return 1;
00230 }
00231 
00232 // ---------------------------------------------------------
00233 int BCMTF::AddChannel(const char * name)
00234 {
00235    // check if channel exists
00236    for (int i = 0; i < fNChannels; ++i) {
00237       // compare names
00238       if (GetChannelIndex(name) >= 0) {
00239       BCLog::OutWarning("BCMultitemplateFitter::AddChannel() : Channel with this name exists already.");
00240          return -1;
00241       }
00242    }
00243 
00244    // create new channel
00245    BCMTFChannel * channel = new BCMTFChannel(name);
00246 
00247    // create new data
00248    BCMTFTemplate * bctemplate = new BCMTFTemplate(channel->GetName().c_str(), "data");
00249 
00250    // add data
00251    channel->SetData(bctemplate);
00252 
00253    // add process templates
00254    for (int i = 0; i < fNProcesses; ++i) {
00255       // get process
00256       BCMTFProcess * process = GetProcess(i);
00257 
00258       // create new template
00259       BCMTFTemplate * bctemplate = new BCMTFTemplate(name, process->GetName().c_str());
00260 
00261       // add template
00262       channel->AddTemplate(bctemplate);
00263    }
00264 
00265    // loop over all systematics
00266    for (int i = 0; i < fNSystematics; ++i) {
00267       // get systematic
00268       BCMTFSystematic * systematic = GetSystematic(i);
00269 
00270       // create new systematic variation
00271       BCMTFSystematicVariation * variation = new BCMTFSystematicVariation(name, systematic->GetName().c_str(), fNProcesses);
00272 
00273       // add systematic variation
00274       channel->AddSystematicVariation(variation);
00275    }
00276 
00277    // add channel
00278    fChannelContainer.push_back(channel);
00279 
00280    // increase number of channels
00281    fNChannels++;
00282 
00283    // no error
00284    return 1;
00285 }
00286 
00287 // ---------------------------------------------------------
00288 int BCMTF::AddProcess(const char * name, double nmin, double nmax)
00289 {
00290    // check if process exists
00291    for (int i = 0; i < fNProcesses; ++i) {
00292       // compare names
00293       if (GetProcessIndex(name) >= 0) {
00294       BCLog::OutWarning("BCMultitemplateFitter::AddProcess() : Process with this name exists already.");
00295          return -1;
00296       }
00297    }
00298 
00299    // create new process
00300    BCMTFProcess * process = new BCMTFProcess(name);
00301 
00302    // add process
00303    fProcessContainer.push_back(process);
00304 
00305    // add process templates
00306    for (int i = 0; i < fNChannels; ++i) {
00307       // get channel
00308       BCMTFChannel * channel = GetChannel(i);
00309 
00310       // create new template
00311       BCMTFTemplate * bctemplate = new BCMTFTemplate(channel->GetName().c_str(), name);
00312 
00313       // add template
00314       channel->AddTemplate(bctemplate);
00315 
00316       // loop over all systematic
00317       for (int j = 0; j < fNSystematics; ++j) {
00318          // get systematic variation
00319          BCMTFSystematicVariation * variation = channel->GetSystematicVariation(j);
00320 
00321          // add histogram
00322          variation->AddHistograms(0, 0);
00323       }
00324    }
00325 
00326    // increase number of processes
00327    fNProcesses++;
00328 
00329    // add parameter index to container
00330    fProcessParIndexContainer.push_back(GetNParameters());
00331 
00332    // add parameter
00333    AddParameter(name, nmin, nmax);
00334 
00335    // add a functional form for the expectation
00336    fExpectationFunctionContainer.push_back(0);
00337 
00338    // no error
00339    return 1;
00340 }
00341 
00342 // ---------------------------------------------------------
00343 int BCMTF::AddSystematic(const char * name, double min, double max)
00344 {
00345    // check if systematic exists
00346    for (int i = 0; i < fNSystematics; ++i) {
00347       // compare names
00348       if (GetSystematicIndex(name) >= 0) {
00349       BCLog::OutWarning("BCMultitemplateFitter::AddSystematic() : Systematic with this name exists already.");
00350          return -1;
00351       }
00352    }
00353 
00354    // create new systematic
00355    BCMTFSystematic * systematic = new BCMTFSystematic(name);
00356 
00357    // add systematic
00358    fSystematicContainer.push_back(systematic);
00359 
00360    // add systematic variations
00361    for (int i = 0; i < fNChannels; ++i) {
00362       // get channel
00363       BCMTFChannel * channel = GetChannel(i);
00364 
00365       // create new systematic variation
00366       BCMTFSystematicVariation * variation = new BCMTFSystematicVariation(channel->GetName().c_str(), name, fNProcesses);
00367 
00368       // add systematic variation
00369       channel->AddSystematicVariation(variation);
00370    }
00371    // ...
00372 
00373    // increase number of systematices
00374    fNSystematics++;
00375 
00376    // add parameter index to container
00377    fSystematicParIndexContainer.push_back(GetNParameters());
00378 
00379    // add parameter
00380    AddParameter(name, min, max);
00381 
00382    // add a functional form for the expectation
00383    fExpectationFunctionContainer.push_back(0);
00384 
00385    // no error
00386    return 1;
00387 }
00388 
00389 // ---------------------------------------------------------
00390 int BCMTF::SetSystematicVariation(const char * channelname, const char * processname,  const char * systematicname, double variation_up, double variation_down)
00391 {
00392 
00393    // get channel index
00394    int channelindex = GetChannelIndex(channelname);
00395 
00396    // check if channel exists
00397    if (channelindex < 0) {
00398       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Channel does not exist.");
00399          return -1;
00400    }
00401 
00402    // get process index
00403    int processindex = GetProcessIndex(processname);
00404 
00405    // check if process exists
00406    if (processindex < 0) {
00407       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Process does not exist.");
00408          return -1;
00409    }
00410 
00411    // get systematic index
00412    int systematicindex = GetSystematicIndex(systematicname);
00413 
00414    // check if systematic exists
00415    if (systematicindex < 0) {
00416       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Systematic does not exist.");
00417          return -1;
00418    }
00419 
00420    // get channel
00421    BCMTFChannel * channel = GetChannel(channelindex);
00422 
00423    BCMTFTemplate * bctemplate = channel->GetTemplate(processindex);
00424 
00425    TH1D * hist_template = bctemplate->GetHistogram();
00426 
00427    TH1D hist_up = TH1D(*hist_template);
00428    TH1D hist_down = TH1D(*hist_template);
00429 
00430    int nbins = hist_up.GetNbinsX();
00431 
00432    // loop over all bins
00433    for (int ibin = 1; ibin <= nbins; ++ibin) {
00434       hist_up.SetBinContent(ibin, variation_up);
00435       hist_down.SetBinContent(ibin, variation_down);
00436    }
00437 
00438    // get systematic variation
00439    BCMTFSystematicVariation * variation = channel->GetSystematicVariation(systematicindex);
00440 
00441    // set histogram
00442    variation->SetHistograms(processindex, new TH1D(hist_up), new TH1D(hist_down));
00443 
00444    // no error
00445    return 1;
00446 }
00447 
00448 // ---------------------------------------------------------
00449 int BCMTF::SetSystematicVariation(const char * channelname, const char * processname,  const char * systematicname, TH1D hist_up, TH1D hist_down)
00450 {
00451    // get channel index
00452    int channelindex = GetChannelIndex(channelname);
00453 
00454    // check if channel exists
00455    if (channelindex < 0) {
00456       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Channel does not exist.");
00457          return -1;
00458    }
00459 
00460    // get process index
00461    int processindex = GetProcessIndex(processname);
00462 
00463    // check if process exists
00464    if (processindex < 0) {
00465       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Process does not exist.");
00466          return -1;
00467    }
00468 
00469    // get systematic index
00470    int systematicindex = GetSystematicIndex(systematicname);
00471 
00472    // check if systematic exists
00473    if (systematicindex < 0) {
00474       BCLog::OutWarning("BCMultitemplateFitter::SetTemplate() : Systematic does not exist.");
00475          return -1;
00476    }
00477 
00478    // get channel
00479    BCMTFChannel * channel = GetChannel(channelindex);
00480 
00481    // get systematic variation
00482    BCMTFSystematicVariation * variation = channel->GetSystematicVariation(systematicindex);
00483 
00484    // set histogram
00485    variation->SetHistograms(processindex, new TH1D(hist_up), new TH1D(hist_down));
00486 
00487    // no error
00488    return 1;
00489 }
00490 
00491 // ---------------------------------------------------------
00492 int BCMTF::SetSystematicVariation(const char * channelname, const char * processname,  const char * systematicname, TH1D hist, TH1D hist_up, TH1D hist_down)
00493 {
00494    // get number of bins
00495    int nbins = hist.GetNbinsX();
00496 
00497    TH1D * hist_up_rel = new TH1D(hist);
00498    TH1D * hist_down_rel = new TH1D(hist);
00499 
00500    // loop over all bins
00501    for (int ibin = 1; ibin <= nbins; ++ibin) {
00502       hist_up_rel->SetBinContent(ibin, (hist_up.GetBinContent(ibin) - hist.GetBinContent(ibin)) / hist.GetBinContent(ibin));
00503       hist_down_rel->SetBinContent(ibin, (hist.GetBinContent(ibin) - hist_down.GetBinContent(ibin)) / hist.GetBinContent(ibin));
00504    }
00505 
00506    // set the systematic variation
00507    return SetSystematicVariation(channelname, processname, systematicname, *hist_up_rel, *hist_down_rel);
00508 }
00509 
00510 // ---------------------------------------------------------
00511 int BCMTF::PrintSummary(const char * filename)
00512 {
00513    // open file
00514    std::ofstream ofi(filename);
00515    ofi.precision(3);
00516 
00517    // check if file is open
00518    if(!ofi.is_open()) {
00519       BCLog::OutWarning(Form("BCMultitemplateFitter::PrintSummary() : Could not open file %s", filename));
00520       return 0;
00521    }
00522 
00523    ofi
00524       << " Multi template fitter summary " << std::endl
00525       << " ----------------------------- " << std::endl
00526       << std::endl
00527       << " Number of channels      : " << fNChannels << std::endl
00528       << " Number of processes     : " << fNProcesses << std::endl
00529       << " Number of systematics   : " << fNSystematics << std::endl
00530       << std::endl;
00531 
00532    ofi
00533       << " Channels :" << std::endl;
00534       for (int i = 0; i < GetNChannels(); ++i) {
00535          ofi
00536             << " " << i
00537             << " : \"" << GetChannel(i)->GetName().c_str() << "\""
00538             << std::endl;
00539       }
00540       ofi
00541          << std::endl;
00542 
00543    ofi
00544       << " Processes :" << std::endl;
00545       for (int i = 0; i < GetNProcesses(); ++i) {
00546          ofi
00547             << " " << i
00548             << " : \"" << GetProcess(i)->GetName().c_str()  << "\""
00549             << " (par index " << GetParIndexProcess(i) << ")"
00550             << std::endl;
00551       }
00552       ofi
00553          << std::endl;
00554 
00555    ofi
00556       << " Systematics :" << std::endl;
00557       for (int i = 0; i < GetNSystematics(); ++i) {
00558          ofi
00559             << " " << i
00560             << " : \"" << GetSystematic(i)->GetName().c_str()  << "\""
00561             << " (par index " << GetParIndexSystematic(i) << ")"
00562             << std::endl;
00563       }
00564       ofi
00565          << std::endl;
00566       if (GetNSystematics() == 0)
00567          ofi
00568             << " - none - " << std::endl;
00569 
00570       ofi
00571          << " Goodness-of-fit: " << std::endl;
00572       for (int i = 0; i < GetNChannels(); ++i) {
00573          ofi
00574             << " i : \"" << GetChannel(i)->GetName().c_str() << "\" : chi2 = "
00575             << CalculateChi2( i, GetBestFitParameters() )
00576             << std::endl;
00577       }
00578       ofi
00579          << std::endl;
00580 
00581 
00582    // close file
00583    ofi.close();
00584 
00585    // no error
00586    return 1;
00587 }
00588 
00589 // ---------------------------------------------------------
00590 double BCMTF::Expectation(int channelindex, int binindex, const std::vector<double> & parameters)
00591 {
00592    double expectation = 0.;
00593     
00594    // loop over all processes
00595    for (int i = 0; i < fNProcesses; ++i) {
00596        // get efficiency
00597        double efficiency = Efficiency(channelindex, i, binindex, parameters);
00598        
00599        // get probability
00600        double probability = Probability(channelindex, i, binindex, parameters);
00601        
00602        // get parameter index
00603        int parindex = fProcessParIndexContainer[i];
00604        
00605        // add to expectation
00606        expectation += ExpectationFunction(parindex, channelindex, i, parameters)
00607           * efficiency
00608           * probability;
00609    }
00610     
00611    // check if expectation is positive
00612    if (expectation < 0)
00613       expectation = 0.;
00614 
00615    return expectation;
00616 }
00617 
00618 // ---------------------------------------------------------
00619 double BCMTF::ExpectationFunction(int parindex, int channelindex, int processindex, const std::vector<double> & parameters)
00620 {
00621   // get function container
00622   std::vector<TF1 *> * funccont = fChannelContainer[channelindex]->GetTemplate(processindex)->GetFunctionContainer();
00623 
00624   if (funccont->size()>0)
00625     return 1.;
00626 
00627   else if (!fExpectationFunctionContainer[parindex])
00628     return parameters[parindex];
00629 
00630   else {
00631     TF1 * func = fExpectationFunctionContainer[parindex];
00632     return func->Eval(parameters[parindex]);
00633   }
00634 }
00635 
00636 // ---------------------------------------------------------
00637 double BCMTF::Efficiency(int channelindex, int processindex, int binindex, const std::vector<double> & parameters)
00638 {
00639    // get channel
00640    BCMTFChannel * channel = fChannelContainer[channelindex];
00641 
00642    double efficiency = channel->GetTemplate(processindex)->GetEfficiency();
00643 
00644    double defficiency = 1.;
00645 
00646    // loop over all systematics
00647    for (int i = 0; i < fNSystematics; ++i) {
00648       if (!(fSystematicContainer[i]->GetFlagSystematicActive()))
00649          continue;
00650 
00651       // get parameter index
00652       int parindex = fSystematicParIndexContainer[i];
00653 
00654       // get histogram
00655       TH1D * hist = 0;
00656 
00657       if (parameters[parindex] > 0)
00658          hist = channel->GetSystematicVariation(i)->GetHistogramUp(processindex);
00659       else
00660          hist = channel->GetSystematicVariation(i)->GetHistogramDown(processindex);
00661 
00662       // check if histogram exists
00663       if (!hist)
00664          continue;
00665 
00666       // multiply efficiency
00667       defficiency += parameters[parindex] * hist->GetBinContent(binindex);
00668    }
00669 
00670    // calculate efficiency
00671    efficiency *= defficiency;
00672 
00673    // make sure efficiency is between 0 and 1
00674    if (fFlagEfficiencyConstraint) {
00675       if (efficiency < 0.)
00676          efficiency = 0.;
00677       if (efficiency > 1.)
00678          efficiency = 1.;
00679    }
00680 
00681    return efficiency;
00682 }
00683 
00684 // ---------------------------------------------------------
00685 double BCMTF::Probability(int channelindex, int processindex, int binindex, const std::vector<double> & parameters)
00686 {
00687    // get histogram
00688    TH1D * hist = fChannelContainer[channelindex]->GetTemplate(processindex)->GetHistogram();
00689 
00690    // get function container
00691    std::vector<TF1 *> * funccont = fChannelContainer[channelindex]->GetTemplate(processindex)->GetFunctionContainer();
00692 
00693    // this needs to be fast
00694    if (!hist && !(funccont->size()>0))
00695       return 0.;
00696 
00697    if (hist)
00698      return hist->GetBinContent(binindex);
00699    else {
00700      int parindex = fProcessParIndexContainer[processindex];
00701      return funccont->at(binindex-1)->Eval(parameters[parindex]);
00702    }
00703 }
00704 
00705 // ---------------------------------------------------------
00706 int BCMTF::PrintStack(const char * channelname, const std::vector<double> & parameters, const char * filename, const char * options)
00707 {
00708    int index = GetChannelIndex(channelname);
00709 
00710    return PrintStack(index, parameters, filename, options);
00711 }
00712 
00713 // ---------------------------------------------------------
00714 int BCMTF::PrintStack(int channelindex, const std::vector<double> & parameters, const char * filename, const char * options)
00715 {
00716    // check if parameters are filled
00717    if (!parameters.size())
00718       return -1;
00719 
00720    // check options
00721    bool flag_logx   = false; // plot x-axis in log-scale
00722    bool flag_logy   = false; // plot y-axis in log-scale
00723    bool flag_bw     = false; // plot in black and white
00724 
00725    if (std::string(options).find("logx") < std::string(options).size())
00726       flag_logx = true;
00727 
00728    if (std::string(options).find("logy") < std::string(options).size())
00729       flag_logy = true;
00730 
00731    if (std::string(options).find("bw") < std::string(options).size())
00732       flag_bw = true;
00733 
00734    // get channel
00735    BCMTFChannel * channel = GetChannel(channelindex);
00736 
00737    // create canvas
00738    TCanvas * c1 = new TCanvas();
00739    c1->cd();
00740 
00741    // set log or linear scale
00742    if (flag_logx)
00743       c1->SetLogx();
00744 
00745    if (flag_logy)
00746       c1->SetLogy();
00747 
00748    // create stack
00749    THStack * stack = new THStack("", "");
00750 
00751    // create a container of temporary histograms
00752    std::vector<TH1D *> histcontainer;
00753 
00754    // get number of templates
00755    unsigned int ntemplates = GetNProcesses();
00756 
00757    // loop over all templates
00758    for (unsigned int i = 0; i < ntemplates; ++i) {
00759 
00760       // get histogram
00761       TH1D * temphist = channel->GetTemplate(i)->GetHistogram();
00762 
00763       // get function container
00764       std::vector<TF1 *> * funccont = channel->GetTemplate(i)->GetFunctionContainer();
00765 
00766       // create new histogram
00767       TH1D * hist(0);
00768 
00769       if (temphist)
00770          hist = new TH1D( *(temphist) );
00771       else if (funccont)
00772          hist = new TH1D( *(channel->GetData()->GetHistogram()));
00773       else
00774          continue;
00775 
00776       if (flag_bw) {
00777          hist->SetFillColor(0);
00778          hist->SetLineWidth(0);
00779          hist->SetLineStyle(int(1+i));
00780       }
00781       else {
00782          hist->SetFillColor(2 + i);
00783          hist->SetFillStyle(1001);
00784       }
00785 
00786       // get number of bins
00787       int nbins = hist->GetNbinsX();
00788 
00789       // scale histogram
00790       for (int ibin = 1; ibin <= nbins; ++ibin) {
00791 
00792          // get efficiency
00793          double efficiency = Efficiency(channelindex, i, ibin, parameters);
00794 
00795          // get probability
00796          double probability = Probability(channelindex, i, ibin, parameters);
00797 
00798          // get parameter index
00799          int parindex = GetParIndexProcess(i);
00800 
00801          // add to expectation
00802          double expectation = parameters[parindex] * efficiency * probability;
00803 
00804          // set bin content
00805          hist->SetBinContent(ibin, expectation);
00806       }
00807 
00808       // add histogram to container (for memory management)
00809       histcontainer.push_back(hist);
00810 
00811       // add histogram to stack
00812       stack->Add(hist);
00813    }
00814 
00815    //draw data
00816    channel->GetData()->GetHistogram()->Draw("PE");
00817 
00818    // set range user
00819    channel->GetData()->GetHistogram()->GetYaxis()->SetRangeUser(0., channel->GetData()->GetHistogram()->GetMaximum() + 2. * sqrt(channel->GetData()->GetHistogram()->GetMaximum()));
00820 
00821    // draw stack
00822    stack->Draw("SAMEHIST");
00823 
00824    //draw data again
00825    channel->GetData()->GetHistogram()->Draw("SAMEPE");
00826 
00827    // redraw the axes
00828    gPad->RedrawAxis();
00829 
00830    // print
00831    c1->Print(filename);
00832 
00833    // free memory
00834    for (unsigned int i = 0; i < histcontainer.size(); ++i) {
00835       TH1D * hist = histcontainer.at(i);
00836       delete hist;
00837    }
00838    delete stack;
00839    delete c1;
00840 
00841    // no error
00842    return 1;
00843 }
00844 
00845 // ---------------------------------------------------------
00846 double BCMTF::CalculateChi2(int channelindex, const std::vector<double> & parameters)
00847 {
00848    if (parameters.size() == 0)
00849       return -1;
00850 
00851    double chi2 = 0;
00852 
00853    // get channel
00854    BCMTFChannel * channel = GetChannel(channelindex);
00855 
00856    // get data
00857    BCMTFTemplate * data = channel->GetData();
00858 
00859    // get histogram
00860    TH1D * hist = data->GetHistogram();
00861 
00862    // check if histogram exists
00863    if (hist) {
00864       // get number of bins in data
00865       int nbins = hist->GetNbinsX();
00866 
00867       // loop over all bins
00868       for (int ibin = 1; ibin <= nbins; ++ibin) {
00869          // get expectation value
00870          double expectation = Expectation(channelindex, ibin, parameters);
00871 
00872          // get observation
00873          double observation = hist->GetBinContent(ibin);
00874 
00875          // add Poisson term
00876          chi2 += (expectation - observation) * (expectation - observation) / expectation;
00877       }
00878    }
00879 
00880    // return chi2;
00881    return chi2;
00882 }
00883 
00884 // ---------------------------------------------------------
00885 double BCMTF::CalculateChi2(const std::vector<double> & parameters)
00886 {
00887    if (parameters.size() == 0)
00888       return -1;
00889 
00890    double chi2 = 0;
00891 
00892    // get number of channels
00893    int nchannels = GetNChannels();
00894 
00895    // loop over all channels
00896    for (int i = 0; i < nchannels; ++i) {
00897       chi2 += CalculateChi2(i, parameters);
00898    }
00899 
00900    // return chi2
00901    return chi2;
00902 }
00903 
00904 // ---------------------------------------------------------
00905 double BCMTF::CalculateCash(int channelindex, const std::vector<double> & parameters)
00906 {
00907    if (parameters.size() == 0)
00908       return -1;
00909 
00910    double cash = 0;
00911 
00912    // get channel
00913    BCMTFChannel * channel = GetChannel(channelindex);
00914 
00915    // get data
00916    BCMTFTemplate * data = channel->GetData();
00917 
00918    // get histogram
00919    TH1D * hist = data->GetHistogram();
00920 
00921    // check if histogram exists
00922    if (hist) {
00923       // get number of bins in data
00924       int nbins = hist->GetNbinsX();
00925 
00926       // loop over all bins
00927       for (int ibin = 1; ibin <= nbins; ++ibin) {
00928          // get expectation value
00929          double expectation = Expectation(channelindex, ibin, parameters);
00930 
00931          // get observation
00932          double observation = hist->GetBinContent(ibin);
00933 
00934          // calculate Cash statistic
00935          cash += 2. * (expectation - observation);
00936 
00937          // check negative log
00938          if (observation > 0)
00939             cash += 2. * observation * log (observation/expectation);
00940       }
00941    }
00942 
00943    // return cash;
00944    return cash;
00945 
00946 }
00947 
00948 // ---------------------------------------------------------
00949 double BCMTF::CalculateCash(const std::vector<double> & parameters)
00950 {
00951    if (parameters.size() == 0)
00952       return -1;
00953 
00954    double cash = 0;
00955 
00956    // get number of channels
00957    int nchannels = GetNChannels();
00958 
00959    // loop over all channels
00960    for (int i = 0; i < nchannels; ++i) {
00961       cash += CalculateCash(i, parameters);
00962    }
00963 
00964    // return cash
00965    return cash;
00966 }
00967 
00968 // ---------------------------------------------------------
00969 double BCMTF::LogLikelihood(const std::vector<double> & parameters)
00970 {
00971    double logprob = 0.;
00972 
00973    // loop over all channels
00974    for (int ichannel = 0; ichannel < fNChannels; ++ichannel) {
00975 
00976       // get channel
00977       BCMTFChannel * channel = fChannelContainer[ichannel];
00978 
00979       // check if channel is active
00980       if (!(channel->GetFlagChannelActive()))
00981          continue;
00982 
00983      // get data
00984       BCMTFTemplate * data = channel->GetData();
00985 
00986       // get histogram
00987       TH1D * hist = data->GetHistogram();
00988 
00989       // check if histogram exists
00990       if (!hist)
00991          continue;
00992 
00993       // get number of bins in data
00994       int nbins = data->GetNBins();
00995 
00996       // loop over all bins
00997       for (int ibin = 1; ibin <= nbins; ++ibin) {
00998 
00999          // get expectation value
01000          double expectation = Expectation(ichannel, ibin, parameters);
01001 
01002          // get observation
01003          double observation = hist->GetBinContent(ibin);
01004 
01005          // add Poisson term
01006          logprob += BCMath::LogPoisson(observation, expectation);
01007       }
01008    }
01009 
01010    return logprob;
01011 }
01012 
01013 // ---------------------------------------------------------

Generated by  doxygen 1.7.1