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
utilities.cxx
Go to the documentation of this file.
1#include "utilities.h"
2#include <algorithm>
3#include <memory>
4/*
5 void utils::buildTrackCollection(std::vector<Track*>& tracks,
6 Event* event,
7 const char* LCTrackCollection)
8 {
9
10 EVENT::LCCollection* lc_tracks event->getLCCollection(LCTrackCollection);
11
12
13 }
14
15*/
16
17
18bool utils::hasCollection(EVENT::LCEvent* lc_event,const std::string& collection) {
19
20 if (!lc_event || collection.empty())
21 return false;
22
23 auto evColls = lc_event->getCollectionNames();
24 auto it = std::find(evColls->begin(),evColls->end(), collection);
25 if (it!=evColls->end())
26 return true;
27 return false;
28}
29
30
31Vertex* utils::buildVertex(EVENT::Vertex* lc_vertex) {
32
33 if (!lc_vertex)
34 return nullptr;
35
36 //TODO move the static cast outside?
37
38 Vertex* vertex = new Vertex();
39 vertex->setChi2 (lc_vertex->getChi2());
40 vertex->setProbability (lc_vertex->getProbability());
41 vertex->setID (lc_vertex->id());
42 vertex->setType (lc_vertex->getAlgorithmType());
43 vertex->setVtxParameters((std::vector<float>)lc_vertex->getParameters());
44
45 //TODO Rotate the covariance matrix!
46 vertex->setCovariance ((std::vector<float>)lc_vertex->getCovMatrix());
47 //std::cout<<lc_vertex->getVertexParameterNames[0]<<std::endl;
48 vertex->setType (lc_vertex->getAlgorithmType());
49
50 vertex->setPos (lc_vertex->getPosition(),false);
51
52
53 return vertex;
54}
55
56Particle* utils::buildParticle(EVENT::ReconstructedParticle* lc_particle,
57 std::string trackstate_location,
58 EVENT::LCCollection* gbl_kink_data,
59 EVENT::LCCollection* track_data)
60
61{
62
63 if (!lc_particle)
64 return nullptr;
65
66 Particle* part = new Particle();
67 // Set the charge of the HpsParticle
68 part->setCharge(lc_particle->getCharge());
69
70 // Set the HpsParticle type
71 part->setType(lc_particle->getType());
72
73 // Set the energy of the HpsParticle
74 part->setEnergy(lc_particle->getEnergy());
75
76 // Set the momentum of the HpsParticle
77 part->setMomentum(lc_particle->getMomentum());
78
79 // Set the mass of the HpsParticle
80 part->setMass(lc_particle->getMass());
81
82 // Set the goodness of PID for the HpsParticle
83 part->setGoodnessOfPID(lc_particle->getGoodnessOfPID());
84
85 // Set the PDG ID for the HpsParticle
86 part->setPDG(lc_particle->getParticleIDUsed()->getPDG());
87
88 // Set the Track for the HpsParticle
89 if (lc_particle->getTracks().size()>0)
90 {
91 Track * trkPtr = utils::buildTrack(lc_particle->getTracks()[0],trackstate_location, gbl_kink_data, track_data);
92 part->setTrack(trkPtr);
93 delete trkPtr;
94 }
95
96 // Set the Track for the HpsParticle
97 if (lc_particle->getClusters().size() > 0)
98 {
99 CalCluster * clusBuf = utils::buildCalCluster(lc_particle->getClusters()[0]);
100 part->setCluster(clusBuf);
101 delete clusBuf;
102 }
103
104 return part;
105}
106
107CalCluster* utils::buildCalCluster(EVENT::Cluster* lc_cluster)
108{
109
110 if (!lc_cluster)
111 return nullptr;
112
113 CalCluster* cluster = new CalCluster();
114 // Set the cluster position
115 cluster->setPosition(lc_cluster->getPosition());
116
117 // Set the cluster energy
118 cluster->setEnergy(lc_cluster->getEnergy());
119
120 // Get the ecal hits used to create the cluster
121 EVENT::CalorimeterHitVec lc_hits = lc_cluster->getCalorimeterHits();
122
123 // Loop over all of the Ecal hits and add them to the Ecal cluster. The
124 // seed hit is set to be the hit with the highest energy. The cluster time
125 // is set to be the hit time of the seed hit.
126 double senergy = 0;
127 double stime = 0;
128 for(int ihit = 0; ihit < (int) lc_hits.size(); ++ihit) {
129 // Get an Ecal hit
130 EVENT::CalorimeterHit* lc_hit = lc_hits[ihit];
131 if (senergy < lc_hit->getEnergy()) {
132 senergy = lc_hit->getEnergy();
133 stime = lc_hit->getTime();
134 }
135 }
136
137 // Set the time of the cluster
138 cluster->setTime(stime);
139
140 return cluster;
141}
142
143bool utils::IsSameTrack(Track* trk1, Track* trk2) {
144 double tol = 1e-6;
145 if (fabs(trk1->getD0() - trk2->getD0()) > tol ||
146 fabs(trk1->getPhi() - trk2->getPhi()) > tol ||
147 fabs(trk1->getOmega() - trk2->getOmega()) > tol ||
148 fabs(trk1->getTanLambda() - trk2->getTanLambda()) > tol ||
149 fabs(trk1->getZ0() - trk2->getZ0()) > tol ||
150 fabs(trk1->getChi2Ndf() - trk2->getChi2Ndf()) > tol
151 )
152 return false;
153
154 return true;
155}
156
157
158Track* utils::buildTrack(EVENT::Track* lc_track,
159 std::string trackstate_location,
160 EVENT::LCCollection* gbl_kink_data,
161 EVENT::LCCollection* track_data) {
162
163 if (!lc_track)
164 return nullptr;
165
166 //TrackState Location map
167 std::map<std::string, int> trackstateLocationMap_ = {
168 {"", EVENT::TrackState::AtIP},
169 {"AtTarget", EVENT::TrackState::LastLocation}
170 };
171
172 int loc;
173 auto it = trackstateLocationMap_.find(trackstate_location);
174 if (it != trackstateLocationMap_.end()){
175 loc = it->second;
176 }
177 else{
178 std::cout << "[utilities]::ERROR Track State Location " << trackstate_location << " Doesn't Exist!" << std::endl;
179 std::cout << "Check map in utilities::buildTrack for defined locations" << std::endl;
180 return nullptr;
181 }
182
183 Track* track = new Track();
184 //If using track AtIP, get params from lc_track
185 if (loc == trackstateLocationMap_[""]){
186 // Set the track parameters
187 track->setTrackParameters(lc_track->getD0(),
188 lc_track->getPhi(),
189 lc_track->getOmega(),
190 lc_track->getTanLambda(),
191 lc_track->getZ0());
192
193 // Set the track covariance matrix
194 track->setCov(static_cast<std::vector<float> > (lc_track->getCovMatrix()));
195 }
196
197 //If other TrackState specified, get track params from track state
198 else {
199 // If track state doesn't exist, no track returned
200 const EVENT::TrackState* ts = lc_track->getTrackState(loc);
201 if (ts == nullptr){
202 return nullptr;
203 }
204 // Set the track parameters using trackstate
205 track->setTrackParameters(ts->getD0(),
206 ts->getPhi(),
207 ts->getOmega(),
208 ts->getTanLambda(),
209 ts->getZ0());
210
211 double position[3] = {
212 ts->getReferencePoint()[1],
213 ts->getReferencePoint()[2],
214 ts->getReferencePoint()[0]
215 };
216
217 track->setCov(static_cast<std::vector<float> > (ts->getCovMatrix()));
218
219 track->setPosition(position);
220 }
221
222 // Set the track id
223 track->setID(lc_track->id());
224
225 // Set the track type
226 track->setType(lc_track->getType());
227
228 // Set the track fit chi^2
229 track->setChi2(lc_track->getChi2());
230
231 // Set the track ndf
232 track->setNdf(lc_track->getNdf());
233
234 // Set the position of the extrapolated track at the ECal face. The
235 // extrapolation uses the full 3D field map.
236 const EVENT::TrackState* track_state
237 = lc_track->getTrackState(EVENT::TrackState::AtCalorimeter);
238
239 if (track_state) {
240 double position_at_ecal[3] = {
241 track_state->getReferencePoint()[1],
242 track_state->getReferencePoint()[2],
243 track_state->getReferencePoint()[0]
244 };
245 track->setPositionAtEcal(position_at_ecal);
246 }
247
248 if (gbl_kink_data) {
249 // Instantiate an LCRelation navigator which will allow faster access
250 // to GBLKinkData object
251 std::shared_ptr<UTIL::LCRelationNavigator> gbl_kink_data_nav = std::make_shared<UTIL::LCRelationNavigator>(gbl_kink_data);
252
253 // Get the list of GBLKinkData associated with the LCIO Track
254 EVENT::LCObjectVec gbl_kink_data_list
255 = gbl_kink_data_nav->getRelatedFromObjects(lc_track);
256
257 // The container of GBLKinkData objects should only contain a
258 // single object. If not, throw an exception
259 if (gbl_kink_data_list.size() == 1) {
260
261 /*
262 std::cout<<"[ Utilities ]: The collection "
263 + std::string(Collections::KINK_DATA)
264 + " has the wrong data structure for this track"<<std::endl;
265 */
266
267 // Get the list GBLKinkData GenericObject associated with the LCIO Track
268 IMPL::LCGenericObjectImpl* gbl_kink_datum
269 = static_cast<IMPL::LCGenericObjectImpl*>(gbl_kink_data_list.at(0));
270
271 // Set the lambda and phi kink values
272 for (int ikink = 0; ikink < gbl_kink_datum->getNDouble(); ++ikink) {
273 track->setLambdaKink(ikink, gbl_kink_datum->getFloatVal(ikink));
274 track->setPhiKink(ikink, gbl_kink_datum->getDoubleVal(ikink));
275 }
276
277 }//gbl_kink_data has right structure
278
279 } // add gbl kink data
280
281 if (track_data) {
282
283 // Instantiate an LCRelation navigator which will allow faster access
284 // to TrackData objects
285 std::shared_ptr<UTIL::LCRelationNavigator> track_data_nav = std::make_shared<UTIL::LCRelationNavigator>(track_data);
286
287 // Get the list of TrackData associated with the LCIO Track
288 EVENT::LCObjectVec track_data_list = track_data_nav->getRelatedFromObjects(lc_track);
289
290 // The container of TrackData objects should only contain a single
291 // object. If not, throw an exception.
292 if (track_data_list.size() == 1) {
293
294 // Get the TrackData GenericObject associated with the LCIO Track
295 IMPL::LCGenericObjectImpl* track_datum = static_cast<IMPL::LCGenericObjectImpl*>(track_data_list.at(0));
296
297 // Check that the TrackData data structure is correct. If it's
298 // not, throw a runtime exception.
299 if (track_datum->getNDouble() > 14 || track_datum->getNFloat() > 8 || track_datum->getNInt() != 1) {
300 throw std::runtime_error("[ TrackingProcessor ]: The collection "
301 + std::string(Collections::TRACK_DATA)
302 + " has the wrong structure.");
303 }
304
305 // Set the SvtTrack isolation values
306 for (int iso_index = 0; iso_index < track_datum->getNDouble(); ++iso_index) {
307 track->setIsolation(iso_index, track_datum->getDoubleVal(iso_index));
308 }
309
310 // Set the SvtTrack time
311 track->setTrackTime(track_datum->getFloatVal(0));
312
313 // Set the Track momentum
314 if (track_datum->getNFloat()>3)
315 track->setMomentum(track_datum->getFloatVal(1),track_datum->getFloatVal(2),track_datum->getFloatVal(3));
316
317 // Set the volume (top/bottom) in which the SvtTrack resides
318 track->setTrackVolume(track_datum->getIntVal(0));
319
320 // Set the BfieldY for track state
321 double bfieldY = -999.9;
322 if(track_datum->getNFloat() > 4){
323 if (loc == trackstateLocationMap_[""])
324 bfieldY = track_datum->getFloatVal(4);
325 if (loc == trackstateLocationMap_["AtTarget"]){
326 bfieldY = track_datum->getFloatVal(5);
327 }
328 if (loc == trackstateLocationMap_["AtCalorimeter"])
329 bfieldY = track_datum->getFloatVal(6);
330 //Bfield needs factor of -1, not sure why... <-TODO investigate
331 track->setMomentum(-bfieldY);
332 }
333 }
334
335 } //add track data
336
337 return track;
338}
339
340RawSvtHit* utils::buildRawHit(EVENT::TrackerRawData* rawTracker_hit,
341 EVENT::LCCollection* raw_svt_hit_fits) {
342
343 EVENT::long64 value =
344 EVENT::long64(rawTracker_hit->getCellID0() & 0xffffffff) |
345 ( EVENT::long64(rawTracker_hit->getCellID1() ) << 32 );
346 decoder.setValue(value);
347
348 RawSvtHit* rawHit = new RawSvtHit();
349 rawHit->setSystem(decoder["system"]);
350 rawHit->setBarrel(decoder["barrel"]);
351 rawHit->setLayer(decoder["layer"]);
352 rawHit->setModule(decoder["module"]);
353 rawHit->setSensor(decoder["sensor"]);
354 rawHit->setSide(decoder["side"]);
355 rawHit->setStrip(decoder["strip"]);
356
357 // Extract ADC values for this hit
358 int hit_adcs[6] = {
359 (int)rawTracker_hit->getADCValues().at(0),
360 (int)rawTracker_hit->getADCValues().at(1),
361 (int)rawTracker_hit->getADCValues().at(2),
362 (int)rawTracker_hit->getADCValues().at(3),
363 (int)rawTracker_hit->getADCValues().at(4),
364 (int)rawTracker_hit->getADCValues().at(5)};
365
366 rawHit->setADCs(hit_adcs);
367 if (raw_svt_hit_fits) {
368 std::shared_ptr<UTIL::LCRelationNavigator> rawTracker_hit_fits_nav = std::make_shared<UTIL::LCRelationNavigator>(raw_svt_hit_fits);
369
370
371 // Get the list of fit params associated with the raw tracker hit
372 EVENT::LCObjectVec rawTracker_hit_fits_list
373 = rawTracker_hit_fits_nav->getRelatedToObjects(rawTracker_hit);
374
375 // Get the list SVTFittedRawTrackerHit GenericObject associated with the SVTRawTrackerHit
376 IMPL::LCGenericObjectImpl* hit_fit_param
377 = static_cast<IMPL::LCGenericObjectImpl*>(rawTracker_hit_fits_list.at(0));
378
379 double fit_params[5] = {
380 (double)hit_fit_param->getDoubleVal(0),
381 (double)hit_fit_param->getDoubleVal(1),
382 (double)hit_fit_param->getDoubleVal(2),
383 (double)hit_fit_param->getDoubleVal(3),
384 (double)hit_fit_param->getDoubleVal(4)
385 };
386 rawHit->setFit(fit_params, 0);
387 if(rawTracker_hit_fits_list.size()>1)
388 {
389 hit_fit_param = static_cast<IMPL::LCGenericObjectImpl*>(rawTracker_hit_fits_list.at(1));
390 fit_params[0] = (double)hit_fit_param->getDoubleVal(0);
391 fit_params[1] = (double)hit_fit_param->getDoubleVal(1);
392 fit_params[2] = (double)hit_fit_param->getDoubleVal(2);
393 fit_params[3] = (double)hit_fit_param->getDoubleVal(3);
394 fit_params[4] = (double)hit_fit_param->getDoubleVal(4);
395
396 rawHit->setFit(fit_params, 1);
397 }
398 rawHit->setFitN(rawTracker_hit_fits_list.size());
399 }//raw svt hit fits
400
401 return rawHit;
402
403}//build raw hit
404
405//type = 0 RotatedHelicalTrackHit type = 1 SiCluster
406TrackerHit* utils::buildTrackerHit(IMPL::TrackerHitImpl* lc_tracker_hit, bool rotate, int type) {
407
408 if (!lc_tracker_hit)
409 return nullptr;
410
411 TrackerHit* tracker_hit = new TrackerHit();
412
413 // Get the position of the LCIO TrackerHit and set the position of
414 // the TrackerHit
415 double hit_position[3] = {
416 lc_tracker_hit->getPosition()[0], //lcio x
417 lc_tracker_hit->getPosition()[1], //lcio y
418 lc_tracker_hit->getPosition()[2] //lcio z
419 };
420 tracker_hit->setPosition(hit_position, rotate, type);
421
422 // Set the covariance matrix of the SvtHit
423 tracker_hit->setCovarianceMatrix(lc_tracker_hit->getCovMatrix());
424
425 // Set the time of the SvtHit
426 tracker_hit->setTime(lc_tracker_hit->getTime());
427
428 // Set the charge of the SvtHit
429 tracker_hit->setCharge(lc_tracker_hit->getEDep());
430
431 // Set the LCIO id
432 tracker_hit->setID(lc_tracker_hit->id());
433
434 return tracker_hit;
435
436
437}
438
439//type 0 rotatedHelicalHit type 1 SiClusterHit
441 IMPL::TrackerHitImpl* lc_tracker_hit,
442 EVENT::LCCollection* raw_svt_fits, std::vector<RawSvtHit*>* rawHits,int type, bool storeRawHit) {
443
444 if (!tracker_hit || !lc_tracker_hit)
445 return false;
446
447 float rawcharge = 0;
448 //0 top 1 bottom
449 int volume = -1;
450 //1-6(7) for rotated 0-13 for SiCluster
451 int layer = -1;
452
453 // Get decoders to read cellids
454 UTIL::BitField64 decoder("system:6,barrel:3,layer:4,module:12,sensor:1,side:32:-2,strip:12");
455
456 //Get the Raw content of the tracker hits
457 EVENT::LCObjectVec lc_rawHits = lc_tracker_hit->getRawHits();
458
459 std::vector<int> rawhit_strips;
460 for (unsigned int irh = 0 ; irh < lc_rawHits.size(); ++irh) {
461 // Get a raw hit from the list of hits
462 EVENT::TrackerRawData* rawTracker_hit
463 = static_cast<EVENT::TrackerRawData*>(lc_rawHits.at(irh));
464
465 //Decode the cellid
466 EVENT::long64 value = EVENT::long64( rawTracker_hit->getCellID0() & 0xffffffff ) |
467 ( EVENT::long64( rawTracker_hit->getCellID1() ) << 32 ) ;
468 decoder.setValue(value);
469
470 //Get rawhit strip number
471 int stripnumber = decoder["strip"];
472 rawhit_strips.push_back(stripnumber);
473
474 //TODO useless to build all of it?
475 RawSvtHit* rawHit = buildRawHit(rawTracker_hit,raw_svt_fits);
476 rawcharge += rawHit->getAmp(0);
477 int currentHitVolume = rawHit->getModule() % 2 ? 1 : 0;
478 int currentHitLayer = (rawHit->getLayer() - 1 ) / 2;
479 if (type == 1)
480 currentHitLayer = rawHit->getLayer() - 1;
481 if (volume == -1 )
482 volume = currentHitVolume;
483 else {
484 if ( currentHitVolume != volume)
485 std::cout<<"[ ERROR ] : utils::addRawInfoTo3dHit raw hits with inconsistent volume found" <<std::endl;
486 }
487
488 if (layer == -1 )
489 layer = currentHitLayer;
490 else {
491 if (currentHitLayer != layer)
492 std::cout<<"[ ERROR ] : utils::addRawInfoTo3dHit raw hits with inconsistent layer found" <<std::endl;
493 }
494
495 if(storeRawHit){
496 tracker_hit->addRawHit(rawHit);
497 if (rawHits)
498 rawHits->push_back(rawHit);
499 }
500 else
501 delete rawHit;
502
503 }
504
505 tracker_hit->setRawCharge(rawcharge);
506 tracker_hit->setVolume(volume);
507 tracker_hit->setLayer(layer);
508 tracker_hit->setRawHitStripNumbers(rawhit_strips);
509
510 return true;
511}
512
513//TODO-improve shared finding algorithm
514
515bool utils::isUsedByTrack(IMPL::TrackerHitImpl* lc_tracker_hit,
516 EVENT::Track* lc_track) {
517
518 EVENT::TrackerHitVec trk_lc_tracker_hits = lc_track->getTrackerHits();
519
520 for (auto trk_lc_tracker_hit : trk_lc_tracker_hits) {
521 //std::cout<<lc_tracker_hit->id()<<" " <<trk_lc_tracker_hit->id()<<std::endl;
522 if (lc_tracker_hit -> id() == trk_lc_tracker_hit -> id())
523 return true;
524 }
525 return false;
526}
527
529 EVENT::Track* lc_track) {
530
531 EVENT::TrackerHitVec trk_lc_tracker_hits = lc_track->getTrackerHits();
532
533 for (auto trk_lc_tracker_hit : trk_lc_tracker_hits) {
534 if (tracker_hit->getID() == trk_lc_tracker_hit->id())
535 return true;
536 }
537 return false;
538}
539
540
542
543 for (int ipart = 0; ipart < vtx->getParticles().GetEntries(); ++ipart) {
544 int pdg_id = ((Particle*)vtx->getParticles().At(ipart))->getPDG();
545 if (pdg_id == 11) {
546 ele = (Particle*)vtx->getParticles().At(ipart);
547 }
548 else if (pdg_id == -11) {
549 pos = (Particle*)vtx->getParticles().At(ipart);
550 }
551
552 else {
553 std::cout<<"Utilities::Wrong particle ID "<< pdg_id <<"associated to vertex. Skip."<<std::endl;
554 return false;
555 }
556 }
557
558 if (!ele || !pos) {
559 std::cout<<"Utilities::Vertex formed without ele/pos. Skip."<<std::endl;
560 return false;
561 }
562
563 return true;
564}
565
566double utils::getKalmanTrackL1Isolations(Track* track, std::vector<TrackerHit*>* siClusters){
567 double L1_axial_iso = 999999.9;
568 double L1_stereo_iso = 999999.9;
569 //Loop over hits on track
570 for(int i = 0; i < track->getSvtHits().GetEntries(); i++){
571 TrackerHit* track_hit = (TrackerHit*)track->getSvtHits().At(i);
572 //Track hit info
573 int trackhit_id = track_hit->getID();
574 int trackhit_layer = track_hit->getLayer();
575 int trackhit_volume = track_hit->getVolume();
576 double trackhit_y = track_hit->getGlobalY();
577 double trackhit_charge = track_hit->getRawCharge();
578 double trackhit_time = track_hit->getTime();
579
580 //Only look at L1
581 if(trackhit_layer > 1)
582 continue;
583
584 //Get rawhit strip information
585 std::vector<int> trackhit_rawhits = track_hit->getRawHitStripNumbers();
586 if(trackhit_rawhits.size() < 1)
587 continue;
588 int trackhit_maxstrip = *max_element(trackhit_rawhits.begin(), trackhit_rawhits.end());
589 int trackhit_minstrip = *min_element(trackhit_rawhits.begin(), trackhit_rawhits.end());
590
591
592 //Isolation only calculated for axial sensors (sensitive to global y)
593 bool isAxial = false;
594 //Axial/Stereo Top/Bot mapping
595 if(trackhit_volume == 1){
596 if(trackhit_layer%2 == 1)
597 isAxial = true;
598 }
599 if(trackhit_volume == 0){
600 if(trackhit_layer%2 == 0)
601 isAxial = true;
602 }
603
604 //Loop over all SiClusters in event
605 //Find closest/best alternative SiCluster to track hit
606 TrackerHit* closest_althit = nullptr;
607 double isohit_dy = 999999.9;
608 double closest_dcharge = 99999.9;
609 double closest_dt = 999999.9;
610 double isohit_y = 999999.9;
611
612 if ( siClusters == nullptr )
613 return isohit_dy;
614 for(int j = 0; j < siClusters->size(); j++){
615 TrackerHit* althit = siClusters->at(j);
616 int althit_id = althit->getID();
617 int althit_volume = althit->getVolume();
618 int althit_layer = althit->getLayer();
619 double althit_y = althit->getGlobalY();
620 double althit_charge = althit->getRawCharge();
621 double althit_time = althit->getTime();
622
623 //Skip if SiCluster not on same layer as track hit
624 if (althit_layer != trackhit_layer)
625 continue;
626
627 //Skip if volume doesn't match
628 if(althit_volume != trackhit_volume)
629 continue;
630
631 //Skip same hit
632 if (althit_id == trackhit_id)
633 continue;
634
635 //Only look at hits that are further from beam-axis in Global Y
636 if ( (trackhit_volume == 0 && althit_y < trackhit_y) ||
637 (trackhit_volume == 1 && althit_y > trackhit_y))
638 continue;
639
640 //Require alternative hits to be within +-30ns (based on SiClustersOnTrack t distr)
641 if(std::abs(althit_time) > 30.0)
642 continue;
643
644 //Skip adjacent rawhits
645 std::vector<int> althit_rawhits = althit->getRawHitStripNumbers();
646 int althit_maxstrip = *max_element(althit_rawhits.begin(), althit_rawhits.end());
647 int althit_minstrip = *min_element(althit_rawhits.begin(), althit_rawhits.end());
648 if(trackhit_minstrip - althit_maxstrip <= 1 && althit_minstrip - trackhit_maxstrip <= 1)
649 continue;
650
651 //Pick closest alternative hit
652 if (std::abs(trackhit_y - althit_y) < isohit_dy){
653 isohit_dy = std::abs(trackhit_y-althit_y);
654 closest_althit = althit;
655 closest_dcharge = trackhit_charge - althit_charge;
656 closest_dt = trackhit_time - althit_time;
657 isohit_y = althit_y;
658 }
659 }
660
661 if(isAxial)
662 L1_axial_iso = isohit_dy;
663 else
664 L1_stereo_iso = isohit_dy;
665
666 }
667
668 if(L1_axial_iso < L1_stereo_iso){
669 return L1_axial_iso;
670 }
671 else
672 return L1_stereo_iso;
673}
674
675void utils::get2016KFMCTruthHitCodes(Track* ele_trk, Track* pos_trk, int& L1L2hitCode, int& L1hitCode, int& L2hitCode){
676 //Count the number of hits per part on the ele track
677 std::map<int, int> nHits4part_ele;
678 for(int i =0; i < ele_trk->getMcpHits().size(); i++)
679 {
680 int partID = ele_trk->getMcpHits().at(i).second;
681 if ( nHits4part_ele.find(partID) == nHits4part_ele.end() )
682 {
683 // not found
684 nHits4part_ele[partID] = 1;
685 }
686 else
687 {
688 // found
689 nHits4part_ele[partID]++;
690 }
691 }
692
693 //Determine the MC part with the most hits on the track
694 int maxNHits_ele = 0;
695 int maxID_ele = 0;
696 for (std::map<int,int>::iterator it=nHits4part_ele.begin(); it!=nHits4part_ele.end(); ++it)
697 {
698 if(it->second > maxNHits_ele)
699 {
700 maxNHits_ele = it->second;
701 maxID_ele = it->first;
702 }
703 }
704
705 //Count the number of hits per part on the pos track
706 std::map<int, int> nHits4part_pos;
707 for(int i =0; i < pos_trk->getMcpHits().size(); i++)
708 {
709 int partID = pos_trk->getMcpHits().at(i).second;
710 if ( nHits4part_pos.find(partID) == nHits4part_pos.end() )
711 {
712 // not found
713 nHits4part_pos[partID] = 1;
714 }
715 else
716 {
717 // found
718 nHits4part_pos[partID]++;
719 }
720 }
721
722 //Determine the MC part with the most hits on the track
723 int maxNHits_pos = 0;
724 int maxID_pos = 0;
725 for (std::map<int,int>::iterator it=nHits4part_pos.begin(); it!=nHits4part_pos.end(); ++it)
726 {
727 if(it->second > maxNHits_pos)
728 {
729 maxNHits_pos = it->second;
730 maxID_pos = it->first;
731 }
732 }
733
734 //Determine Ele L1 and L2 truth information
735 bool ele_trueAxialL1 = false;
736 bool ele_trueStereoL1 = false;
737 bool ele_trueAxialL2 = false;
738 bool ele_trueStereoL2 = false;
739 bool pos_trueAxialL1 = false;
740 bool pos_trueStereoL1 = false;
741 bool pos_trueAxialL2 = false;
742 bool pos_trueStereoL2 = false;
743
744 if(ele_trk->isKalmanTrack()){
745 for(int i = 0; i < ele_trk->getMcpHits().size(); i++){
746 int mcpid = ele_trk->getMcpHits().at(i).second;
747 int layer = ele_trk->getMcpHits().at(i).first;
748 int volume = -1;
749 if(ele_trk->isTopTrack() == 1)
750 volume = 0;
751 else
752 volume = 1;
753
754 bool isAxial = false;
755 bool isStereo = false;
756 bool isL1 = false;
757 bool isL2 = false;
758 bool isGood = false;
759
760 //L1 and L2 only
761 if(layer < 2)
762 isL1 = true;
763 else if(layer > 1 && layer < 4)
764 isL2 = true;
765 else
766 continue;
767
768 if(volume == 1){
769 if(layer%2 == 1)
770 isAxial = true;
771 else
772 isStereo = true;
773 }
774 if(volume == 0){
775 if(layer%2 == 0)
776 isAxial = true;
777 else
778 isStereo = true;
779 }
780
781 if(mcpid == maxID_ele)
782 isGood = true;
783
784 if(isGood){
785 if(isAxial){
786 if(isL1)
787 ele_trueAxialL1 = true;
788 if(isL2)
789 ele_trueAxialL2 = true;
790 }
791 if(isStereo){
792 if(isL1)
793 ele_trueStereoL1 = true;
794 if(isL2)
795 ele_trueStereoL2 = true;
796 }
797 }
798 }
799 }
800
801 if(pos_trk->isKalmanTrack()){
802 for(int i = 0; i < pos_trk->getMcpHits().size(); i++){
803 int mcpid = pos_trk->getMcpHits().at(i).second;
804 int layer = pos_trk->getMcpHits().at(i).first;
805 int volume = -1;
806 if(pos_trk->isTopTrack() == 1)
807 volume = 0;
808 else
809 volume = 1;
810
811 bool isAxial = false;
812 bool isStereo = false;
813 bool isL1 = false;
814 bool isL2 = false;
815 bool isGood = false;
816
817 //L1 and L2 only
818 if(layer < 2)
819 isL1 = true;
820 else if(layer > 1 && layer < 4)
821 isL2 = true;
822 else
823 continue;
824
825 if(volume == 1){
826 if(layer%2 == 1)
827 isAxial = true;
828 else
829 isStereo = true;
830 }
831 if(volume == 0){
832 if(layer%2 == 0)
833 isAxial = true;
834 else
835 isStereo = true;
836 }
837
838 if(mcpid == maxID_pos)
839 isGood = true;
840
841 if(isGood){
842 if(isAxial){
843 if(isL1)
844 pos_trueAxialL1 = true;
845 if(isL2)
846 pos_trueAxialL2 = true;
847 }
848 if(isStereo){
849 if(isL1)
850 pos_trueStereoL1 = true;
851 if(isL2)
852 pos_trueStereoL2 = true;
853 }
854 }
855 }
856 }
857
858
859 //Require both Axial and Stereo truth hits to be 'Good' hit
860 if(ele_trueAxialL1 && ele_trueStereoL1) L1L2hitCode = L1L2hitCode | (0x1 << 3);
861 if(pos_trueAxialL1 && pos_trueStereoL1) L1L2hitCode = L1L2hitCode | (0x1 << 2);
862 if(ele_trueAxialL2 && ele_trueStereoL2) L1L2hitCode = L1L2hitCode | (0x1 << 1);
863 if(pos_trueAxialL2 && pos_trueStereoL2) L1L2hitCode = L1L2hitCode | (0x1 << 0);
864
865
866 //Set L1 axial/stereo hit code for ele and positron
867 if(ele_trueAxialL1) L1hitCode = L1hitCode | (0x1 << 3);
868 if(ele_trueStereoL1) L1hitCode = L1hitCode | (0x1 << 2);
869 if(pos_trueAxialL1) L1hitCode = L1hitCode | (0x1 << 1);
870 if(pos_trueStereoL1) L1hitCode = L1hitCode | (0x1 << 0);
871
872 //Set L2 axial/stereo hit code for ele and positron
873 if(ele_trueAxialL2) L2hitCode = L2hitCode | (0x1 << 3);
874 if(ele_trueStereoL2) L2hitCode = L2hitCode | (0x1 << 2);
875 if(pos_trueAxialL2) L2hitCode = L2hitCode | (0x1 << 1);
876 if(pos_trueStereoL2) L2hitCode = L2hitCode | (0x1 << 0);
877}
878
879double utils::v0_projection_to_target_significance(json v0proj_fits, int run, double &vtx_proj_x, double &vtx_proj_y,
880 double &vtx_proj_x_signif, double &vtx_proj_y_signif, double vtx_x, double vtx_y, double vtx_z,
881 double vtx_px, double vtx_py, double vtx_pz){
882 //V0 Projection fit parameters are calculated externally by projecting vertices to the target z position,
883 //and then fitting the 2D distribution vtx_x vs vtx_y with a rotated 2D Gaussian.
884 //The fit parameters are defined along the rotated coordinate system.
885 //Therefore, the vertex position must be rotated into this coordinate system before calculating significance.
886 //The rotation angle corresponding to the fit is provided in the json file containing the rotated fit values.
887
888 //Read v0 projection fits from json file
889 int closest_run;
890 for(auto entry : v0proj_fits.items()){
891 int check_run = std::stoi(entry.key());
892 if(check_run > run)
893 break;
894 else{
895 closest_run = check_run;
896 }
897 }
898 double target_pos = v0proj_fits[std::to_string(closest_run)]["target_position"];
899 double rot_mean_x = v0proj_fits[std::to_string(closest_run)]["rotated_mean_x"];
900 double rot_mean_y = v0proj_fits[std::to_string(closest_run)]["rotated_mean_y"];
901 double rot_sigma_x = v0proj_fits[std::to_string(closest_run)]["rotated_sigma_x"];
902 double rot_sigma_y = v0proj_fits[std::to_string(closest_run)]["rotated_sigma_y"];
903 double rotation_angle = (double)v0proj_fits[std::to_string(closest_run)]["rotation_angle_mrad"]/1000.0;
904
905 //project vertex to target position
906 vtx_proj_x = vtx_x - ((vtx_z - target_pos)*(vtx_px/vtx_pz));
907 vtx_proj_y = vtx_y - ((vtx_z - target_pos)*(vtx_py/vtx_pz));
908
909 //Rotate projected vertex by angle corresponding to run number
910 double rot_vtx_proj_x = vtx_proj_x*std::cos(rotation_angle) - vtx_proj_y*std::sin(rotation_angle);
911 double rot_vtx_proj_y = vtx_proj_x*std::sin(rotation_angle) + vtx_proj_y*std::cos(rotation_angle);
912
913 //Calculate significance
914 vtx_proj_x_signif = (rot_vtx_proj_x - rot_mean_x)/rot_sigma_x;
915 vtx_proj_y_signif = (rot_vtx_proj_y - rot_mean_y)/rot_sigma_y;
916
917 //
918 double significance = std::sqrt( vtx_proj_x_signif*vtx_proj_x_signif + vtx_proj_y_signif*vtx_proj_y_signif );
919
920 return significance;
921}
nlohmann::json json
void setPosition(const float *position)
void setEnergy(const double energy)
Definition CalCluster.h:70
void setTime(const double time)
Definition CalCluster.h:80
void setGoodnessOfPID(const double goodness_pid)
Definition Particle.h:84
void setPDG(const int pdg)
Definition Particle.h:98
void setMass(const double mass)
Definition Particle.h:112
void setCharge(const int charge)
Definition Particle.h:77
void setType(const int type)
Definition Particle.h:91
void setMomentum(const double *momentum)
Definition Particle.cxx:23
void setTrack(Track *track)
Definition Particle.h:42
void setEnergy(const double energy)
Definition Particle.h:105
void setCluster(CalCluster *cluster)
Definition Particle.h:56
void setBarrel(int barrel)
Definition RawSvtHit.cxx:48
void setSystem(int system)
Definition RawSvtHit.cxx:44
void setFit(double fit[5], int fitI)
Definition RawSvtHit.cxx:27
int getLayer()
Definition RawSvtHit.cxx:88
void setStrip(int strip)
Definition RawSvtHit.cxx:68
void setSensor(int sensor)
Definition RawSvtHit.cxx:60
void setSide(int side)
Definition RawSvtHit.cxx:64
double getAmp(int fitI)
Definition RawSvtHit.h:101
void setADCs(int adcs[6])
Definition RawSvtHit.cxx:35
int getModule()
Definition RawSvtHit.cxx:92
void setFitN(int fitN)
Definition RawSvtHit.cxx:23
void setLayer(int layer)
Definition RawSvtHit.cxx:52
void setModule(int module)
Definition RawSvtHit.cxx:56
Definition Track.h:32
bool isKalmanTrack() const
Definition Track.h:231
bool isTopTrack() const
Definition Track.h:330
double getPhi() const
Definition Track.h:88
double getChi2Ndf() const
Definition Track.h:126
double getOmega() const
Definition Track.h:90
double getD0() const
Definition Track.h:87
double getTanLambda() const
Definition Track.h:91
std::vector< std::pair< int, int > > getMcpHits()
Definition Track.h:103
double getZ0() const
Definition Track.h:92
double getGlobalY() const
Definition TrackerHit.h:53
void setLayer(const int layer)
Definition TrackerHit.h:106
int getLayer() const
Definition TrackerHit.h:109
void setVolume(const int volume)
Definition TrackerHit.h:100
void addRawHit(TObject *rawhit)
Definition TrackerHit.h:115
void setCovarianceMatrix(std::vector< float > cov)
void setID(const int id)
Definition TrackerHit.h:122
std::vector< int > getRawHitStripNumbers()
Definition TrackerHit.h:130
void setRawHitStripNumbers(std::vector< int > rawhit_strips)
Definition TrackerHit.h:133
void setRawCharge(const double rawcharge)
Definition TrackerHit.h:95
int getID() const
Definition TrackerHit.h:124
void setCharge(const double charge)
Definition TrackerHit.h:90
int getVolume()
Definition TrackerHit.h:103
double getRawCharge() const
Definition TrackerHit.h:98
void setTime(const double time)
Definition TrackerHit.h:73
void setPosition(const double *position, bool rotate=false, int type=0)
double getTime() const
Definition TrackerHit.h:83
void setPos(const float *pos, bool rotate=false)
Definition Vertex.cxx:39
void setVtxParameters(const std::vector< float > &parameters)
Definition Vertex.cxx:81
void setID(const int id)
Definition Vertex.h:83
void setType(const std::string &type)
Definition Vertex.h:80
void setProbability(const float probability)
Definition Vertex.h:92
void setCovariance(const std::vector< float > &vec)
Definition Vertex.cxx:35
void setChi2(const double chi2)
Definition Vertex.h:45
constexpr const char * TRACK_DATA
Definition Collections.h:41
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
CalCluster * buildCalCluster(EVENT::Cluster *lc_cluster)
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
bool hasCollection(EVENT::LCEvent *lc_event, const std::string &collection)
description
Definition utilities.cxx:18
bool isUsedByTrack(IMPL::TrackerHitImpl *lc_tracker_hit, EVENT::Track *lc_track)
description
TrackerHit * buildTrackerHit(IMPL::TrackerHitImpl *lc_trackerHit, bool rotate=true, int type=0)
description
double getKalmanTrackL1Isolations(Track *track, std::vector< TrackerHit * > *siClusters)
description
bool getParticlesFromVertex(Vertex *vtx, Particle *ele, Particle *pos)
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
double v0_projection_to_target_significance(json v0proj_fits, int run, double &vtx_proj_x, double &vtx_proj_y, double &vtx_proj_x_signif, double &vtx_proj_y_signif, double vtx_x, double vtx_y, double vtx_z, double vtx_px, double vtx_py, double vtx_pz)
description
void get2016KFMCTruthHitCodes(Track *ele_trk, Track *pos_trk, int &L1L2hitCode, int &L1hitCode, int &L2hitCode)
description
bool IsSameTrack(Track *trk1, Track *trk2)
description
RawSvtHit * buildRawHit(EVENT::TrackerRawData *rawTracker_hit, EVENT::LCCollection *raw_svt_hit_fits)