DD4hep - The AIDA detector description toolkit for high energy physics experiments
DD4hep  Rev:Unversioneddirectory
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
GenericEventHandler.cpp
Go to the documentation of this file.
1 // $Id: $
2 //==========================================================================
3 // AIDA Detector description implementation for LCD
4 //--------------------------------------------------------------------------
5 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
6 // All rights reserved.
7 //
8 // For the licensing terms see $DD4hepINSTALL/LICENSE.
9 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
10 //
11 // Author : M.Frank
12 //
13 //==========================================================================
14 
15 // Framework include files
17 #include "DD4hep/Primitives.h"
18 #include "DD4hep/Factories.h"
19 #include "DD4hep/Plugins.h"
20 #include <stdexcept>
21 
23 #include "TGMsgBox.h"
24 #include "TSystem.h"
25 #include <climits>
26 
27 using namespace std;
28 using namespace DD4hep;
29 
31 
32 GenericEventHandler::GenericEventHandler() : m_current(0) {
34 }
35 
37 GenericEventHandler::~GenericEventHandler() {
38  m_subscriptions.clear();
39  deletePtr(m_current);
40 }
41 
42 EventHandler* GenericEventHandler::current() const {
43  if ( m_current ) {
44  return m_current;
45  }
46  throw runtime_error("Invalid event handler");
47 }
48 
50 void GenericEventHandler::NotifySubscribers(void (EventConsumer::*pmf)(EventHandler*)) {
51  for(Subscriptions::iterator i=m_subscriptions.begin(); i!=m_subscriptions.end();++i)
52  ((*i)->*pmf)(this);
53 }
54 
56 void GenericEventHandler::Subscribe(EventConsumer* consumer) {
57  m_subscriptions.insert(consumer);
58 }
59 
61 void GenericEventHandler::Unsubscribe(EventConsumer* consumer) {
62  Subscriptions::iterator i=m_subscriptions.find(consumer);
63  if ( i != m_subscriptions.end() ) m_subscriptions.erase(i);
64 }
65 
67 long GenericEventHandler::numEvents() const {
68  return current()->numEvents();
69 }
70 
72 string GenericEventHandler::datasourceName() const {
73  return current()->datasourceName();
74 }
75 
77 EventHandler::CollectionType GenericEventHandler::collectionType(const string& collection) const {
78  if ( m_current && m_current->hasEvent() ) {
79  return m_current->collectionType(collection);
80  }
81  return NO_COLLECTION;
82 }
83 
85 size_t GenericEventHandler::collectionLoop(const string& collection, DDEveHitActor& actor) {
86  if ( m_current && m_current->hasEvent() ) {
87  return m_current->collectionLoop(collection,actor);
88  }
89  return 0;
90 }
91 
93 size_t GenericEventHandler::collectionLoop(const string& collection, DDEveParticleActor& actor) {
94  if ( m_current && m_current->hasEvent() ) {
95  return m_current->collectionLoop(collection,actor);
96  }
97  return 0;
98 }
99 
101 bool GenericEventHandler::Open(const string& file_type, const string& file_name) {
102  size_t idx = file_name.find("lcio");
103  size_t idr = file_name.find("root");
104  string err;
105  m_hasFile = false;
106  m_hasEvent = false;
107  try {
108  deletePtr(m_current);
109  // prefer event handler configured in xml
110  if ( file_type.find("FCC") != string::npos ) {
111  m_current = (EventHandler*)PluginService::Create<void*>("DDEve_FCCEventHandler",(const char*)0);
112  }
113  // fall back to defaults according to file ending
114  else if ( idx != string::npos ) {
115  m_current = (EventHandler*)PluginService::Create<void*>("DDEve_LCIOEventHandler",(const char*)0);
116  }
117  else if ( idr != string::npos ) {
118  m_current = (EventHandler*)PluginService::Create<void*>("DDEve_DDG4EventHandler",(const char*)0);
119  }
120  else {
121  throw runtime_error("Attempt to open file:"+file_name+" of unknown type:"+file_type);
122  }
123  if ( m_current ) {
124  if ( m_current->Open(file_type, file_name) ) {
125  m_hasFile = true;
126  NotifySubscribers(&EventConsumer::OnFileOpen);
127  return true;
128  }
129  err = "+++ Failed to open the data file:"+file_name;
130  deletePtr(m_current);
131  }
132  else {
133  err = "+++ Failed to create fikle reader for file '"+file_name+"' of type '"+file_type+"'";
134  }
135  }
136  catch(const exception& e) {
137  err = "\nAn exception occurred \n"
138  "while opening event data:\n" + string(e.what()) + "\n\n";
139  }
140  string path = TString::Format("%s/icons/stop_t.xpm", gSystem->Getenv("ROOTSYS")).Data();
141  const TGPicture* pic = gClient->GetPicture(path.c_str());
142  new TGMsgBox(gClient->GetRoot(),0,"Failed to open event data",err.c_str(),pic,
143  kMBDismiss,0,kVerticalFrame,kTextLeft|kTextCenterY);
144  return false;
145 }
146 
148 bool GenericEventHandler::NextEvent() {
149  m_hasEvent = false;
150  try {
151  if ( m_hasFile ) {
152  if ( current()->NextEvent() > 0 ) {
153  m_hasEvent = true;
154  NotifySubscribers(&EventConsumer::OnNewEvent);
155  return 1;
156  }
157  }
158  throw runtime_error("+++ EventHandler::readEvent: No file open!");
159  }
160  catch(const exception& e) {
161  string path = TString::Format("%s/icons/stop_t.xpm", gSystem->Getenv("ROOTSYS")).Data();
162  string err = "\nAn exception occurred \n"
163  "while reading a new event:\n" + string(e.what()) + "\n\n";
164  const TGPicture* pic = gClient->GetPicture(path.c_str());
165  new TGMsgBox(gClient->GetRoot(),0,"Failed to read event", err.c_str(),pic,
166  kMBDismiss,0,kVerticalFrame,kTextLeft|kTextCenterY);
167  }
168  return -1;
169 }
170 
172 bool GenericEventHandler::PreviousEvent() {
173  m_hasEvent = false;
174  if ( m_hasFile && current()->PreviousEvent() > 0 ) {
175  m_hasEvent = true;
176  NotifySubscribers(&EventConsumer::OnNewEvent);
177  return 1;
178  }
179  return -1;
180 }
181 
183 bool GenericEventHandler::GotoEvent(long event_number) {
184  m_hasEvent = false;
185  if ( m_hasFile && current()->GotoEvent(event_number) > 0 ) {
186  m_hasEvent = true;
187  NotifySubscribers(&EventConsumer::OnNewEvent);
188  return 1;
189  }
190  return -1;
191 }
Event data actor base class for hits. Used to extract data from concrete classes. ...
Definition: EventHandler.h:44
void deletePtr(T *&p)
Helper to delete objects from heap and reset the pointer. Saves many many lines of code...
Definition: Primitives.h:234
Event handler base class. Interface to all DDEve I/O actions.
Definition: EventHandler.h:66
return e
Definition: Volumes.cpp:297
Event handler base class. Interface to all DDEve I/O actions.
Event data actor base class for particles. Used to extract data from concrete classes.
Definition: EventHandler.h:55
ClassImp(GenericEventHandler) GenericEventHandler
Standard constructor.