EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ReadEICFiles.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file ReadEICFiles.cc
1 #include "ReadEICFiles.h"
2 
3 #include "EicEventHeader.h"
4 #include "EicEventHeaderv1.h"
5 
7 
9 #include <phool/PHIODataNode.h>
10 #include <phool/PHNode.h> // for PHNode
11 #include <phool/PHNodeIterator.h>
12 #include <phool/PHObject.h> // for PHObject
13 #include <phool/getClass.h>
14 #include <phool/phool.h> // for PHWHERE
15 
16 #include <HepMC/GenEvent.h>
17 #include <HepMC/GenParticle.h> // for GenParticle
18 #include <HepMC/GenVertex.h>
19 #include <HepMC/PdfInfo.h> // for PdfInfo
20 #include <HepMC/SimpleVector.h> // for FourVector
21 #include <HepMC/Units.h> // for GEV, MM
22 
23 // eicsmear classes
24 #include <eicsmear/erhic/EventMC.h>
26 #include <eicsmear/erhic/ParticleMC.h> // for ParticleMC
27 #include <eicsmear/erhic/Pid.h> // for Pid
28 
30 
31 // General Root and C++ classes
32 #include <TBranch.h> // for TBranch
33 #include <TChain.h>
34 #include <TVector3.h> // for TVector3
35 
36 #include <iostream> // for operator<<, endl, basic_ostream
37 #include <vector> // for vector
38 
39 class PHHepMCGenEvent;
40 
41 using namespace std;
42 
44 
46  : SubsysReco(name)
47  , filename("")
48  , Tin(nullptr)
49  , nEntries(0)
50  , entry(0)
51  , m_EvtGenId(EvtGen::Unknown)
52  , GenEvent(nullptr)
53  , _node_name("PHHepMCGenEvent")
54 {
55  PHHepMCGenHelper::set_embedding_id(1); // default embedding ID to 1
56  return;
57 }
58 
60 
62 {
63  delete Tin;
64 }
65 
67 
68 bool ReadEICFiles::OpenInputFile(const string &name)
69 {
70  filename = name;
71  Tin = new TChain("EICTree", "EICTree");
72  Tin->Add(name.c_str());
73  GetTree();
74  return true;
75 }
76 
78 
80 {
81  /* Print the actual class of the event branch record,
82  i.e. erhic::EventMilou or other */
83  string EventClass(Tin->GetBranch("event")->GetClassName());
84  if (Verbosity() > 0)
85  {
86  cout << "ReadEICFiles: Input Branch Event Class = "
87  << EventClass << endl;
88  }
89  if (EventClass.find("Milou") != string::npos)
90  {
91  m_EvtGenId = EvtGen::Milou;
92  }
93 
94  if (EventClass.find("DEMP") != string::npos)
95  {
96  m_EvtGenId = EvtGen::DEMP;
97  }
98 
99  Tin->SetBranchAddress("event", &GenEvent);
100  nEntries = Tin->GetEntries();
101 }
102 
104 {
105  /* Create node tree */
106  CreateNodeTree(topNode);
107  return 0;
108 }
109 
111 {
112  /* Check if there is an unused event left in input file */
113  if (entry >= nEntries)
114  {
115  if (filename.size() > 0)
116  {
117  cout << "input file " << filename << " exhausted" << endl;
118  }
119  else
120  {
121  cout << "no file opened" << endl;
122  }
124  }
125 
126  /* Get event record from input file */
127  Tin->GetEntry(entry);
128 
129  EicEventHeader *evthead = findNode::getClass<EicEventHeader>(topNode, "EicEventHeader");
130  switch (m_EvtGenId)
131  {
132  case EvtGen::Milou:
133  {
134  erhic::EventMilou *gen = dynamic_cast<erhic::EventMilou *>(GenEvent);
135  evthead->set_eventgenerator_type(EicEventHeader::EvtGen::Milou);
136  evthead->set_milou_weight(gen->weight);
137  evthead->set_milou_trueX(gen->trueX);
138  evthead->set_milou_trueQ2(gen->trueQ2);
139  }
140  break;
141  case EvtGen::DEMP:
142  {
143  erhic::EventDEMP *gen = dynamic_cast<erhic::EventDEMP *>(GenEvent);
144  evthead->set_eventgenerator_type(EicEventHeader::EvtGen::DEMP);
145  evthead->set_demp_weight(gen->weight);
146  }
147  break;
148  case EvtGen::Unknown:
149  cout << "unknown event generator" << endl;
150  break;
151  default:
152  cout << "what is this " << m_EvtGenId << " ????" << endl;
153  break;
154  }
155 
156  /* Create GenEvent */
157  HepMC::GenEvent *evt = new HepMC::GenEvent();
158 
159  /* define the units (Pythia uses GeV and mm) */
160  evt->use_units(HepMC::Units::GEV, HepMC::Units::MM);
161 
162  /* add global information to the event */
163  evt->set_event_number(entry);
164 
165  /* process ID from pythia */
166  evt->set_signal_process_id(GenEvent->GetProcess());
167 
168  /* Set the PDF information */
169  HepMC::PdfInfo pdfinfo;
170  pdfinfo.set_x1(1);
171  pdfinfo.set_x2(GenEvent->GetX());
172  pdfinfo.set_scalePDF(GenEvent->GetQ2());
173  evt->set_pdf_info(pdfinfo);
174 
175  /* Loop over all particles for this event in input file and fill
176  * vector with HepMC particles */
177  vector<HepMC::GenParticle *> hepmc_particles;
178  vector<unsigned> origin_index;
179 
180  /* save pointers to beam particles */
181  HepMC::GenParticle *hepmc_beam1 = nullptr;
182  HepMC::GenParticle *hepmc_beam2 = nullptr;
183 
184  for (unsigned ii = 0; ii < GenEvent->GetNTracks(); ii++)
185  {
186  /* Get particle / track from event records.
187  * Class documentation for erhic::VirtualParticle at
188  * http://www.star.bnl.gov/~tpb/eic-smear/classerhic_1_1_virtual_particle.html */
189  erhic::ParticleMC *track_ii = GenEvent->GetTrack(ii);
190 
191  if (Verbosity() > 1)
192  {
193  cout << __PRETTY_FUNCTION__ << " : " << __LINE__ << endl;
194  cout << "\ttrack " << ii << endl;
195  cout << "\t4mom\t" << track_ii->GetPx() << "\t" << track_ii->GetPy() << "\t" << track_ii->GetPz() << "\t" << track_ii->GetE() << endl;
196  cout << "\tstatus= " << track_ii->GetStatus() << "\tindex= " << track_ii->GetIndex() << "\tmass= " << track_ii->GetM() << endl;
197  }
198 
199  /* Create HepMC particle record */
200  HepMC::GenParticle *hepmcpart = new HepMC::GenParticle(HepMC::FourVector(track_ii->GetPx(),
201  track_ii->GetPy(),
202  track_ii->GetPz(),
203  track_ii->GetE()),
204  track_ii->Id());
205 
206  /* translate eic-smear status codes to HepMC status codes */
207  switch (track_ii->GetStatus())
208  {
209  case 1:
210  hepmcpart->set_status(1);
211  break; // 'stable particle'
212 
213  case 21:
214  hepmcpart->set_status(3);
215  break; // 'documentation line'
216 
217  default:
218  hepmcpart->set_status(0);
219  break; // 'null entry'
220  }
221 
222  /* assume the first two particles are the beam particles (which getsHepMC status 4)*/
223  if (ii < 2)
224  hepmcpart->set_status(4);
225 
226  /* add particle information */
227  hepmcpart->setGeneratedMass(track_ii->GetM());
228 
229  /* append particle to vector */
230  hepmc_particles.push_back(hepmcpart);
231  origin_index.push_back(track_ii->GetIndex());
232 
233  /* if first particle, call this the first beam particle */
234  if (ii == 0)
235  hepmc_beam1 = hepmcpart;
236 
237  /* if second particle, call this the second beam particle */
238  if (ii == 1)
239  hepmc_beam2 = hepmcpart;
240  }
241 
242  /* Check if hepmc_particles and origin_index vectors are the same size */
243  if (hepmc_particles.size() != origin_index.size())
244  {
245  cout << "ReadEICFiles::process_event - Lengths of HepMC particles and Origin index vectors do not match!" << endl;
246 
247  delete evt;
249  }
250 
251  /* add HepMC particles to Hep MC vertices; skip first two particles
252  * in loop, assuming that they are the beam particles */
253  vector<HepMC::GenVertex *> hepmc_vertices;
254 
255  for (unsigned p = 2; p < hepmc_particles.size(); p++)
256  {
257  HepMC::GenParticle *pp = hepmc_particles.at(p);
258 
259  /* continue if vertices for particle are already set */
260  if (pp->production_vertex() && pp->end_vertex())
261  continue;
262 
263  /* access mother particle vertex */
264  erhic::ParticleMC *track_pp = GenEvent->GetTrack(p);
265 
266  unsigned parent_index = track_pp->GetParentIndex();
267 
268  HepMC::GenParticle *pmother = nullptr;
269  for (unsigned m = 0; m < hepmc_particles.size(); m++)
270  {
271  if (origin_index.at(m) == parent_index)
272  {
273  pmother = hepmc_particles.at(m);
274  break;
275  }
276  }
277 
278  /* if mother does not exist: create new vertex and add this particle as outgoing to vertex */
279  if (!pmother)
280  {
281  HepMC::GenVertex *hepmcvtx = new HepMC::GenVertex(HepMC::FourVector(track_pp->GetVertex()[0],
282  track_pp->GetVertex()[1],
283  track_pp->GetVertex()[2],
284  0));
285  hepmc_vertices.push_back(hepmcvtx);
286  hepmcvtx->add_particle_out(pp);
287  continue;
288  }
289  /* if mother exists and has end vertex: add this particle as outgoing to the mother's end vertex */
290  else if (pmother->end_vertex())
291  {
292  pmother->end_vertex()->add_particle_out(pp);
293  }
294  /* if mother exists and has no end vertex: create new vertex */
295  else
296  {
297  HepMC::GenVertex *hepmcvtx = new HepMC::GenVertex(HepMC::FourVector(track_pp->GetVertex()[0],
298  track_pp->GetVertex()[1],
299  track_pp->GetVertex()[2],
300  0));
301  hepmc_vertices.push_back(hepmcvtx);
302  hepmcvtx->add_particle_in(pmother);
303  pmother->end_vertex()->add_particle_out(pp);
304  }
305  }
306 
307  /* Add end vertex to beam particles if they don't have one yet */
308  for (unsigned p = 0; p < 2; p++)
309  {
310  HepMC::GenParticle *pp = hepmc_particles.at(p);
311 
312  if (!pp->end_vertex())
313  {
314  /* create collision vertex */
315  HepMC::GenVertex *hepmcvtx = new HepMC::GenVertex(HepMC::FourVector(0,
316  0,
317  0,
318  0));
319  hepmc_vertices.push_back(hepmcvtx);
320  hepmcvtx->add_particle_in(pp);
321  }
322  }
323 
324  /* Check that all particles (except beam particles) have a production vertex */
325  for (unsigned p = 2; p < hepmc_particles.size(); p++)
326  {
327  if (!hepmc_particles.at(p)->production_vertex())
328  {
329  cout << "ReadEICFiles::process_event - Missing production vertex for one or more non-beam particles!" << endl;
330  delete evt;
332  }
333  }
334 
335  /* Add HepMC vertices to event */
336  for (unsigned v = 0; v < hepmc_vertices.size(); v++)
337  {
338  evt->add_vertex(hepmc_vertices.at(v));
339  }
340 
341  /* set beam particles */
342  evt->set_beam_particles(hepmc_beam1, hepmc_beam2);
343 
344  /* pass HepMC to PHNode*/
346  if (Verbosity() > 1)
347  {
348  cout << __PRETTY_FUNCTION__ << " : " << __LINE__ << endl;
349  evt->print();
350  cout << endl << endl;
351  }
352 
353  if (!success)
354  {
355  cout << "ReadEICFiles::process_event - Failed to add event to HepMC record!" << endl;
357  }
358  /* Count up number of 'used' events from input file */
359  entry++;
360 
361  /* Done */
362  return 0;
363 }
364 
366 {
368 
369  PHNodeIterator iter(topNode);
370  // Looking for the DST node
371  PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
372  if (!dstNode)
373  {
374  cout << PHWHERE << "DST Node missing, doing nothing." << endl;
376  }
377  EicEventHeader *evthead = findNode::getClass<EicEventHeader>(topNode, "EicEventHeader");
378  if (!evthead)
379  {
380  evthead = new EicEventHeaderv1();
381  PHIODataNode<PHObject> *HeaderNode = new PHIODataNode<PHObject>(evthead, "EicEventHeader", "PHObject");
382  dstNode->addNode(HeaderNode);
383  }
385 }