More information on this package is here.
//=================================================================
//  TDSTReader is a calss implements tools for Reading ROOT DST's
//  without any run group specific details. 
//  It will contain tools only to read EVENT and Scaller trees.
// 
// ________________________________________________________________
//
//   Author : G.Gavalian (UNH) 04/02/2002
//  
//

#ifndef __TDST_READER__
#define __TDST_READER__

#define WRITE_MCHD
#define WRITE_GSIM
#define WRITE_PHOTON
#define WRITE_VERTEX
#define WRITE_TBER
// #define WRITE_TGBI

// Standart Includes 

#include <iostream>
using namespace std;
#include <stdlib.h>
#include <unistd.h>

// ROOT includes 

#include "TROOT.h"
#include "TSystem.h"
#include "TObject.h"
#include "TString.h"
#include "TRegexp.h"
#include "TChain.h"
#include "TFile.h"
#include "TGFileDialog.h"
#include "TTree.h"
#include "TClonesArray.h"
#include "TMatrix.h"
#include "TVector2.h"
#include "TVector3.h"
#include "TRotation.h"
#include "TLorentzVector.h"
#include "TLorentzRotation.h"
#include "TMath.h"
#include "TH1.h"
#include "TRandom.h"
#include "TString.h"
#include "TObjString.h"

#ifdef __USE_MYSQL__
#include "TSQLServer.h"
#include "TSQLRow.h"
#include "TSQLResult.h"
#endif

//
// Base class.
//
#include "TVirtualReader.h"
//
// CLAS Banks Includes
//

#include "THEADERClass.h"
#include "TEVNTClass.h"
#include "TECPBClass.h"
#include "TSCPBClass.h"
#include "TDCPBClass.h"
#include "TCCPBClass.h"
#include "TLCPBClass.h"
#include "THEADClass.h"
#include "TTGBIClass.h"
#include "TEPICClass.h"

// Photon Banks
#ifdef WRITE_PHOTON
#include "TSTPBClass.h"
#include "TTGPBClass.h"
#include "TTAGRClass.h"
#endif

// GSIM Bank 
#ifdef WRITE_MCHD
#include "TMCHDClass.h"
#endif

#ifdef WRITE_GSIM
#include "TGSIMClass.h"
#endif

#ifdef WRITE_VERTEX
#include "TVERTClass.h"
#include "TMVRTClass.h"
#endif

#ifdef WRITE_TBER
#include "TTBERClass.h"
#endif

//
// Scaler Banks
//
#include "THLSClass.h"
#include "TS1STClass.h"
#include "TTRGSClass.h"


class TDSTReader : public TVirtualReader{

 public:

  // Data specific to a ROOT DST file.

  THEADERClass *fEvntHeader;   // Event Header pointer
  THEADClass   *fScalerHead;   // Scaler Event Header pointer
  TEPICClass   *fEPIC;         // EPICClass pointer

  TClonesArray *fcaDataPointer;  // Points to fcaEVNTStore or fcaGSIMStore
  TClonesArray *fcaEVNTStore;  // Clones Array for EVNT bank
  TClonesArray *fcaECPBStore;  // Clones Array for ECPB bank
  TClonesArray *fcaSCPBStore;  // Clones Array for SCPB bank
  TClonesArray *fcaDCPBStore;  // Clones Array for DCPB bank
  TClonesArray *fcaCCPBStore;  // Clones Array for CCPB bank
  TClonesArray *fcaLCPBStore;  // Clones Array for LCPB bank
#ifdef WRITE_PHOTON
  TClonesArray *fcaSTPBStore;  // Clones Array for STPB bank
  TClonesArray *fcaTGPBStore;  // Clones Array for TGPB bank
  TClonesArray *fcaTAGRStore;  // Clones Array for TAGR bank
#endif

#ifdef WRITE_VERTEX
  TClonesArray *fcaVERTStore;  // Clones Array for VERT bank
  TClonesArray *fcaMVRTStore;  // Clones Array for MVRT bank
#endif

#ifdef WRITE_TBER
  TClonesArray *fcaTBERStore; // Clones Array for the TBER bank.
#endif

#ifdef WRITE_TGBI
  TClonesArray *fcaTGBIStore;  // Clones Array for TGBI bank
#endif

  //---------------------------------- For MCHD Simulation header bank
#ifdef WRITE_MCHD
  TClonesArray *fcaMCHDStore;
#endif

  //---------------------------------- For GSIM Simulation bank
#ifdef WRITE_GSIM
  TClonesArray *fcaGSIMStore;
#endif
  //---------------------------------- For Scaler Events.
  TClonesArray *fcaHLSStore;  // Clones Array for HLS bank
  TClonesArray *fcaS1STStore; // Clones Array for S1ST bank
  TClonesArray *fcaTRGSStore; // Clones Array for TRGS bank

  TChain *fEventChain;    // Chain of files with data. TTree points to CLASEVENT 
  TChain *fScalerChain;   // Chain of files with data. TTree points to SCALER
  TChain *fEPICChain;     // Chain of files with data. TTree points to CLASEPIC

  Int_t fCurrentEPIC;  // Keeps track of position in EPIC bank.Perhaps should be in TVirtualReader?
  TArrayI *EPICTable;   // Table storing the last physics event before each EPIC. Used in GetEPICValue

 public:

  TDSTReader();   // Default Constructor.
  ~TDSTReader(); // Default Destructor.
  TDSTReader( const TDSTReader& reader );  // Copy Constructor

  Int_t   Add( const Char_t *filename );     // Add data files to the chain, using wildcards.
  Int_t   AddFile(const Char_t *filename);   // Add a data file to the chain
  Int_t   AddFile(const Char_t *filename,Int_t nentries); // Add a data file to the chain
  Int_t   AddEPICFile(const Char_t *filename,Int_t nentries); // Add an Epics file to the EPICChain
  void    CleanEventClones();   // Clean the TClonesArray's before getting next event.
  void    CleanScalerClones();  // Clean the scalar TClonesArray's before getting next scalar event.
  void    InitClones();         // Initialize the TClonesArray's for retreiving data.
  void    InitBranches();       // Initialize the branch addresses for all the data members.
  void    InitEPICBranch();     // Initialize the branch address for the EPIC branch.

  void     MakeEPICTable();      //Builds the table used to reference which EPIC belongs with which Physics
  Float_t  GetEPICValue(const char* fName); //Return the EPIC value corresponding to fName
  TEPICClass *GetEPICClass(){return(fEPIC);}; //Return the current EPIC. This does not check to see if the current EPIC corresponds to the current physics event, to make sure it is the "right" EPIC call GetEPICValue before calling this.
  TObject  *GetBankRow(const Char_t *bankname,Int_t nrow);
  TChain  *GetChain(char *which="Event"); // Return the Chain requested ("Event"=default, "Scaler")
  Int_t    GetCurrentEvent();
  Int_t    GetEntry(Int_t entry=-1, Int_t getall=0); // Return entry numbered "entry".
  GetEntries_return_t GetEntries() const;                // Get the number of events in the chain.
  GetEntries_return_t GetEntriesFast() const;            // Same as GetEntries, but use previously stored value (much faster.)
  Int_t    GetEventNumber();
  THEADERClass *GetHEADER(){return(fEvntHeader);};            // Return the HEADER.
  Int_t    GetHelicity();               // Return the event helicity bit.
  Int_t    GetNFiles();                 // Return the number of files in the chain.
  Int_t    GetNPart(){ return ( fEvntHeader->GetNPGP() % 100 ); };     // Returns number of particles in event

  Int_t    GetNRows(const Char_t *bankname="EVNT");
  TVector3 GetPart3Vector(Int_t indx);    // Get 3 Vector of Particle in indx-th row
  Int_t    GetPartID(Int_t indx);         // Get Particle ID from indx-th row
  Int_t    GetPartCharge(Int_t indx);     // Get Particle ID from indx-th row
  TVector3 GetPartVertex(Int_t indx);     // Get the vertex vector for particle at indx-th row.
  TString  GetReaderType();             // Returns the reader type
  Int_t    GetRunNumber();
  TChain   *GetTree(Char_t *treename="EVNT"); // Return the TTree.
  Int_t   GetTreeNumber(void){return(fEventChain->GetTreeNumber());}; // Return the number of the current tree.

  Int_t   Next();                          // Read next event from chain
  Bool_t  Notify();                    // Hook that will be notified when a new file is opened by the TChain.
  void    Print(Int_t full=0);  // Print information about the event. 
  void    PrintEventStatus();   // Print the status of the current event.
  void    PrintEPICPosition();    // Print the position of the current EPIC in terms of run and evt numbers
  void    Reset();                     // Reset everything

  ClassDef(TDSTReader,1) // A reader interface for ROOT files.
};

inline Int_t TDSTReader::GetNFiles(){
 return(((TObjArray *)fEventChain->GetListOfFiles())->GetEntries());
}

inline void TDSTReader::Reset(){
  // Set the fCurrentEvent and the fCurrentScaler to 0.
  // Thus the next Next() call will start reading at the beginning of the 
  // current Chain of data files.
  fCurrentEvent=0;
  fCurrentScaler=0;
  fCurrentEPIC=0;
}

inline TChain *TDSTReader::GetChain(char *which){
  if(which[0]=='S'){
    return(fScalerChain);
  }else{
    return(fEventChain);
  }
}
#endif


Author: Maurik Holtrop
Last update: ClasTool/TDSTReader:$Name: $:$Id: TDSTReader.cc,v 1.23 2007/04/05 19:55:25 holtrop Exp $
CopyLeft - This code is freely available.

ROOT Home page - ClasTool Class index - Class Hierarchy - Top of the page