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
SvtDataProcessor.cxx
Go to the documentation of this file.
1
2#include "SvtDataProcessor.h"
3
4SvtDataProcessor::SvtDataProcessor(const std::string& name, Process& process)
5 : Processor(name, process) {
6}
7
10
12 // Add branches to tree
13 tree->Branch(Collections::GBL_TRACKS, &tracks_);
14 tree->Branch(Collections::TRACKER_HITS, &hits_);
15}
16
18
19 for(int i = 0; i < tracks_.size(); i++) delete tracks_.at(i);
20 for(int i = 0; i < hits_.size(); i++) delete hits_.at(i);
21 tracks_.clear();
22 hits_.clear();
23
24 Event* event = static_cast<Event*> (ievent);
25
26 // Get the collection of 3D hits from the LCIO event. If no such collection
27 // exist, a DataNotAvailableException is thrown
28 EVENT::LCCollection* tracker_hits = event->getLCCollection(Collections::TRACKER_HITS);
29
30 // Create a map from an LCIO TrackerHit to a SvtHit. This will be used when
31 // assigning references to a track
32 // TODO: Use an unordered map for faster access
33 std::map<EVENT::TrackerHit*, TrackerHit*> hit_map;
34
35 // Loop over all of the 3D hits in the LCIO event and add them to the
36 // HPS event
37 for (int ihit = 0; ihit < tracker_hits->getNumberOfElements(); ++ihit) {
38
39 // Get a 3D hit from the list of hits
40 IMPL::TrackerHitImpl* lc_tracker_hit = static_cast<IMPL::TrackerHitImpl*>(tracker_hits->getElementAt(ihit));
41
42 // Add a tracker hit to the event
43 TrackerHit* tracker_hit = new TrackerHit();
44
45 // Rotate the position of the LCIO TrackerHit and set the position of
46 // the TrackerHit
47 double hit_position[3] = {
48 lc_tracker_hit->getPosition()[1],
49 lc_tracker_hit->getPosition()[2],
50 lc_tracker_hit->getPosition()[0]
51 };
52 tracker_hit->setPosition(hit_position);
53
54 // Set the covariance matrix of the SvtHit
55 tracker_hit->setCovarianceMatrix(lc_tracker_hit->getCovMatrix());
56
57 // Set the time of the SvtHit
58 tracker_hit->setTime(lc_tracker_hit->getTime());
59
60 // Set the charge of the SvtHit
61 tracker_hit->setCharge(lc_tracker_hit->getEDep());
62
63 // Map the TrackerHit object to the corresponding SvtHit object. This
64 // will be used later when setting references for hits on tracks.
65 hit_map[lc_tracker_hit] = tracker_hit;
66 hits_.push_back(tracker_hit);
67
68 }
69
70 // Add the hit collection to the event
71 //event->addCollection(Collections::TRACKER_HITS, hits_);
72
73 // Get all track collections from the event
74 EVENT::LCCollection* tracks = event->getLCCollection(Collections::GBL_TRACKS);
75
76
77 // Loop over all the LCIO Tracks and add them to the HPS event.
78 for (int itrack = 0; itrack < tracks->getNumberOfElements(); ++itrack) {
79
80 // Get a LCIO Track from the LCIO event
81 EVENT::Track* lc_track = static_cast<EVENT::Track*>(tracks->getElementAt(itrack));
82
83 // Add a track to the event
84 Track* track = new Track();
85
86 // Set the track parameters
87 track->setTrackParameters(lc_track->getD0(),
88 lc_track->getPhi(),
89 lc_track->getOmega(),
90 lc_track->getTanLambda(),
91 lc_track->getZ0());
92
93 // Set the track type
94 track->setType(lc_track->getType());
95
96 // Set the track fit chi^2
97 track->setChi2(lc_track->getChi2());
98
99 // Set the position of the extrapolated track at the ECal face. The
100 // extrapolation uses the full 3D field map.
101 const EVENT::TrackState* track_state
102 = lc_track->getTrackState(EVENT::TrackState::AtCalorimeter);
103 double position_at_ecal[3] = {
104 track_state->getReferencePoint()[1],
105 track_state->getReferencePoint()[2],
106 track_state->getReferencePoint()[0]
107 };
108 track->setPositionAtEcal(position_at_ecal);
109
110 // Get the collection of LCRelations between GBL kink data variables
111 // (GBLKinkData) and the corresponding track.
112 EVENT::LCCollection* gbl_kink_data =
113 static_cast<EVENT::LCCollection*>(event->getLCCollection(Collections::KINK_DATA_REL));
114
115 // Instantiate an LCRelation navigator which will allow faster access
116 // to GBLKinkData object
117 UTIL::LCRelationNavigator* gbl_kink_data_nav
118 = new UTIL::LCRelationNavigator(gbl_kink_data);
119
120 // Get the list of GBLKinkData associated with the LCIO Track
121 EVENT::LCObjectVec gbl_kink_data_list
122 = gbl_kink_data_nav->getRelatedFromObjects(lc_track);
123
124 // The container of GBLKinkData objects should only contain a
125 // single object. If not, throw an exception
126 if (gbl_kink_data_list.size() != 1) {
127 throw std::runtime_error("[ SvtDataProcessor ]: The collection "
128 + std::string(Collections::TRACK_DATA_REL)
129 + " has the wrong data structure.");
130 }
131
132 // Get the list GBLKinkData GenericObject associated with the LCIO Track
133 IMPL::LCGenericObjectImpl* gbl_kink_datum
134 = static_cast<IMPL::LCGenericObjectImpl*>(gbl_kink_data_list.at(0));
135
136 // Set the lambda and phi kink values
137 for (int ikink = 0; ikink < gbl_kink_datum->getNDouble(); ++ikink) {
138 track->setLambdaKink(ikink, gbl_kink_datum->getFloatVal(ikink));
139 track->setPhiKink(ikink, gbl_kink_datum->getDoubleVal(ikink));
140 }
141
142 delete gbl_kink_data_nav;
143
144 // Get the collection of LCRelations between track data variables
145 // (TrackData) and the corresponding track.
146 EVENT::LCCollection* track_data = static_cast<EVENT::LCCollection*>(
147 event->getLCCollection(Collections::TRACK_DATA_REL));
148
149 // Instantiate an LCRelation navigator which will allow faster access
150 // to TrackData objects
151 UTIL::LCRelationNavigator* track_data_nav
152 = new UTIL::LCRelationNavigator(track_data);
153
154 // Get the list of TrackData associated with the LCIO Track
155 EVENT::LCObjectVec track_data_list = track_data_nav->getRelatedFromObjects(lc_track);
156
157 // The container of TrackData objects should only contain a single
158 // object. If not, throw an exception.
159 if (track_data_list.size() == 1) {
160
161 // Get the TrackData GenericObject associated with the LCIO Track
162 IMPL::LCGenericObjectImpl* track_datum = static_cast<IMPL::LCGenericObjectImpl*>(track_data_list.at(0));
163
164 // Check that the TrackData data structure is correct. If it's
165 // not, throw a runtime exception.
166 if (track_datum->getNDouble() > 14 || track_datum->getNFloat() > 7
167 || track_datum->getNInt() != 1) {
168 throw std::runtime_error("[ SvtDataProcessor ]: The collection "
169 + std::string(Collections::TRACK_DATA)
170 + " has the wrong structure.");
171 }
172
173 // Set the SvtTrack isolation values
174 for (int iso_index = 0; iso_index < track_datum->getNDouble(); ++iso_index) {
175 track->setIsolation(iso_index, track_datum->getDoubleVal(iso_index));
176 }
177
178 // Set the SvtTrack time
179 track->setTrackTime(track_datum->getFloatVal(0));
180
181 // Set the volume (top/bottom) in which the SvtTrack resides
182 track->setTrackVolume(track_datum->getIntVal(0));
183 }
184
185 delete track_data_nav;
186
187 // Get the collection of 3D hits associated with a LCIO Track
188 EVENT::TrackerHitVec lc_tracker_hits = lc_track->getTrackerHits();
189
190 // Iterate through the collection of 3D hits (TrackerHit objects)
191 // associated with a track, find the corresponding hits in the HPS
192 // event and add references to the track
193 for (auto hit : lc_tracker_hits) {
194
195 // Add a reference to the hit
196 track->addHit(hit_map[hit]);
197 }
198 tracks_.push_back(track);
199 }
200
201 //event->addCollection(Collections::GBL_TRACKS, tracks_);
202
203 return true;
204}
205
208
#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
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
Insert description here. more details.
virtual void finalize()
Callback for the Processor to take any necessary action when the processing of events finishes.
std::vector< Track * > tracks_
virtual void initialize(TTree *tree)
Callback for the Processor to take any necessary action when the processing of events starts.
SvtDataProcessor(const std::string &name, Process &process)
Class constructor.
std::vector< TrackerHit * > hits_
Definition Track.h:32
void setCovarianceMatrix(std::vector< float > cov)
void setCharge(const double charge)
Definition TrackerHit.h:90
void setTime(const double time)
Definition TrackerHit.h:73
void setPosition(const double *position, bool rotate=false, int type=0)
constexpr const char * GBL_TRACKS
Definition Collections.h:14
constexpr const char * TRACK_DATA_REL
Definition Collections.h:44
constexpr const char * KINK_DATA_REL
Definition Collections.h:38
constexpr const char * TRACK_DATA
Definition Collections.h:41
constexpr const char * TRACKER_HITS
Definition Collections.h:29