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
FinalStateParticleProcessor.cxx
Go to the documentation of this file.
1
7#include "utilities.h"
8
10 : Processor(name, process) {
11 }
12
15
17
18 std::cout << "Configuring FinalStateParticleProcessor" << std::endl;
19 try
20 {
21 debug_ = parameters.getInteger("debug", debug_);
22 fspCollLcio_ = parameters.getString("fspCollLcio", fspCollLcio_);
23 fspCollRoot_ = parameters.getString("fspCollRoot", fspCollRoot_);
24 kinkRelCollLcio_ = parameters.getString("kinkRelCollLcio", kinkRelCollLcio_);
25 trkRelCollLcio_ = parameters.getString("trkRelCollLcio", trkRelCollLcio_);
26 hitFitsCollLcio_ = parameters.getString("hitFitsCollLcio", hitFitsCollLcio_);
27 trkhitCollRoot_ = parameters.getString("trkhitCollRoot",trkhitCollRoot_);
28 rawhitCollRoot_ = parameters.getString("rawhitCollRoot",rawhitCollRoot_);
29 bfield_ = parameters.getDouble("bfield",bfield_);
30 }
31 catch (std::runtime_error& error)
32 {
33 std::cout << error.what() << std::endl;
34 }
35}
36
38 // Add branches to tree
39
40 if (!trkhitCollRoot_.empty())
41 tree->Branch(trkhitCollRoot_.c_str(), &hits_);
42
43 if (!rawhitCollRoot_.empty())
44 tree->Branch(rawhitCollRoot_.c_str(), &rawhits_);
45 tree->Branch(fspCollRoot_.c_str(), &fsps_);
46}
47
49
50 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Clear output vector" << std::endl;
51
52 //Clean up
53 if (hits_.size() > 0) {
54 for (std::vector<TrackerHit *>::iterator it = hits_.begin(); it != hits_.end(); ++it) {
55 delete *it;
56 }
57 hits_.clear();
58 }
59
60 if (rawhits_.size() > 0) {
61 for (std::vector<RawSvtHit *>::iterator it = rawhits_.begin(); it != rawhits_.end(); ++it) {
62 delete *it;
63 }
64 rawhits_.clear();
65 }
66
67 for(int i = 0; i < fsps_.size(); i++) delete fsps_.at(i);
68 fsps_.clear();
69
70 Event* event = static_cast<Event*> (ievent);
71
72 // Get the collection of vertices from the LCIO event. If no such collection
73 // exist, a DataNotAvailableException is thrown
74 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Get LCIO Collection " << fspCollLcio_ << std::endl;
75 EVENT::LCCollection* lc_fsps= nullptr;
76 try
77 {
78 lc_fsps = event->getLCCollection(fspCollLcio_.c_str());
79 }
80 catch (EVENT::DataNotAvailableException e)
81 {
82 std::cout << e.what() << std::endl;
83 return false;
84 }
85
86 // Get the collection of LCRelations between GBL tracks and kink data and track data variables.
87 EVENT::LCCollection* gbl_kink_data{nullptr};
88 EVENT::LCCollection* track_data{nullptr};
89
90 try
91 {
92 if (!kinkRelCollLcio_.empty())
93 gbl_kink_data = static_cast<EVENT::LCCollection*>(event->getLCCollection(kinkRelCollLcio_.c_str()));
94 if (!trkRelCollLcio_.empty())
95 track_data = static_cast<EVENT::LCCollection*>(event->getLCCollection(trkRelCollLcio_.c_str()));
96 }
97 catch (EVENT::DataNotAvailableException e)
98 {
99 std::cout << e.what() << std::endl;
100 if (!gbl_kink_data)
101 std::cout<<"Failed retrieving " << kinkRelCollLcio_ <<std::endl;
102 if (!track_data)
103 std::cout<<"Failed retrieving " << trkRelCollLcio_ <<std::endl;
104
105 }
106
107
108 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Converting"<< std::endl;
109
110 bool rotateHits = true;
111 int hitType = 0;
112 EVENT::LCCollection* raw_svt_hit_fits = nullptr;
113
114 auto evColls = event->getLCEvent()->getCollectionNames();
115 auto it = std::find (evColls->begin(), evColls->end(), hitFitsCollLcio_.c_str());
116 bool hasFits = true;
117 if(it == evColls->end()) hasFits = false;
118 if(hasFits)
119 {
120 raw_svt_hit_fits = event->getLCCollection(hitFitsCollLcio_.c_str());
121 }
122 for (int ifsp = 0 ; ifsp < lc_fsps->getNumberOfElements(); ++ifsp)
123 {
124 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Converting FinalStateParticle " << ifsp << std::endl;
125 EVENT::ReconstructedParticle* lc_fsp{nullptr};
126 lc_fsp = static_cast<EVENT::ReconstructedParticle*>(lc_fsps->getElementAt(ifsp));
127 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Build Particle" << std::endl;
128
129 Particle * fsp = utils::buildParticle(lc_fsp,"", gbl_kink_data, track_data);
130 if (lc_fsp->getTracks().size()>0){
131 EVENT::Track* lc_track = static_cast<EVENT::Track*>(lc_fsp->getTracks()[0]);
132 Track* track = utils::buildTrack(lc_track,"",gbl_kink_data,track_data);
133 if (bfield_ > 0.0) track->setMomentum(bfield_);
134 if (track->isKalmanTrack()) hitType = 1; //SiClusters
135 EVENT::TrackerHitVec lc_tracker_hits = lc_track->getTrackerHits();
136 for (auto lc_tracker_hit : lc_tracker_hits) {
137 TrackerHit* tracker_hit = utils::buildTrackerHit(static_cast<IMPL::TrackerHitImpl*>(lc_tracker_hit),rotateHits,hitType);
138 std::vector<RawSvtHit*> rawSvthitsOn3d;
139 utils::addRawInfoTo3dHit(tracker_hit,static_cast<IMPL::TrackerHitImpl*>(lc_tracker_hit),
140 raw_svt_hit_fits,&rawSvthitsOn3d,hitType);
141 for (auto rhit : rawSvthitsOn3d)
142 rawhits_.push_back(rhit);
143 //rawhits_->addHit(rhit);
144
145 track->addHit(tracker_hit);
146 hits_.push_back(tracker_hit);
147 rawSvthitsOn3d.clear();
148 // loop on j>i tracks
149 }
150 fsp->setTrack(track);
151 delete track;
152 }
153
154 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: Add Particle" << std::endl;
155 fsps_.push_back(fsp);
156 }
157
158 if (debug_ > 0) std::cout << "FinalStateParticleProcessor: End process" << std::endl;
159 return true;
160}
161
164
#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
Insert description here. 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.
FinalStateParticleProcessor(const std::string &name, Process &process)
Class constructor.
virtual void initialize(TTree *tree)
Callback for the Processor to take any necessary action when the processing of events starts.
std::vector< RawSvtHit * > rawhits_
std::vector< TrackerHit * > hits_
Definition IEvent.h:7
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
Definition Track.h:32
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
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