hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
MCParticleProcessor.cxx
Go to the documentation of this file.
1
9
10MCParticleProcessor::MCParticleProcessor(const std::string& name, Process& process)
11 : Processor(name, process) {
12 }
13
16
18
19 std::cout << "Configuring MCParticleProcessor" << std::endl;
20 try
21 {
22 debug_ = parameters.getInteger("debug", debug_ );
23 mcPartCollLcio_ = parameters.getString("mcPartCollLcio", mcPartCollLcio_);
24 mcPartCollRoot_ = parameters.getString("mcPartCollRoot", mcPartCollRoot_);
25 }
26 catch (std::runtime_error& error)
27 {
28 std::cout << error.what() << std::endl;
29 }
30
31
32}
33
35 // Add branch to tree
36 tree->Branch(mcPartCollRoot_.c_str(),&mc_particles_);
37
38}
39
41
42 Event* event = static_cast<Event*> (ievent);
43
44 // Get the collection from the event
45 EVENT::LCCollection* lc_particles{nullptr};
46 try
47 {
48 lc_particles = event->getLCCollection(mcPartCollLcio_.c_str());
49 }
50 catch (EVENT::DataNotAvailableException e)
51 {
52 std::cout << e.what() << std::endl;
53 }
54
55
56 //Clean up
57 if (mc_particles_.size() > 0 )
58 {
59 for (std::vector<MCParticle*>::iterator it = mc_particles_.begin(); it != mc_particles_.end(); ++it)
60 {
61 delete *it;
62 }
63 mc_particles_.clear();
64 }
65
66
67 // Loop through all of the particles in the event
68 for (int iparticle = 0; iparticle < lc_particles->getNumberOfElements(); ++iparticle) {
69
70 // Get a particle from the LCEvent
71 IMPL::MCParticleImpl* lc_particle
72 = static_cast<IMPL::MCParticleImpl*>(lc_particles->getElementAt(iparticle));
73
74 // Make an MCParticle to build and add to vector
75 MCParticle* particle = new MCParticle();
76
77 // Set the charge of the HpsMCParticle
78 particle->setCharge(lc_particle->getCharge());
79
80 // Set the HpsMCParticle type
81 particle->setTime(lc_particle->getTime());
82
83 // Set the energy of the HpsMCParticle
84 particle->setEnergy(lc_particle->getEnergy());
85
86 // Set the momentum of the HpsMCParticle
87 particle->setMomentum(lc_particle->getMomentum());
88
89 // Set the momentum of HpsMCParticle at Endpoint
90 particle->setEndpointMomentum(lc_particle->getMomentumAtEndpoint());
91
92 // Set the mass of the HpsMCParticle
93 particle->setMass(lc_particle->getMass());
94
95 // Set the PDG of the particle
96 particle->setPDG(lc_particle->getPDG());
97
98 // Set the LCIO id of the particle
99 particle->setID(lc_particle->id());
100
101 // Set the PDG of the particle
102 std::vector<EVENT::MCParticle*> parentVec = lc_particle->getParents();
103 if(parentVec.size() > 0) particle->setMomPDG(parentVec.at(parentVec.size()-1)->getPDG());
104
105 // Set the generator status of the particle
106 particle->setGenStatus(lc_particle->getGeneratorStatus());
107
108 // Set the generator status of the particle
109 particle->setSimStatus(lc_particle->getSimulatorStatus());
110
111 // Loop through all of the tracks associated with the particle
112 // and add references to the MCParticle object.
113 /*for (auto const &lc_track : lc_particle->getTracks()) {
114
115 TClonesArray* tracks = event->getCollection(Collections::GBL_TRACKS);
116
117 // Loop through all of the tracks in the HpsEvent and find the one
118 // that matches the track associated with the particle
119 for (int itrack = 0; itrack < tracks->GetEntriesFast(); ++itrack) {
120 Track* track = static_cast<Track*>(tracks->At(itrack));
121
122 // Use the track chi^2 to find the match
123 // TODO: Verify that the chi^2 is unique enough to find the match
124 if (lc_track->getChi2() == track->getChi2()) {
125
126 // Add a reference to the track
127 particle->addTrack(track);
128
129 // If the particle is a final state particle, add a
130 // reference from the corresponding track to the particle
131 if ((collections.first.compare(Collections::FINAL_STATE_PARTICLES) == 0)
132 || (collections.first.compare(Collections::OTHER_ELECTRONS) == 0) ) {
133 track->setMCParticle(particle);
134 track->setMomentum(particle->getMomentum());
135 track->setCharge(particle->getCharge());
136 }
137 break;
138 }
139 }
140
141 }*/
142
143 // Only add vertex information if the particle is not a final state particle
144 //if ((collections.first.compare(Collections::FINAL_STATE_PARTICLES) == 0) ||
145 // (collections.first.compare(Collections::OTHER_ELECTRONS) == 0)) {
146 // // Set the PDG ID of the particle
147 // particle->setPDG(lc_particle->getParticleIDUsed()->getPDG());
148 // continue;
149 //}
150
151 // Set the vertex position of the particle
152 particle->setVertexPosition(lc_particle->getVertex());
153
154 // Set the vertex position of the particle
155 particle->setEndPoint(lc_particle->getEndpoint());
156
157 mc_particles_.push_back(particle);
158 // If the particle has daughter particles, add the daughters to the
159 // MCParticle.
160 //
161
162 // Loop through all of the daughter particles associated with the particle
163 /*for (auto const &daughter : lc_particle->getParticles()) {
164
165 // Loop through all of the final state particles in the event
166 // and find the one that matches the daughters associated with
167 // the particles.
168 for (int iparticle = 0;
169 iparticle < fs_particles->GetEntriesFast(); ++iparticle) {
170
171 MCParticle* dparticle = static_cast<MCParticle*>(fs_particles->At(iparticle));
172
173 // Try to find the match between a final state particle
174 // and ReconstructedParticle daughter. For now, use the
175 // momentum as the matching criterion.
176 // TODO: Verify that the track momentum is always unique in an event.
177 if ((dparticle->getMomentum()[0] == lc_particle->getMomentum()[0])
178 && (dparticle->getMomentum()[1] == lc_particle->getMomentum()[1])
179 && (dparticle->getMomentum()[2] == lc_particle->getMomentum()[2])) {
180
181 particle->addDaughter(dparticle);
182
183 if (dparticle->getTracks()->GetEntriesFast() != 0)
184 particle->addTrack(dparticle->getTracks()->At(0));
185
186 if (dparticle->getClusters()->GetEntriesFast() != 0)
187 particle->addCluster(dparticle->getClusters()->At(0));
188
189 break;
190 }
191 }
192 }*/
193 }
194
195 return true;
196}
197
200
Processor used to translate LCIO MCParticles to DST MCParticle objects.
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
Definition Event.h:35
EVENT::LCCollection * getLCCollection(std::string name)
Definition Event.h:102
Definition IEvent.h:7
Processor used to translate LCIO MCParticles to DST MCParticle objects. more details.
virtual void configure(const ParameterSet &parameters)
Callback for the Processor to configure itself from the given set of parameters.
virtual void finalize()
Callback for the Processor to take any necessary action when the processing of events finishes.
virtual void initialize(TTree *tree)
Callback for the Processor to take any necessary action when the processing of events starts.
std::vector< MCParticle * > mc_particles_
std::string mcPartCollLcio_
description
MCParticleProcessor(const std::string &name, Process &process)
Class constructor.
std::string mcPartCollRoot_
description
void setVertexPosition(const double *vtx_pos)
void setEndpointMomentum(const double *momentum_ep)
void setPDG(const int pdg)
Definition MCParticle.h:60
void setMass(const double mass)
Definition MCParticle.h:102
void setID(const int id)
Definition MCParticle.h:67
void setGenStatus(const int gen)
Definition MCParticle.h:81
void setMomPDG(const int momPDG)
Definition MCParticle.h:74
void setCharge(const int charge)
Definition MCParticle.h:53
void setMomentum(const double *momentum)
void setEnergy(const double energy)
Definition MCParticle.h:95
void setTime(const double time)
Definition MCParticle.h:109
void setEndPoint(const double *ep_pos)
void setSimStatus(const int sim)
Definition MCParticle.h:88
description
Base class for all event processing components.
Definition Processor.h:34
virtual bool process()
Process the histograms and generate analysis output.
Definition Processor.h:95