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

BCTemplateEnsembleTest.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 <TFile.h>
00011 #include <TTree.h>
00012 #include <TROOT.h>
00013 #include <TRandom3.h>
00014 #include <TPostScript.h>
00015 #include <TCanvas.h>
00016 
00017 #include <math.h>
00018 #include <iostream>
00019 
00020 #include "../../BAT/BCLog.h"
00021 #include "../../BAT/BCH1D.h"
00022 
00023 #include "BCTemplateFitter.h"
00024 #include "BCTemplateEnsembleTest.h"
00025 
00026 using namespace std;
00027 
00028 // ---------------------------------------------------------
00029 BCTemplateEnsembleTest::BCTemplateEnsembleTest()
00030    : fTemplateFitter(0)
00031    , fFile(0)
00032    , fTree(0)
00033    , fEnsembleCounter(0)
00034    , fEnsembleExpectation(0)
00035    , fNEnsembles(0)
00036    , fFlagMCMC(false)
00037 {
00038    fRandom = new TRandom3(0);
00039 
00040    BCLog::OutWarning("Class BCTemplateFitter is depreceted and it will be removed");
00041    BCLog::OutWarning("in the future version of BAT. The same is true for the class");
00042    BCLog::OutWarning("BCTemplateEnsembleTest. Use the new class BCMTF and its ensemble");
00043    BCLog::OutWarning("testing capabilities instead.");
00044 }
00045 
00046 // ---------------------------------------------------------
00047 BCTemplateEnsembleTest::~BCTemplateEnsembleTest()
00048 {
00049    if (fRandom)
00050       delete fRandom;
00051 }
00052 
00053 // ---------------------------------------------------------
00054 int BCTemplateEnsembleTest::SetEnsembleTemplate(TH1D hist)
00055 {
00056    // calculate integral
00057    double integral = hist.Integral();
00058 
00059    // check if integral is ok
00060    if (integral <= 0) {
00061       std::cout << "Template not valid. Integral is lower or equal to 0." << std::endl;
00062       return 0;
00063    }
00064 
00065    // set template
00066    fEnsembleTemplate = hist;
00067 
00068    // scale template
00069    fEnsembleTemplate.Scale(1.0/integral);
00070 
00071    // no error
00072    return 1;
00073 }
00074 
00075 // ---------------------------------------------------------
00076 int BCTemplateEnsembleTest::PerformEnsembleTest(TTree* tree)
00077 {
00078    // set log level to nothing
00079    BCLog::LogLevel ll = BCLog::GetLogLevelScreen();
00080    BCLog::SetLogLevel(BCLog::nothing);
00081 
00082    // initialize template fitter
00083    fTemplateFitter->Initialize();
00084 
00085    // Prepare tree
00086    PrepareTree();
00087 
00088    // get number of parameters
00089    int npar = fTemplateFitter->GetNParameters();
00090 
00091    // connect tree
00092    if (tree) {
00093       fTemplateParameters = std::vector<double>(npar);
00094       for (int i = 0; i < npar; ++i) {
00095          tree->SetBranchAddress(Form("Parameter%i", i), &(fTemplateParameters[i]));
00096       }
00097    }
00098 
00099    // loop over ensembles
00100    for(int j = 0; j < fNEnsembles; j++){
00101 
00102       // print status
00103       if ((j+1) % 100 == 0 && j > 0)
00104          std::cout << "Fraction of ensembles analyzed: " << double(j+1) / double(fNEnsembles) * 100 << "%" << std::endl;
00105 
00106       // get parameters from tree
00107       if (tree) {
00108         int index = (int) fRandom->Uniform(tree->GetEntries());
00109         tree->GetEntry(index);
00110       }
00111 
00112       // create new ensemble
00113       TH1D* ensemble = BuildEnsemble();
00114 
00115       // set ensemble as new data set
00116       fTemplateFitter->SetData(*ensemble);
00117 
00118       // find mode
00119       fTemplateFitter->FindMode();
00120 
00121       // perform MCMC
00122       if(fFlagMCMC) {
00123          fTemplateFitter->MCMCInitialize();
00124 
00125          fTemplateFitter->MarginalizeAll();
00126 
00127          // find mode with MCMC best fit
00128          fTemplateFitter->FindMode(fTemplateFitter->GetBestFitParameters());
00129 
00130          // loop over parameters and set tree variables
00131          for (int i = 0; i < npar; ++i) {
00132             BCH1D* hist = fTemplateFitter->GetMarginalized(fTemplateFitter->GetParameter(i));
00133             fOutParModeMarg[i]       = hist->GetMode();
00134             fOutParMedianMarg[i]     = hist->GetMedian();
00135             fOutParMeanMarg[i]       = hist->GetMean();
00136             fOutParRMSMarg[i]        = hist->GetRMS();
00137             fOutParErrorUpMarg[i]    = hist->GetQuantile(0.84)-hist->GetMode();
00138             fOutParErrorDownMarg[i]  = hist->GetMode()-hist->GetQuantile(0.16);
00139             fOutParQuantile5Marg[i]  = hist->GetQuantile(0.05);
00140             fOutParQuantile10Marg[i] = hist->GetQuantile(0.10);
00141             fOutParQuantile90Marg[i] = hist->GetQuantile(0.90);
00142             fOutParQuantile95Marg[i] = hist->GetQuantile(0.95);
00143             fOutParPullMarg[i]       = -1;
00144             double error = 0.5 * (fOutParErrorUpMarg.at(i) + fOutParErrorDownMarg.at(i));
00145             if (error > 0)
00146                fOutParPullMarg[i]     = (fOutParModeMarg.at(i) - fTemplateParameters.at(i)) / ( error );
00147          }
00148          fOutChi2Marg             = fTemplateFitter->CalculateChi2( fTemplateFitter->GetBestFitParametersMarginalized() );
00149          fOutChi2ProbMarg         = fTemplateFitter->CalculateChi2Prob(fTemplateFitter->GetBestFitParametersMarginalized());
00150       }
00151 
00152       if (fFlagMCMC) {
00153          int nratios = fTemplateFitter->GetNRatios();
00154          for (int i = 0; i < nratios; ++i) {
00155             TH1D histtemp = fTemplateFitter->GetHistRatio1D(i);
00156              BCH1D * hist = new BCH1D( &histtemp );
00157              fOutRatioModeMarg[i]       = hist->GetMode();
00158              fOutRatioMedianMarg[i]     = hist->GetMedian();
00159              fOutRatioMeanMarg[i]       = hist->GetMean();
00160              fOutRatioRMSMarg[i]        = hist->GetRMS();
00161              fOutRatioErrorUpMarg[i]    = hist->GetQuantile(0.84)-hist->GetMode();
00162             fOutRatioErrorDownMarg[i]  = hist->GetMode()-hist->GetQuantile(0.16);
00163              fOutRatioQuantile5Marg[i]  = hist->GetQuantile(0.05);
00164              fOutRatioQuantile10Marg[i] = hist->GetQuantile(0.10);
00165              fOutRatioQuantile90Marg[i] = hist->GetQuantile(0.90);
00166              fOutRatioQuantile95Marg[i] = hist->GetQuantile(0.95);
00167          }
00168       }
00169 
00170       // set tree variables
00171       fOutParValue           = fTemplateParameters;
00172       fOutParModeGlobal      = fTemplateFitter->GetBestFitParameters();
00173       fOutParErrorUpGlobal   = fTemplateFitter->GetBestFitParameterErrors();
00174       fOutParErrorDownGlobal = fTemplateFitter->GetBestFitParameterErrors();
00175       fOutChi2Global         = fTemplateFitter->CalculateChi2( fTemplateFitter->GetBestFitParameters() );
00176       fOutNDF                = fTemplateFitter->GetNDF();
00177       fOutChi2ProbGlobal     = fTemplateFitter->CalculateChi2Prob(fTemplateFitter->GetBestFitParameters());
00178       fOutKSProb             = fTemplateFitter->CalculateKSProb();
00179       fOutPValue             = fTemplateFitter->CalculatePValue();
00180       fOutNEvents            = int(fTemplateFitter->GetData().Integral());
00181 
00182       for (int i = 0; i < npar; ++i) {
00183          fOutParPullGlobal[i] = -1;
00184          double error = 0.5 * (fOutParErrorUpGlobal.at(i) + fOutParErrorDownGlobal.at(i));
00185          if (error > 0)
00186             fOutParPullGlobal[i] = (fOutParModeGlobal.at(i) - fTemplateParameters.at(i)) / ( error );
00187       }
00188 
00189       // fill the tree
00190       fTree->Fill();
00191    }
00192 
00193    // reset log level
00194    BCLog::SetLogLevel(ll);
00195 
00196    // no error
00197    return 1;
00198 }
00199 
00200 //---------------------------------------------------------------------------------------------------------
00201 TH1D* BCTemplateEnsembleTest::BuildEnsemble()
00202 {
00203    // get histogram parameters
00204   int nbins   = fTemplateFitter->GetData().GetNbinsX();
00205 
00206    // create new ensemble
00207    TH1D* ensemble = new TH1D(fTemplateFitter->GetData());
00208 
00209    // increase ensemble counter
00210    fEnsembleCounter++;
00211 
00212    // get new parameter if needed
00213    std::vector<double> parameters = fTemplateParameters;
00214 
00215    // loop over bins and fill them
00216   for(int ibin = 1; ibin <= nbins; ++ibin){
00217       double nexp = fTemplateFitter->Expectation(ibin, parameters);
00218       double nobs = gRandom->Poisson(nexp);
00219 
00220       // set the bin content
00221     ensemble->SetBinContent(ibin, nobs);
00222   }
00223 
00224    // return the ensemble histogram
00225   return ensemble;
00226 }
00227 
00228 //---------------------------------------------------------------------------------------------------------
00229 int BCTemplateEnsembleTest::Write(const char * filename)
00230 {
00231    // open file
00232    fFile = new TFile(filename, "RECREATE");
00233 
00234    // write tree
00235    fTree->Write();
00236 
00237    // close file
00238    fFile->Close();
00239 
00240    // free memory
00241    delete fFile;
00242 
00243    // no error
00244    return 1;
00245 }
00246 
00247 //---------------------------------------------------------------------------------------------------------
00248 int BCTemplateEnsembleTest::PrepareTree()
00249 {
00250    // delete old tree if necessary
00251    if (fTree)
00252       delete fTree;
00253 
00254    // create new tree
00255    fTree = new TTree("fTree", "fTree");
00256 
00257    // get number of parameters and ratios
00258    int npar = fTemplateFitter->GetNParameters();
00259    int nratios = fTemplateFitter->GetNRatios();
00260 
00261    // initialize variables
00262    fOutParValue.assign(npar, 0);
00263    fOutParModeGlobal.assign(npar, 0);
00264    fOutParErrorUpGlobal.assign(npar, 0);
00265    fOutParErrorDownGlobal.assign(npar, 0);
00266    fOutParPullGlobal.assign(npar, 0);
00267    fOutParModeMarg.assign(npar, 0);
00268    fOutParMeanMarg.assign(npar, 0);
00269    fOutParMedianMarg.assign(npar, 0);
00270    fOutParRMSMarg.assign(npar, 0);
00271    fOutParErrorUpMarg.assign(npar, 0);
00272    fOutParErrorDownMarg.assign(npar, 0);
00273    fOutParPullMarg.assign(npar, 0);
00274    fOutParQuantile5Marg.assign(npar, 0);
00275    fOutParQuantile10Marg.assign(npar, 0);
00276    fOutParQuantile90Marg.assign(npar, 0);
00277    fOutParQuantile95Marg.assign(npar, 0);
00278    fOutRatioModeMarg.assign(nratios, 0);
00279    fOutRatioMeanMarg.assign(nratios, 0);
00280    fOutRatioMedianMarg.assign(nratios, 0);
00281    fOutRatioRMSMarg.assign(nratios, 0);
00282    fOutRatioErrorUpMarg.assign(nratios, 0);
00283    fOutRatioErrorDownMarg.assign(nratios, 0);
00284    fOutRatioQuantile5Marg.assign(nratios, 0);
00285    fOutRatioQuantile10Marg.assign(nratios, 0);
00286    fOutRatioQuantile90Marg.assign(nratios, 0);
00287    fOutRatioQuantile95Marg.assign(nratios, 0);
00288 
00289    fTree->Branch("chi2_global", &fOutChi2Global,     "chi2 (global)/D");
00290    fTree->Branch("ndf",      &fOutNDF,      "ndf/I");
00291    fTree->Branch("chi2prob_global", &fOutChi2ProbGlobal, "chi2 prob probability (global)/D");
00292    fTree->Branch("KSprob",   &fOutKSProb,   "KS probability/D");
00293    fTree->Branch("pvalue",   &fOutPValue,   "p-value/D");
00294    fTree->Branch("nevents",  &fOutNEvents,  "n events/I");
00295 
00296    for (int i = 0; i < npar; ++i) {
00297       // add branches
00298       fTree->Branch(Form("parameter_%i", i),                 &fOutParValue[i],      Form("parameter_%i/D", i));
00299       fTree->Branch(Form("par_global_mode_par_%i", i),       &fOutParModeGlobal[i],      Form("par_global_Mode_par_%i/D", i));
00300       fTree->Branch(Form("par_global_error_up_par_%i", i),   &fOutParErrorUpGlobal[i],   Form("par_global_error_up_par_%i/D", i));
00301       fTree->Branch(Form("par_global_error_down_par_%i", i), &fOutParErrorDownGlobal[i], Form("par_global_error_down_par_%i/D", i));
00302       fTree->Branch(Form("par_global_pull_par_%i", i),       &fOutParPullGlobal[i],      Form("par_global_pull_par_%i/D", i));
00303 
00304       if(fFlagMCMC) {
00305          fTree->Branch(Form("par_marg_mode_par_%i", i),       &fOutParModeMarg[i],       Form("par_marg_mode_par_%i/D", i));
00306          fTree->Branch(Form("par_marg_mean_par_%i", i),       &fOutParMeanMarg[i],       Form("par_marg_mean_par_%i/D", i));
00307          fTree->Branch(Form("par_marg_median_par_%i", i),     &fOutParMedianMarg[i],     Form("par_marg_median_par_%i/D", i));
00308          fTree->Branch(Form("par_marg_rms_par_%i", i),        &fOutParRMSMarg[i],        Form("par_marg_rms_par_%i/D", i));
00309          fTree->Branch(Form("par_marg_error_up_par_%i", i),   &fOutParErrorUpMarg[i],    Form("par_marg_ErrorUp_par_%i/D", i));
00310          fTree->Branch(Form("par_marg_error_down_par_%i", i), &fOutParErrorDownMarg[i],  Form("par_marg_error_down_par_%i/D", i));
00311          fTree->Branch(Form("par_marg_pull_par_%i", i),       &fOutParPullMarg[i],       Form("par_marg_pull_par_%i/D", i));
00312          fTree->Branch(Form("par_marg_quantile5_par_%i", i),  &fOutParQuantile5Marg[i],  Form("par_marg_Quantile5_par_%i/D", i));
00313          fTree->Branch(Form("par_marg_quantile10_par_%i", i), &fOutParQuantile10Marg[i], Form("par_marg_Quantile10_par_%i/D", i));
00314          fTree->Branch(Form("par_marg_quantile90_par_%i", i), &fOutParQuantile90Marg[i], Form("par_marg_Quantile90_par_%i/D", i));
00315          fTree->Branch(Form("par_marg_quantile95_par_%i", i), &fOutParQuantile95Marg[i], Form("par_marg_Quantile95_par_%i/D", i));
00316       }
00317    }
00318 
00319    if (fFlagMCMC) {
00320       fTree->Branch("chi2_marg", &fOutChi2Marg,     "chi2 (marginalized)/D");
00321       fTree->Branch("chi2prob_marg", &fOutChi2ProbMarg, "chi2 prob probability (marginalized)/D");
00322 
00323       for (int i = 0; i < nratios; ++i) {
00324          fTree->Branch(Form("ratio_marg_mode_ratio_%i", i),       &fOutRatioModeMarg[i],       Form("ratio_marg_mode_ratio_%i/D", i));
00325          fTree->Branch(Form("ratio_marg_mean_ratio_%i", i),       &fOutRatioMeanMarg[i],       Form("ratio_marg_mean_ratio_%i/D", i));
00326          fTree->Branch(Form("ratio_marg_median_ratio_%i", i),     &fOutRatioMedianMarg[i],     Form("ratio_marg_median_ratio_%i/D", i));
00327          fTree->Branch(Form("ratio_marg_rms_ratio_%i", i),        &fOutRatioRMSMarg[i],        Form("ratio_marg_rms_ratio_%i/D", i));
00328          fTree->Branch(Form("ratio_marg_error_up_ratio_%i", i),   &fOutRatioErrorUpMarg[i],    Form("ratio_marg_ErrorUp_ratio_%i/D", i));
00329          fTree->Branch(Form("ratio_marg_error_down_ratio_%i", i), &fOutRatioErrorDownMarg[i],  Form("ratio_marg_error_down_ratio_%i/D", i));
00330          fTree->Branch(Form("ratio_marg_quantile5_ratio_%i", i),  &fOutRatioQuantile5Marg[i],  Form("ratio_marg_Quantile5_ratio_%i/D", i));
00331          fTree->Branch(Form("ratio_marg_quantile10_ratio_%i", i), &fOutRatioQuantile10Marg[i], Form("ratio_marg_Quantile10_ratio_%i/D", i));
00332          fTree->Branch(Form("ratio_marg_quantile90_ratio_%i", i), &fOutRatioQuantile90Marg[i], Form("ratio_marg_Quantile90_ratio_%i/D", i));
00333          fTree->Branch(Form("ratio_marg_quantile95_ratio_%i", i), &fOutRatioQuantile95Marg[i], Form("ratio_marg_Quantile95_ratio_%i/D", i));
00334       }
00335    }
00336 
00337    // no error
00338    return 1;
00339 }
00340 
00341 //---------------------------------------------------------------------------------------------------------
00342 void BCTemplateEnsembleTest::PrintPulls(const char* filename)
00343 {
00344    // create postscript
00345    TPostScript * ps = new TPostScript(filename);
00346 
00347    // create canvas and prepare postscript
00348    TCanvas * canvas = new TCanvas();
00349 
00350    canvas->Update();
00351    ps->NewPage();
00352    canvas->cd();
00353 
00354    // get number of parameters
00355    int npar = fTemplateFitter->GetNParameters();
00356 
00357    // create histogram container
00358    std::vector<TH1D*> histcontainer = std::vector<TH1D*>(0);
00359 
00360    // loop over all parameters
00361    for (int j = 0; j < npar; ++j) {
00362       TH1D* hist = new TH1D("", "Pull;N", 11, -5.5, 5.5);
00363       histcontainer.push_back(hist);
00364    }
00365 
00366    // loop over all ensembles
00367    for (int i = 0; i < fNEnsembles; ++i) {
00368 
00369       // get ensemble
00370       fTree->GetEntry(i);
00371 
00372       // loop over all parameters
00373       for (int j = 0; j < npar; ++j) {
00374          histcontainer[j]->Fill(fOutParPullGlobal[j]);
00375       }
00376    }
00377 
00378    // loop over all parameters
00379    for (int j = 0; j < npar; ++j) {
00380       // update post script
00381       canvas->Update();
00382       ps->NewPage();
00383       canvas->cd();
00384 
00385       canvas->cd();
00386       histcontainer.at(j)->Draw();
00387    }
00388    // close ps
00389    canvas->Update();
00390    ps->Close();
00391 
00392    // free memory
00393    delete ps;
00394    delete canvas;
00395 
00396 }
00397 
00398 //---------------------------------------------------------------------------------------------------------

Generated by  doxygen 1.7.1