00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <sys/stat.h>
00011 #include <unistd.h>
00012 #include <iostream>
00013
00014 #include <TROOT.h>
00015 #include <TCanvas.h>
00016 #include <TTree.h>
00017 #include <TH1D.h>
00018 #include <TRandom3.h>
00019 #include <TFile.h>
00020
00021 #include "../../BAT/BCLog.h"
00022 #include "../../BAT/BCH1D.h"
00023
00024 #include "BCMTF.h"
00025 #include "BCMTFChannel.h"
00026 #include "BCMTFTemplate.h"
00027 #include "BCMTFSystematic.h"
00028 #include "BCMTFComparisonTool.h"
00029
00030 #include "BCMTFAnalysisFacility.h"
00031
00032
00033 BCMTFAnalysisFacility::BCMTFAnalysisFacility(BCMTF * mtf)
00034 : fRandom(new TRandom3(0))
00035 , fFlagMCMC(false)
00036 , fLogLevel(BCLog::nothing)
00037 {
00038 fMTF = mtf;
00039 BCLog::OutDetail(Form("Prepared Analysis Facility for MTF model \'%s\'",mtf->GetName().c_str()));
00040 }
00041
00042
00043 BCMTFAnalysisFacility::~BCMTFAnalysisFacility()
00044 {
00045 }
00046
00047
00048 std::vector<TH1D> BCMTFAnalysisFacility::BuildEnsemble(const std::vector<double> & parameters)
00049 {
00050
00051 int nchannels = fMTF->GetNChannels();
00052
00053
00054 std::vector<TH1D> histograms;
00055
00056
00057 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00058
00059
00060 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00061
00062
00063 TH1D hist( *(channel->GetData()->GetHistogram()) );
00064
00065
00066 int nbins = hist.GetNbinsX();
00067
00068
00069 for (int ibin = 1; ibin <= nbins; ++ibin) {
00070 double expectation = fMTF->Expectation(ichannel, ibin, parameters);
00071
00072 double observation = gRandom->Poisson(expectation);
00073
00074 hist.SetBinContent(ibin, observation);
00075 }
00076
00077
00078 histograms.push_back(hist);
00079 }
00080
00081
00082 return histograms;
00083 }
00084
00085
00086 TTree * BCMTFAnalysisFacility::BuildEnsembles(TTree * tree, int nensembles)
00087 {
00088
00089 int nchannels = fMTF->GetNChannels();
00090
00091 BCLog::OutDetail(Form("MTF Building %d ensambles for %d channels.",nensembles,nchannels));
00092
00093
00094 int nparameters = fMTF->GetNParameters();
00095
00096
00097 std::vector<double> parameters(nparameters);
00098
00099
00100 for (int i = 0; i < nparameters; ++i) {
00101 tree->SetBranchAddress(Form("Parameter%i", i), &(parameters[i]));
00102 }
00103
00104
00105 TTree * tree_out = new TTree("ensembles", "ensembles");
00106
00107
00108 std::vector< std::vector<double> > nbins_matrix;
00109
00110
00111
00112 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00113
00114 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00115
00116
00117 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00118
00119
00120 std::vector<double> nbins_column(nbins);
00121
00122
00123 nbins_matrix.push_back(nbins_column);
00124 }
00125
00126 std::vector<double> in_parameters(nparameters);
00127
00128
00129
00130 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00131
00132 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00133
00134
00135 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00136
00137
00138 for (int ibin = 1; ibin <= nbins; ++ibin) {
00139
00140 tree_out->Branch(Form("channel_%i_bin_%i", ichannel, ibin),
00141 &(nbins_matrix[ichannel])[ibin-1], "n/D");
00142 }
00143 }
00144
00145 for (int i = 0; i < nparameters; ++i) {
00146 tree_out->Branch(Form("parameter_%i", i), &in_parameters[i], Form("parameter_%i/D", i));
00147 }
00148
00149
00150 std::vector<TH1D> histograms;
00151
00152
00153 for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
00154
00155 int index = (int) fRandom->Uniform(tree->GetEntries());
00156 tree->GetEntry(index);
00157
00158
00159 histograms = BuildEnsemble(parameters);
00160
00161
00162
00163 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00164
00165 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00166
00167
00168 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00169
00170
00171 for (int ibin = 1; ibin <= nbins; ++ibin) {
00172
00173 (nbins_matrix[ichannel])[ibin-1] = histograms.at(ichannel).GetBinContent(ibin);
00174 }
00175 }
00176
00177
00178 for (int i = 0; i < nparameters; ++i) {
00179 in_parameters[i] = parameters.at(i);
00180 }
00181
00182
00183 tree_out->Fill();
00184 }
00185
00186
00187 return tree_out;
00188 }
00189
00190
00191 TTree * BCMTFAnalysisFacility::BuildEnsembles(const std::vector<double> & parameters, int nensembles)
00192 {
00193
00194 int nchannels = fMTF->GetNChannels();
00195
00196 BCLog::OutDetail(Form("MTF Building %d ensambles for %d channels.",nensembles,nchannels));
00197
00198
00199 TTree * tree = new TTree("ensembles", "ensembles");
00200
00201
00202 std::vector< std::vector<double> > nbins_matrix;
00203
00204
00205
00206 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00207
00208 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00209
00210
00211 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00212
00213
00214 std::vector<double> nbins_column(nbins);
00215
00216
00217 nbins_matrix.push_back(nbins_column);
00218 }
00219
00220
00221 int nparameters = fMTF->GetNParameters();
00222
00223 std::vector<double> in_parameters(nparameters);
00224
00225
00226
00227 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00228
00229 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00230
00231
00232 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00233
00234
00235 for (int ibin = 1; ibin <= nbins; ++ibin) {
00236
00237 tree->Branch(Form("channel_%i_bin_%i", ichannel, ibin),
00238 &(nbins_matrix[ichannel])[ibin-1], "n/D");
00239 }
00240 }
00241
00242 for (int i = 0; i < nparameters; ++i) {
00243 tree->Branch(Form("parameter_%i", i), &in_parameters[i], Form("parameter_%i/D", i));
00244 }
00245
00246
00247 std::vector<TH1D> histograms;
00248
00249
00250 for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
00251
00252 histograms = BuildEnsemble(parameters);
00253
00254
00255
00256 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00257
00258 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00259
00260
00261 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00262
00263
00264 for (int ibin = 1; ibin <= nbins; ++ibin) {
00265
00266 (nbins_matrix[ichannel])[ibin-1] = histograms.at(ichannel).GetBinContent(ibin);
00267 }
00268 }
00269
00270
00271 for (int i = 0; i < nparameters; ++i) {
00272 in_parameters[i] = parameters.at(i);
00273 }
00274
00275
00276 tree->Fill();
00277 }
00278
00279
00280 return tree;
00281 }
00282
00283
00284 TTree * BCMTFAnalysisFacility::PerformEnsembleTest(const std::vector<double> & parameters, int nensembles)
00285 {
00286
00287 TTree * tree = 0;
00288
00289
00290 tree = BuildEnsembles(parameters, nensembles);
00291
00292
00293 return PerformEnsembleTest(tree, nensembles);
00294 }
00295
00296
00297 TTree * BCMTFAnalysisFacility::PerformEnsembleTest(TTree * tree, int nensembles, int start)
00298 {
00299 BCLog::OutSummary("Running ensemble test.");
00300 if (fFlagMCMC) {
00301 BCLog::OutSummary("Fit for each ensemble is going to be run using MCMC. It can take a while.");
00302 }
00303
00304
00305
00306
00307 BCLog::LogLevel lls = BCLog::GetLogLevelScreen();
00308 BCLog::LogLevel llf = BCLog::GetLogLevelFile();
00309 if(fLogLevel==BCLog::nothing) {
00310 BCLog::OutSummary("No log messages for the ensemble fits are going to be printed.");
00311 BCLog::SetLogLevel(fLogLevel);
00312 }
00313 else if(fLogLevel!=lls) {
00314 BCLog::OutSummary(Form("The log level for the ensemble test is set to \'%s\'.",BCLog::ToString(fLogLevel)));
00315 BCLog::SetLogLevel(fLogLevel);
00316 }
00317
00318
00319 int nchannels = fMTF->GetNChannels();
00320
00321
00322 std::vector<TH1D *> histograms_data(nchannels);
00323
00324
00325 std::vector< std::vector<double> > nbins_matrix;
00326
00327
00328
00329 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00330
00331 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00332
00333
00334 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00335
00336
00337 std::vector<double> nbins_column(nbins);
00338
00339
00340 nbins_matrix.push_back(nbins_column);
00341 }
00342
00343
00344 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00345
00346 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00347
00348
00349 int nbins = channel->GetData()->GetHistogram()->GetNbinsX();
00350
00351
00352 for (int ibin = 1; ibin <= nbins; ++ibin) {
00353
00354 tree->SetBranchAddress(Form("channel_%i_bin_%i", ichannel, ibin),
00355 &(nbins_matrix[ichannel])[ibin-1]);
00356 }
00357 }
00358
00359
00360 int nparameters = fMTF->GetNParameters();
00361
00362
00363 std::vector<double> out_parameters(nparameters);
00364
00365 for (int i = 0; i < nparameters; ++i) {
00366 tree->SetBranchAddress(Form("parameter_%i", i), &out_parameters[i]);
00367 }
00368
00369
00370
00371 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00372
00373 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00374
00375
00376 histograms_data[ichannel] = channel->GetData()->GetHistogram();
00377 }
00378
00379
00380 TTree * tree_out = new TTree("ensemble_test", "ensemble test");
00381
00382
00383 std::vector<double> out_mode_global(nparameters);
00384 std::vector<double> out_std_global(nparameters);
00385 std::vector<double> out_mode_marginalized(nparameters);
00386 std::vector<double> out_mean_marginalized(nparameters);
00387 std::vector<double> out_median_marginalized(nparameters);
00388 std::vector<double> out_5quantile_marginalized(nparameters);
00389 std::vector<double> out_10quantile_marginalized(nparameters);
00390 std::vector<double> out_16quantile_marginalized(nparameters);
00391 std::vector<double> out_84quantile_marginalized(nparameters);
00392 std::vector<double> out_90quantile_marginalized(nparameters);
00393 std::vector<double> out_95quantile_marginalized(nparameters);
00394 std::vector<double> out_std_marginalized(nparameters);
00395 std::vector<double> out_chi2_generated(nchannels);
00396 std::vector<double> out_chi2_mode(nchannels);
00397 std::vector<double> out_cash_generated(nchannels);
00398 std::vector<double> out_cash_mode(nchannels);
00399 std::vector<int> out_nevents(nchannels);
00400 double out_chi2_generated_total;
00401 double out_chi2_mode_total;
00402 double out_cash_generated_total;
00403 double out_cash_mode_total;
00404 int out_nevents_total;
00405
00406
00407 for (int i = 0; i < nparameters; ++i) {
00408 tree_out->Branch(Form("parameter_%i", i), &out_parameters[i], Form("parameter %i/D", i));
00409 tree_out->Branch(Form("mode_global_%i", i), &out_mode_global[i], Form("global mode of par. %i/D", i));
00410 tree_out->Branch(Form("std_global_%i", i), &out_std_global[i], Form("global std of par. %i/D", i));
00411 if (fFlagMCMC) {
00412 tree_out->Branch(Form("mode_marginalized_%i", i), &out_mode_marginalized[i], Form("marginalized mode of par. %i/D", i));
00413 tree_out->Branch(Form("mean_marginalized_%i", i), &out_mean_marginalized[i], Form("marginalized mean of par. %i/D", i));
00414 tree_out->Branch(Form("median_marginalized_%i", i), &out_median_marginalized[i], Form("median of par. %i/D", i));
00415 tree_out->Branch(Form("5quantile_marginalized_%i", i), &out_5quantile_marginalized[i], Form("marginalized 5 per cent quantile of par. %i/D", i));
00416 tree_out->Branch(Form("10quantile_marginalized_%i", i), &out_10quantile_marginalized[i], Form("marginalized 10 per cent quantile of par. %i/D", i));
00417 tree_out->Branch(Form("16quantile_marginalized_%i", i), &out_16quantile_marginalized[i], Form("marginalized 16 per cent quantile of par. %i/D", i));
00418 tree_out->Branch(Form("84quantile_marginalized_%i", i), &out_84quantile_marginalized[i], Form("marginalized 84 per cent quantile of par. %i/D", i));
00419 tree_out->Branch(Form("90quantile_marginalized_%i", i), &out_90quantile_marginalized[i], Form("marginalized 90 per cent quantile of par. %i/D", i));
00420 tree_out->Branch(Form("95quantile_marginalized_%i", i), &out_95quantile_marginalized[i], Form("marginalized 95 per cent quantile of par. %i/D", i));
00421 tree_out->Branch(Form("std_marginalized_%i", i), &out_std_marginalized[i], Form("marginalized std of par. %i/D", i));
00422 }
00423 }
00424 for (int i = 0; i < nchannels; ++i) {
00425 tree_out->Branch(Form("chi2_generated_%i", i), &out_chi2_generated[i], Form("chi2 (generated par.) in channel %i/D", i));
00426 tree_out->Branch(Form("chi2_mode_%i", i), &out_chi2_mode[i], Form("chi2 (mode of par.) in channel %i/D", i));
00427 tree_out->Branch(Form("cash_generated_%i", i), &out_cash_generated[i], Form("cash statistic (generated par.) in channel %i/D", i));
00428 tree_out->Branch(Form("cash_mode_%i", i), &out_cash_mode[i], Form("cash statistic (mode of par.) in channel %i/D", i));
00429 tree_out->Branch(Form("nevents_%i", i), &out_nevents[i], Form("nevents in channel %i/I", i));
00430 }
00431 tree_out->Branch("chi2_generated_total", &out_chi2_generated_total, "chi2 (generated par.) in all channels/D");
00432 tree_out->Branch("chi2_mode_total", &out_chi2_mode_total, "chi2 (mode of par.) in all channels/D");
00433 tree_out->Branch("cash_generated_total", &out_cash_generated_total, "cash statistic (generated par.) in all channels/D");
00434 tree_out->Branch("cash_mode_total", &out_cash_mode_total, "cash statistic (mode of par.) in all channels/D");
00435 tree_out->Branch("nevents_total", &out_nevents_total, "total number of events/I");
00436
00437
00438 for (int iensemble = 0; iensemble < nensembles; ++iensemble) {
00439
00440 if ((iensemble+1)%100 == 0 && iensemble > 0) {
00441 BCLog::SetLogLevel(lls,llf);
00442 int frac = double(iensemble+1) / double(nensembles) * 100.;
00443 BCLog::OutDetail(Form("Fraction of ensembles analyzed: %i%%",frac));
00444 BCLog::SetLogLevel(fLogLevel);
00445 }
00446
00447
00448
00449 int index = iensemble + start;
00450 tree->GetEntry(index);
00451
00452
00453 std::vector<TH1D> histograms;
00454
00455
00456 histograms = MatrixToHistograms(nbins_matrix);
00457
00458
00459 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00460
00461 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00462
00463
00464 if (iensemble == 0)
00465 channel->GetData()->SetHistogram(0);
00466
00467
00468 fMTF->SetData(channel->GetName().c_str(), histograms.at(ichannel));
00469 }
00470
00471
00472 if (fFlagMCMC) {
00473 BCLog::SetLogLevel(lls,llf);
00474 BCLog::OutDetail(Form("Running MCMC for ensemble %i",iensemble));
00475 BCLog::SetLogLevel(fLogLevel);
00476
00477
00478 fMTF->MCMCResetResults();
00479
00480
00481 fMTF->MarginalizeAll();
00482
00483
00484 fMTF->FindMode( fMTF->GetBestFitParameters() );
00485 }
00486 else {
00487
00488 fMTF->FindMode();
00489 }
00490
00491
00492 out_mode_global = fMTF->GetBestFitParameters();
00493 out_std_global = fMTF->GetBestFitParameterErrors();
00494 out_chi2_generated_total = 0;
00495 out_chi2_mode_total = 0;
00496 out_cash_generated_total = 0;
00497 out_cash_mode_total = 0;
00498 out_nevents_total = 0;
00499
00500 for (int i = 0; i < nparameters; ++i) {
00501 if (fFlagMCMC) {
00502 BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
00503 out_mode_marginalized[i] = hist->GetMode();
00504 out_mean_marginalized[i] = hist->GetMean();
00505 out_median_marginalized[i] = hist->GetMedian();
00506 out_5quantile_marginalized[i] = hist->GetQuantile(0.05);
00507 out_10quantile_marginalized[i] = hist->GetQuantile(0.10);
00508 out_16quantile_marginalized[i] = hist->GetQuantile(0.16);
00509 out_84quantile_marginalized[i] = hist->GetQuantile(0.84);
00510 out_90quantile_marginalized[i] = hist->GetQuantile(0.90);
00511 out_95quantile_marginalized[i] = hist->GetQuantile(0.95);
00512 out_std_marginalized[i]=hist->GetRMS();
00513 }
00514 }
00515
00516 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00517
00518 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00519
00520
00521 out_nevents[ichannel] = (int) channel->GetData()->GetHistogram()->Integral();
00522
00523
00524 out_chi2_generated[ichannel] = fMTF->CalculateChi2( ichannel, out_parameters );
00525 out_chi2_mode[ichannel] = fMTF->CalculateChi2( ichannel, fMTF->GetBestFitParameters() );
00526
00527
00528 out_cash_generated[ichannel] = fMTF->CalculateCash( ichannel, out_parameters );
00529 out_cash_mode[ichannel] = fMTF->CalculateCash( ichannel, fMTF->GetBestFitParameters() );
00530
00531
00532 out_nevents_total += out_nevents[ichannel];
00533
00534
00535 out_chi2_generated_total += out_chi2_generated[ichannel];
00536 out_chi2_mode_total += out_chi2_mode[ichannel];
00537
00538
00539 out_cash_generated_total += out_cash_generated[ichannel];
00540 out_cash_mode_total += out_cash_mode[ichannel];
00541 }
00542
00543
00544 tree_out->Fill();
00545 }
00546
00547
00548
00549 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00550
00551 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00552
00553
00554 channel->GetData()->SetHistogram(histograms_data.at(ichannel));
00555 }
00556
00557
00558 BCLog::SetLogLevel(lls,llf);
00559
00560
00561 fMTF->MCMCResetResults();
00562
00563 BCLog::OutSummary("Ensemble test ran successfully.");
00564
00565
00566 return tree_out;
00567 }
00568
00569
00570 std::vector<TH1D> BCMTFAnalysisFacility::MatrixToHistograms(const std::vector< std::vector<double> > & matrix)
00571 {
00572
00573 std::vector<TH1D> histograms;
00574
00575
00576 int nchannels = matrix.size();;
00577
00578
00579 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00580
00581 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00582
00583
00584 std::vector<double> nbins_column = matrix[ichannel];
00585
00586
00587 TH1D hist( *(channel->GetData()->GetHistogram()) );
00588
00589
00590 int nbins = hist.GetNbinsX();
00591
00592
00593 for (int ibin = 1; ibin <= nbins; ++ibin) {
00594 hist.SetBinContent(ibin, nbins_column.at(ibin-1));
00595 }
00596
00597
00598 histograms.push_back(hist);
00599 }
00600
00601
00602 return histograms;
00603
00604 }
00605
00606
00607 int BCMTFAnalysisFacility::PerformSingleChannelAnalyses(const char * dirname, const char * options)
00608 {
00609 BCLog::OutSummary(Form("Running single channel analysis in directory \'%s\'.",dirname));
00610
00611
00612
00613 mkdir(dirname, 0777);
00614 chdir(dirname);
00615
00616
00617
00618 bool flag_syst = true;
00619 bool flag_mcmc = true;
00620
00621 if (std::string(options).find("nosyst") < std::string(options).size())
00622 flag_syst = false;
00623
00624 if (std::string(options).find("mcmc") < std::string(options).size())
00625 flag_mcmc = true;
00626
00627
00628 int nchannels = fMTF->GetNChannels();
00629
00630
00631 int nsystematics = fMTF->GetNSystematics();
00632
00633
00634 std::vector<bool> flag_channel(nchannels);
00635
00636
00637 std::vector<bool> flag_systematic(nsystematics);
00638
00639
00640 std::vector<BCMTFComparisonTool *> ctc;
00641 std::vector<BCMTFComparisonTool *> ctc_mcmc;
00642
00643
00644 int nparameters = fMTF->GetNParameters();
00645
00646
00647
00648
00649 for (int i = 0; i < nparameters; ++ i) {
00650
00651 BCMTFComparisonTool * ct = new BCMTFComparisonTool(fMTF->GetParameter(i)->GetName().c_str());
00652 BCMTFComparisonTool * ct_mcmc = new BCMTFComparisonTool(fMTF->GetParameter(i)->GetName().c_str());
00653
00654
00655 ctc.push_back(ct);
00656 ctc_mcmc.push_back(ct_mcmc);
00657 }
00658
00659
00660
00661
00662 for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
00663
00664 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
00665
00666
00667 flag_systematic[isystematic] = systematic->GetFlagSystematicActive();
00668
00669
00670 if (flag_systematic[isystematic])
00671 systematic->SetFlagSystematicActive(flag_syst);
00672 }
00673
00674
00675
00676
00677 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00678
00679 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00680
00681
00682 flag_channel[ichannel] = channel->GetFlagChannelActive();
00683
00684
00685 channel->SetFlagChannelActive(true);
00686 }
00687
00688
00689 if (flag_mcmc) {
00690
00691 fMTF->MCMCResetResults();
00692
00693
00694 fMTF->MarginalizeAll();
00695
00696
00697 fMTF->FindMode( fMTF->GetBestFitParameters() );
00698 }
00699 else {
00700
00701 fMTF->FindMode();
00702 }
00703
00704
00705 if (flag_mcmc)
00706 fMTF->PrintAllMarginalized("marginalized_combined.ps");
00707 fMTF->PrintResults("results_combined.txt");
00708
00709
00710 for (int i = 0; i < nparameters; ++ i) {
00711
00712 BCMTFComparisonTool * ct = ctc.at(i);
00713 BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
00714
00715 ct->AddContribution("all channels",
00716 fMTF->GetBestFitParameters().at(i),
00717 fMTF->GetBestFitParameterErrors().at(i));
00718 if (flag_mcmc) {
00719 BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
00720
00721 ct_mcmc->AddContribution("all channels",
00722 hist->GetMean(),
00723 hist->GetRMS());
00724 }
00725 }
00726
00727
00728
00729
00730 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00731
00732 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00733
00734
00735 channel->SetFlagChannelActive(false);
00736 }
00737
00738
00739
00740
00741 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00742
00743 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00744
00745
00746 channel->SetFlagChannelActive(true);
00747
00748
00749
00750 if (flag_mcmc) {
00751
00752 fMTF->MCMCResetResults();
00753
00754
00755 fMTF->MarginalizeAll();
00756
00757
00758 fMTF->FindMode( fMTF->GetBestFitParameters() );
00759 }
00760 else {
00761
00762 fMTF->FindMode();
00763 }
00764
00765
00766 if (flag_mcmc)
00767 fMTF->PrintAllMarginalized(Form("marginalized_channel_%i.ps", ichannel));
00768 fMTF->PrintResults(Form("results_channel_%i.txt", ichannel));
00769
00770
00771
00772
00773 for (int i = 0; i < nparameters; ++ i) {
00774
00775 BCMTFComparisonTool * ct = ctc.at(i);
00776 BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
00777
00778 ct->AddContribution(channel->GetName().c_str(),
00779 fMTF->GetBestFitParameters().at(i),
00780 fMTF->GetBestFitParameterErrors().at(i));
00781 if (flag_mcmc) {
00782 BCH1D * hist = fMTF->GetMarginalized( fMTF->GetParameter(i) );
00783
00784 ct_mcmc->AddContribution(channel->GetName().c_str(),
00785 hist->GetMean(),
00786 hist->GetRMS());
00787 }
00788
00789
00790 channel->SetFlagChannelActive(false);
00791 }
00792 }
00793
00794
00795
00796 for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
00797
00798 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
00799
00800
00801 if (flag_systematic[isystematic])
00802 systematic->SetFlagSystematicActive(flag_systematic[isystematic]);
00803 }
00804
00805
00806
00807
00808 for (int ichannel = 0; ichannel < nchannels; ++ichannel) {
00809
00810 BCMTFChannel * channel = fMTF->GetChannel(ichannel);
00811
00812
00813 channel->SetFlagChannelActive(flag_channel[ichannel]);
00814 }
00815
00816
00817 fMTF->MCMCResetResults();
00818
00819
00820 TCanvas * c1 = new TCanvas();
00821 c1->cd();
00822
00823
00824 BCMTFComparisonTool * ct = ctc.at(0);
00825 ct->DrawOverview();
00826
00827 c1->Print((std::string("overview.ps")+std::string("(")).c_str());
00828
00829
00830 for (int i = 1; i < nparameters-1; ++ i) {
00831
00832 BCMTFComparisonTool * ct = ctc.at(i);
00833
00834 ct->DrawOverview();
00835 c1->Print("overview.ps");
00836 }
00837
00838
00839 ct = ctc.at(nparameters-1);
00840 ct->DrawOverview();
00841 c1->Print((std::string("overview.ps")+std::string(")")).c_str());
00842
00843
00844 if (flag_mcmc) {
00845 TCanvas * c2 = new TCanvas();
00846 c2->cd();
00847
00848
00849 BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(0);
00850 ct_mcmc->DrawOverview();
00851
00852 c2->Print((std::string("overview_mcmc.ps")+std::string("(")).c_str());
00853
00854
00855 for (int i = 1; i < nparameters-1; ++ i) {
00856
00857 BCMTFComparisonTool * ct_mcmc = ctc_mcmc.at(i);
00858
00859 ct_mcmc->DrawOverview();
00860 c2->Print("overview_mcmc.ps");
00861 }
00862
00863
00864 ct_mcmc = ctc_mcmc.at(nparameters-1);
00865 ct_mcmc->DrawOverview();
00866 c2->Print((std::string("overview_mcmc.ps")+std::string(")")).c_str());
00867 delete c2;
00868 }
00869
00870
00871 for (int i = 0; i < nparameters; ++i) {
00872 BCMTFComparisonTool * ct = ctc[i];
00873 BCMTFComparisonTool * ct_mcmc = ctc_mcmc[i];
00874 delete ct;
00875 delete ct_mcmc;
00876 }
00877 ctc.clear();
00878 ctc_mcmc.clear();
00879
00880 delete c1;
00881
00882
00883
00884 chdir("../");
00885
00886 BCLog::OutSummary("Single channel analysis ran successfully");
00887
00888
00889 return 1;
00890 }
00891
00892
00893 int BCMTFAnalysisFacility::PerformSingleSystematicAnalyses(const char * dirname, const char * options)
00894 {
00895 BCLog::OutSummary(Form("Running single channel systematic analysis in directory \'%s\'.",dirname));
00896
00897
00898
00899 mkdir(dirname, 0777);
00900 chdir(dirname);
00901
00902
00903
00904 bool flag_mcmc = true;
00905
00906 if (std::string(options).find("mcmc") < std::string(options).size())
00907 flag_mcmc = true;
00908
00909
00910 int nchannels = fMTF->GetNChannels();
00911
00912
00913 int nsystematics = fMTF->GetNSystematics();
00914
00915
00916 std::vector<bool> flag_channel(nchannels);
00917
00918
00919 std::vector<bool> flag_systematic(nsystematics);
00920
00921
00922 std::vector<BCMTFComparisonTool *> ctc;
00923
00924
00925 int nparameters = fMTF->GetNParameters();
00926
00927
00928
00929
00930 for (int i = 0; i < nparameters; ++ i) {
00931
00932 BCMTFComparisonTool * ct = new BCMTFComparisonTool(fMTF->GetParameter(i)->GetName().c_str());
00933
00934
00935 ctc.push_back(ct);
00936 }
00937
00938
00939
00940 for (int isystematic = 0; isystematic < nsystematics; ++ isystematic) {
00941
00942 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
00943
00944
00945 flag_systematic[isystematic] = systematic->GetFlagSystematicActive();
00946
00947
00948 systematic->SetFlagSystematicActive(true);
00949 }
00950
00951 if (flag_mcmc) {
00952
00953 fMTF->MCMCResetResults();
00954
00955
00956 fMTF->MarginalizeAll();
00957
00958
00959 fMTF->FindMode( fMTF->GetBestFitParameters() );
00960 }
00961 else {
00962
00963 fMTF->FindMode();
00964 }
00965
00966
00967 if (flag_mcmc)
00968 fMTF->PrintAllMarginalized("marginalized_all.ps");
00969 fMTF->PrintResults("results_all.txt");
00970
00971
00972 for (int i = 0; i < nparameters; ++ i) {
00973
00974 BCMTFComparisonTool * ct = ctc.at(i);
00975
00976 ct->AddContribution("all systematics",
00977 fMTF->GetBestFitParameters().at(i),
00978 fMTF->GetBestFitParameterErrors().at(i));
00979 }
00980
00981
00982
00983
00984 for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
00985
00986 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
00987
00988
00989 systematic->SetFlagSystematicActive(false);
00990 }
00991
00992
00993
00994
00995 for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
00996
00997 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
00998
00999
01000 systematic->SetFlagSystematicActive(true);
01001
01002
01003 if (flag_mcmc) {
01004
01005 fMTF->MCMCResetResults();
01006
01007
01008 fMTF->MarginalizeAll();
01009
01010
01011 fMTF->FindMode( fMTF->GetBestFitParameters() );
01012 }
01013 else {
01014
01015 fMTF->FindMode();
01016 }
01017
01018
01019 if (flag_mcmc)
01020 fMTF->PrintAllMarginalized(Form("marginalized_systematic_%i.ps", isystematic));
01021 fMTF->PrintResults(Form("results_systematic_%i.txt", isystematic));
01022
01023
01024
01025
01026 for (int i = 0; i < nparameters; ++ i) {
01027
01028 BCMTFComparisonTool * ct = ctc.at(i);
01029
01030 ct->AddContribution(systematic->GetName().c_str(),
01031 fMTF->GetBestFitParameters().at(i),
01032 fMTF->GetBestFitParameterErrors().at(i));
01033 }
01034
01035
01036 systematic->SetFlagSystematicActive(false);
01037 }
01038
01039
01040
01041 if (flag_mcmc) {
01042
01043 fMTF->MCMCResetResults();
01044
01045
01046 fMTF->MarginalizeAll();
01047
01048
01049 fMTF->FindMode( fMTF->GetBestFitParameters() );
01050 }
01051 else {
01052
01053 fMTF->FindMode();
01054 }
01055
01056
01057 if (flag_mcmc)
01058 fMTF->PrintAllMarginalized("marginalized_none.ps");
01059 fMTF->PrintResults("results_none.txt");
01060
01061
01062 for (int i = 0; i < nparameters; ++ i) {
01063
01064 BCMTFComparisonTool * ct = ctc.at(i);
01065
01066 ct->AddContribution("no systematics",
01067 fMTF->GetBestFitParameters().at(i),
01068 fMTF->GetBestFitParameterErrors().at(i));
01069 }
01070
01071
01072
01073
01074
01075
01076 for (int isystematic = 0; isystematic < nsystematics; ++isystematic) {
01077
01078 BCMTFSystematic * systematic = fMTF->GetSystematic(isystematic);
01079
01080
01081 if (flag_systematic[isystematic])
01082 systematic->SetFlagSystematicActive(flag_systematic[isystematic]);
01083 }
01084
01085
01086 fMTF->MCMCResetResults();
01087
01088
01089 TCanvas * c1 = new TCanvas();
01090 c1->cd();
01091
01092
01093 BCMTFComparisonTool * ct = ctc.at(0);
01094 ct->DrawOverview();
01095
01096 c1->Print((std::string("overview.ps")+std::string("(")).c_str());
01097
01098
01099 for (int i = 1; i < nparameters-1; ++i) {
01100
01101 BCMTFComparisonTool * ct = ctc.at(i);
01102
01103 ct->DrawOverview();
01104 c1->Print("overview.ps");
01105 }
01106
01107
01108 ct = ctc.at(nparameters-1);
01109 ct->DrawOverview();
01110 c1->Print((std::string("overview.ps")+std::string(")")).c_str());
01111
01112
01113 for (int i = 0; i < nparameters; ++i) {
01114 BCMTFComparisonTool * ct = ctc[i];
01115 delete ct;
01116 }
01117 ctc.clear();
01118
01119 delete c1;
01120
01121
01122
01123 chdir("../");
01124
01125 BCLog::OutSummary("Single channel analysis ran successfully");
01126
01127
01128 return 1;
01129 }
01130
01131
01132 int BCMTFAnalysisFacility::PerformCalibrationAnalysis(const char * dirname, const std::vector<double> & default_parameters, int index, const std::vector<double> & parametervalues, int nensembles)
01133 {
01134 BCLog::OutSummary(Form("Running calibration analysis in directory \'%s\'.",dirname));
01135
01136
01137
01138 mkdir(dirname, 0777);
01139 chdir(dirname);
01140
01141
01142
01143 int nvalues = int(parametervalues.size());
01144 for (int ivalue = 0; ivalue < nvalues; ++ivalue) {
01145
01146
01147 TFile * file = new TFile(Form("ensemble_%i.root", ivalue), "RECREATE");
01148 file->cd();
01149
01150
01151 std::vector<double> parameters = default_parameters;
01152 parameters[index] = parametervalues.at(ivalue);
01153
01154
01155 TTree * tree = PerformEnsembleTest(parameters, nensembles);
01156
01157
01158 tree->Write();
01159
01160
01161 file->Close();
01162
01163
01164 delete file;
01165 }
01166
01167
01168
01169 chdir("../");
01170
01171 BCLog::OutSummary("Calibration analysis ran successfully");
01172
01173
01174 return 1;
01175 }
01176
01177