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