BayesianAnalysisToolkit  0.9.3
BCMTFAnalysisFacility.cxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2007-2013, the BAT core developer team
3  * All rights reserved.
4  *
5  * For the licensing terms see doc/COPYING.
6  * For documentation see http://mpp.mpg.de/bat
7  */
8 
9 // ---------------------------------------------------------
10 
11 #include "BCMTFAnalysisFacility.h"
12 
13 #include "BCMTF.h"
14 #include "BCMTFChannel.h"
15 #include "BCMTFComparisonTool.h"
16 #include "BCMTFSystematic.h"
17 #include "BCMTFTemplate.h"
18 
19 #include "../../BAT/BCLog.h"
20 #include "../../BAT/BCH1D.h"
21 #include "../../BAT/BCParameter.h"
22 
23 #include <TCanvas.h>
24 #include <TFile.h>
25 #include <TH1D.h>
26 #include <TRandom3.h>
27 #include <TROOT.h>
28 #include <TTree.h>
29 
30 #include <sys/stat.h>
31 #include <unistd.h>
32 #include <iostream>
33 
34 // ---------------------------------------------------------
36  : fRandom(new TRandom3(0))
37  , fFlagMarginalize(false)
38  , fLogLevel(BCLog::nothing)
39 {
40  fMTF = mtf;
41  BCLog::OutDetail(Form("Prepared Analysis Facility for MTF model \'%s\'",mtf->GetName().c_str()));
42 }
43 
44 // ---------------------------------------------------------
46 {
47 }
48 
49 // ---------------------------------------------------------
50 std::vector<TH1D> BCMTFAnalysisFacility::BuildEnsemble(const std::vector<double> & parameters, std::string options)
51 {
52  // option flags
53  bool flag_data = false;
54 
55  // check content of options string
56  if (options.find("data") < options.size()) {
57  flag_data = true;
58  }
59 
60  // get number of channels
61  int nchannels = fMTF->GetNChannels();
62 
63  // create vector of histograms
64  std::vector<TH1D> histograms;
65 
66  // loop over channels
67  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
68 
69  // get channel
70  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
71 
72  // create new histogram
73  TH1D hist( *(channel->GetData()->GetHistogram()) );
74 
75  // get number of bins
76  int nbins = hist.GetNbinsX();
77 
78  // loop over all bins
79  for (int ibin = 1; ibin <= nbins; ++ibin) {
80  if (!flag_data) {
81  double expectation = fMTF->Expectation(ichannel, ibin, parameters);
82  double observation = fRandom->Poisson(expectation);
83 
84  hist.SetBinContent(ibin, observation);
85  }
86  }
87 
88  // add histogram
89  histograms.push_back(hist);
90  }
91 
92  // return histograms
93  return histograms;
94 }
95 
96 // ---------------------------------------------------------
97 TTree * BCMTFAnalysisFacility::BuildEnsembles(TTree * tree, int nensembles, std::string options)
98 {
99  // get number of channels
100  int nchannels = fMTF->GetNChannels();
101 
102  BCLog::OutDetail(Form("MTF Building %d ensembles for %d channels.",nensembles,nchannels));
103 
104  // get number of parameters
105  int nparameters = fMTF->GetNParameters();
106 
107  // create tree variables for input tree
108  std::vector<double> parameters(nparameters);
109 
110  // set branch addresses
111  for (int i = 0; i < nparameters; ++i) {
112  tree->SetBranchAddress(Form("Parameter%i", i), &(parameters[i]));
113  }
114 
115  // create tree
116  TTree * tree_out = new TTree("ensembles", "ensembles");
117 
118  // create matrix of number of bins
119  std::vector< std::vector<double> > nbins_matrix;
120 
121  // prepare the tree variables
122  // loop over channels
123  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
124  // get channel
125  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
126 
127  // get number of bins
128  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
129 
130  // create new matrix row
131  std::vector<double> nbins_column(nbins);
132 
133  // add matrix
134  nbins_matrix.push_back(nbins_column);
135  }
136 
137  std::vector<double> in_parameters(nparameters);
138 
139  // create branches
140  // loop over channels
141  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
142  // get channel
143  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
144 
145  // get number of bins
146  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
147 
148  // loop over bins
149  for (int ibin = 1; ibin <= nbins; ++ibin) {
150  // create branches
151  tree_out->Branch(Form("channel_%i_bin_%i", ichannel, ibin),
152  &(nbins_matrix[ichannel])[ibin-1], "n/D");
153  }
154  }
155 
156  for (int i = 0; i < nparameters; ++i) {
157  tree_out->Branch(Form("parameter_%i", i), &in_parameters[i], Form("parameter_%i/D", i));
158  }
159 
160  // create vector of histograms
161  std::vector<TH1D> histograms;
162 
163  // loop over ensembles
164  for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
165  // get random event from tree
166  int index = (int) fRandom->Uniform(tree->GetEntries());
167  tree->GetEntry(index);
168 
169  // create ensembles
170  histograms = BuildEnsemble(parameters, options);
171 
172  // copy information from histograms into tree variables
173  // loop over channels
174  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
175  // get channel
176  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
177 
178  // get number of bins
179  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
180 
181  // loop over bins
182  for (int ibin = 1; ibin <= nbins; ++ibin) {
183  // fill tree variable
184  (nbins_matrix[ichannel])[ibin-1] = histograms.at(ichannel).GetBinContent(ibin);
185  }
186  }
187 
188  // copy parameter information
189  for (int i = 0; i < nparameters; ++i) {
190  in_parameters[i] = parameters.at(i);
191  }
192 
193  // fill tree
194  tree_out->Fill();
195  }
196 
197  // return tree
198  return tree_out;
199 }
200 
201 // ---------------------------------------------------------
202 TTree * BCMTFAnalysisFacility::BuildEnsembles(const std::vector<double> & parameters, int nensembles, std::string options)
203 {
204  // get number of channels
205  int nchannels = fMTF->GetNChannels();
206 
207  BCLog::OutDetail(Form("MTF Building %d ensambles for %d channels.",nensembles,nchannels));
208 
209  // create tree
210  TTree * tree = new TTree("ensembles", "ensembles");
211 
212  // create matrix of number of bins
213  std::vector< std::vector<double> > nbins_matrix;
214 
215  // prepare the tree variables
216  // loop over channels
217  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
218  // get channel
219  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
220 
221  // get number of bins
222  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
223 
224  // create new matrix row
225  std::vector<double> nbins_column(nbins);
226 
227  // add matrix
228  nbins_matrix.push_back(nbins_column);
229  }
230 
231  // get number of parameters
232  int nparameters = fMTF->GetNParameters();
233 
234  std::vector<double> in_parameters(nparameters);
235 
236  // create branches
237  // loop over channels
238  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
239  // get channel
240  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
241 
242  // get number of bins
243  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
244 
245  // loop over bins
246  for (int ibin = 1; ibin <= nbins; ++ibin) {
247  // create branches
248  tree->Branch(Form("channel_%i_bin_%i", ichannel, ibin),
249  &(nbins_matrix[ichannel])[ibin-1], "n/D");
250  }
251  }
252 
253  for (int i = 0; i < nparameters; ++i) {
254  tree->Branch(Form("parameter_%i", i), &in_parameters[i], Form("parameter_%i/D", i));
255  }
256 
257  // create vector of histograms
258  std::vector<TH1D> histograms;
259 
260  // loop over ensembles
261  for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
262  // create ensembles
263  histograms = BuildEnsemble(parameters, options);
264 
265  // copy information from histograms into tree variables
266  // loop over channels
267  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
268  // get channel
269  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
270 
271  // get number of bins
272  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
273 
274  // loop over bins
275  for (int ibin = 1; ibin <= nbins; ++ibin) {
276  // fill tree variable
277  (nbins_matrix[ichannel])[ibin-1] = histograms.at(ichannel).GetBinContent(ibin);
278  }
279  }
280 
281  // copy parameter information
282  for (int i = 0; i < nparameters; ++i) {
283  if (parameters.size() > 0)
284  in_parameters[i] = parameters.at(i);
285  else
286  in_parameters[i] = 0;
287  }
288 
289  // fill tree
290  tree->Fill();
291  }
292 
293  // return tree
294  return tree;
295 }
296 
297 // ---------------------------------------------------------
298 TTree * BCMTFAnalysisFacility::PerformEnsembleTest(const std::vector<double> & parameters, int nensembles, std::string options)
299 {
300  // create new tree
301  TTree * tree = 0;
302 
303  // create ensembles
304  tree = BuildEnsembles(parameters, nensembles, options);
305 
306  // perform ensemble test
307  return PerformEnsembleTest(tree, nensembles, 0, options);
308 }
309 
310 // ---------------------------------------------------------
311 TTree * BCMTFAnalysisFacility::PerformEnsembleTest(TTree * tree, int nensembles, int start, std::string options)
312 {
313  BCLog::OutSummary("Running ensemble test.");
314  if (fFlagMarginalize) {
315  BCLog::OutSummary("Fit for each ensemble is going to be run using MCMC. It can take a while.");
316  }
317 
318  // option flags
319  bool flag_mc = false;
320 
321  // check content of options string
322  if (options.find("MC") < options.size()) {
323  flag_mc = true;
324  }
325 
326  // set log level
327  // It would be better to set the level for the screen and for the file
328  // separately. Perhaps in the future.
332  BCLog::OutSummary("No log messages for the ensemble fits are going to be printed.");
334  }
335  else if(fLogLevel!=lls) {
336  BCLog::OutSummary(Form("The log level for the ensemble test is set to \'%s\'.",BCLog::ToString(fLogLevel)));
338  }
339 
340  // get number of channels
341  int nchannels = fMTF->GetNChannels();
342 
343  // define set of histograms for the original data sets
344  std::vector<TH1D *> histograms_data(nchannels);
345 
346  // create matrix of number of bins
347  std::vector< std::vector<double> > nbins_matrix;
348 
349  // prepare the tree
350  // loop over channels
351  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
352  // get channel
353  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
354 
355  // get number of bins
356  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
357 
358  // create new matrix row
359  std::vector<double> nbins_column(nbins);
360 
361  // add matrix
362  nbins_matrix.push_back(nbins_column);
363  }
364 
365  // loop over channels
366  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
367  // get channel
368  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
369 
370  // get number of bins
371  int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
372 
373  // loop over bins
374  for (int ibin = 1; ibin <= nbins; ++ibin) {
375  // create branches
376  tree->SetBranchAddress(Form("channel_%i_bin_%i", ichannel, ibin),
377  &(nbins_matrix[ichannel])[ibin-1]);
378  }
379  }
380 
381  // get number of parameters
382  int nparameters = fMTF->GetNParameters();
383 
384  // define tree variables
385  std::vector<double> out_parameters(nparameters);
386 
387  for (int i = 0; i < nparameters; ++i) {
388  tree->SetBranchAddress(Form("parameter_%i", i), &out_parameters[i]);
389  }
390 
391  // copy the original data sets
392  // loop over channels
393  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
394  // get channel
395  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
396 
397  // set data pointer
398  histograms_data[ichannel] = channel->GetData()->GetHistogram();
399  }
400 
401  // create output tree
402  TTree * tree_out = new TTree("ensemble_test", "ensemble test");
403 
404  // define tree variables
405  std::vector<double> out_mode_global(nparameters);
406  std::vector<double> out_std_global(nparameters);
407  std::vector<double> out_mode_marginalized(nparameters);
408  std::vector<double> out_mean_marginalized(nparameters);
409  std::vector<double> out_median_marginalized(nparameters);
410  std::vector<double> out_5quantile_marginalized(nparameters);
411  std::vector<double> out_10quantile_marginalized(nparameters);
412  std::vector<double> out_16quantile_marginalized(nparameters);
413  std::vector<double> out_84quantile_marginalized(nparameters);
414  std::vector<double> out_90quantile_marginalized(nparameters);
415  std::vector<double> out_95quantile_marginalized(nparameters);
416  std::vector<double> out_std_marginalized(nparameters);
417  std::vector<double> out_chi2_generated(nchannels);
418  std::vector<double> out_chi2_mode(nchannels);
419  std::vector<double> out_cash_generated(nchannels);
420  std::vector<double> out_cash_mode(nchannels);
421  std::vector<int> out_nevents(nchannels);
422  double out_chi2_generated_total;
423  double out_chi2_mode_total;
424  double out_cash_generated_total;
425  double out_cash_mode_total;
426  int out_nevents_total;
427 
428  // create branches
429  for (int i = 0; i < nparameters; ++i) {
430  tree_out->Branch(Form("parameter_%i", i), &out_parameters[i], Form("parameter %i/D", i));
431  tree_out->Branch(Form("mode_global_%i", i), &out_mode_global[i], Form("global mode of par. %i/D", i));
432  tree_out->Branch(Form("std_global_%i", i), &out_std_global[i], Form("global std of par. %i/D", i));
433  if (fFlagMarginalize) {
434  tree_out->Branch(Form("mode_marginalized_%i", i), &out_mode_marginalized[i], Form("marginalized mode of par. %i/D", i));
435  tree_out->Branch(Form("mean_marginalized_%i", i), &out_mean_marginalized[i], Form("marginalized mean of par. %i/D", i));
436  tree_out->Branch(Form("median_marginalized_%i", i), &out_median_marginalized[i], Form("median of par. %i/D", i));
437  tree_out->Branch(Form("5quantile_marginalized_%i", i), &out_5quantile_marginalized[i], Form("marginalized 5 per cent quantile of par. %i/D", i));
438  tree_out->Branch(Form("10quantile_marginalized_%i", i), &out_10quantile_marginalized[i], Form("marginalized 10 per cent quantile of par. %i/D", i));
439  tree_out->Branch(Form("16quantile_marginalized_%i", i), &out_16quantile_marginalized[i], Form("marginalized 16 per cent quantile of par. %i/D", i));
440  tree_out->Branch(Form("84quantile_marginalized_%i", i), &out_84quantile_marginalized[i], Form("marginalized 84 per cent quantile of par. %i/D", i));
441  tree_out->Branch(Form("90quantile_marginalized_%i", i), &out_90quantile_marginalized[i], Form("marginalized 90 per cent quantile of par. %i/D", i));
442  tree_out->Branch(Form("95quantile_marginalized_%i", i), &out_95quantile_marginalized[i], Form("marginalized 95 per cent quantile of par. %i/D", i));
443  tree_out->Branch(Form("std_marginalized_%i", i), &out_std_marginalized[i], Form("marginalized std of par. %i/D", i));
444  }
445  }
446  for (int i = 0; i < nchannels; ++i) {
447  tree_out->Branch(Form("chi2_generated_%i", i), &out_chi2_generated[i], Form("chi2 (generated par.) in channel %i/D", i));
448  tree_out->Branch(Form("chi2_mode_%i", i), &out_chi2_mode[i], Form("chi2 (mode of par.) in channel %i/D", i));
449  tree_out->Branch(Form("cash_generated_%i", i), &out_cash_generated[i], Form("cash statistic (generated par.) in channel %i/D", i));
450  tree_out->Branch(Form("cash_mode_%i", i), &out_cash_mode[i], Form("cash statistic (mode of par.) in channel %i/D", i));
451  tree_out->Branch(Form("nevents_%i", i), &out_nevents[i], Form("nevents in channel %i/I", i));
452  }
453  tree_out->Branch("chi2_generated_total", &out_chi2_generated_total, "chi2 (generated par.) in all channels/D");
454  tree_out->Branch("chi2_mode_total", &out_chi2_mode_total, "chi2 (mode of par.) in all channels/D");
455  tree_out->Branch("cash_generated_total", &out_cash_generated_total, "cash statistic (generated par.) in all channels/D");
456  tree_out->Branch("cash_mode_total", &out_cash_mode_total, "cash statistic (mode of par.) in all channels/D");
457  tree_out->Branch("nevents_total", &out_nevents_total, "total number of events/I");
458 
459  // define temporary vector of histogram for fluctated templates
460  std::vector<TH1D*> histlist(0);
461 
462  // loop over ensembles
463  for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
464  // print status
465  if ((iensemble+1)%100 == 0 && iensemble > 0) {
466  BCLog::SetLogLevel(lls,llf);
467  int frac = int (double(iensemble+1) / double(nensembles) * 100.);
468  BCLog::OutDetail(Form("Fraction of ensembles analyzed: %i%%",frac));
470  }
471 
472  // get next (commented out: random) event from tree
473  // int index = (int) fRandom->Uniform(tree->GetEntries());
474  int index = iensemble + start;
475  tree->GetEntry(index);
476 
477  // define histograms
478  std::vector<TH1D> histograms;
479 
480  // transform matrix into histograms
481  histograms = MatrixToHistograms(nbins_matrix);
482 
483  // loop over channels and set data
484  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
485  // get channel
486  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
487 
488  // overwrite data
489  if (iensemble == 0)
490  channel->GetData()->SetHistogram(0);
491 
492  // set data histogram
493  fMTF->SetData(channel->GetName().c_str(), histograms.at(ichannel));
494  }
495 
496  // fluctuate templates if option "MC" is chosen
497  if (flag_mc) {
498  // get number of templates
499  unsigned int ntemplates = fMTF->GetNProcesses();
500 
501  // loop over channels
502  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
503  // get channel
504  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
505 
506  // loop over all templates
507  for (unsigned int i = 0; i < ntemplates; ++i) {
508 
509  // get histogram
510  TH1D * temphist = channel->GetTemplate(i)->GetHistogram();
511  histlist.push_back(temphist);
512 
513  // replace by fluctuated histogram
514  if (temphist) {
515  TH1D* temphistfluc = new TH1D(channel->GetTemplate(i)->FluctuateHistogram(options, channel->GetTemplate(i)->GetOriginalNorm()));
516  channel->GetTemplate(i)->SetHistogram(temphistfluc, channel->GetTemplate(i)->GetNorm());
517  }
518  }
519  }
520  }
521 
522  // work-around: force initialization
523  fMTF->ResetResults();
524 
525  // check if marginalization should be run and perform analysis
526  if (fFlagMarginalize) {
527  BCLog::SetLogLevel(lls,llf);
528  BCLog::OutDetail(Form("Running MCMC for ensemble %i",iensemble));
530 
531  // run mcmc
532  fMTF->MarginalizeAll();
533 
534  // find mode
536  }
537  else {
538  // find mode
539  fMTF->FindMode();
540  }
541 
542  // fill tree variables
543  out_mode_global = fMTF->GetBestFitParameters();
544  out_std_global = fMTF->GetBestFitParameterErrors();
545  out_chi2_generated_total = 0;
546  out_chi2_mode_total = 0;
547  out_cash_generated_total = 0;
548  out_cash_mode_total = 0;
549  out_nevents_total = 0;
550 
551  for (int i = 0; i < nparameters; ++i) {
552  if (fFlagMarginalize) {
553  BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
554  out_mode_marginalized[i] = hist->GetMode();
555  out_mean_marginalized[i] = hist->GetMean();
556  out_median_marginalized[i] = hist->GetMedian();
557  out_5quantile_marginalized[i] = hist->GetQuantile(0.05);
558  out_10quantile_marginalized[i] = hist->GetQuantile(0.10);
559  out_16quantile_marginalized[i] = hist->GetQuantile(0.16);
560  out_84quantile_marginalized[i] = hist->GetQuantile(0.84);
561  out_90quantile_marginalized[i] = hist->GetQuantile(0.90);
562  out_95quantile_marginalized[i] = hist->GetQuantile(0.95);
563  out_std_marginalized[i]=hist->GetRMS();
564  }
565  }
566 
567  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
568  // get channel
569  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
570 
571  // get number of events
572  out_nevents[ichannel] = (int) channel->GetData()->GetHistogram()->Integral();
573 
574  // calculate chi2
575  out_chi2_generated[ichannel] = fMTF->CalculateChi2( ichannel, out_parameters );
576  out_chi2_mode[ichannel] = fMTF->CalculateChi2( ichannel, fMTF->GetBestFitParameters() );
577 
578  // calculate cash statistic
579  out_cash_generated[ichannel] = fMTF->CalculateCash( ichannel, out_parameters );
580  out_cash_mode[ichannel] = fMTF->CalculateCash( ichannel, fMTF->GetBestFitParameters() );
581 
582  // increase the total number of events
583  out_nevents_total += out_nevents[ichannel];
584 
585  // increase the total chi2
586  out_chi2_generated_total += out_chi2_generated[ichannel];
587  out_chi2_mode_total += out_chi2_mode[ichannel];
588 
589  // increase the total cash statistic
590  out_cash_generated_total += out_cash_generated[ichannel];
591  out_cash_mode_total += out_cash_mode[ichannel];
592  }
593 
594  // fill tree
595  tree_out->Fill();
596 
597  // but original template back if options "MC" is chosen
598  if (flag_mc) {
599  // get number of templates
600  unsigned int ntemplates = fMTF->GetNProcesses();
601 
602  // loop over channels
603  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
604  // get channel
605  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
606 
607  // loop over all templates
608  for (unsigned int i = 0; i < ntemplates; ++i) {
609 
610  // get histogram
611  TH1D * temphist = histlist.at(ichannel * ntemplates + i);
612  temphist->Scale(channel->GetTemplate(i)->GetOriginalNorm()/ temphist->Integral());
613 
614  // replace by fluctuated histogram
615  TH1D* temphistfluc = channel->GetTemplate(i)->GetHistogram();
616  delete temphistfluc;
617  channel->GetTemplate(i)->SetHistogram(temphist, channel->GetTemplate(i)->GetNorm());
618  }
619  }
620  histlist.clear();
621  }
622  }
623 
624  // put the original data back in place
625  // loop over channels
626  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
627  // get channel
628  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
629 
630  // set data pointer
631  channel->GetData()->SetHistogram(histograms_data.at(ichannel));
632  }
633 
634  // reset log level
635  BCLog::SetLogLevel(lls,llf);
636 
637  // work-around: force initialization
638  fMTF->ResetResults();
639 
640  BCLog::OutSummary("Ensemble test ran successfully.");
641 
642  // return output tree
643  return tree_out;
644 }
645 
646 // ---------------------------------------------------------
647 std::vector<TH1D> BCMTFAnalysisFacility::MatrixToHistograms(const std::vector< std::vector<double> > & matrix)
648 {
649  // create vector of histograms
650  std::vector<TH1D> histograms;
651 
652  // get number of channels
653  int nchannels = matrix.size();;
654 
655  // loop over channels
656  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
657  // get channel
658  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
659 
660  // get column
661  std::vector<double> nbins_column = matrix[ichannel];
662 
663  // create new histogram
664  TH1D hist( *(channel->GetData()->GetHistogram()) );
665 
666  // get number of bins
667  int nbins = hist.GetNbinsX();
668 
669  // fill bin content
670  for (int ibin = 1; ibin <= nbins; ++ibin) {
671  hist.SetBinContent(ibin, nbins_column.at(ibin-1));
672  }
673 
674  // add histogram to container
675  histograms.push_back(hist);
676  }
677 
678  // return histograms
679  return histograms;
680 
681 }
682 
683 // ---------------------------------------------------------
684 int BCMTFAnalysisFacility::PerformSingleChannelAnalyses(const char * dirname, const char * options)
685 {
686  BCLog::OutSummary(Form("Running single channel analysis in directory \'%s\'.",dirname));
687 
688  // ---- create new directory ---- //
689 
690  mkdir(dirname, 0777);
691  chdir(dirname);
692 
693  // ---- check options ---- //
694 
695  bool flag_syst = true;
696  bool flag_mcmc = true;
697 
698  if (std::string(options).find("nosyst") < std::string(options).size())
699  flag_syst = false;
700 
701  if (std::string(options).find("mcmc") < std::string(options).size())
702  flag_mcmc = true;
703 
704  // get number of channels
705  int nchannels = fMTF->GetNChannels();
706 
707  // get number of systematics
708  int nsystematics = fMTF->GetNSystematics();
709 
710  // container of flags for channels
711  std::vector<bool> flag_channel(nchannels);
712 
713  // container of flags for systematics
714  std::vector<bool> flag_systematic(nsystematics);
715 
716  // create new container of comparison tools
717  std::vector<BCMTFComparisonTool *> ctc;
718  std::vector<BCMTFComparisonTool *> ctc_mcmc;
719 
720  // get number of parameters
721  int nparameters = fMTF->GetNParameters();
722 
723  // ---- add one comparison tool for each parameter ---- //
724 
725  // loop over all parameters
726  for (int i = 0; i < nparameters; ++ i) {
727  // create new comparison tool
729  BCMTFComparisonTool * ct_mcmc = new BCMTFComparisonTool(fMTF->GetParameter(i)->GetName().c_str());
730 
731  // add comparison tool to container
732  ctc.push_back(ct);
733  ctc_mcmc.push_back(ct_mcmc);
734  }
735 
736  // ---- switch on/off all systematics ---- //
737 
738  // loop over all systematics
739  for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
740  // get systematic
741  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
742 
743  // remember old setting
744  flag_systematic[isystematic] = systematic->GetFlagSystematicActive();
745 
746  // switch off systematic
747  if (flag_systematic[isystematic])
748  systematic->SetFlagSystematicActive(flag_syst);
749  }
750 
751  // ---- switch on all channels ---- //
752 
753  // loop over all channels
754  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
755  // get channel
756  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
757 
758  // remember old setting
759  flag_channel[ichannel] = channel->GetFlagChannelActive();
760 
761  // switch channel on/off
762  channel->SetFlagChannelActive(true);
763  }
764 
765  // ---- perform analysis for combination ---- //
766  if (flag_mcmc) {
767  // work-around: force initialization
768  fMTF->ResetResults();
769 
770  // run mcmc
771  fMTF->MarginalizeAll();
772 
773  // find mode
775  }
776  else {
777  // find mode
778  fMTF->FindMode();
779  }
780 
781  // print results
782  if (flag_mcmc)
783  fMTF->PrintAllMarginalized("marginalized_combined.pdf");
784  fMTF->PrintResults("results_combined.txt");
785 
786  // loop over all parameters
787  for (int i = 0; i < nparameters; ++ i) {
788  // get comparison tool
789  BCMTFComparisonTool * ct = ctc.at(i);
790  BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
791 
792  ct->AddContribution("all channels",
793  fMTF->GetBestFitParameters().at(i),
795  if (flag_mcmc) {
796  BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
797 
798  ct_mcmc->AddContribution("all channels",
799  hist->GetMean(),
800  hist->GetRMS());
801  }
802  }
803 
804  // ---- switch off all channels ----//
805 
806  // loop over all channels
807  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
808  // get channel
809  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
810 
811  // switch off channel
812  channel->SetFlagChannelActive(false);
813  }
814 
815  // ---- perform analysis on all channels separately ---- //
816 
817  // loop over all channels
818  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
819  // get channel
820  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
821 
822  // switch on one channel
823  channel->SetFlagChannelActive(true);
824 
825  // perform analysis
826 
827  if (flag_mcmc) {
828  // work-around: force initialization
829  fMTF->ResetResults();
830 
831  // run mcmc
832  fMTF->MarginalizeAll();
833 
834  // find mode
836  }
837  else {
838  // find mode
839  fMTF->FindMode();
840  }
841 
842  // print results
843  if (flag_mcmc)
844  fMTF->PrintAllMarginalized(Form("marginalized_channel_%i.pdf", ichannel));
845  fMTF->PrintResults(Form("results_channel_%i.txt", ichannel));
846 
847  // ---- update comparison tools ---- //
848 
849  // loop over all parameters
850  for (int i = 0; i < nparameters; ++ i) {
851  // get comparison tool
852  BCMTFComparisonTool * ct = ctc.at(i);
853  BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
854 
855  ct->AddContribution(channel->GetName().c_str(),
856  fMTF->GetBestFitParameters().at(i),
858  if (flag_mcmc) {
859  BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
860 
861  ct_mcmc->AddContribution(channel->GetName().c_str(),
862  hist->GetMean(),
863  hist->GetRMS());
864  }
865 
866  // switch off channel
867  channel->SetFlagChannelActive(false);
868  }
869  }
870 
871  // ---- reset all systematics ---- //
872  // loop over all systematics
873  for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
874  // get systematic
875  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
876 
877  // switch off systematic
878  if (flag_systematic[isystematic])
879  systematic->SetFlagSystematicActive(flag_systematic[isystematic]);
880  }
881 
882  // ---- reset all channels ---- //
883 
884  // loop over all channels
885  for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
886  // get channel
887  BCMTFChannel * channel = fMTF->GetChannel(ichannel);
888 
889  // switch channel on/off
890  channel->SetFlagChannelActive(flag_channel[ichannel]);
891  }
892 
893  // ---- workaround: reset MCMC ---- //
894  fMTF->ResetResults();
895 
896  // ---- print everything ---- //
897  TCanvas * c1 = new TCanvas();
898  c1->cd();
899 
900  // draw first one
901  BCMTFComparisonTool * ct = ctc.at(0);
902  ct->DrawOverview();
903  // c1->Print((std::string(filename)+std::string("(")).c_str());
904  c1->Print((std::string("overview.pdf")+std::string("(")).c_str());
905 
906  // loop over all parameters
907  for (int i = 1; i < nparameters-1; ++ i) {
908  // create new comparison tool
909  BCMTFComparisonTool * ct = ctc.at(i);
910 
911  ct->DrawOverview();
912  c1->Print("overview.pdf");
913  }
914 
915  // draw last one
916  ct = ctc.at(nparameters-1);
917  ct->DrawOverview();
918  c1->Print((std::string("overview.pdf")+std::string(")")).c_str());
919 
920  // ---- print everything (mcmc) ---- //
921  if (flag_mcmc) {
922  TCanvas * c2 = new TCanvas();
923  c2->cd();
924 
925  // draw first one
926  BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(0);
927  ct_mcmc->DrawOverview();
928  // c2->Print((std::string(filename)+std::string("(")).c_str());
929  c2->Print((std::string("overview_mcmc.pdf")+std::string("(")).c_str());
930 
931  // loop over all parameters
932  for (int i = 1; i < nparameters-1; ++ i) {
933  // create new comparison tool
934  BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
935 
936  ct_mcmc->DrawOverview();
937  c2->Print("overview_mcmc.pdf");
938  }
939 
940  // draw last one
941  ct_mcmc = ctc_mcmc.at(nparameters-1);
942  ct_mcmc->DrawOverview();
943  c2->Print((std::string("overview_mcmc.pdf")+std::string(")")).c_str());
944  delete c2;
945  }
946 
947  // ---- free memory ---- //
948  for (int i = 0; i < nparameters; ++i) {
949  BCMTFComparisonTool * ct = ctc[i];
950  BCMTFComparisonTool * ct_mcmc = ctc_mcmc[i];
951  delete ct;
952  delete ct_mcmc;
953  }
954  ctc.clear();
955  ctc_mcmc.clear();
956 
957  delete c1;
958 
959  // ---- change directory ---- //
960 
961  chdir("../");
962 
963  BCLog::OutSummary("Single channel analysis ran successfully");
964 
965  // no error
966  return 1;
967 }
968 
969 // ---------------------------------------------------------
970 int BCMTFAnalysisFacility::PerformSingleSystematicAnalyses(const char * dirname, const char * options)
971 {
972  BCLog::OutSummary(Form("Running single channel systematic analysis in directory \'%s\'.",dirname));
973 
974  // ---- create new directory ---- //
975 
976  mkdir(dirname, 0777);
977  chdir(dirname);
978 
979  // ---- check options ---- //
980 
981  bool flag_mcmc = true;
982 
983  if (std::string(options).find("mcmc") < std::string(options).size())
984  flag_mcmc = true;
985 
986  // get number of channels
987  int nchannels = fMTF->GetNChannels();
988 
989  // get number of systematics
990  int nsystematics = fMTF->GetNSystematics();
991 
992  // container of flags for channels
993  std::vector<bool> flag_channel(nchannels);
994 
995  // container of flags for systematics
996  std::vector<bool> flag_systematic(nsystematics);
997 
998  // create new container of comparison tools
999  std::vector<BCMTFComparisonTool *> ctc;
1000 
1001  // get number of parameters
1002  int nparameters = fMTF->GetNParameters();
1003 
1004  // ---- add one comparison tool for each systematic ---- //
1005 
1006  // loop over all parameters
1007  for (int i = 0; i < nparameters; ++ i) {
1008  // create new comparison tool
1010 
1011  // add comparison tool to container
1012  ctc.push_back(ct);
1013  }
1014 
1015  // ---- switch on all systematics ---- //
1016 
1017  for (int isystematic = 0; isystematic < nsystematics; ++ isystematic) {
1018  // get systematic
1019  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
1020 
1021  // remember old setting
1022  flag_systematic[isystematic] = systematic->GetFlagSystematicActive();
1023 
1024  // switch on
1025  systematic->SetFlagSystematicActive(true);
1026  }
1027 
1028  if (flag_mcmc) {
1029  // work-around: force initialization
1030  fMTF->ResetResults();
1031 
1032  // run mcmc
1033  fMTF->MarginalizeAll();
1034 
1035  // find mode
1037  }
1038  else {
1039  // find mode
1040  fMTF->FindMode();
1041  }
1042 
1043  // print results
1044  if (flag_mcmc)
1045  fMTF->PrintAllMarginalized("marginalized_all.pdf");
1046  fMTF->PrintResults("results_all.txt");
1047 
1048  // loop over all parameters
1049  for (int i = 0; i < nparameters; ++ i) {
1050  // get comparison tool
1051  BCMTFComparisonTool * ct = ctc.at(i);
1052 
1053  ct->AddContribution("all systematics",
1054  fMTF->GetBestFitParameters().at(i),
1055  fMTF->GetBestFitParameterErrors().at(i));
1056  }
1057 
1058  // ---- switch off all systematics ---- //
1059 
1060  // loop over all systematics
1061  for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
1062  // get systematic
1063  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
1064 
1065  // switch off
1066  systematic->SetFlagSystematicActive(false);
1067  }
1068 
1069  // ---- perform analysis with all systematics separately ---- //
1070 
1071  // loop over all channels
1072  for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
1073  // get systematic
1074  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
1075 
1076  // switch on systematic
1077  systematic->SetFlagSystematicActive(true);
1078 
1079  // perform analysis
1080  if (flag_mcmc) {
1081  // work-around: force initialization
1082  fMTF->ResetResults();
1083 
1084  // run mcmc
1085  fMTF->MarginalizeAll();
1086 
1087  // find mode
1089  }
1090  else {
1091  // find mode
1092  fMTF->FindMode();
1093  }
1094 
1095  // print results
1096  if (flag_mcmc)
1097  fMTF->PrintAllMarginalized(Form("marginalized_systematic_%i.pdf", isystematic));
1098  fMTF->PrintResults(Form("results_systematic_%i.txt", isystematic));
1099 
1100  // ---- update comparison tools ---- //
1101 
1102  // loop over all parameters
1103  for (int i = 0; i < nparameters; ++ i) {
1104  // get comparison tool
1105  BCMTFComparisonTool * ct = ctc.at(i);
1106 
1107  ct->AddContribution(systematic->GetName().c_str(),
1108  fMTF->GetBestFitParameters().at(i),
1109  fMTF->GetBestFitParameterErrors().at(i));
1110  }
1111 
1112  // switch off systematic
1113  systematic->SetFlagSystematicActive(false);
1114  }
1115 
1116  // ---- analysis without any systematic uncertainty ---- //
1117 
1118  if (flag_mcmc) {
1119  // work-around: force initialization
1120  fMTF->ResetResults();
1121 
1122  // run mcmc
1123  fMTF->MarginalizeAll();
1124 
1125  // find mode
1127  }
1128  else {
1129  // find mode
1130  fMTF->FindMode();
1131  }
1132 
1133  // print results
1134  if (flag_mcmc)
1135  fMTF->PrintAllMarginalized("marginalized_none.pdf");
1136  fMTF->PrintResults("results_none.txt");
1137 
1138  // loop over all parameters
1139  for (int i = 0; i < nparameters; ++ i) {
1140  // get comparison tool
1141  BCMTFComparisonTool * ct = ctc.at(i);
1142 
1143  ct->AddContribution("no systematics",
1144  fMTF->GetBestFitParameters().at(i),
1145  fMTF->GetBestFitParameterErrors().at(i));
1146  }
1147 
1148 
1149 
1150  // ---- reset all systematics ---- //
1151 
1152  // loop over all systematics
1153  for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
1154  // get systematic
1155  BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
1156 
1157  // switch off systematic
1158  if (flag_systematic[isystematic])
1159  systematic->SetFlagSystematicActive(flag_systematic[isystematic]);
1160  }
1161 
1162  // ---- workaround: reset MCMC ---- //
1163  fMTF->ResetResults();
1164 
1165  // ---- print everything ---- //
1166  TCanvas * c1 = new TCanvas();
1167  c1->cd();
1168 
1169  // draw first one
1170  BCMTFComparisonTool * ct = ctc.at(0);
1171  ct->DrawOverview();
1172  // c1->Print((std::string(filename)+std::string("(")).c_str());
1173  c1->Print((std::string("overview.pdf")+std::string("(")).c_str());
1174 
1175  // loop over all parameters
1176  for (int i = 1; i < nparameters-1; ++i) {
1177  // create new comparison tool
1178  BCMTFComparisonTool * ct = ctc.at(i);
1179 
1180  ct->DrawOverview();
1181  c1->Print("overview.pdf");
1182  }
1183 
1184  // draw last one
1185  ct = ctc.at(nparameters-1);
1186  ct->DrawOverview();
1187  c1->Print((std::string("overview.pdf")+std::string(")")).c_str());
1188 
1189  // ---- free memory ---- //
1190  for (int i = 0; i < nparameters; ++i) {
1191  BCMTFComparisonTool * ct = ctc[i];
1192  delete ct;
1193  }
1194  ctc.clear();
1195 
1196  delete c1;
1197 
1198  // ---- change directory ---- //
1199 
1200  chdir("../");
1201 
1202  BCLog::OutSummary("Single channel analysis ran successfully");
1203 
1204  // no error
1205  return 1;
1206 }
1207 
1208 // ---------------------------------------------------------
1209 int BCMTFAnalysisFacility::PerformCalibrationAnalysis(const char * dirname, const std::vector<double> & default_parameters, int index, const std::vector<double> & parametervalues, int nensembles)
1210 {
1211  BCLog::OutSummary(Form("Running calibration analysis in directory \'%s\'.",dirname));
1212 
1213  // ---- create new directory ---- //
1214 
1215  mkdir(dirname, 0777);
1216  chdir(dirname);
1217 
1218  // ---- loop over parameter values and perform analysis ---- //
1219 
1220  int nvalues = int(parametervalues.size());
1221  for (int ivalue = 0; ivalue < nvalues; ++ivalue) {
1222 
1223  // open file
1224  TFile * file = TFile::Open(Form("ensemble_%i.root", ivalue), "RECREATE");
1225  file->cd();
1226 
1227  // set parameters
1228  std::vector<double> parameters = default_parameters;
1229  parameters[index] = parametervalues.at(ivalue);
1230 
1231  // create ensemble
1232  TTree * tree = PerformEnsembleTest(parameters, nensembles);
1233 
1234  // write tree
1235  tree->Write();
1236 
1237  // close file
1238  file->Close();
1239 
1240  // free memory
1241  delete file;
1242  }
1243 
1244  // ---- change directory ---- //
1245 
1246  chdir("../");
1247 
1248  BCLog::OutSummary("Calibration analysis ran successfully");
1249 
1250  // no error
1251  return 1;
1252 }
1253 
1254 // ---------------------------------------------------------