DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
DDG4EventHandler.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
14 // Framework include files
15 #include <DDEve/DDG4EventHandler.h>
16 #include <DD4hep/Printout.h>
17 #include <DD4hep/Objects.h>
18 #include <DD4hep/Factories.h>
19 
20 #include <TFile.h>
21 #include <TTree.h>
22 #include <TBranch.h>
23 
24 // C/C++ include files
25 #include <stdexcept>
26 
27 using namespace dd4hep;
28 
30 namespace {
31  union FCN {
32  FCN(void* p) { ptr = p; }
35  void* ptr;
36  };
38  void* _create(const char*) {
39  EventHandler* h = new DDG4EventHandler();
40  return h;
41  }
42 }
43 using namespace dd4hep::detail;
44 DECLARE_CONSTRUCTOR(DD4hep_DDEve_DDG4EventHandler,_create)
45 
46 DDG4EventHandler::DDG4EventHandler() : EventHandler(), m_file(0,0), m_entry(-1) {
48  void* ptr = PluginService::Create<void*>("DD4hep_DDEve_DDG4HitAccess",(const char*)"");
49  if ( 0 == ptr ) {
50  throw std::runtime_error("FATAL: Failed to access function pointer from factory DD4hep_DDEve_DDG4HitAccess");
51  }
52  m_simhitConverter = FCN(ptr).hits;
53  ptr = PluginService::Create<void*>("DD4hep_DDEve_DDG4ParticleAccess",(const char*)"");
54  if ( 0 == ptr ) {
55  throw std::runtime_error("FATAL: Failed to access function pointer from factory DD4hep_DDEve_DDG4ParticleAccess");
56  }
57  m_particleConverter = FCN(ptr).particles;
58 }
59 
61 DDG4EventHandler::~DDG4EventHandler() {
62  if ( m_file.first ) {
63  m_file.first->Close();
64  delete m_file.first;
65  }
66  m_file.first = 0;
67  m_file.second = 0;
68 }
69 
72  return ReadEvent(++m_entry) > 0;
73 }
74 
77  return ReadEvent(--m_entry) > 0;
78 }
79 
81 bool DDG4EventHandler::GotoEvent(long event_number) {
82  return ReadEvent(m_entry = event_number) > 0;
83 }
84 
87  if ( hasFile() ) {
88  return m_file.second->GetEntries();
89  }
90  return -1;
91 }
92 
94 std::string DDG4EventHandler::datasourceName() const {
95  if ( hasFile() ) {
96  return m_file.first->GetName();
97  }
98  return "UNKNOWN";
99 }
100 
103  Branches::const_iterator i = m_branches.find(collection);
104  if ( i != m_branches.end() ) {
105  const char* cl = (*i).second.first->GetClassName();
106  if ( ::strstr(cl,"sim::Geant4Calorimeter::Hit") ) return CALO_HIT_COLLECTION;
107  else if ( ::strstr(cl,"sim::Geant4Tracker::Hit") ) return TRACKER_HIT_COLLECTION;
108  else if ( ::strstr(cl,"sim::Geant4Particle") ) return PARTICLE_COLLECTION;
109  // These are OLD types. Eventually remove these lines.....
110  else if ( ::strstr(cl,"sim::SimpleCalorimeter::Hit") ) return CALO_HIT_COLLECTION;
111  else if ( ::strstr(cl,"sim::SimpleTracker::Hit") ) return TRACKER_HIT_COLLECTION;
112  else if ( ::strstr(cl,"sim::Particle") ) return PARTICLE_COLLECTION;
113  }
114  return NO_COLLECTION;
115 }
116 
118 size_t DDG4EventHandler::collectionLoop(const std::string& collection, DDEveHitActor& actor) {
119  typedef std::vector<void*> _P;
120  Branches::const_iterator i = m_branches.find(collection);
121  if ( i != m_branches.end() ) {
122  const _P* data_ptr = (_P*)(*i).second.second;
123  if ( data_ptr ) {
124  DDEveHit hit;
125  actor.setSize(data_ptr->size());
126  for(_P::const_iterator j=data_ptr->begin(); j!=data_ptr->end(); ++j) {
127  if ( (*m_simhitConverter)(*j,&hit) ) {
128  actor(hit);
129  }
130  }
131  return data_ptr->size();
132  }
133  }
134  return 0;
135 }
136 
138 size_t DDG4EventHandler::collectionLoop(const std::string& collection, DDEveParticleActor& actor) {
139  typedef std::vector<void*> _P;
140  Branches::const_iterator i = m_branches.find(collection);
141  if ( i != m_branches.end() ) {
142  const _P* data_ptr = (_P*)(*i).second.second;
143  if ( data_ptr ) {
144  DDEveParticle part;
145  actor.setSize(data_ptr->size());
146  for(_P::const_iterator j=data_ptr->begin(); j!=data_ptr->end(); ++j) {
147  if ( (*m_particleConverter)(*j,&part) ) {
148  actor(part);
149  }
150  }
151  return data_ptr->size();
152  }
153  }
154  return 0;
155 }
156 
158 Int_t DDG4EventHandler::ReadEvent(Long64_t event_number) {
159  m_data.clear();
160  m_hasEvent = false;
161  if ( hasFile() ) {
162  if ( event_number >= m_file.second->GetEntries() ) {
163  event_number = m_file.second->GetEntries()-1;
164  printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Cannot read across End-of-file! Reading last event:%d.",event_number);
165  }
166  else if ( event_number < 0 ) {
167  event_number = 0;
168  printout(ERROR,"DDG4EventHandler","+++ nextEvent: Cannot read across Start-of-file! Reading first event:%d.",event_number);
169  }
170 
171  Int_t nbytes = m_file.second->GetEntry(event_number);
172  if ( nbytes >= 0 ) {
173  printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Read %d bytes of event data for entry:%d",nbytes,event_number);
174  for(Branches::const_iterator i=m_branches.begin(); i != m_branches.end(); ++i) {
175  TBranch* b = (*i).second.first;
176  std::vector<void*>* ptr_data = *(std::vector<void*>**)b->GetAddress();
177  m_data[b->GetClassName()].emplace_back(b->GetName(),ptr_data->size());
178  }
179  m_hasEvent = true;
180  return nbytes;
181  }
182  printout(ERROR,"DDG4EventHandler","+++ ReadEvent: Cannot read event data for entry:%d",event_number);
183  throw std::runtime_error("+++ EventHandler::readEvent: Failed to read event");
184  }
185  throw std::runtime_error("+++ EventHandler::readEvent: No file open!");
186 }
187 
189 bool DDG4EventHandler::Open(const std::string&, const std::string& name) {
190  if ( m_file.first ) m_file.first->Close();
191  m_hasFile = false;
192  m_hasEvent = false;
193  TFile* f = TFile::Open(name.c_str());
194  if ( f && !f->IsZombie() ) {
195  m_file.first = f;
196  TTree* t = (TTree*)f->Get("EVENT");
197  if ( t ) {
198  TObjArray* br = t->GetListOfBranches();
199  m_file.second = t;
200  m_entry = -1;
201  m_branches.clear();
202  for(Int_t i=0; i<br->GetSize(); ++i) {
203  TBranch* b = (TBranch*)br->At(i);
204  if ( !b ) continue;
205  m_branches[b->GetName()] = std::make_pair(b,(void*)0);
206  printout(INFO,"DDG4EventHandler::open","+++ Branch %s has %ld entries.",b->GetName(),b->GetEntries());
207  }
208  for(Int_t i=0; i<br->GetSize(); ++i) {
209  TBranch* b = (TBranch*)br->At(i);
210  if ( !b ) continue;
211  b->SetAddress(&m_branches[b->GetName()].second);
212  }
213  m_hasFile = true;
214  return true;
215  }
216  throw std::runtime_error("+++ Failed to access tree EVENT in ROOT file:"+name);
217  }
218  throw std::runtime_error("+++ Failed to open ROOT file:"+name);
219 }
dd4hep::DDEveParticleActor
Event data actor base class for particles. Used to extract data from concrete classes.
Definition: EventHandler.h:55
dd4hep::EventHandler::NO_COLLECTION
@ NO_COLLECTION
Definition: EventHandler.h:70
dd4hep::DDEveHit
DDEve event classes: Basic hit.
Definition: DDEveEventData.h:29
dd4hep::DDG4EventHandler::Open
virtual bool Open(const std::string &type, const std::string &file_name) override
Open new data file.
Definition: DDG4EventHandler.cpp:189
Objects.h
dd4hep::EventHandler::m_hasEvent
bool m_hasEvent
Flag to indicate that an event is loaded.
Definition: EventHandler.h:86
dd4hep::DDG4EventHandler::m_file
std::pair< TFile *, TTree * > m_file
Reference to data file.
Definition: DDG4EventHandler.h:45
dd4hep::EventHandler
Event handler base class: Interface to all DDEve I/O actions.
Definition: EventHandler.h:67
DDG4EventHandler.h
dd4hep::EventHandler::TRACKER_HIT_COLLECTION
@ TRACKER_HIT_COLLECTION
Definition: EventHandler.h:74
dd4hep::DDG4EventHandler
Event I/O handler class for the dd4hep event display.
Definition: DDG4EventHandler.h:38
dd4hep::DDG4EventHandler::ReadEvent
Int_t ReadEvent(Long64_t n)
Load the specified event.
Definition: DDG4EventHandler.cpp:158
dd4hep::DDG4EventHandler::GotoEvent
virtual bool GotoEvent(long event_number) override
Goto a specified event in the file.
Definition: DDG4EventHandler.cpp:81
dd4hep::DDEveHitActor
Event data actor base class for hits. Used to extract data from concrete classes.
Definition: EventHandler.h:43
Factories.h
dd4hep::EventHandler::m_hasFile
bool m_hasFile
Flag to indicate that a file is opened.
Definition: EventHandler.h:84
dd4hep::DDG4EventHandler::PreviousEvent
virtual bool PreviousEvent() override
User overloadable function: Load the previous event.
Definition: DDG4EventHandler.cpp:76
DECLARE_CONSTRUCTOR
#define DECLARE_CONSTRUCTOR(name, func)
Definition: Factories.h:286
dd4hep::EventHandler::hasFile
virtual bool hasFile() const
Check if a data file is connected to the handler.
Definition: EventHandler.h:95
dd4hep::DDG4EventHandler::collectionLoop
virtual size_t collectionLoop(const std::string &collection, DDEveHitActor &actor) override
Call functor on hit collection.
Definition: DDG4EventHandler.cpp:118
dd4hep::DDG4EventHandler::m_particleConverter
ParticleAccessor_t m_particleConverter
Function pointer to interprete particles.
Definition: DDG4EventHandler.h:53
dd4hep::detail
DD4hep internal namespace.
Definition: Alignments.h:32
dd4hep::EventHandler::CALO_HIT_COLLECTION
@ CALO_HIT_COLLECTION
Definition: EventHandler.h:73
dd4hep::DDG4EventHandler::m_simhitConverter
HitAccessor_t m_simhitConverter
Function pointer to interprete hits.
Definition: DDG4EventHandler.h:51
dd4hep::DDG4EventHandler::m_branches
Branches m_branches
Branch map.
Definition: DDG4EventHandler.h:47
dd4hep::EventHandler::CollectionType
CollectionType
Definition: EventHandler.h:69
dd4hep::DDG4EventHandler::numEvents
virtual long numEvents() const override
Access the number of events on the current input data source (-1 if no data source connected)
Definition: DDG4EventHandler.cpp:86
dd4hep::DDEveParticleActor::setSize
virtual void setSize(size_t)
Definition: EventHandler.h:58
dd4hep::DDG4EventHandler::datasourceName
std::string datasourceName() const override
Access the data source name.
Definition: DDG4EventHandler.cpp:94
dd4hep::DDG4EventHandler::ParticleAccessor_t
void *(* ParticleAccessor_t)(void *, DDEveParticle *)
Definition: DDG4EventHandler.h:41
dd4hep::DDG4EventHandler::m_data
TypedEventCollections m_data
Data collection map.
Definition: DDG4EventHandler.h:55
dd4hep::EventHandler::PARTICLE_COLLECTION
@ PARTICLE_COLLECTION
Definition: EventHandler.h:71
dd4hep::DDG4EventHandler::HitAccessor_t
void *(* HitAccessor_t)(void *, DDEveHit *)
Definition: DDG4EventHandler.h:42
ClassImp
ClassImp(DDG4EventHandler) namespace
Definition: DDG4EventHandler.cpp:29
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::DDG4EventHandler::NextEvent
virtual bool NextEvent() override
User overloadable function: Load the next event.
Definition: DDG4EventHandler.cpp:71
dd4hep::DDG4EventHandler::collectionType
virtual CollectionType collectionType(const std::string &collection) const override
Access to the collection type by name.
Definition: DDG4EventHandler.cpp:102
dd4hep::DDEveParticle
Data structure to store the MC particle information.
Definition: DDEveEventData.h:56
dd4hep::DDG4EventHandler::DDG4EventHandler
DDG4EventHandler()
Standard constructor.
Definition: DDG4EventHandler.cpp:47
dd4hep::DDG4EventHandler::m_entry
Long64_t m_entry
File entry number.
Definition: DDG4EventHandler.h:49
Printout.h
dd4hep::DDEveHitActor::setSize
virtual void setSize(size_t)
Definition: EventHandler.h:46