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
VertexProcessor.cxx
Go to the documentation of this file.
1
6#include "VertexProcessor.h"
7#include "utilities.h"
8
9VertexProcessor::VertexProcessor(const std::string& name, Process& process)
10 : Processor(name, process) {
11 }
12
15
17
18 std::cout << "Configuring VertexProcessor" << std::endl;
19 try
20 {
21 debug_ = parameters.getInteger("debug", debug_);
22 vtxCollLcio_ = parameters.getString("vtxCollLcio", vtxCollLcio_);
23 vtxCollRoot_ = parameters.getString("vtxCollRoot", vtxCollRoot_);
24 partCollRoot_ = parameters.getString("partCollRoot", partCollRoot_);
25 kinkRelCollLcio_ = parameters.getString("kinkRelCollLcio", kinkRelCollLcio_);
26 trkRelCollLcio_ = parameters.getString("trkRelCollLcio", trkRelCollLcio_);
27 trackStateLocation_= parameters.getString("trackStateLocation", trackStateLocation_);
28 hitFitsCollLcio_ = parameters.getString("hitFitsCollLcio", hitFitsCollLcio_);
29 trkhitCollRoot_ = parameters.getString("trkhitCollRoot",trkhitCollRoot_);
30 rawhitCollRoot_ = parameters.getString("rawhitCollRoot",rawhitCollRoot_);
31 mcPartRelLcio_ = parameters.getString("mcPartRelLcio", mcPartRelLcio_);
32 useTrackerHits_ = parameters.getInteger("useTrackerHits",useTrackerHits_);
33 bfield_ = parameters.getDouble("bfield",bfield_);
34 }
35 catch (std::runtime_error& error)
36 {
37 std::cout << error.what() << std::endl;
38 }
39}
40
41void VertexProcessor::initialize(TTree* tree) {
42 // Add branches to tree
43 tree->Branch(vtxCollRoot_.c_str(), &vtxs_);
44 tree->Branch(partCollRoot_.c_str(), &parts_);
45 if (!trkhitCollRoot_.empty())
46 tree->Branch(trkhitCollRoot_.c_str(), &hits_);
47
48 if (!rawhitCollRoot_.empty())
49 tree->Branch(rawhitCollRoot_.c_str(), &rawhits_);
50}
51
53
54 if (debug_ > 0) std::cout << "VertexProcessor: Clear output vector" << std::endl;
55
56 if (hits_.size() > 0) {
57 for (std::vector<TrackerHit *>::iterator it = hits_.begin(); it != hits_.end(); ++it) {
58 delete *it;
59 }
60 hits_.clear();
61 }
62
63 if (rawhits_.size() > 0) {
64 for (std::vector<RawSvtHit *>::iterator it = rawhits_.begin(); it != rawhits_.end(); ++it) {
65 delete *it;
66 }
67 rawhits_.clear();
68 }
69
70 for(int i = 0; i < vtxs_.size(); i++) delete vtxs_.at(i);
71 vtxs_.clear();
72 for(int i = 0; i < parts_.size(); i++) delete parts_.at(i);
73 parts_.clear();
74
75 Event* event = static_cast<Event*> (ievent);
76 UTIL::LCRelationNavigator* mcPartRel_nav;
77
78 // Get the collection of vertices from the LCIO event. If no such collection
79 // exist, a DataNotAvailableException is thrown
80 if (debug_ > 0) std::cout << "VertexProcessor: Get LCIO Collection " << vtxCollLcio_ << std::endl;
81 EVENT::LCCollection* lc_vtxs = nullptr;
82 try
83 {
84 lc_vtxs = event->getLCCollection(vtxCollLcio_.c_str());
85 }
86 catch (EVENT::DataNotAvailableException e)
87 {
88 std::cout << e.what() << std::endl;
89 return false;
90 }
91
92 // Get the collection of LCRelations between GBL tracks and kink data and track data variables.
93 EVENT::LCCollection* gbl_kink_data{nullptr};
94 EVENT::LCCollection* track_data{nullptr};
95
96 try
97 {
98 if (!kinkRelCollLcio_.empty())
99 gbl_kink_data = static_cast<EVENT::LCCollection*>(event->getLCCollection(kinkRelCollLcio_.c_str()));
100 if (!trkRelCollLcio_.empty())
101 track_data = static_cast<EVENT::LCCollection*>(event->getLCCollection(trkRelCollLcio_.c_str()));
102 }
103 catch (EVENT::DataNotAvailableException e)
104 {
105 std::cout << e.what() << std::endl;
106 if (!gbl_kink_data)
107 std::cout<<"Failed retrieving " << kinkRelCollLcio_ <<std::endl;
108 if (!track_data)
109 std::cout<<"Failed retrieving " << trkRelCollLcio_ <<std::endl;
110
111 }
112
113 bool rotateHits = true;
114 int hitType = 0;
115 EVENT::LCCollection* raw_svt_hit_fits = nullptr;
116
117 auto evColls = event->getLCEvent()->getCollectionNames();
118 auto it = std::find (evColls->begin(), evColls->end(), hitFitsCollLcio_.c_str());
119 bool hasFits = true;
120 if(it == evColls->end()) hasFits = false;
121 if(hasFits)
122 {
123 raw_svt_hit_fits = event->getLCCollection(hitFitsCollLcio_.c_str());
124 }
125 //Check to see if MC Particles are in the file
126 it = std::find (evColls->begin(), evColls->end(), mcPartRelLcio_.c_str());
127 bool hasMCParts = true;
128 EVENT::LCCollection* mcPartRel;
129 if(it == evColls->end()) hasMCParts = false;
130 if(hasMCParts)
131 {
132 mcPartRel = event->getLCCollection(mcPartRelLcio_.c_str());
133 // Heap an LCRelation navigator which will allow faster access
134 mcPartRel_nav = new UTIL::LCRelationNavigator(mcPartRel);
135
136 }
137
138 if (debug_ > 0) std::cout << "VertexProcessor: Converting Verteces" << std::endl;
139 for (int ivtx = 0 ; ivtx < lc_vtxs->getNumberOfElements(); ++ivtx)
140 {
141 if (debug_ > 0) std::cout << "VertexProcessor: Converting Vertex " << ivtx << std::endl;
142 EVENT::Vertex * lc_vtx{nullptr};
143 lc_vtx = static_cast<EVENT::Vertex*>(lc_vtxs->getElementAt(ivtx));
144
145 if (debug_ > 0) std::cout << "VertexProcessor: Build Vertex" << std::endl;
146 Vertex* vtx = utils::buildVertex(lc_vtx);
147
148 if (debug_ > 0) std::cout << "VertexProcessor: Get Particles" << std::endl;
149 std::vector<EVENT::ReconstructedParticle*> lc_parts = lc_vtx->getAssociatedParticle()->getParticles();
150 for(auto lc_part : lc_parts)
151 {
152 if (debug_ > 0) std::cout << "VertexProcessor: Build particle" << std::endl;
153 Particle * part = utils::buildParticle(lc_part,trackStateLocation_, gbl_kink_data, track_data);
154 if (lc_part->getTracks().size()>0){
155 EVENT::Track* lc_track = static_cast<EVENT::Track*>(lc_part->getTracks()[0]);
156 Track* track = utils::buildTrack(lc_track,trackStateLocation_,gbl_kink_data,track_data);
157 int nHits = 0;
158 if (bfield_ > 0.0) track->setMomentum(bfield_);
159 if (track->isKalmanTrack()) hitType = 1; //SiClusters
160
161 if(!useTrackerHits_){
162 auto hitPattern = lc_track->getSubdetectorHitNumbers();
163 for( int h=0; h<hitPattern.size(); h++){
164 if( hitPattern[h]>0 ){
165 track->addHitLayer(h);
166 nHits++;
167 }
168 }
169 if(nHits==0){
170 std::cout << "[VertexProcessor::process] WARNING Track with no hits -- likely getSubdetectorHitNumbers isn't filled!" << std::endl;
171 std::cout << "[VertexProcessor::process] WARNING If hit collections are saved, you may want to turn on useTrackerHits" << std::endl;
172 }
173 track->setTrackerHitCount(nHits);
174 }
175 else{
176 EVENT::TrackerHitVec lc_tracker_hits = lc_track->getTrackerHits();
177 for (auto lc_tracker_hit : lc_tracker_hits) {
178 TrackerHit* tracker_hit = utils::buildTrackerHit(static_cast<IMPL::TrackerHitImpl*>(lc_tracker_hit),rotateHits,hitType);
179 std::vector<RawSvtHit*> rawSvthitsOn3d;
180 utils::addRawInfoTo3dHit(tracker_hit,static_cast<IMPL::TrackerHitImpl*>(lc_tracker_hit),
181 raw_svt_hit_fits,&rawSvthitsOn3d,hitType);
182
183 int hitLayer = tracker_hit->getLayer();
184 nHits++;
185 EVENT::LCObjectVec rawHits = lc_tracker_hit->getRawHits();
186 for(int irawhit = 0; irawhit < rawHits.size(); ++irawhit){
187 IMPL::TrackerHitImpl* rawhit = static_cast<IMPL::TrackerHitImpl*>(rawHits.at(irawhit));
188 if(debug_ > 0)
189 std::cout << "rawhit on track has lcio id: " << rawhit->id() << std::endl;
190
191 if (hasMCParts)
192 {
193 // Get the list of fit params associated with the raw tracker hit
194 EVENT::LCObjectVec lc_simtrackerhits = mcPartRel_nav->getRelatedToObjects(rawhit);
195
196 //Loop over SimTrackerHits to get MCParticles
197 for(int isimhit = 0; isimhit < lc_simtrackerhits.size(); isimhit++){
198 IMPL::SimTrackerHitImpl* lc_simhit = static_cast<IMPL::SimTrackerHitImpl*>(lc_simtrackerhits.at(isimhit));
199 IMPL::MCParticleImpl* lc_mcp = static_cast<IMPL::MCParticleImpl*>(lc_simhit->getMCParticle());
200 if(lc_mcp == nullptr)
201 std::cout << "mcp is null" << std::endl;
202 track->addMcpHit(hitLayer, lc_mcp->id());
203 if(debug_ > 0) {
204 std::cout << "simtrackerhit lcio id: " << lc_simhit->id() << std::endl;
205 std::cout << "mcp lcio id: " << lc_mcp->id() << std::endl;
206 }
207 }
208 }
209 } // Loop over raw tracker hits
210
211 track->addHitLayer(hitLayer);
212 hits_.push_back(tracker_hit);
213 for (std::vector<RawSvtHit *>::iterator it = rawSvthitsOn3d.begin(); it != rawSvthitsOn3d.end(); ++it) {
214 delete *it;
215 }
216 rawSvthitsOn3d.clear();
217
218 } // Loop over hits on track
219
220 track->setTrackerHitCount(nHits);
221
222 } // !useTrackerHits
223
224 part->setTrack(track);
225 delete track;
226
227 }
228 // =========================================
229 if (debug_ > 0) std::cout << "VertexProcessor: Add particle" << std::endl;
230 parts_.push_back(part);
231 vtx->addParticle(part);
232 }
233
234 if (debug_ > 0) std::cout << "VertexProcessor: Add Vertex" << std::endl;
235 vtxs_.push_back(vtx);
236 }
237
238 if(hasMCParts) delete mcPartRel_nav;
239 if (debug_ > 0) std::cout << "VertexProcessor: End process" << std::endl;
240 return true;
241}
242
245
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
Class used to convert LCIO Vertex collections into ROOT collections.
Definition Event.h:35
EVENT::LCCollection * getLCCollection(std::string name)
Definition Event.h:102
Definition IEvent.h:7
description
void setTrack(Track *track)
Definition Particle.h:42
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
Definition Track.h:32
int getLayer() const
Definition TrackerHit.h:109
Insert description here. more details.
std::string mcPartRelLcio_
std::string partCollRoot_
description
std::string trkRelCollLcio_
description
std::string trackStateLocation_
select track state for tracks DEFAULT AtIP
int debug_
Debug Level.
std::string trkhitCollRoot_
double bfield_
magnetic field
int useTrackerHits_
Load hit collections, otherwise get from getSubdetectorHitNumbers.
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::string hitFitsCollLcio_
std::vector< RawSvtHit * > rawhits_
std::vector< Vertex * > vtxs_
description
std::vector< TrackerHit * > hits_
std::vector< Particle * > parts_
description
std::string vtxCollRoot_
description
std::string kinkRelCollLcio_
description
std::string vtxCollLcio_
description
std::string rawhitCollRoot_
VertexProcessor(const std::string &name, Process &process)
Class constructor.
bool addRawInfoTo3dHit(TrackerHit *tracker_hit, IMPL::TrackerHitImpl *lc_tracker_hit, EVENT::LCCollection *raw_svt_fits, std::vector< RawSvtHit * > *rawHits=nullptr, int type=0, bool storeRawHit=true)
description
Track * buildTrack(EVENT::Track *lc_track, std::string trackstate_location, EVENT::LCCollection *gbl_kink_data, EVENT::LCCollection *track_data)
description
Vertex * buildVertex(EVENT::Vertex *lc_vertex)
description
Definition utilities.cxx:31
TrackerHit * buildTrackerHit(IMPL::TrackerHitImpl *lc_trackerHit, bool rotate=true, int type=0)
description
Particle * buildParticle(EVENT::ReconstructedParticle *lc_particle, std::string trackstate_location, EVENT::LCCollection *gbl_kink_data, EVENT::LCCollection *track_data)
description
Definition utilities.cxx:56