#ifndef _sdet_T2LifeROOTFileManager_h_ #define _sdet_T2LifeROOTFileManager_h_ /** \file \author Ralf Ulrich \version $Id: T2LifeROOTFileManager.h 19041 2011-04-15 18:49:23Z javierg $ \date Wed Nov 22 00:58:48 CET 2006 */ static const char CVSId_sdet_T2LifeROOTFileManager[] = "$Id: T2LifeROOTFileManager.h 19041 2011-04-15 18:49:23Z javierg $"; #include class TTree; namespace sdet { /** \brief Class to manage Sd life time ROOT files Identical functionality as T2LifeFileManager, but much more user-friendly and better performance. The input data file is a single relatively small ROOT file containing a TTree, which contains an identical image of the corresponding ASCII T2-information. The needed ROOT data file can be obtained from http://augerobserver.fzk.de/fileadmin/user_upload/AugerData/UpTime/v1r0_ICRC07/T2Life.root (~64MB)
A small and simple program, which can be used for the conversion from T2 ASCII data to ROOT data looks like: \code // ************************************************** // // T2toRoot.cc // // converts T2 ASCII files to root tree // // v1.0 04/02/06 Michael.Unger@ik.fzk.de // // ************************************************** #include #include #include #include #include #include #include using namespace std; int main(int argc, char *argv[] ) { if(argc<2) { cout << "usage: T2toRoot " << endl; return 1; } int j=1; string line; ifstream * inFile; const int kMaxStation=1600; TBits * stationArray = new TBits(kMaxStation); unsigned int gpsStart; unsigned int gpsStop; TFile rootFile("/userdata/Fabian/AugerUpTime/NewFiles/T2Life.root","RECREATE"); TTree * theTree = new TTree("T2Life","T2 lifetimes"); theTree->Branch("gpsStart",&gpsStart,"gpsStart/i"); theTree->Branch("gpsStop",&gpsStop,"gpsStop/i"); theTree->Branch("stationArray","TBits",&stationArray,16000,99); while ( argv[j] != NULL ) { cout << " ---> " << argv[j] << endl; inFile = new ifstream(argv[j]); while (getline(*inFile, line)) { // find next valid line stationArray->ResetAllBits(); istringstream buffer; buffer.str(line); buffer >> gpsStart >> gpsStop; int iStation=100; int flag; while ( buffer >> flag) { if(flag) stationArray->SetBitNumber(iStation,true); iStation++; } theTree->Fill(); } delete inFile; j++; } rootFile.cd(); rootFile.Write(); rootFile.Close(); } \endcode \author Ralf Ulrich \date Wed Nov 22 00:59:29 CET 2006 \ingroup T2LifeROOTFileManager \ingroup managers */ class T2LifeROOTFileManager : public det::VManager { public: T2LifeROOTFileManager() {} virtual ~T2LifeROOTFileManager() { } VMANAGER_GETDATA_CALL(GetOkFlag, int) VMANAGER_GETDATA_NOTFOUND(double) VMANAGER_GETDATA_NOTFOUND(std::string) VMANAGER_GETDATA_NOTFOUND(std::list) VMANAGER_GETDATA_NOTFOUND(std::list) VMANAGER_GETDATA_NOTFOUND(std::list) VMANAGER_GETDATA_NOTFOUND(std::list >) VMANAGER_GETDATA_NOTFOUND(std::vector) VMANAGER_GETDATA_NOTFOUND(std::vector) VMANAGER_GETDATA_NOTFOUND(std::vector) VMANAGER_GETDATA_NOTFOUND(std::vector) VMANAGER_GETDATA_NOTFOUND(utl::TabulatedFunction) VMANAGER_GETDATA_NOTFOUND(utl::TabulatedFunctionComplexLgAmpPhase) VMANAGER_GETDATA_NOTFOUND(std::map) VMANAGER_GETDATA_NOTFOUND(std::map) void Init(const std::string& configLink); private: Status GetOkFlag(int& returnData, const std::string& componentProperty, const std::string& componentName, const IndexMap& componentIndex) const; int GetStationIndex(const IndexMap& componentIndex) const; bool ReadData() const; private: int fVerbosity; std::string fFileName; TTree* fDataTree; // number of pixels per telescope (as used in AugerUpTime-tree) //static const unsigned int kNPixelPerTel = 440; // data cache mutable det::ValidityStamp fDataValidity; mutable std::vector fStationList; // index map static const double fgSearchMapBinning; typedef std::pair IndexEntry; typedef std::map Index; typedef std::map::iterator IndexIterator; typedef std::map::const_iterator ConstIndexIterator; Index fIndex; }; } // sdet #endif // _sdet_T2LifeROOTFileManager_h_ // Configure (x)emacs for this file ... // Local Variables: // mode: c++ // compile-command: "make -C .. -k" // End: