BCModelOutput Class Reference

#include <BCModelOutput.h>

Collaboration diagram for BCModelOutput:

Collaboration graph
[legend]

List of all members.


Detailed Description

A class for creating an (ROOT) output file.

Author:
Daniel Kollar

Kevin Kröninger

Version:
1.0
Date:
08.2008 This class defines an output interface for the analysis. It creates a ROOT file which can contain summary information, histograms and Markov chains.

Definition at line 39 of file BCModelOutput.h.


Public Member Functions

Constructors and destructors


 BCModelOutput (const BCModelOutput &modeloutput)
 BCModelOutput (BCModel *model, const char *filenname)
 BCModelOutput ()
virtual ~BCModelOutput ()
Member functions (miscellaneous methods)


void Close ()
void FillAnalysisTree ()
void Write (TObject *o)
void WriteErrorBand ()
void WriteMarginalizedDistributions ()
void WriteMarkovChain (bool flag)
Getters


TTree * GetAnalysisTree ()
TFile * GetFile ()
Assignment operators


BCModelOutputoperator= (const BCModelOutput &modeloutput)
Setters


void SetFile (const char *filename)
void SetModel (BCModel *model)

Private Member Functions

void Copy (BCModelOutput &modeloutput) const
void InitializeAnalysisTree ()
void InitializeMarkovChainTrees ()

Private Attributes

TTree * fAnalysisTree
const char * fFilename
int fIndex
std::vector< int > * fIteration
std::vector< double > * fLogLikelihood
TTree * fMarkovChainTree
std::vector< TTree * > fMarkovChainTrees
double fMean_marginalized [MAXNPARAMETERS]
double fMedian_marginalized [MAXNPARAMETERS]
double fMode_global [MAXNPARAMETERS]
double fMode_marginalized [MAXNPARAMETERS]
BCModelfModel
unsigned int fNParameters
TFile * fOutputFile
std::vector< double > * fParameters
double fProbability_aposteriori
double fProbability_apriori
double fQuantile_05 [MAXNPARAMETERS]
double fQuantile_10 [MAXNPARAMETERS]
double fQuantile_16 [MAXNPARAMETERS]
double fQuantile_84 [MAXNPARAMETERS]
double fQuantile_90 [MAXNPARAMETERS]
double fQuantile_95 [MAXNPARAMETERS]

Constructor & Destructor Documentation

BCModelOutput::BCModelOutput (  ) 

The default constructor.

Definition at line 27 of file BCModelOutput.cxx.

00028 {
00029    fIndex = 0;
00030 
00031    fOutputFile = 0;
00032    fMarkovChainTree = 0;
00033    fAnalysisTree = 0;
00034    fModel = 0;
00035    fOutputFile = 0;
00036 }

BCModelOutput::BCModelOutput ( BCModel model,
const char *  filenname 
)

A constructor.

Parameters:
model The model to which this output class is assigned.
filename Name of the output file.

Definition at line 40 of file BCModelOutput.cxx.

00041 {
00042    // remeber current directory
00043    TDirectory * dir = gDirectory;
00044    fFilename = filename;
00045 
00046    // create a new file
00047    fOutputFile = new TFile(fFilename, "RECREATE");
00048 
00049    BCModelOutput();
00050    fModel = model;
00051 
00052    fModel -> MCMCInitialize();
00053 
00054    // initialize trees
00055    this -> InitializeAnalysisTree();
00056    this -> InitializeMarkovChainTrees();
00057 
00058    // change again to the old directory
00059    gDirectory = dir;
00060 }

BCModelOutput::BCModelOutput ( const BCModelOutput modeloutput  ) 

The default copy constructor.

Definition at line 86 of file BCModelOutput.cxx.

00087 {
00088    modeloutput.Copy(* this);
00089 }

BCModelOutput::~BCModelOutput (  )  [virtual]

The default destructor.

Definition at line 64 of file BCModelOutput.cxx.

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 }


Member Function Documentation

void BCModelOutput::Close (  ) 

Closes the TFile.

Definition at line 245 of file BCModelOutput.cxx.

00246 {
00247    // remember current directory
00248    TDirectory * dir = gDirectory;
00249 
00250    // change to file
00251    fOutputFile -> cd();
00252 
00253    // write analysis tree to file
00254    if (fAnalysisTree -> GetEntries() > 0)
00255       fAnalysisTree -> Write();
00256 
00257    // write markov chain tree to file
00258    for (int i = 0; i < fModel -> MCMCGetNChains(); ++i)
00259       if (fMarkovChainTrees[i] -> GetEntries() > 0)
00260          fMarkovChainTrees[i] -> Write();
00261 
00262    // write control plots to file
00263 
00264    // if (fModel -> MCMCGetH1RValue())
00265    //   fModel -> MCMCGetH1RValue() -> Write();
00266 
00267    // if (fModel -> MCMCGetH1Efficiency())
00268    //   fModel -> MCMCGetH1Efficiency() -> Write();
00269 
00270    // close file
00271    fOutputFile -> Close();
00272 
00273    // return to old directory
00274    gDirectory = dir;
00275 }

void BCModelOutput::Copy ( BCModelOutput modeloutput  )  const [private]

Definition at line 340 of file BCModelOutput.cxx.

00341 {
00342    // don't copy the content
00343    modeloutput.fModel           = this -> fModel;
00344    modeloutput.fAnalysisTree    = this -> fAnalysisTree;
00345    modeloutput.fMarkovChainTree = this -> fMarkovChainTree;
00346 }

void BCModelOutput::FillAnalysisTree (  ) 

Fill the output TTree with the current information.

Definition at line 135 of file BCModelOutput.cxx.

00136 {
00137    // get output values from model
00138    fNParameters = fModel -> GetNParameters();
00139    fProbability_apriori   = fModel -> GetModelAPrioriProbability();
00140    fProbability_aposteriori = fModel -> GetModelAPosterioriProbability();
00141 
00142    // loop over parameters
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    // fill tree
00164    fAnalysisTree -> Fill();
00165 
00166    // increase index
00167    fIndex++;
00168 }

TTree* BCModelOutput::GetAnalysisTree (  )  [inline]

Returns the output TTree tree.

Returns:
The pointer to the output TTree

Definition at line 81 of file BCModelOutput.h.

00082          { return fAnalysisTree; };

TFile* BCModelOutput::GetFile (  )  [inline]

Returns the output TFile.

Returns:
The pointer to the output TFile

Definition at line 87 of file BCModelOutput.h.

00088          { return fOutputFile; };

void BCModelOutput::InitializeAnalysisTree (  )  [private]

Initialize the output TTree.

Definition at line 279 of file BCModelOutput.cxx.

00280 {
00281    // create new tree
00282    fAnalysisTree = new TTree("AnalysisTree", "AnalysisTree");
00283 
00284    // set branch addresses
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 }

void BCModelOutput::InitializeMarkovChainTrees (  )  [private]

Initialize the Markov Chain TTree.

Definition at line 303 of file BCModelOutput.cxx.

00304 {
00305    // create new tree
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    // connect pointer to parameter vectors
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    // loop over all parameters
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 }

BCModelOutput & BCModelOutput::operator= ( const BCModelOutput modeloutput  ) 

The defaut assignment operator

Definition at line 93 of file BCModelOutput.cxx.

00094 {
00095    if (this != &modeloutput)
00096       modeloutput.Copy(* this);
00097 
00098    return * this;
00099 }

void BCModelOutput::SetFile ( const char *  filename  ) 

Sets the output filename.

Parameters:
filename The filename

Definition at line 103 of file BCModelOutput.cxx.

00104 {
00105    // delete the old file
00106    if (fOutputFile)
00107    {
00108       fOutputFile -> Close();
00109       delete fOutputFile;
00110    }
00111 
00112    // remember current directory
00113    TDirectory * dir = gDirectory;
00114 
00115    // create a new file
00116    fOutputFile = new TFile(fFilename, "RECREATE");
00117 
00118    // change back to the old directory
00119    gDirectory = dir;
00120 
00121    // initialize tree
00122    this -> InitializeAnalysisTree();
00123 }

void BCModelOutput::SetModel ( BCModel model  )  [inline]

Assign a BCModel to this output class.

Parameters:
model A pointer to the BCModel

Definition at line 98 of file BCModelOutput.h.

00099          { fModel = model; };

void BCModelOutput::Write ( TObject *  o  ) 

Writes any object derived from TObject to TFile.

Definition at line 229 of file BCModelOutput.cxx.

00230 {
00231    // remember current directory
00232    TDirectory * dir = gDirectory;
00233 
00234    // change to file
00235    fOutputFile -> cd();
00236 
00237    o -> Write();
00238 
00239    // return to old directory
00240    gDirectory = dir;
00241 }

void BCModelOutput::WriteErrorBand (  ) 

Writes the error band histogram into the TFile.

Definition at line 196 of file BCModelOutput.cxx.

00197 {
00198    // remember current directory
00199    TDirectory * dir = gDirectory;
00200 
00201    // change to file
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    // return to old directory
00224    gDirectory = dir;
00225 }

void BCModelOutput::WriteMarginalizedDistributions (  ) 

Writes the marginalized histograms to the TFile.

Definition at line 172 of file BCModelOutput.cxx.

00173 {
00174    // remember current directory
00175    TDirectory * dir = gDirectory;
00176 
00177    // change to file
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    // return to old directory
00191    gDirectory = dir;
00192 }

void BCModelOutput::WriteMarkovChain ( bool  flag  ) 

Flag for writing Markov chain to file

Parameters:
flag Writes (true) or does not write (false) the Markov chain

Definition at line 127 of file BCModelOutput.cxx.

00128 {
00129    if (fModel)
00130       fModel -> WriteMarkovChain(flag);
00131 }


Member Data Documentation

TTree* BCModelOutput::fAnalysisTree [private]

Pointer to the TTree containing the summary output information.

Definition at line 154 of file BCModelOutput.h.

const char* BCModelOutput::fFilename [private]

The output filename

Definition at line 167 of file BCModelOutput.h.

int BCModelOutput::fIndex [private]

The analysis tree variables

Definition at line 179 of file BCModelOutput.h.

std::vector<int>* BCModelOutput::fIteration [private]

Definition at line 198 of file BCModelOutput.h.

std::vector<double>* BCModelOutput::fLogLikelihood [private]

Definition at line 197 of file BCModelOutput.h.

Pointer to the TTree containing the Markov chain.

Definition at line 158 of file BCModelOutput.h.

std::vector<TTree *> BCModelOutput::fMarkovChainTrees [private]

Definition at line 163 of file BCModelOutput.h.

Definition at line 185 of file BCModelOutput.h.

Definition at line 186 of file BCModelOutput.h.

Definition at line 183 of file BCModelOutput.h.

Definition at line 184 of file BCModelOutput.h.

Pointer to the model this output class is assigned to

Definition at line 175 of file BCModelOutput.h.

unsigned int BCModelOutput::fNParameters [private]

Definition at line 180 of file BCModelOutput.h.

TFile* BCModelOutput::fOutputFile [private]

Pointer to the output TFile.

Definition at line 171 of file BCModelOutput.h.

std::vector<double>* BCModelOutput::fParameters [private]

The markov chain tree variables

Definition at line 196 of file BCModelOutput.h.

Definition at line 182 of file BCModelOutput.h.

Definition at line 181 of file BCModelOutput.h.

Definition at line 187 of file BCModelOutput.h.

Definition at line 188 of file BCModelOutput.h.

Definition at line 189 of file BCModelOutput.h.

Definition at line 190 of file BCModelOutput.h.

Definition at line 191 of file BCModelOutput.h.

Definition at line 192 of file BCModelOutput.h.


Generated on Fri Jan 16 10:24:11 2009 for BAT - Bayesian Analysis Toolkit by  doxygen 1.5.6