Go to the documentation of this file.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 #include "BAT/BCLog.h"
00017
00018 #include <TDirectory.h>
00019 #include <TFile.h>
00020 #include <TTree.h>
00021 #include <TObject.h>
00022 #include <TH2D.h>
00023 #include <TH1D.h>
00024 #include <TString.h>
00025
00026 #include <iostream>
00027
00028
00029 BCModelOutput::BCModelOutput()
00030 {
00031 Init();
00032 }
00033
00034
00035 BCModelOutput::BCModelOutput(BCModel * model, const char * filename)
00036 {
00037 Init();
00038 SetModel(model);
00039 SetFile(filename);
00040 }
00041
00042
00043 BCModelOutput::~BCModelOutput()
00044 {
00045 if (fOutputFile) {
00046 fOutputFile -> Close();
00047 delete fOutputFile;
00048 }
00049 }
00050
00051
00052 BCModelOutput::BCModelOutput(const BCModelOutput & modeloutput)
00053 {
00054 modeloutput.Copy(* this);
00055 }
00056
00057
00058 BCModelOutput & BCModelOutput::operator = (const BCModelOutput & modeloutput)
00059 {
00060 if (this != &modeloutput)
00061 modeloutput.Copy(* this);
00062
00063 return * this;
00064 }
00065
00066
00067 void BCModelOutput::Init()
00068 {
00069 fIndex = 0;
00070 fOutputFile = 0;
00071 fAnalysisTree = 0;
00072 fTreeSA = 0;
00073 fModel = 0;
00074 fFileName = 0;
00075 }
00076
00077
00078 void BCModelOutput::SetModel(BCModel * model)
00079 {
00080 fModel = model;
00081 fModel->MCMCInitialize();
00082 fModel->SAInitialize();
00083 }
00084
00085
00086 void BCModelOutput::SetFile(const char * filename)
00087 {
00088 if(!fModel) {
00089 BCLog::OutError("BCModelOutput::SetFile : Cannot set file if model is not set.");
00090 return;
00091 }
00092
00093
00094 if (fOutputFile)
00095 {
00096 fOutputFile->Close();
00097 delete fOutputFile;
00098 }
00099
00100
00101 TDirectory * dir = gDirectory;
00102
00103
00104 fFileName = const_cast<char *>(filename);
00105 fOutputFile = new TFile(fFileName, "RECREATE");
00106
00107
00108 InitializeAnalysisTree();
00109 fModel->MCMCInitializeMarkovChainTrees();
00110 fModel->InitializeSATree();
00111
00112
00113 gDirectory = dir;
00114 }
00115
00116
00117 void BCModelOutput::WriteMarkovChain(bool flag)
00118 {
00119 if (fModel)
00120 fModel -> WriteMarkovChain(flag);
00121 }
00122
00123
00124 void BCModelOutput::FillAnalysisTree()
00125 {
00126 if(!fOutputFile)
00127 {
00128 BCLog::OutError("BCModelOutput::FillAnalysisTree : No file to write to.");
00129 return;
00130 }
00131
00132
00133 fNParameters = fModel -> GetNParameters();
00134 fProbability_apriori = fModel -> GetModelAPrioriProbability();
00135 fProbability_aposteriori = fModel -> GetModelAPosterioriProbability();
00136
00137
00138 int nparameters = fModel -> GetNParameters();
00139 for (int i = 0; i < nparameters; ++i)
00140 {
00141 BCParameter * parameter = fModel -> GetParameter(i);
00142 if (fModel -> GetBestFitParameters().size() > 0)
00143 fMode_global[i] = fModel -> GetBestFitParameters().at(i);
00144 if (fModel -> GetMarginalized(parameter -> GetName().data()))
00145 {
00146 fMode_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMode();
00147 fMean_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMean();
00148 fMedian_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMedian();
00149 fQuantile_05[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.05);
00150 fQuantile_10[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.10);
00151 fQuantile_16[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.16);
00152 fQuantile_84[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.84);
00153 fQuantile_90[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.90);
00154 fQuantile_95[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.95);
00155 }
00156 }
00157
00158
00159 fAnalysisTree -> Fill();
00160
00161
00162 fIndex++;
00163 }
00164
00165
00166 void BCModelOutput::WriteMarginalizedDistributions()
00167 {
00168 if(!fOutputFile)
00169 {
00170 BCLog::OutError("BCModelOutput::WriteMarginalizedDistributions : No file to write to.");
00171 return;
00172 }
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 void BCModelOutput::WriteErrorBand()
00196 {
00197 if(!fOutputFile)
00198 {
00199 BCLog::OutError("BCModelOutput::WriteErrorBand : No file to write to.");
00200 return;
00201 }
00202
00203
00204 TDirectory * dir = gDirectory;
00205
00206
00207 fOutputFile -> cd();
00208
00209 TH2D * h0 = fModel -> GetErrorBandXY();
00210 if (h0)
00211 {
00212 TH2D * h1 = (TH2D*)h0 -> Clone("errorbandxy");
00213 h1 -> Write();
00214
00215 double levels[] = { .68, .90, .95 };
00216 int nlevels = sizeof(levels)/sizeof(double);
00217 for (int i=0;i<nlevels;i++)
00218 {
00219 TH2D * htmp = fModel -> GetErrorBandXY_yellow(levels[i]);
00220 htmp -> SetName(TString::Format("%s_sub_%f.2",h1->GetName(),levels[i]));
00221 htmp -> Write();
00222 delete htmp;
00223 }
00224
00225 delete h1;
00226 }
00227
00228
00229 gDirectory = dir;
00230 }
00231
00232
00233 void BCModelOutput::Write(TObject * o)
00234 {
00235 if(!fOutputFile)
00236 {
00237 BCLog::OutError("BCModelOutput::Write : No file to write to.");
00238 return;
00239 }
00240
00241
00242 TDirectory * dir = gDirectory;
00243
00244
00245 fOutputFile -> cd();
00246
00247 o -> Write();
00248
00249
00250 gDirectory = dir;
00251 }
00252
00253
00254 void BCModelOutput::Close()
00255 {
00256
00257 TDirectory * dir = gDirectory;
00258
00259
00260 fOutputFile -> cd();
00261
00262
00263 if (fAnalysisTree -> GetEntries() > 0)
00264 fAnalysisTree -> Write();
00265
00266
00267 for (int i = 0; i < fModel->MCMCGetNChains(); ++i)
00268 if (fModel->MCMCGetMarkovChainTree(i)->GetEntries() > 0)
00269 fModel->MCMCGetMarkovChainTree(i)->Write();
00270
00271
00272 if (fModel->GetSATree()->GetEntries() > 0)
00273 fModel->GetSATree()->Write();
00274
00275
00276 fOutputFile -> Close();
00277
00278
00279 gDirectory = dir;
00280 }
00281
00282
00283 void BCModelOutput::InitializeAnalysisTree()
00284 {
00285
00286 fAnalysisTree = new TTree("AnalysisTree", "AnalysisTree");
00287
00288
00289 fAnalysisTree -> Branch("fIndex", &fIndex, "index/I");
00290 fAnalysisTree -> Branch("fNParameters", &fNParameters, "parameters/I");
00291 fAnalysisTree -> Branch("fProbability_apriori" , &fProbability_apriori, "apriori probability/D");
00292 fAnalysisTree -> Branch("fProbability_aposteriori", &fProbability_aposteriori, "aposteriori probability/D");
00293 fAnalysisTree -> Branch("fMode_global", fMode_global, "mode (global) [parameters]/D");
00294 fAnalysisTree -> Branch("fMode_marginalized", fMode_marginalized, "mode (marginalized) [parameters]/D");
00295 fAnalysisTree -> Branch("fMean_marginalized", fMean_marginalized, "mean (marginalized)[parameters]/D");
00296 fAnalysisTree -> Branch("fMedian_marginalized", fMedian_marginalized, "median (marginalized)[parameters]/D");
00297 fAnalysisTree -> Branch("fQuantile_05" , fQuantile_05, "quantile 5% [parameters]/D");
00298 fAnalysisTree -> Branch("fQuantile_10" , fQuantile_10, "quantile 10% [parameters]/D");
00299 fAnalysisTree -> Branch("fQuantile_16" , fQuantile_16, "quantile 16% [parameters]/D");
00300 fAnalysisTree -> Branch("fQuantile_84" , fQuantile_84, "quantile 84% [parameters]/D");
00301 fAnalysisTree -> Branch("fQuantile_90" , fQuantile_90, "quantile 90% [parameters]/D");
00302 fAnalysisTree -> Branch("fQuantile_95" , fQuantile_95, "quantile 95% [parameters]/D");
00303 }
00304
00305
00306 void BCModelOutput::Copy(BCModelOutput & modeloutput) const
00307 {
00308
00309 modeloutput.fModel = this -> fModel;
00310 modeloutput.fAnalysisTree = this -> fAnalysisTree;
00311
00312 }
00313
00314