Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PHHepMCGenEvent.cc
Go to the documentation of this file.
1 #include "PHHepMCGenEvent.h"
2 #include <HepMC/GenEvent.h>
3 
4 #include <TBuffer.h>
5 #include <TClass.h>
6 
7 #include <RVersion.h> // root version
8 
9 #include <boost/foreach.hpp>
10 
11 #include <algorithm>
12 #include <cassert>
13 #include <cstdlib>
14 #include <iomanip>
15 #include <sstream>
16 #include <stdexcept>
17 #include <vector>
18 
19 //ClassImp(PHHepMCGenEvent)
20 
21 using namespace std;
22 
24  : _embedding_id(0)
25  , _isSimulated(false)
26  , _collisionVertex(0, 0, 0, 0)
27  , _theEvt(nullptr)
28 {
29 }
30 
32  : _embedding_id(event.get_embedding_id())
33  , _isSimulated(event.is_simulated())
34  , _collisionVertex(event.get_collision_vertex())
35  , _theEvt(nullptr)
36 {
37  _theEvt = new HepMC::GenEvent(*event.getEvent());
38  return;
39 }
40 
42 {
43  if (&event == this) return *this;
44 
45  Reset();
46 
47  _embedding_id = event.get_embedding_id();
48  _isSimulated = event.is_simulated();
49 
50  const HepMC::GenEvent* hepmc = event.getEvent();
51  _theEvt = new HepMC::GenEvent(*(hepmc));
52 
53  return *this;
54 }
55 
57 {
58  Reset();
59 }
60 
62 {
63  _embedding_id = 0;
64  _isSimulated = false;
65  _collisionVertex.set(0, 0, 0, 0);
66  if (_theEvt)
67  {
68  delete _theEvt;
69  _theEvt = NULL;
70  }
71 }
72 
73 HepMC::GenEvent* PHHepMCGenEvent::getEvent()
74 {
75  return _theEvt;
76 }
77 
78 const HepMC::GenEvent* PHHepMCGenEvent::getEvent() const
79 {
80  return _theEvt;
81 }
82 
83 bool PHHepMCGenEvent::addEvent(HepMC::GenEvent* evt)
84 {
85  if (_theEvt) delete _theEvt;
86 
87  _theEvt = evt;
88  if (!_theEvt) return false;
89  return true;
90 }
91 
92 bool PHHepMCGenEvent::swapEvent(HepMC::GenEvent*& evt)
93 {
94  swap(_theEvt, evt);
95 
96  if (!_theEvt) return false;
97  return true;
98 }
99 
100 bool PHHepMCGenEvent::addEvent(HepMC::GenEvent& evt)
101 {
102  return addEvent(new HepMC::GenEvent(evt));
103 }
104 
106 {
107  if (_theEvt) _theEvt->clear();
108 }
109 
110 void PHHepMCGenEvent::moveVertex(double x, double y, double z, double t)
111 {
112  _collisionVertex.setX(_collisionVertex.x() + x);
113  _collisionVertex.setY(_collisionVertex.y() + y);
114  _collisionVertex.setZ(_collisionVertex.z() + z);
115  _collisionVertex.setT(_collisionVertex.t() + t);
116 }
117 
118 int PHHepMCGenEvent::size(void) const
119 {
120  if (_theEvt)
121  return _theEvt->particles_size();
122  else
123  return 0;
124 }
125 
127 {
128  if (_theEvt)
129  return _theEvt->vertices_size();
130  else
131  return 0;
132 }
133 
134 //_____________________________________________________________________________
135 void PHHepMCGenEvent::identify(std::ostream& os) const
136 {
137  os << "identify yourself: PHHepMCGenEvent Object, ";
138  os << " embedding_id = " << _embedding_id;
139  os << " isSimulated = " << _isSimulated;
140  os << " collisionVertex = (" << _collisionVertex.x()<<","<< _collisionVertex.y()<<","<< _collisionVertex.z()<<") cm, "<< _collisionVertex.t()<<" ns";
141  os << ", No of Particles: " << size();
142  os << ", No of Vertices: " << vertexSize();
143  os << endl;
144  return;
145 }
146 
147 void PHHepMCGenEvent::print(std::ostream& out) const
148 {
149  identify(out);
150 }
151 
153 {
154  _theEvt->print();
155 }
bool swapEvent(HepMC::GenEvent *&evt)
bool addEvent(HepMC::GenEvent *evt)
host an HepMC event
HepMC::GenEvent * _theEvt
The HEP MC record from event generator. Note the units are recorded in GenEvent.
virtual void print(std::ostream &os=std::cout) const
virtual void moveVertex(double x, double y, double z, double t=0)
move the collision vertex position in the Hall coordinate system, use PHENIX units of cm...
virtual HepMC::GenEvent * getEvent()
#define NULL
Definition: Pdb.h:9
virtual int size(void) const
void swap(array< T, N > &x, array< T, N > &y)
Definition: array.hpp:366
virtual ~PHHepMCGenEvent()
virtual void Reset()
Clear Event.
bool _isSimulated
whether this event has been processed in Geant4 simulation
virtual int vertexSize(void) const
int _embedding_id
Embedding ID for this generated event positive ID is the embedded event of interest, e.g. jetty event from pythia negative IDs are backgrounds, .e.g out of time pile up collisions Usually, ID = 0 means the primary Au+Au collision background.
HepMC::FourVector _collisionVertex
collision vertex position in the Hall coordinate system, use PHENIX units of cm, ns ...
PHHepMCGenEvent & operator=(const PHHepMCGenEvent &event)
virtual void identify(std::ostream &os=std::cout) const