A class for creating an (ROOT) output file. More...
#include <BCModelOutput.h>
A class for creating an (ROOT) output file.
Definition at line 39 of file BCModelOutput.h.
BCModelOutput::BCModelOutput | ( | ) |
BCModelOutput::BCModelOutput | ( | BCModel * | model, | |
const char * | filenname | |||
) |
A constructor.
model | The model to which this output class is assigned. | |
filename | Name of the output file. |
Definition at line 35 of file BCModelOutput.cxx.
BCModelOutput::BCModelOutput | ( | const BCModelOutput & | modeloutput | ) |
The default copy constructor.
Definition at line 52 of file BCModelOutput.cxx.
{ modeloutput.Copy(* this); }
BCModelOutput::~BCModelOutput | ( | ) | [virtual] |
The default destructor.
Definition at line 43 of file BCModelOutput.cxx.
{ if (fOutputFile) { fOutputFile -> Close(); delete fOutputFile; } }
void BCModelOutput::Close | ( | ) |
Closes the TFile.
Definition at line 254 of file BCModelOutput.cxx.
{ // remember current directory TDirectory * dir = gDirectory; // change to file fOutputFile -> cd(); // write analysis tree to file if (fAnalysisTree -> GetEntries() > 0) fAnalysisTree -> Write(); // write markov chain tree to file for (int i = 0; i < fModel->MCMCGetNChains(); ++i) if (fModel->MCMCGetMarkovChainTree(i)->GetEntries() > 0) fModel->MCMCGetMarkovChainTree(i)->Write(); // write SA tree to file if (fModel->GetSATree()->GetEntries() > 0) fModel->GetSATree()->Write(); // close file fOutputFile -> Close(); // return to old directory gDirectory = dir; }
void BCModelOutput::Copy | ( | BCModelOutput & | modeloutput | ) | const [private] |
Definition at line 306 of file BCModelOutput.cxx.
{ // don't copy the content modeloutput.fModel = this -> fModel; modeloutput.fAnalysisTree = this -> fAnalysisTree; // modeloutput.fMarkovChainTrees = this -> fMarkovChainTrees; }
void BCModelOutput::FillAnalysisTree | ( | ) |
Fill the output TTree with the current information.
Definition at line 124 of file BCModelOutput.cxx.
{ if(!fOutputFile) { BCLog::OutError("BCModelOutput::FillAnalysisTree : No file to write to."); return; } // get output values from model fNParameters = fModel -> GetNParameters(); fProbability_apriori = fModel -> GetModelAPrioriProbability(); fProbability_aposteriori = fModel -> GetModelAPosterioriProbability(); // loop over parameters int nparameters = fModel -> GetNParameters(); for (int i = 0; i < nparameters; ++i) { BCParameter * parameter = fModel -> GetParameter(i); if (fModel -> GetBestFitParameters().size() > 0) fMode_global[i] = fModel -> GetBestFitParameters().at(i); if (fModel -> GetMarginalized(parameter -> GetName().data())) { fMode_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMode(); fMean_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMean(); fMedian_marginalized[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetMedian(); fQuantile_05[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.05); fQuantile_10[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.10); fQuantile_16[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.16); fQuantile_84[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.84); fQuantile_90[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.90); fQuantile_95[i] = fModel -> GetMarginalized(parameter -> GetName().data()) -> GetQuantile(0.95); } } // fill tree fAnalysisTree -> Fill(); // increase index fIndex++; }
TTree* BCModelOutput::GetAnalysisTree | ( | ) | [inline] |
Returns the output TTree tree.
Definition at line 81 of file BCModelOutput.h.
{ return fAnalysisTree; };
TFile* BCModelOutput::GetFile | ( | ) | [inline] |
Returns the output TFile.
Definition at line 87 of file BCModelOutput.h.
{ return fOutputFile; };
void BCModelOutput::Init | ( | ) | [private] |
Initialize the variables
Definition at line 67 of file BCModelOutput.cxx.
{ fIndex = 0; fOutputFile = 0; fAnalysisTree = 0; fTreeSA = 0; fModel = 0; fFileName = 0; }
void BCModelOutput::InitializeAnalysisTree | ( | ) | [private] |
Initialize the output TTree.
Definition at line 283 of file BCModelOutput.cxx.
{ // create new tree fAnalysisTree = new TTree("AnalysisTree", "AnalysisTree"); // set branch addresses fAnalysisTree -> Branch("fIndex", &fIndex, "index/I"); fAnalysisTree -> Branch("fNParameters", &fNParameters, "parameters/I"); fAnalysisTree -> Branch("fProbability_apriori" , &fProbability_apriori, "apriori probability/D"); fAnalysisTree -> Branch("fProbability_aposteriori", &fProbability_aposteriori, "aposteriori probability/D"); fAnalysisTree -> Branch("fMode_global", fMode_global, "mode (global) [parameters]/D"); fAnalysisTree -> Branch("fMode_marginalized", fMode_marginalized, "mode (marginalized) [parameters]/D"); fAnalysisTree -> Branch("fMean_marginalized", fMean_marginalized, "mean (marginalized)[parameters]/D"); fAnalysisTree -> Branch("fMedian_marginalized", fMedian_marginalized, "median (marginalized)[parameters]/D"); fAnalysisTree -> Branch("fQuantile_05" , fQuantile_05, "quantile 5% [parameters]/D"); fAnalysisTree -> Branch("fQuantile_10" , fQuantile_10, "quantile 10% [parameters]/D"); fAnalysisTree -> Branch("fQuantile_16" , fQuantile_16, "quantile 16% [parameters]/D"); fAnalysisTree -> Branch("fQuantile_84" , fQuantile_84, "quantile 84% [parameters]/D"); fAnalysisTree -> Branch("fQuantile_90" , fQuantile_90, "quantile 90% [parameters]/D"); fAnalysisTree -> Branch("fQuantile_95" , fQuantile_95, "quantile 95% [parameters]/D"); }
void BCModelOutput::InitializeSATree | ( | ) | [private] |
Initialize SA TTree.
BCModelOutput & BCModelOutput::operator= | ( | const BCModelOutput & | modeloutput | ) |
The defaut assignment operator
Definition at line 58 of file BCModelOutput.cxx.
{ if (this != &modeloutput) modeloutput.Copy(* this); return * this; }
void BCModelOutput::SetFile | ( | const char * | filename | ) |
Sets the output filename.
filename | The filename |
Definition at line 86 of file BCModelOutput.cxx.
{ if(!fModel) { BCLog::OutError("BCModelOutput::SetFile : Cannot set file if model is not set."); return; } // delete the old file if (fOutputFile) { fOutputFile->Close(); delete fOutputFile; } // remember current directory TDirectory * dir = gDirectory; // create a new file fFileName = const_cast<char *>(filename); fOutputFile = new TFile(fFileName, "RECREATE"); // initialize trees InitializeAnalysisTree(); fModel->MCMCInitializeMarkovChainTrees(); fModel->InitializeSATree(); // change back to the old directory gDirectory = dir; }
void BCModelOutput::SetModel | ( | BCModel * | model | ) |
Assign a BCModel to this output class.
model | A pointer to the BCModel |
Definition at line 78 of file BCModelOutput.cxx.
{ fModel = model; fModel->MCMCInitialize(); fModel->SAInitialize(); }
void BCModelOutput::Write | ( | TObject * | o | ) |
Writes any object derived from TObject to TFile.
Definition at line 233 of file BCModelOutput.cxx.
{ if(!fOutputFile) { BCLog::OutError("BCModelOutput::Write : No file to write to."); return; } // remember current directory TDirectory * dir = gDirectory; // change to file fOutputFile -> cd(); o -> Write(); // return to old directory gDirectory = dir; }
void BCModelOutput::WriteErrorBand | ( | ) |
Writes the error band histogram into the TFile.
Definition at line 195 of file BCModelOutput.cxx.
{ if(!fOutputFile) { BCLog::OutError("BCModelOutput::WriteErrorBand : No file to write to."); return; } // remember current directory TDirectory * dir = gDirectory; // change to file fOutputFile -> cd(); TH2D * h0 = fModel -> GetErrorBandXY(); if (h0) { TH2D * h1 = (TH2D*)h0 -> Clone("errorbandxy"); h1 -> Write(); double levels[] = { .68, .90, .95 }; int nlevels = sizeof(levels)/sizeof(double); for (int i=0;i<nlevels;i++) { TH2D * htmp = fModel -> GetErrorBandXY_yellow(levels[i]); htmp -> SetName(TString::Format("%s_sub_%f.2",h1->GetName(),levels[i])); htmp -> Write(); delete htmp; } delete h1; } // return to old directory gDirectory = dir; }
void BCModelOutput::WriteMarginalizedDistributions | ( | ) |
Writes the marginalized histograms to the TFile.
Definition at line 166 of file BCModelOutput.cxx.
{ if(!fOutputFile) { BCLog::OutError("BCModelOutput::WriteMarginalizedDistributions : No file to write to."); return; } // remember current directory TDirectory * dir = gDirectory; // change to file fOutputFile -> cd(); int nparameters = fModel -> GetNParameters(); for (int i = 0; i < nparameters; ++i) fModel -> GetMarginalized(fModel -> GetParameter(i)) -> GetHistogram() -> Write(); if (nparameters > 1) for (int i = 0; i < nparameters - 1; ++i) for (int j = i + 1; j < nparameters; ++j) fModel -> GetMarginalized(fModel -> GetParameter(i), fModel -> GetParameter(j)) -> GetHistogram() -> Write(); // return to old directory gDirectory = dir; }
void BCModelOutput::WriteMarkovChain | ( | bool | flag = true |
) |
Flag for writing Markov chain to file
flag | Writes (true) or does not write (false) the Markov chain |
Definition at line 117 of file BCModelOutput.cxx.
{ if (fModel) fModel -> WriteMarkovChain(flag); }
TTree* BCModelOutput::fAnalysisTree [private] |
Pointer to the TTree containing the summary output information.
Definition at line 158 of file BCModelOutput.h.
char* BCModelOutput::fFileName [private] |
The output filename
Definition at line 171 of file BCModelOutput.h.
int BCModelOutput::fIndex [private] |
The analysis tree variables
Definition at line 183 of file BCModelOutput.h.
std::vector<TTree *> BCModelOutput::fMarkovChainTrees [private] |
Definition at line 163 of file BCModelOutput.h.
double BCModelOutput::fMean_marginalized[MAXNPARAMETERS] [private] |
Definition at line 189 of file BCModelOutput.h.
double BCModelOutput::fMedian_marginalized[MAXNPARAMETERS] [private] |
Definition at line 190 of file BCModelOutput.h.
double BCModelOutput::fMode_global[MAXNPARAMETERS] [private] |
Definition at line 187 of file BCModelOutput.h.
double BCModelOutput::fMode_marginalized[MAXNPARAMETERS] [private] |
Definition at line 188 of file BCModelOutput.h.
BCModel* BCModelOutput::fModel [private] |
Pointer to the model this output class is assigned to
Definition at line 179 of file BCModelOutput.h.
unsigned int BCModelOutput::fNParameters [private] |
Definition at line 184 of file BCModelOutput.h.
TFile* BCModelOutput::fOutputFile [private] |
Pointer to the output TFile.
Definition at line 175 of file BCModelOutput.h.
double BCModelOutput::fProbability_aposteriori [private] |
Definition at line 186 of file BCModelOutput.h.
double BCModelOutput::fProbability_apriori [private] |
Definition at line 185 of file BCModelOutput.h.
double BCModelOutput::fQuantile_05[MAXNPARAMETERS] [private] |
Definition at line 191 of file BCModelOutput.h.
double BCModelOutput::fQuantile_10[MAXNPARAMETERS] [private] |
Definition at line 192 of file BCModelOutput.h.
double BCModelOutput::fQuantile_16[MAXNPARAMETERS] [private] |
Definition at line 193 of file BCModelOutput.h.
double BCModelOutput::fQuantile_84[MAXNPARAMETERS] [private] |
Definition at line 194 of file BCModelOutput.h.
double BCModelOutput::fQuantile_90[MAXNPARAMETERS] [private] |
Definition at line 195 of file BCModelOutput.h.
double BCModelOutput::fQuantile_95[MAXNPARAMETERS] [private] |
Definition at line 196 of file BCModelOutput.h.
TTree* BCModelOutput::fTreeSA [private] |
Definition at line 167 of file BCModelOutput.h.