Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:14:11

0001 //==========================================================================
0002 //  AIDA Detector description implementation 
0003 //--------------------------------------------------------------------------
0004 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
0005 // All rights reserved.
0006 //
0007 // For the licensing terms see $DD4hepINSTALL/LICENSE.
0008 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
0009 //
0010 // Author     : M.Frank
0011 //
0012 //==========================================================================
0013 
0014 // Framework include files
0015 #include <DDEve/DDG4EventHandler.h>
0016 #include <DD4hep/Printout.h>
0017 #include <DD4hep/Objects.h>
0018 #include <DD4hep/Factories.h>
0019 
0020 #include <TFile.h>
0021 #include <TTree.h>
0022 #include <TBranch.h>
0023 
0024 // C/C++ include files
0025 #include <stdexcept>
0026 
0027 using namespace dd4hep;
0028 
0029 ClassImp(DDG4EventHandler)
0030 namespace {
0031   union FCN  {
0032     FCN(void* p) { ptr = p; }
0033     DDG4EventHandler::HitAccessor_t hits;
0034     DDG4EventHandler::ParticleAccessor_t particles;
0035     void* ptr;
0036   };
0037   /// Factory entry point
0038   void* _create(const char*)  {
0039     EventHandler* h = new DDG4EventHandler();
0040     return h;
0041   }
0042 }
0043 using namespace dd4hep::detail;
0044 DECLARE_CONSTRUCTOR(DD4hep_DDEve_DDG4EventHandler,_create)
0045 
0046 /// Standard constructor
0047 DDG4EventHandler::DDG4EventHandler() : EventHandler(), m_file(0,0), m_entry(-1) {
0048   void* ptr = PluginService::Create<void*>("DD4hep_DDEve_DDG4HitAccess",(const char*)"");
0049   if ( 0 == ptr )   {
0050     throw std::runtime_error("FATAL: Failed to access function pointer from factory DD4hep_DDEve_DDG4HitAccess");
0051   }
0052   m_simhitConverter = FCN(ptr).hits;
0053   ptr = PluginService::Create<void*>("DD4hep_DDEve_DDG4ParticleAccess",(const char*)"");
0054   if ( 0 == ptr )   {
0055     throw std::runtime_error("FATAL: Failed to access function pointer from factory DD4hep_DDEve_DDG4ParticleAccess");
0056   }
0057   m_particleConverter = FCN(ptr).particles;
0058 }
0059 
0060 /// Default destructor
0061 DDG4EventHandler::~DDG4EventHandler()   {
0062   if ( m_file.first )  {
0063     m_file.first->Close();
0064     delete m_file.first;
0065   }
0066   m_file.first = 0;
0067   m_file.second = 0;
0068 }
0069 
0070 /// Load the next event
0071 bool DDG4EventHandler::NextEvent()   {
0072   return ReadEvent(++m_entry) > 0;
0073 }
0074 
0075 /// Load the previous event
0076 bool DDG4EventHandler::PreviousEvent()   {
0077   return ReadEvent(--m_entry) > 0;
0078 }
0079 
0080 /// Goto a specified event in the file
0081 bool DDG4EventHandler::GotoEvent(long event_number)   {
0082   return ReadEvent(m_entry = event_number) > 0;
0083 }
0084 
0085 /// Access the number of events on the current input data source (-1 if no data source connected)
0086 long DDG4EventHandler::numEvents() const   {
0087   if ( hasFile() )  {
0088     return m_file.second->GetEntries();
0089   }
0090   return -1;
0091 }
0092 
0093 /// Access the data source name
0094 std::string DDG4EventHandler::datasourceName() const   {
0095   if ( hasFile() )  {
0096     return m_file.first->GetName();
0097   }
0098   return "UNKNOWN";
0099 }
0100 
0101 /// Access to the collection type by name
0102 EventHandler::CollectionType DDG4EventHandler::collectionType(const std::string& collection) const {
0103   Branches::const_iterator i = m_branches.find(collection);
0104   if ( i != m_branches.end() )   {
0105     const char* cl = (*i).second.first->GetClassName();
0106     if ( ::strstr(cl,"sim::Geant4Calorimeter::Hit") )  return CALO_HIT_COLLECTION;
0107     else if ( ::strstr(cl,"sim::Geant4Tracker::Hit") ) return TRACKER_HIT_COLLECTION;
0108     else if ( ::strstr(cl,"sim::Geant4Particle") )     return PARTICLE_COLLECTION;
0109     // These are OLD types. Eventually remove these lines.....
0110     else if ( ::strstr(cl,"sim::SimpleCalorimeter::Hit") ) return CALO_HIT_COLLECTION;
0111     else if ( ::strstr(cl,"sim::SimpleTracker::Hit") )     return TRACKER_HIT_COLLECTION;
0112     else if ( ::strstr(cl,"sim::Particle") )               return PARTICLE_COLLECTION;
0113   }
0114   return NO_COLLECTION;
0115 }
0116 
0117 /// Call functor on hit collection
0118 size_t DDG4EventHandler::collectionLoop(const std::string& collection, DDEveHitActor& actor)   {
0119   typedef std::vector<void*> _P;
0120   Branches::const_iterator i = m_branches.find(collection);
0121   if ( i != m_branches.end() )   {
0122     const _P* data_ptr = (_P*)(*i).second.second;
0123     if ( data_ptr )  {
0124       DDEveHit hit;
0125       actor.setSize(data_ptr->size());
0126       for(_P::const_iterator j=data_ptr->begin(); j!=data_ptr->end(); ++j)   {
0127         if ( (*m_simhitConverter)(*j,&hit) )    {
0128           actor(hit);
0129         }
0130       }
0131       return data_ptr->size();
0132     }
0133   }
0134   return 0;
0135 }
0136 
0137 /// Loop over collection and extract particle data
0138 size_t DDG4EventHandler::collectionLoop(const std::string& collection, DDEveParticleActor& actor)    {
0139   typedef std::vector<void*> _P;
0140   Branches::const_iterator i = m_branches.find(collection);
0141   if ( i != m_branches.end() )   {
0142     const _P* data_ptr = (_P*)(*i).second.second;
0143     if ( data_ptr )  {
0144       DDEveParticle part;
0145       actor.setSize(data_ptr->size());
0146       for(_P::const_iterator j=data_ptr->begin(); j!=data_ptr->end(); ++j)   {
0147         if ( (*m_particleConverter)(*j,&part) )    {
0148           actor(part);
0149         }
0150       }
0151       return data_ptr->size();
0152     }
0153   }
0154   return 0;
0155 }
0156 
0157 /// Load the specified event
0158 Int_t DDG4EventHandler::ReadEvent(Long64_t event_number)   {
0159   m_data.clear();
0160   m_hasEvent = false;
0161   if ( hasFile() )  {
0162     if ( event_number >= m_file.second->GetEntries() )  {
0163       event_number = m_file.second->GetEntries()-1;
0164       printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Cannot read across End-of-file! Reading last event:%d.",event_number);
0165     }
0166     else if ( event_number < 0 )  {
0167       event_number = 0;
0168       printout(ERROR,"DDG4EventHandler","+++ nextEvent: Cannot read across Start-of-file! Reading first event:%d.",event_number);
0169     }
0170 
0171     Int_t nbytes = m_file.second->GetEntry(event_number);
0172     if ( nbytes >= 0 )   {
0173       printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Read %d bytes of event data for entry:%d",nbytes,event_number);
0174       for(Branches::const_iterator i=m_branches.begin(); i != m_branches.end(); ++i)  {
0175         TBranch* b = (*i).second.first;
0176         std::vector<void*>* ptr_data = *(std::vector<void*>**)b->GetAddress();
0177         m_data[b->GetClassName()].emplace_back(b->GetName(),ptr_data->size());
0178       }
0179       m_hasEvent = true;
0180       return nbytes;
0181     }
0182     printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Cannot read event data for entry:%d",event_number);
0183     throw std::runtime_error("+++ EventHandler::readEvent: Failed to read event");
0184   }
0185   throw std::runtime_error("+++ EventHandler::readEvent: No file open!");
0186 }
0187 
0188 /// Open new data file
0189 bool DDG4EventHandler::Open(const std::string&, const std::string& name)   {
0190   if ( m_file.first ) m_file.first->Close();
0191   m_hasFile = false;
0192   m_hasEvent = false;
0193   TFile* f = TFile::Open(name.c_str());
0194   if ( f && !f->IsZombie() )  {
0195     m_file.first = f;
0196     TTree* t = (TTree*)f->Get("EVENT");
0197     if ( t )   {
0198       TObjArray* br = t->GetListOfBranches();
0199       m_file.second = t;
0200       m_entry = -1;
0201       m_branches.clear();
0202       for(Int_t i=0; i<br->GetSize(); ++i)  {
0203         TBranch* b = (TBranch*)br->At(i);
0204         if ( !b ) continue;
0205         m_branches[b->GetName()] = std::make_pair(b,(void*)0);
0206         printout(INFO,"DDG4EventHandler::open","+++ Branch %s has %ld entries.",b->GetName(),b->GetEntries());
0207       }
0208       for(Int_t i=0; i<br->GetSize(); ++i)  {
0209         TBranch* b = (TBranch*)br->At(i);
0210         if ( !b ) continue;
0211         b->SetAddress(&m_branches[b->GetName()].second);
0212       }
0213       m_hasFile = true;
0214       return true;
0215     }
0216     throw std::runtime_error("+++ Failed to access tree EVENT in ROOT file:"+name);
0217   }
0218   throw std::runtime_error("+++ Failed to open ROOT file:"+name);
0219 }