Class Reference for E1039 Core & Analysis Software
Fun4AllHepMCInputManager.cc
Go to the documentation of this file.
2 #include "PHHepMCGenEvent.h"
3 #include "PHHepMCGenEventMap.h"
4 
8 #include <phool/getClass.h>
9 #include <phool/recoConsts.h>
10 
11 #include <ffaobjects/RunHeader.h>
12 
13 //#include <frog/FROG.h>
14 #include <phool/PHCompositeNode.h>
15 #include <phool/PHDataNode.h>
16 
17 #include <HepMC/GenEvent.h>
18 #include <HepMC/IO_GenEvent.h>
19 
20 #include <TPRegexp.h>
21 #include <TString.h>
22 
23 #include <fstream>
24 #include <iostream>
25 #include <istream>
26 #include <memory>
27 #include <sstream>
28 
29 #include <boost/iostreams/filter/bzip2.hpp>
30 #include <boost/iostreams/filter/gzip.hpp>
31 
32 #include <cstdlib>
33 #include <memory>
34 
35 #include <phool/PHRandomSeed.h>
36 
37 #include <gsl/gsl_const.h>
38 
39 using namespace std;
40 
41 static const double toMM = 1.e-12;
42 
43 Fun4AllHepMCInputManager::Fun4AllHepMCInputManager(const string &name, const string &nodename, const string &topnodename)
44  : Fun4AllInputManager(name, nodename, topnodename)
45  , isopen(0)
46  , events_total(0)
47  , events_thisfile(0)
48  , readoscar(0)
49  , topNodeName(topnodename)
50  , ascii_in(NULL)
51  , evt(NULL)
52  , save_evt(NULL)
53  , filestream(NULL)
54  , unzipstream(NULL)
55 {
56  set_embedding_id(0); // default embedding ID. Welcome to change via macro
57 
59  topNode = se->topNode(topNodeName.c_str());
60  PHNodeIterator iter(topNode);
61  PHCompositeNode *dstNode = se->getNode(InputNode.c_str(), topNodeName.c_str());
62 
63  PHHepMCGenEventMap *geneventmap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
64  if (!geneventmap)
65  {
66  geneventmap = new PHHepMCGenEventMap();
67  PHIODataNode<PHObject> *newmapnode = new PHIODataNode<PHObject>(geneventmap, "PHHepMCGenEventMap", "PHObject");
68  dstNode->addNode(newmapnode);
69  }
70 
71  hepmc_helper.set_geneventmap(geneventmap);
72 
73  return;
74 }
75 
77 {
78  fileclose();
79 
80  delete ascii_in;
81  delete filestream;
82  delete unzipstream;
83 }
84 
85 int Fun4AllHepMCInputManager::fileopen(const string &filenam)
86 {
87  if (!mySyncManager)
88  {
89  cout << "Call fileopen only after you registered your Input Manager " << Name() << " with the Fun4AllServer" << endl;
90  exit(1);
91  }
92  if (isopen)
93  {
94  cout << "Closing currently open file "
95  << filename
96  << " and opening " << filenam << endl;
97  fileclose();
98  }
99  filename = filenam;
100  //FROG frog;
101  //string fname(frog.location(filename.c_str()));
102  string fname(filename);
103  if (verbosity > 0)
104  {
105  cout << ThisName << ": opening file " << fname << endl;
106  }
107 
108  if (readoscar)
109  {
110  theOscarFile.open(fname.c_str());
111  }
112  else
113  {
114  TString tstr(fname);
115  TPRegexp bzip_ext(".bz2$");
116  TPRegexp gzip_ext(".gz$");
117  if (tstr.Contains(bzip_ext))
118  {
119  // use boost iosteam library to decompress bz2 on the fly
120  filestream = new ifstream(fname.c_str(), std::ios::in | std::ios::binary);
121  zinbuffer.push(boost::iostreams::bzip2_decompressor());
122  zinbuffer.push(*filestream);
123  unzipstream = new istream(&zinbuffer);
124  ascii_in = new HepMC::IO_GenEvent(*unzipstream);
125  }
126  else if (tstr.Contains(gzip_ext))
127  {
128  // use boost iosream to decompress the gzip file on the fly
129  filestream = new ifstream(fname.c_str(), std::ios::in | std::ios::binary);
130  zinbuffer.push(boost::iostreams::gzip_decompressor());
131  zinbuffer.push(*filestream);
132  unzipstream = new istream(&zinbuffer);
133  ascii_in = new HepMC::IO_GenEvent(*unzipstream);
134  }
135  else
136  {
137  // expects normal ascii hepmc file
138  ascii_in = new HepMC::IO_GenEvent(fname, std::ios::in);
139  }
140  }
141 
143  static bool run_number_forced = rc->FlagExist("RUNNUMBER");
144  if (run_number_forced)
145  {
146  mySyncManager->CurrentRun(rc->get_IntFlag("RUNNUMBER"));
147  }
148  else
149  {
151  }
152  events_thisfile = 0;
153  isopen = 1;
154  AddToFileOpened(fname); // add file to the list of files which were opened
155  return 0;
156 }
157 
159 {
160  // attempt to retrieve a valid event from inputs
161  while (true)
162  {
163  if (!isopen)
164  {
165  if (!filelist.size())
166 
167  {
168  if (verbosity > 0)
169  {
170  cout << "Fun4AllHepMCInputManager::run::" << Name() << ": No Input file open" << endl;
171  }
172  return -1;
173  }
174  else
175  {
176  if (OpenNextFile())
177  {
178  cout << "Fun4AllHepMCInputManager::run::" << Name() << ": No Input file from filelist opened" << endl;
179  return -1;
180  }
181  }
182  }
183 
184  if (save_evt) // if an event was pushed back, copy saved pointer and reset save_evt pointer
185  {
186  evt = save_evt;
187  save_evt = NULL;
188  }
189  else
190  {
191  if (readoscar)
192  {
193  evt = ConvertFromOscar();
194  }
195  else
196  {
197  evt = ascii_in->read_next_event();
198  }
199  }
200 
201  if (!evt)
202  {
203  if (verbosity > 1)
204  {
205  cout << "Fun4AllHepMCInputManager::run::" << Name()
206  << ": error type: " << ascii_in->error_type()
207  << ", rdstate: " << ascii_in->rdstate() << endl;
208  }
209  fileclose();
210  }
211  else
212  {
213  mySyncManager->CurrentEvent(evt->event_number());
214  if (verbosity > 0)
215  {
216  cout << "Fun4AllHepMCInputManager::run::" << Name()
217  << ": hepmc evt no: " << evt->event_number() << endl;
218  }
219 
222  if (ievt != hepmc_helper.get_geneventmap()->end())
223  {
224  // override existing event
225  ievt->second->addEvent(evt);
226  }
227  else
229 
230  events_total++;
231  events_thisfile++;
232 
233  // check if the local SubsysReco discards this event
235  {
236  ResetEvent();
237  }
238  else
239  break; // have a good event, move on
240  }
241  } // attempt to retrieve a valid event from inputs
242 
243  return 0;
244 }
245 
247 {
248  if (!isopen)
249  {
250  cout << Name() << ": fileclose: No Input file open" << endl;
251  return -1;
252  }
253  if (readoscar)
254  {
255  theOscarFile.close();
256  }
257  else
258  {
259  delete ascii_in;
260  ascii_in = NULL;
261  }
262  isopen = 0;
263  // if we have a file list, move next entry to top of the list
264  // or repeat the same entry again
265  if (filelist.size() > 0)
266  {
267  if (repeat)
268  {
269  filelist.push_back(*(filelist.begin()));
270  if (repeat > 0)
271  {
272  repeat--;
273  }
274  }
275  filelist.pop_front();
276  }
277  return 0;
278 }
279 
280 void Fun4AllHepMCInputManager::Print(const string &what) const
281 {
283  return;
284 }
285 
287 {
288  while (filelist.size() > 0)
289  {
290  list<string>::const_iterator iter = filelist.begin();
291  if (verbosity)
292  {
293  cout << PHWHERE << " opening next file: " << *iter << endl;
294  }
295  if (fileopen((*iter).c_str()))
296  {
297  cout << PHWHERE << " could not open file: " << *iter << endl;
298  filelist.pop_front();
299  }
300  else
301  {
302  return 0;
303  }
304  }
305  return -1;
306 }
307 
309 {
310  // PushBackEvents is supposedly pushing events back on the stack which works
311  // easily with root trees (just grab a different entry) but hard in these HepMC ASCII files.
312  // A special case is when the synchronization fails and we need to only push back a single
313  // event. In this case we save the evt pointer as save_evt which is used in the run method
314  // instead of getting the next event.
315  if (i > 0)
316  {
317  if (i == 1 && evt) // check on evt pointer makes sure it is not done from the cmd line
318  {
319  if (verbosity > 3)
320  {
321  cout << ThisName << ": pushing back evt no " << evt->event_number() << endl;
322  }
323  save_evt = evt;
324  return 0;
325  }
326  cout << PHWHERE << ThisName
327  << " Fun4AllHepMCInputManager cannot pop back events into file"
328  << endl;
329  return -1;
330  }
331  if (!isopen)
332  {
333  cout << PHWHERE << ThisName
334  << " no file opened yet" << endl;
335  return -1;
336  }
337  // Skipping events is implemented as
338  // pushing a negative number of events on the stack, so in order to implement
339  // the skipping of events we read -i events.
340  int nevents = -i; // negative number of events to push back -> skip num events
341  int errorflag = 0;
342  while (nevents > 0 && !errorflag)
343  {
344  evt = ascii_in->read_next_event();
345  if (!evt)
346  {
347  cout << "Error after skipping " << i - nevents << endl;
348  cout << "error type: " << ascii_in->error_type()
349  << ", rdstate: " << ascii_in->rdstate() << endl;
350  errorflag = -1;
351  fileclose();
352  }
353  else
354  {
355  if (verbosity > 3)
356  {
357  cout << "Skipping evt no: " << evt->event_number() << endl;
358  }
359  }
360  delete evt;
361  nevents--;
362  }
363  return errorflag;
364 }
365 
366 HepMC::GenEvent *
368 {
369  delete evt;
370  evt = NULL;
371  if (theOscarFile.eof()) // if the file is exhausted bail out during this next read
372  {
373  return evt;
374  }
375 
376  //use PHENIX unit
377  evt = new HepMC::GenEvent(HepMC::Units::GEV, HepMC::Units::CM);
378 
379  if (verbosity > 1) cout << "Reading Oscar Event " << events_total << endl;
380  //Grab New Event From Oscar
381  string theLine;
382  vector<vector<double> > theEventVec;
383  vector<HepMC::FourVector> theVtxVec;
384  while (getline(theOscarFile, theLine))
385  {
386  if (theLine.find("#") == 0) continue;
387  vector<double> theInfo; //format: N,pid,px,py,pz,E,mass,xvtx,yvtx,zvtx,?
388  double number;
389  for (istringstream numbers_iss(theLine); numbers_iss >> number;)
390  {
391  theInfo.push_back(number);
392  }
393 
394  if (theInfo.size() == 2 && theInfo[0] == 0 && theInfo[1] == 0)
395  {
396  break;
397  }
398  else if (theInfo.size() == 2 && theInfo[0] == 0 && theInfo[1] > 0)
399  {
400  continue;
401  }
402  else
403  {
404  theEventVec.push_back(theInfo);
405  HepMC::FourVector vert(theInfo[8] * toMM, theInfo[9] * toMM, theInfo[10] * toMM, theInfo[11]);
406  theVtxVec.push_back(vert);
407  }
408 
409  } //while(getline)
410 
411  //Set Event Number
412  evt->set_event_number(events_total);
413 
414  //Loop Over One Event, Fill HepMC
415  for (unsigned int i = 0; i < theEventVec.size(); i++)
416  {
417  //int N = (int)theEventVec[i][0];
418  int pid = (int) theEventVec[i][1];
419  double px = theEventVec[i][3];
420  double py = theEventVec[i][4];
421  double pz = theEventVec[i][5];
422  double E = theEventVec[i][6];
423  double m = theEventVec[i][7];
424  int status = 1; //oscar only writes final state particles
425 
426  HepMC::GenVertex *v = new HepMC::GenVertex(theVtxVec[i]);
427  evt->add_vertex(v);
428 
429  HepMC::GenParticle *p = new HepMC::GenParticle(HepMC::FourVector(px, py, pz, E), pid, status);
430  p->setGeneratedMass(m);
431  p->suggest_barcode(i + 1);
432  v->add_particle_out(p);
433  }
434  if (verbosity > 3)
435  {
436  evt->print();
437  }
438  return evt;
439 }
static const double toMM
#define NULL
Definition: Pdb.h:9
int verbosity
The verbosity level. 0 means not verbose at all.
Definition: Fun4AllBase.h:75
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
std::string ThisName
Definition: Fun4AllBase.h:72
virtual int PushBackEvents(const int i)
boost::iostreams::filtering_streambuf< boost::iostreams::input > zinbuffer
virtual void Print(const std::string &what="ALL") const
virtual int fileopen(const std::string &filenam)
Fun4AllHepMCInputManager(const std::string &name="DUMMY", const std::string &nodename="DST", const std::string &topnodename="TOP")
virtual int run(const int nevents=0)
PHHepMCGenHelper hepmc_helper
helper for insert HepMC event to DST node and add vertex smearing
virtual void Print(const std::string &what="ALL") const
Fun4AllSyncManager * mySyncManager
std::list< std::string > filelist
void AddToFileOpened(const std::string &filename)
static Fun4AllServer * instance()
PHCompositeNode * topNode() const
Definition: Fun4AllServer.h:59
PHCompositeNode * getNode(const char *name, const char *topnodename="TOP")
void CurrentEvent(const int evt)
PHBoolean addNode(PHNode *)
virtual int FlagExist(const std::string &name) const
Definition: PHFlag.cc:255
virtual int get_IntFlag(const std::string &name) const
Definition: PHFlag.cc:117
PHHepMCGenEventMap is collection of HEPMC events input into this simulation map of embedding ID -> PH...
std::map< int, PHHepMCGenEvent * >::iterator Iter
ConstIter end() const
ConstIter find(unsigned int idkey) const
find
void set_geneventmap(PHHepMCGenEventMap *geneventmap)
const PHHepMCGenEventMap * get_geneventmap() const
int get_embedding_id() const
PHHepMCGenEvent * insert_event(HepMC::GenEvent *evt)
send HepMC::GenEvent to DST tree. This function takes ownership of evt
static recoConsts * instance()
Definition: recoConsts.cc:7
#define PHWHERE
Definition: phool.h:23