00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "BAT/BCModelOutput.h"
00011
00012 #include "BAT/BCModel.h"
00013 #include "BAT/BCParameter.h"
00014 #include "BAT/BCH1D.h"
00015 #include "BAT/BCH2D.h"
00016
00017 #include <TDirectory.h>
00018 #include <TFile.h>
00019 #include <TTree.h>
00020 #include <TObject.h>
00021 #include <TH2D.h>
00022 #include <TH1D.h>
00023 #include <TString.h>
00024
00025
00026
00027 BCModelOutput::BCModelOutput()
00028 {
00029 fIndex = 0;
00030
00031 fOutputFile = 0;
00032 fMarkovChainTree = 0;
00033 fAnalysisTree = 0;
00034 fModel = 0;
00035 fOutputFile = 0;
00036 }
00037
00038
00039
00040 BCModelOutput::BCModelOutput(BCModel * model, const char * filename)
00041 {
00042
00043 TDirectory * dir = gDirectory;
00044 fFilename = filename;
00045
00046
00047 fOutputFile = new TFile(fFilename, "RECREATE");
00048
00049 BCModelOutput();
00050 fModel = model;
00051
00052 fModel -> MCMCInitialize();
00053
00054
00055 this -> InitializeAnalysisTree();
00056 this -> InitializeMarkovChainTrees();
00057
00058
00059 gDirectory = dir;
00060 }
00061
00062
00063
00064 BCModelOutput::~BCModelOutput()
00065 {
00066 if (fOutputFile)
00067 {
00068 fOutputFile -> Close();
00069 delete fOutputFile;
00070 }
00071
00072 if (fAnalysisTree)
00073 delete fAnalysisTree;
00074
00075 if (fMarkovChainTree)
00076 delete fMarkovChainTree;
00077
00078 if (fOutputFile)
00079 delete fOutputFile;
00080
00081 fModel = 0;
00082 }
00083
00084
00085
00086 BCModelOutput::BCModelOutput(const BCModelOutput & modeloutput)
00087 {
00088 modeloutput.Copy(* this);
00089 }
00090
00091
00092
00093 BCModelOutput & BCModelOutput::operator = (const BCModelOutput & modeloutput)
00094 {
00095 if (this != &modeloutput)
00096 modeloutput.Copy(* this);
00097
00098 return * this;
00099 }
00100
00101
00102
00103 void BCModelOutput::SetFile(const char * filename)
00104 {
00105
00106 if (fOutputFile)
00107 {
00108 fOutputFile -> Close();
00109 delete fOutputFile;
00110 }
00111
00112
00113 TDirectory * dir = gDirectory;
00114
00115
00116 fOutputFile = new TFile(fFilename, "RECREATE");
00117
00118
00119 gDirectory = dir;
00120
00121
00122 this -> InitializeAnalysisTree();
00123 }
00124
00125
00126
00127 void BCModelOutput::WriteMarkovChain(bool flag)
00128 {
00129 if (fModel)
00130 fModel -> WriteMarkovChain(flag);
00131 }
00132
00133
00134
00135 void BCModelOutput::FillAnalysisTree()
00136 {
00137
00138 fNParameters = fModel -> GetNParameters();
00139 fProbability_apriori = fModel -> GetModelAPrioriProbability();
00140 fProbability_aposteriori = fModel -> GetModelAPosterioriProbability();
00141
00142
00143 int nparameters = fModel -> GetNParameters();
00144 for (int i = 0; i < nparameters; ++i)
00145 {
00146 BCParameter * parameter = fModel -> GetParameter(i);
00147 if (fModel -> GetBestFitParameters().size() > 0)
00148 fMode_global[i] = fModel -> GetBestFitParameters().at(i);
00149 if (fModel -> GetMarginalized(parameter -> GetName().data()))
00150 {
00151 fMode_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMode();
00152 fMean_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMean();
00153 fMedian_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMedian();
00154 fQuantile_05[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.05);
00155 fQuantile_10[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.10);
00156 fQuantile_16[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.16);
00157 fQuantile_84[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.84);
00158 fQuantile_90[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.90);
00159 fQuantile_95[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.95);
00160 }
00161 }
00162
00163
00164 fAnalysisTree -> Fill();
00165
00166
00167 fIndex++;
00168 }
00169
00170
00171
00172 void BCModelOutput::WriteMarginalizedDistributions()
00173 {
00174
00175 TDirectory * dir = gDirectory;
00176
00177
00178 fOutputFile -> cd();
00179
00180 int nparameters = fModel -> GetNParameters();
00181 for (int i = 0; i < nparameters; ++i)
00182 fModel -> GetMarginalized(fModel -> GetParameter(i)) -> GetHistogram() -> Write();
00183
00184 if (nparameters > 1)
00185 for (int i = 0; i < nparameters - 1; ++i)
00186 for (int j = i + 1; j < nparameters; ++j)
00187 fModel -> GetMarginalized(fModel -> GetParameter(i),
00188 fModel -> GetParameter(j)) -> GetHistogram() -> Write();
00189
00190
00191 gDirectory = dir;
00192 }
00193
00194
00195
00196 void BCModelOutput::WriteErrorBand()
00197 {
00198
00199 TDirectory * dir = gDirectory;
00200
00201
00202 fOutputFile -> cd();
00203
00204 TH2D * h0 = fModel -> GetErrorBandXY();
00205 if (h0)
00206 {
00207 TH2D * h1 = (TH2D*)h0 -> Clone("errorbandxy");
00208 h1 -> Write();
00209
00210 double levels[] = { .68, .90, .95 };
00211 int nlevels = sizeof(levels)/sizeof(double);
00212 for (int i=0;i<nlevels;i++)
00213 {
00214 TH2D * htmp = fModel -> GetErrorBandXY_yellow(levels[i]);
00215 htmp -> SetName(TString::Format("%s_sub_%f.2",h1->GetName(),levels[i]));
00216 htmp -> Write();
00217 delete htmp;
00218 }
00219
00220 delete h1;
00221 }
00222
00223
00224 gDirectory = dir;
00225 }
00226
00227
00228
00229 void BCModelOutput::Write(TObject * o)
00230 {
00231
00232 TDirectory * dir = gDirectory;
00233
00234
00235 fOutputFile -> cd();
00236
00237 o -> Write();
00238
00239
00240 gDirectory = dir;
00241 }
00242
00243
00244
00245 void BCModelOutput::Close()
00246 {
00247
00248 TDirectory * dir = gDirectory;
00249
00250
00251 fOutputFile -> cd();
00252
00253
00254 if (fAnalysisTree -> GetEntries() > 0)
00255 fAnalysisTree -> Write();
00256
00257
00258 for (int i = 0; i < fModel -> MCMCGetNChains(); ++i)
00259 if (fMarkovChainTrees[i] -> GetEntries() > 0)
00260 fMarkovChainTrees[i] -> Write();
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271 fOutputFile -> Close();
00272
00273
00274 gDirectory = dir;
00275 }
00276
00277
00278
00279 void BCModelOutput::InitializeAnalysisTree()
00280 {
00281
00282 fAnalysisTree = new TTree("AnalysisTree", "AnalysisTree");
00283
00284
00285 fAnalysisTree -> Branch("fIndex", &fIndex, "index/I");
00286 fAnalysisTree -> Branch("fNParameters", &fNParameters, "parameters/I");
00287 fAnalysisTree -> Branch("fProbability_apriori" , &fProbability_apriori, "apriori probability/D");
00288 fAnalysisTree -> Branch("fProbability_aposteriori", &fProbability_aposteriori, "aposteriori probability/D");
00289 fAnalysisTree -> Branch("fMode_global", fMode_global, "mode (global) [parameters]/D");
00290 fAnalysisTree -> Branch("fMode_marginalized", fMode_marginalized, "mode (marginalized) [parameters]/D");
00291 fAnalysisTree -> Branch("fMean_marginalized", fMean_marginalized, "mean (marginalized)[parameters]/D");
00292 fAnalysisTree -> Branch("fMedian_marginalized", fMedian_marginalized, "median (marginalized)[parameters]/D");
00293 fAnalysisTree -> Branch("fQuantile_05" , fQuantile_05, "quantile 5% [parameters]/D");
00294 fAnalysisTree -> Branch("fQuantile_10" , fQuantile_10, "quantile 10% [parameters]/D");
00295 fAnalysisTree -> Branch("fQuantile_16" , fQuantile_16, "quantile 16% [parameters]/D");
00296 fAnalysisTree -> Branch("fQuantile_84" , fQuantile_84, "quantile 84% [parameters]/D");
00297 fAnalysisTree -> Branch("fQuantile_90" , fQuantile_90, "quantile 90% [parameters]/D");
00298 fAnalysisTree -> Branch("fQuantile_95" , fQuantile_95, "quantile 95% [parameters]/D");
00299 }
00300
00301
00302
00303 void BCModelOutput::InitializeMarkovChainTrees()
00304 {
00305
00306 fMarkovChainTrees.clear();
00307 for (int i = 0; i < fModel -> MCMCGetNChains(); ++i)
00308 {
00309 TTree * tree = new TTree(TString::Format("MarkovChainTree_%i", i), "MarkovChainTree");
00310 fMarkovChainTrees.push_back(tree);
00311 }
00312
00313
00314 fParameters = fModel -> MCMCGetP2x();
00315 fIteration = fModel -> MCMCGetP2NIterations();
00316 fLogLikelihood = fModel -> MCMCGetP2LogProbx();
00317 fNParameters = fModel -> MCMCGetNParameters();
00318
00319 for (int i = 0; i < fModel -> MCMCGetNChains(); ++i)
00320 {
00321 fMarkovChainTrees[i] -> Branch("fIteration", &(*fIteration)[i], "iteration/I");
00322 fMarkovChainTrees[i] -> Branch("fNParameters", &fNParameters, "parameters/I");
00323 fMarkovChainTrees[i] -> Branch("fLogLikelihood", &(*fLogLikelihood)[i], "log (likelihood)/D");
00324 }
00325
00326
00327 for (int i = 0; i < fModel -> MCMCGetNChains(); ++i)
00328 for (int j = 0; j < fModel -> MCMCGetNParameters(); ++j)
00329 {
00330 fMarkovChainTrees[i] -> Branch(TString::Format("fParameter%i", j),
00331 &(*fParameters)[i * fModel -> MCMCGetNParameters() + j],
00332 TString::Format("parameter %i/D", j));
00333 }
00334
00335 fModel -> MCMCSetMarkovChainTrees(fMarkovChainTrees);
00336 }
00337
00338
00339
00340 void BCModelOutput::Copy(BCModelOutput & modeloutput) const
00341 {
00342
00343 modeloutput.fModel = this -> fModel;
00344 modeloutput.fAnalysisTree = this -> fAnalysisTree;
00345 modeloutput.fMarkovChainTree = this -> fMarkovChainTree;
00346 }
00347
00348