CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ConversionTrackMerger.cc
Go to the documentation of this file.
1 //
2 // Package: RecoTracker/FinalTrackSelectors
3 // Class: ConversionTrackMerger
4 //
5 // Description: Merger for ConversionTracks, adapted from SimpleTrackListMerger
6 //
7 // Original Author: J.Bendavid
8 //
9 //
10 
11 #include <memory>
12 #include <string>
13 #include <iostream>
14 #include <cmath>
15 #include <vector>
16 
18 
24 
26 
28 
33 
34 //#include "DataFormats/TrackReco/src/classes.h"
35 
37 
38 
40  conf_(conf)
41  {
42  // retrieve producer name of input TrackCollection(s)
43  trackProducer1 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer1"));
44  trackProducer2 = consumes<reco::ConversionTrackCollection>(conf_.getParameter<edm::InputTag>("TrackProducer2"));
45 
46  produces<reco::ConversionTrackCollection>();
47 
48  }
49 
50 
51  // Virtual destructor needed.
53 
54  // Functions that gets called by framework every event
56  {
57 
58 
59  double shareFrac = conf_.getParameter<double>("ShareFrac");
60  bool allowFirstHitShare = conf_.getParameter<bool>("allowFirstHitShare");
61  bool checkCharge = conf_.getParameter<bool>("checkCharge");
62  double minProb = conf_.getParameter<double>("minProb");
63 
64  int outputPreferCollection = conf_.getParameter<int>("outputPreferCollection");
65  int trackerOnlyPreferCollection = conf_.getParameter<int>("trackerOnlyPreferCollection");
66  int arbitratedEcalSeededPreferCollection = conf_.getParameter<int>("arbitratedEcalSeededPreferCollection");
67  int arbitratedMergedPreferCollection = conf_.getParameter<int>("arbitratedMergedPreferCollection");
68  int arbitratedMergedEcalGeneralPreferCollection = conf_.getParameter<int>("arbitratedMergedEcalGeneralPreferCollection");
69 
70  // get Inputs
71  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
72  // this allows ConversionTrackMerger to be used as a cleaner only if handed just one list
73  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
74  //
75  const reco::ConversionTrackCollection *TC1 = 0;
76  static const reco::ConversionTrackCollection s_empty1, s_empty2;
78  e.getByToken(trackProducer1, trackCollection1);
79  if (trackCollection1.isValid()) {
80  TC1 = trackCollection1.product();
81  //std::cout << "1st collection " << trackProducer1 << " has "<< TC1->size() << " tracks" << std::endl ;
82  } else {
83  TC1 = &s_empty1;
84  edm::LogWarning("ConversionTrackMerger")
85  << "1st TrackCollection not found;"
86  << " will only clean 2nd TrackCollection ";
87  }
89 
90  const reco::ConversionTrackCollection *TC2 = 0;
92  e.getByToken(trackProducer2, trackCollection2);
93  if (trackCollection2.isValid()) {
94  TC2 = trackCollection2.product();
95  //std::cout << "2nd collection " << trackProducer2 << " has "<< TC2->size() << " tracks" << std::endl ;
96  } else {
97  TC2 = &s_empty2;
98  edm::LogWarning("ConversionTrackMerger")
99  << "2nd TrackCollection not found;"
100  <<" will only clean 1st TrackCollection ";
101  }
103 
104  // Step B: create empty output collection
105  outputTrks = std::auto_ptr<reco::ConversionTrackCollection>(new reco::ConversionTrackCollection);
106  int i;
107 
108  std::vector<int> selected1; for (unsigned int i=0; i<tC1.size(); ++i){selected1.push_back(1);}
109  std::vector<int> selected2; for (unsigned int i=0; i<tC2.size(); ++i){selected2.push_back(1);}
110 
111 
112  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh1;
113  std::map<reco::ConversionTrackCollection::const_iterator, std::vector<const TrackingRecHit*> > rh2;
114  for (reco::ConversionTrackCollection::const_iterator track=tC1.begin(); track!=tC1.end(); ++track){
115  trackingRecHit_iterator itB = track->track()->recHitsBegin();
116  trackingRecHit_iterator itE = track->track()->recHitsEnd();
117  for (trackingRecHit_iterator it = itB; it != itE; ++it) {
118  const TrackingRecHit* hit = &(**it);
119  rh1[track].push_back(hit);
120  }
121  }
122  for (reco::ConversionTrackCollection::const_iterator track=tC2.begin(); track!=tC2.end(); ++track){
123  trackingRecHit_iterator jtB = track->track()->recHitsBegin();
124  trackingRecHit_iterator jtE = track->track()->recHitsEnd();
125  for (trackingRecHit_iterator jt = jtB; jt != jtE; ++jt) {
126  const TrackingRecHit* hit = &(**jt);
127  rh2[track].push_back(hit);
128  }
129  }
130 
131  if ( (0<tC1.size())&&(0<tC2.size()) ){
132  i=-1;
133  for (reco::ConversionTrackCollection::iterator track=tC1.begin(); track!=tC1.end(); ++track){
134  i++;
135 
136  //clear flags if preferCollection was set to 0
137  selected1[i] = selected1[i] && outputPreferCollection!=0;
138  track->setIsTrackerOnly ( track->isTrackerOnly() && trackerOnlyPreferCollection!=0 );
139  track->setIsArbitratedEcalSeeded( track->isArbitratedEcalSeeded() && arbitratedEcalSeededPreferCollection!=0 );
140  track->setIsArbitratedMerged( track->isArbitratedMerged() && arbitratedMergedPreferCollection!=0 );
141  track->setIsArbitratedMergedEcalGeneral( track->isArbitratedMergedEcalGeneral() && arbitratedMergedEcalGeneralPreferCollection!=0 );
142 
143 
144  std::vector<const TrackingRecHit*>& iHits = rh1[track];
145  unsigned nh1 = iHits.size();
146  int j=-1;
147  for (reco::ConversionTrackCollection::iterator track2=tC2.begin(); track2!=tC2.end(); ++track2){
148  j++;
149 
150  //clear flags if preferCollection was set to 0
151  selected2[j] = selected2[j] && outputPreferCollection!=0;
152  track2->setIsTrackerOnly ( track2->isTrackerOnly() && trackerOnlyPreferCollection!=0 );
153  track2->setIsArbitratedEcalSeeded( track2->isArbitratedEcalSeeded() && arbitratedEcalSeededPreferCollection!=0 );
154  track2->setIsArbitratedMerged( track2->isArbitratedMerged() && arbitratedMergedPreferCollection!=0 );
155  track2->setIsArbitratedMergedEcalGeneral( track2->isArbitratedMergedEcalGeneral() && arbitratedMergedEcalGeneralPreferCollection!=0 );
156 
157  std::vector<const TrackingRecHit*>& jHits = rh2[track2];
158  unsigned nh2 = jHits.size();
159  int noverlap=0;
160  int firstoverlap=0;
161  for ( unsigned ih=0; ih<nh1; ++ih ) {
162  const TrackingRecHit* it = iHits[ih];
163  if (it->isValid()){
164  int jj=-1;
165  for ( unsigned jh=0; jh<nh2; ++jh ) {
166  const TrackingRecHit* jt = jHits[jh];
167  jj++;
168  if (jt->isValid()){
169  if ( it->sharesInput(jt,TrackingRecHit::some) ) {
170  noverlap++;
171  if ( allowFirstHitShare && ( ih == 0 ) && ( jh == 0 ) ) firstoverlap=1;
172  }
173  }
174  }
175  }
176  }
177  int nhit1 = track->track()->numberOfValidHits();
178  int nhit2 = track2->track()->numberOfValidHits();
179  //if (noverlap>0) printf("noverlap = %i, firstoverlap = %i, nhit1 = %i, nhit2 = %i, algo1 = %i, algo2 = %i, q1 = %i, q2 = %i\n",noverlap,firstoverlap,nhit1,nhit2,track->track()->algo(),track2->track()->algo(),track->track()->charge(),track2->track()->charge());
180  //std::cout << " trk1 trk2 nhits1 nhits2 nover " << i << " " << j << " " << track->numberOfValidHits() << " " << track2->numberOfValidHits() << " " << noverlap << " " << fi << " " << fj <<std::endl;
181  if ( (noverlap-firstoverlap) > (std::min(nhit1,nhit2)-firstoverlap)*shareFrac && (!checkCharge || track->track()->charge()*track2->track()->charge()>0) ) {
182  //printf("overlapping tracks\n");
183  //printf ("ndof1 = %5f, chisq1 = %5f, ndof2 = %5f, chisq2 = %5f\n",track->track()->ndof(),track->track()->chi2(),track2->track()->ndof(),track2->track()->chi2());
184 
185  double probFirst = ChiSquaredProbability(track->track()->chi2(), track->track()->ndof());
186  double probSecond = ChiSquaredProbability(track2->track()->chi2(), track2->track()->ndof());
187 
188  //arbitrate by number of hits and reduced chisq
189  bool keepFirst = ( nhit1>nhit2 || (nhit1==nhit2 && track->track()->normalizedChi2()<track2->track()->normalizedChi2()) );
190 
191  //override decision in case one track is radically worse quality than the other
192  keepFirst |= (probFirst>minProb && probSecond<=minProb);
193  keepFirst &= !(probFirst<=minProb && probSecond>minProb);
194 
195  bool keepSecond = !keepFirst;
196 
197  //set flags based on arbitration decision and precedence settings
198 
199  selected1[i] = selected1[i] && ( (keepFirst && outputPreferCollection==3) || outputPreferCollection==-1 || outputPreferCollection==1 );
200  track->setIsTrackerOnly ( track->isTrackerOnly() && ( (keepFirst && trackerOnlyPreferCollection==3) || trackerOnlyPreferCollection==-1 || trackerOnlyPreferCollection==1 ) );
201  track->setIsArbitratedEcalSeeded( track->isArbitratedEcalSeeded() && ( (keepFirst && arbitratedEcalSeededPreferCollection==3) || arbitratedEcalSeededPreferCollection==-1 || arbitratedEcalSeededPreferCollection==1 ) );
202  track->setIsArbitratedMerged( track->isArbitratedMerged() && ( (keepFirst && arbitratedMergedPreferCollection==3) || arbitratedMergedPreferCollection==-1 || arbitratedMergedPreferCollection==1 ) );
203  track->setIsArbitratedMergedEcalGeneral( track->isArbitratedMergedEcalGeneral() && ( (keepFirst && arbitratedMergedEcalGeneralPreferCollection==3) || arbitratedMergedEcalGeneralPreferCollection==-1 || arbitratedMergedEcalGeneralPreferCollection==1 ) );
204 
205  selected2[j] = selected2[j] && ( (keepSecond && outputPreferCollection==3) || outputPreferCollection==-1 || outputPreferCollection==2 );
206  track2->setIsTrackerOnly ( track2->isTrackerOnly() && ( (keepSecond && trackerOnlyPreferCollection==3) || trackerOnlyPreferCollection==-1 || trackerOnlyPreferCollection==2 ) );
207  track2->setIsArbitratedEcalSeeded( track2->isArbitratedEcalSeeded() && ( (keepSecond && arbitratedEcalSeededPreferCollection==3) || arbitratedEcalSeededPreferCollection==-1 || arbitratedEcalSeededPreferCollection==2 ) );
208  track2->setIsArbitratedMerged( track2->isArbitratedMerged() && ( (keepSecond && arbitratedMergedPreferCollection==3) || arbitratedMergedPreferCollection==-1 || arbitratedMergedPreferCollection==2 ) );
209  track2->setIsArbitratedMergedEcalGeneral( track2->isArbitratedMergedEcalGeneral() && ( (keepSecond && arbitratedMergedEcalGeneralPreferCollection==3) || arbitratedMergedEcalGeneralPreferCollection==-1 || arbitratedMergedEcalGeneralPreferCollection==2 ) );
210 
211  }//end got a duplicate
212  }//end track2 loop
213  }//end track loop
214  }//end more than 1 track
215 
216  //
217  // output selected tracks - if any
218  //
219 
220  if ( 0<tC1.size() ){
221  i=0;
222  for (reco::ConversionTrackCollection::const_iterator track=tC1.begin(); track!=tC1.end();
223  ++track, ++i){
224  //don't store tracks rejected as duplicates
225  if (!selected1[i]){
226  continue;
227  }
228  //fill the TrackCollection
229  outputTrks->push_back(*track);
230  }//end faux loop over tracks
231  }//end more than 0 track
232 
233  //Fill the trajectories, etc. for 1st collection
234 
235 
236  if ( 0<tC2.size() ){
237  i=0;
238  for (reco::ConversionTrackCollection::const_iterator track=tC2.begin(); track!=tC2.end();
239  ++track, ++i){
240  //don't store tracks rejected as duplicates
241  if (!selected2[i]){
242  continue;
243  }
244  //fill the TrackCollection
245  outputTrks->push_back(*track);
246  }//end faux loop over tracks
247  }//end more than 0 track
248 
249  e.put(outputTrks);
250  return;
251 
252  }//end produce
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::auto_ptr< reco::ConversionTrackCollection > outputTrks
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer1
std::vector< ConversionTrack > ConversionTrackCollection
collection of ConversionTracks
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
int j
Definition: DBlmapReader.cc:9
float ChiSquaredProbability(double chiSquared, double nrDOF)
ConversionTrackMerger(const edm::ParameterSet &conf)
T min(T a, T b)
Definition: MathUtil.h:58
bool isValid() const
Definition: HandleBase.h:75
tuple conf
Definition: dbtoconf.py:185
T const * product() const
Definition: Handle.h:81
bool isValid() const
virtual void produce(edm::Event &e, const edm::EventSetup &c)
edm::EDGetTokenT< reco::ConversionTrackCollection > trackProducer2
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection