• Main Page
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

BCDataSet.cxx

Go to the documentation of this file.
00001 /*
00002  * Copyright (C) 2008-2010, Daniel Kollar and Kevin Kroeninger.
00003  * All rights reserved.
00004  *
00005  * For the licensing terms see doc/COPYING.
00006  */
00007 
00008 // ---------------------------------------------------------
00009 
00010 #include "BAT/BCDataSet.h"
00011 
00012 #include "BAT/BCDataPoint.h"
00013 #include "BAT/BCLog.h"
00014 #include "BAT/BCErrorCodes.h"
00015 
00016 #include <TFile.h>
00017 #include <TTree.h>
00018 #include <TString.h>
00019 
00020 #include <iostream>
00021 #include <fstream>
00022 
00023 // ---------------------------------------------------------
00024 
00025 BCDataSet::BCDataSet()
00026 {
00027    fBCDataVector = 0;
00028 }
00029 
00030 // ---------------------------------------------------------
00031 
00032 BCDataSet::~BCDataSet()
00033 {
00034    if (fBCDataVector)
00035       delete fBCDataVector;
00036 }
00037 
00038 // ---------------------------------------------------------
00039 
00040 unsigned int BCDataSet::GetNDataPoints()
00041 {
00042    // check if vector exists. Return number of data points if true ...
00043    if (fBCDataVector)
00044       return fBCDataVector -> size();
00045 
00046    // ... or give out warning and return 0 if not.
00047    BCLog::Out(BCLog::warning, BCLog::warning,"BCDataSet::GetNDataPoints : DataSet not yet created.");
00048    return 0;
00049 }
00050 
00051 // ---------------------------------------------------------
00052 
00053 unsigned int BCDataSet::GetNValuesPerPoint()
00054 {
00055    // check if vector exists and contains datapoints
00056    if (fBCDataVector && fBCDataVector -> size() > 0)
00057       return this -> GetDataPoint(0) -> GetNValues();
00058 
00059    BCLog::Out(BCLog::error, BCLog::error,
00060          "BCDataSet::GetNValuesPerPoint : Data set doesn't exist yet");
00061    return 0;
00062 }
00063 
00064 // ---------------------------------------------------------
00065 
00066 BCDataPoint * BCDataSet::GetDataPoint(unsigned int index)
00067 {
00068    if (!fBCDataVector || this -> GetNDataPoints()==0 )
00069    {
00070       BCLog::Out(BCLog::error, BCLog::error,"BCDataSet::GetDataPoint : Dataset is empty.");
00071       return 0;
00072    }
00073 
00074    // check if index is within range. Return the data point if true ...
00075    if(index >= 0 && index < this -> GetNDataPoints())
00076       return fBCDataVector -> at(index);
00077 
00078    // ... or give out warning and return 0 if not.
00079    BCLog::Out(BCLog::error, BCLog::error,"BCDataSet::GetDataPoint : Index out of range. Return 0.");
00080    return 0;
00081 }
00082 
00083 // ---------------------------------------------------------
00084 std::vector<double> BCDataSet::GetDataComponents( int index)
00085 {
00086    unsigned int N = this->GetNDataPoints();
00087    std::vector<double> components( N , 0.0 );
00088 
00089    BCDataPoint* point=0;
00090    for (unsigned int i = 0; i < N; ++i) {
00091       //rely on index checking in Get... methods
00092       point = this -> GetDataPoint(i);
00093       components[i] = point -> GetValue(index);
00094    }
00095 
00096    return components;
00097 }
00098 
00099 
00100 
00101 // ---------------------------------------------------------
00102 
00103 int BCDataSet::ReadDataFromFileTree(const char * filename, const char * treename, const char * branchnames)
00104 {
00105    // open root file
00106    TFile * file = new TFile(filename, "READ");
00107 
00108    // check if file is open and warn if not.
00109    if (!file -> IsOpen())
00110    {
00111       BCLog::Out(BCLog::error, BCLog::error,
00112             Form("BCDataSet::ReadDataFromFileTree : Could not open file %s.", filename));
00113       return ERROR_FILENOTFOUND;
00114    }
00115 
00116    // get tree
00117    TTree * tree = (TTree*) file -> Get(treename);
00118 
00119    // check if tree is there and warn if not.
00120    if (!tree)
00121    {
00122       BCLog::Out(BCLog::error, BCLog::error,
00123             Form("BCDataSet::ReadDataFromFileTree : Could not find TTree %s.", treename));
00124 
00125       // close file
00126       file -> Close();
00127 
00128       return ERROR_TREENOTFOUND;
00129    }
00130 
00131    // if data set contains data, clear data object container ...
00132    if (fBCDataVector != 0)
00133    {
00134       fBCDataVector -> clear();
00135 
00136       BCLog::Out(BCLog::detail, BCLog::detail,"BCDataSet::ReadDataFromFileTree : Overwrite existing data.");
00137    }
00138 
00139    // ... or allocate memory for the vector if not.
00140    else
00141       fBCDataVector = new BCDataVector();
00142 
00143    // get branch names.
00144 
00145    // first, copy the branchnames into a std::string.
00146    std::string branches(branchnames);
00147 
00148    // define a vector of std::strings which contain the tree names.
00149    std::vector<std::string> * branchnamevector = new std::vector<std::string>;
00150 
00151    // the names are supposed to be separated by commas. find first comma
00152    // entry in the string.
00153    int temp_index = branches.find_first_of(",");
00154 
00155    // reset number of branches
00156    int nbranches = 0;
00157 
00158    // repeat until the is nothing left in the string.
00159    while(branches.size() > 0)
00160    {
00161       // temporary string which contains the name of the current branch
00162       std::string branchname;
00163 
00164       // get current branch name
00165 
00166       // if there is no comma the current branchname corresponds to the whole string, ...
00167       if (temp_index == -1)
00168          branchname = branches;
00169 
00170       // ... if there is a comma, copy that part of the string into the current branchname.
00171       else
00172          branchname.assign(branches, 0, temp_index);
00173 
00174       // write branch name to a vector
00175       branchnamevector -> push_back(branchname);
00176 
00177       // increase the number of branches found
00178       nbranches++;
00179 
00180       // cut remaining string with branchnames
00181 
00182       // if there is no comma left empty the string, ...
00183       if (temp_index == -1)
00184             branches = "";
00185 
00186       // ... if there is a comma remove the current branchname from the string.
00187       else
00188          branches.erase(0, temp_index + 1);
00189 
00190       // find the next comma
00191       temp_index = branches.find_first_of(",");
00192    }
00193 
00194    // create temporary vector with data and assign some zeros.
00195    std::vector<double> data;
00196    data.assign(nbranches, 0.0);
00197 
00198    // set the branch address.
00199    for (int i = 0; i < nbranches; i++)
00200       tree -> SetBranchAddress(branchnamevector -> at(i).data(), &data.at(i));
00201 
00202    // calculate maximum number of entries
00203    int nentries = tree -> GetEntries();
00204 
00205    // check if there are any events in the tree and close file if not.
00206    if (nentries <= 0)
00207    {
00208       BCLog::Out(BCLog::error, BCLog::error,
00209             Form("BCDataSet::ReadDataFromFileTree : No events in TTree %s.", treename));
00210 
00211       // close file
00212       file -> Close();
00213 
00214       return ERROR_NOEVENTS;
00215    }
00216 
00217    // loop over entries
00218    for (int ientry = 0; ientry < nentries; ientry++)
00219    {
00220       // get entry
00221       tree -> GetEntry(ientry);
00222 
00223       // create data object
00224       BCDataPoint * datapoint = new BCDataPoint(nbranches);
00225 
00226       // copy data
00227 
00228       for (int i = 0; i < nbranches; i++)
00229          datapoint -> SetValue(i, data.at(i));
00230 
00231       // add data point to this data set.
00232       this -> AddDataPoint(datapoint);
00233    }
00234 
00235    // close file
00236    file -> Close();
00237 
00238    // remove file pointer.
00239    if (file)
00240       delete file;
00241 
00242    return 0;
00243 
00244 }
00245 
00246 // ---------------------------------------------------------
00247 
00248 int BCDataSet::ReadDataFromFileTxt(const char * filename, int nbranches)
00249 {
00250    // open text file.
00251    std::fstream file;
00252    file.open(filename, std::fstream::in);
00253 
00254    // check if file is open and warn if not.
00255    if (!file.is_open())
00256    {
00257       BCLog::Out(BCLog::error, BCLog::error,
00258             Form("BCDataSet::ReadDataFromFileText : Could not open file %s.", filename));
00259 
00260       return ERROR_FILENOTFOUND;
00261    }
00262 
00263    // if data set contains data, clear data object container ...
00264    if (fBCDataVector != 0)
00265    {
00266       fBCDataVector -> clear();
00267 
00268       BCLog::Out(BCLog::detail, BCLog::detail,"BCDataSet::ReadDataFromFileTxt : Overwrite existing data.");
00269    }
00270 
00271    // ... or allocate memory for the vector if not.
00272    else
00273       fBCDataVector = new BCDataVector();
00274 
00275    // create temporary vector with data and assign some zeros.
00276    std::vector<double> data;
00277    data.assign(nbranches, 0.0);
00278 
00279    // reset counter
00280    int nentries = 0;
00281 
00282    // read data and create data points.
00283    while (!file.eof())
00284    {
00285       // read data from file
00286       int i=0;
00287       while(file >> data[i])
00288       {
00289          if (i==nbranches-1)
00290             break;
00291          i++;
00292       }
00293 
00294       // create data point.
00295       if(i == nbranches-1)
00296       {
00297          BCDataPoint * datapoint = new BCDataPoint(nbranches);
00298 
00299          // copy data into data point
00300          for (int i = 0; i < nbranches; i++)
00301             datapoint -> SetValue(i, data.at(i));
00302 
00303          // add data point to this data set.
00304          this -> AddDataPoint(datapoint);
00305 
00306          // increase counter
00307          nentries++;
00308       }
00309    }
00310 
00311    // check if there are any events in the tree and close file if not.
00312    if (nentries <= 0)
00313    {
00314       BCLog::Out(BCLog::error, BCLog::error,
00315             Form("BCDataSet::ReadDataFromFileText : No events in the file %s.", filename));
00316 
00317       // close file
00318       file.close();
00319 
00320       return ERROR_NOEVENTS;
00321    }
00322 
00323    // close file
00324    file.close();
00325 
00326    return 0;
00327 
00328 }
00329 
00330 // ---------------------------------------------------------
00331 
00332 void BCDataSet::AddDataPoint(BCDataPoint * datapoint)
00333 {
00334 
00335    // check if memory for the vector has been allocated and
00336    // allocate if not.
00337    if (fBCDataVector == 0)
00338       fBCDataVector = new BCDataVector();
00339 
00340    // add data point to the data set.
00341    fBCDataVector -> push_back(datapoint);
00342 
00343 }
00344 
00345 // ---------------------------------------------------------
00346 
00347 void BCDataSet::Reset()
00348 {
00349 
00350    // if memory has been allocated to the data set
00351    // clear the content.
00352    if (fBCDataVector != 0)
00353       fBCDataVector -> clear();
00354 
00355 }
00356 
00357 // ---------------------------------------------------------
00358 
00359 void BCDataSet::Dump()
00360 {
00361    if (!fBCDataVector)
00362    {
00363       BCLog::Out(BCLog::error, BCLog::error, "BCDataSet::Dump : Data set is empty. Nothing to dump.");
00364       return;
00365    }
00366 
00367    std::cout << std::endl
00368       << "Dumping dataset:" << std::endl
00369       << "----------------" << std::endl
00370       << " - number of points:            " << fBCDataVector -> size() << std::endl
00371       << " - number of values per point:  " << this -> GetDataPoint(0) -> GetNValues() << std::endl
00372       << " - values:" << std::endl;
00373    unsigned int n = this -> GetDataPoint(0) -> GetNValues();
00374    for (unsigned int i=0; i< fBCDataVector -> size(); i++)
00375    {
00376       std::cout << Form("%5d :  ", i);
00377       for (unsigned int j=0; j<n; j++)
00378          std::cout << Form("%12.5g", this -> GetDataPoint(i) -> GetValue(j));
00379       std::cout << std::endl;
00380    }
00381    std::cout << std::endl;
00382 
00383 }
00384 
00385 
00386 // ---------------------------------------------------------

Generated on Mon Aug 30 2010 22:14:54 for Bayesian Analysis Toolkit by  doxygen 1.7.1