CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HICaloCompatibleTrackSelector.cc
Go to the documentation of this file.
1 /*
2 
3  Based on analytical track selector
4 
5  - This track selector assigns a quality bit to tracks deemed compatible with matching calo info
6  - The default mode is to use the matching provided by particle flow,
7  but a delta R matching to calorimeter towers is also supported
8  - No selection is done other then selecting calo-compatible tracks.
9  - The code always keeps all tracks in the input collection (should make configurable)
10  - Note that matching by PF candidate only works on the same track collection used as input to PF
11  - Tower code not up to data
12 
13  Authors: Matthew Nguyen, Andre Yoon, Frank Ma (November 4th, 2011)
14 
15  */
16 
17 
18 // Basic inclusion
20 
24 
27 
29 #include <Math/DistFunc.h>
30 #include "TMath.h"
31 
32 
34 
35 HICaloCompatibleTrackSelector::HICaloCompatibleTrackSelector( const edm::ParameterSet & cfg ) :
36  srcTracks_(consumes<TrackCollection>(cfg.getParameter<edm::InputTag>("srcTracks"))),
37  srcPFCands_(consumes<PFCandidateCollection>(cfg.getParameter<edm::InputTag>("srcPFCands"))),
38  srcTower_(consumes<CaloTowerCollection>(cfg.getParameter<edm::InputTag>("srcTower"))),
39  usePFCandMatching_(cfg.getUntrackedParameter<bool>("usePFCandMatching", true)),
40  trkMatchPtMin_(cfg.getUntrackedParameter<double>("trkMatchPtMin",10.0)),
41  trkCompPtMin_(cfg.getUntrackedParameter<double>("trkCompPtMin",35.0)),
42  trkEtaMax_(cfg.getUntrackedParameter<double>("trkEtaMax",2.4)),
43  towerPtMin_(cfg.getUntrackedParameter<double>("towerPtMin",5.0)),
44  matchConeRadius_(cfg.getUntrackedParameter<double>("matchConeRadius",0.087)),
45  keepAllTracks_(cfg.getUntrackedParameter<bool>("keepAllTracks", true)),
46  copyExtras_(cfg.getUntrackedParameter<bool>("copyExtras", true)),
47  copyTrajectories_(cfg.getUntrackedParameter<bool>("copyTrajectories", true)),
48  qualityToSet_(cfg.getParameter<std::string>("qualityToSet")),
49  qualityToSkip_(cfg.getParameter<std::string>("qualityToSkip")),
50  qualityToMatch_(cfg.getParameter<std::string>("qualityToMatch")),
51  minimumQuality_(cfg.getParameter<std::string>("minimumQuality")),
52  resetQuality_(cfg.getUntrackedParameter<bool>("resetQuality", true)),
53  passMuons_(cfg.getUntrackedParameter<bool>("passMuons", true)),
54  passElectrons_(cfg.getUntrackedParameter<bool>("passElectrons", false)),
55  funcDeltaRTowerMatch_(cfg.getParameter<std::string>("funcDeltaRTowerMatch")),
56  funcCaloComp_(cfg.getParameter<std::string>("funcCaloComp"))
57 {
58  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
59  produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
60  if (copyExtras_) {
61  produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
62  produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
63  }
64  if (copyTrajectories_) {
65  produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
66  produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
67  srcTrackTrajs_ = (consumes<std::vector<Trajectory> >(cfg.getParameter<edm::InputTag>("srcTracks")));
68  srcTrackTrajAssoc_ = (consumes<TrajTrackAssociationCollection>(cfg.getParameter<edm::InputTag>("srcTracks")));
69  }
70 
71  // pt dependence of delta R matching requirement
72  fDeltaRTowerMatch = new TF1("fDeltaRTowerMatch",funcDeltaRTowerMatch_.c_str(),0,200);
73  // pt dependance of calo compatibility, i.e., minimum sum Calo Et vs. track pT
74  fCaloComp = new TF1("fCaloComp",funcCaloComp_.c_str(),0,200); // a parameterization of pt dependent cut
75 
76 
77 }
78 
80 }
81 
82 
84 {
85  using namespace std;
86  using namespace edm;
87  using namespace reco;
88 
89  LogDebug("HICaloCompatibleTrackSelector")<<"min pt for selection = "<< trkMatchPtMin_<<endl;
90 
91 
92  Handle<TrackCollection> hSrcTrack;
96 
97  evt.getByToken(srcTracks_,hSrcTrack);
98 
99  selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
101  if (copyExtras_) {
102  selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
103  selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
106  }
107 
108 
109  if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
110 
111 
114 
115  bool isPFThere = false;
116  bool isTowerThere = false;
117 
118  if(usePFCandMatching_) isPFThere = evt.getByToken(srcPFCands_, pfCandidates);
119  else isTowerThere = evt.getByToken(srcTower_, towers);
120 
121  size_t current = 0;
122  for (TI ti = hSrcTrack->begin(), ed = hSrcTrack->end(); ti != ed; ++ti, ++current) {
123 
124 
125  const reco::Track& trk = *ti;
126 
127  bool isSelected;
128  if(usePFCandMatching_) isSelected = selectByPFCands(ti, hSrcTrack, pfCandidates, isPFThere);
129  else isSelected = selectByTowers(ti, hSrcTrack, towers, isTowerThere);
130 
131  if(!keepAllTracks_ && !isSelected) continue;
132 
133  // Add all tracks to output collection, the rest of the code only sets the quality bit
134  selTracks_->push_back( reco::Track( trk ) ); // clone and store
135 
136 
137  if(isSelected) selTracks_->back().setQuality(reco::TrackBase::qualityByName(qualityToSet_.c_str()));
138 
139 
140  if (copyExtras_) {
141  // TrackExtras
142  selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
143  trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
144  trk.outerStateCovariance(), trk.outerDetId(),
145  trk.innerStateCovariance(), trk.innerDetId(),
146  trk.seedDirection(), trk.seedRef() ) );
147  selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
148  TrackExtra & tx = selTrackExtras_->back();
149  tx.setResiduals(trk.residuals());
150  // TrackingRecHits
151  for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
152  selHits_->push_back( (*hit)->clone() );
153  tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
154  }
155  }
156  if (copyTrajectories_) {
157  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
158  }
159 
160 
161  } // close track loop
162 
163  if ( copyTrajectories_ ) {
166  evt.getByToken(srcTrackTrajs_, hTTAss);
167  evt.getByToken(srcTrackTrajAssoc_, hTraj);
168  selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>());
169  rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
170  selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
171  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
172  Ref< vector<Trajectory> > trajRef(hTraj, i);
174  if (match != hTTAss->end()) {
175  const Ref<TrackCollection> &trkRef = match->val;
176  short oldKey = static_cast<short>(trkRef.key());
177  if (trackRefs_[oldKey].isNonnull()) {
178  selTrajs_->push_back( Trajectory(*trajRef) );
179  selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
180  }
181  }
182  }
183  }
184 
185  static const string emptyString;
186  evt.put(selTracks_);
187  if (copyExtras_ ) {
188  evt.put(selTrackExtras_);
189  evt.put(selHits_);
190  }
191  if ( copyTrajectories_ ) {
192  evt.put(selTrajs_);
193  evt.put(selTTAss_);
194  }
195 }
196 
197 
198 
199 bool
201 {
202 
203  const reco::Track& trk = *ti;
204 
205  // If it passes this quality threshold or is under the minimum match pT, automatically save it
207  return true;
208  }
210  return false;
211  }
212  else
213  {
214 
215  double trkPt = trk.pt();
216  //if(trkPt < trkMatchPtMin_ ) return false;
217 
218  double caloEt = 0.0;
219  if(usePFCandMatching_){
220  if(isPFThere){
221  unsigned int trackKey = ti - hSrcTrack->begin();
222  caloEt = matchPFCandToTrack(pfCandidates, trackKey, trkPt);
223  }
224  }
225 
226  // Set quality bit based on calo matching
227  if(!(caloEt>0.)) return false;
228 
229 
230  if(trkPt<=trkCompPtMin_){
232  else return false;
233  }
234  else{
235  // loose cuts are implied in selectors, make configurable?
236  float compPt = (fCaloComp->Eval(trkPt)!=fCaloComp->Eval(trkPt)) ? 0 : fCaloComp->Eval(trkPt); // protect agains NaN
237  if(caloEt>compPt) return true;
238  else return false;
239  }
240  } // else above trkMatchPtMin_
241 
242  throw cms::Exception("Undefined case in HICaloCompatibleTrackSelector") << "Undefined case in HICaloCompatibleTrackSelector";
243 }
244 
245 
246 bool
248 {
249 
250  // Not up to date! use PF towers instead
251 
252  const reco::Track& trk = *ti;
253 
254  // If it passes the high purity cuts, then consider it confirmed
256  return true;
257  }
258  else{
259  if(trk.pt() < trkMatchPtMin_ || !trk.quality(reco::TrackBase::qualityByName(qualityToMatch_))) return false;
260 
261  double caloEt = 0.0;
262  if(isTowerThere){
263  double matchDr;
264  matchByDrAllowReuse(trk,towers,matchDr,caloEt);
265  float matchConeRadius_pt = (fDeltaRTowerMatch->Eval(trk.pt())!=fDeltaRTowerMatch->Eval(trk.pt())) ? 0 : fDeltaRTowerMatch->Eval(trk.pt()); // protect agains NaN
266  if (caloEt>0 && matchDr>matchConeRadius_pt) caloEt=0.;
267  }
268 
269  if(trk.pt()<trkCompPtMin_||caloEt>0.75*(trk.pt()-trkCompPtMin_)) return true;
270  else return false;
271  }
272 }
273 
274 double
276 {
277 
278  // This function returns the sum of the calo energy in the block containing the track
279  // There is another piece of information which could be useful which is the pT assigned to the charged hadron by PF
280 
281 
282  double sum_ecal=0.0, sum_hcal=0.0;
283 
284  int candType = 0;
285  reco::PFCandidate matchedCand;
286 
287  // loop over the PFCandidates until you find the one whose trackRef points to the track
288 
289  for( CI ci = pfCandidates->begin(); ci!=pfCandidates->end(); ++ci) {
290 
291  const reco::PFCandidate& cand = *ci;
292 
293  int type = cand.particleId();
294 
295  // only charged hadrons and leptons can be asscociated with a track
296  if(!(type == PFCandidate::h || //type1
297  type == PFCandidate::e || //type2
298  type == PFCandidate::mu //type3
299  )
300  ) continue;
301 
302 
303  unsigned candTrackRefKey = cand.trackRef().key();
304 
305  if(it==candTrackRefKey) {
306  matchedCand = cand;
307  candType = type;
308  break;
309  }
310  }
311  // take all muons as compatible, extend to electrons when validataed
312  if(passMuons_ && candType==3) return 9999.;
313  if(passElectrons_ && candType==2) return 9999.;
314 
315  if(trkPt < trkMatchPtMin_ ) return 0.;
316 
317  if(candType>0){
318 
319  // Now that we found the matched PF candidate, loop over the elements in the block summing the calo Et
320  for(unsigned ib=0; ib<matchedCand.elementsInBlocks().size(); ib++) {
321 
322  PFBlockRef blockRef = matchedCand.elementsInBlocks()[ib].first;
323 
324  unsigned indexInBlock = matchedCand.elementsInBlocks()[ib].second;
325  const edm::OwnVector< reco::PFBlockElement>& elements = (*blockRef).elements();
326 
327  //This tells you what type of element it is:
328  //cout<<" block type"<<elements[indexInBlock].type()<<endl;
329 
330  switch (elements[indexInBlock].type()) {
331 
332  case PFBlockElement::ECAL: {
333  reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
334  sum_ecal += clusterRef->energy()/cosh(clusterRef->eta());
335  break;
336  }
337 
338  case PFBlockElement::HCAL: {
339  reco::PFClusterRef clusterRef = elements[indexInBlock].clusterRef();
340  sum_hcal += clusterRef->energy()/cosh(clusterRef->eta());
341  break;
342  }
343  case PFBlockElement::TRACK: {
344  //Do nothing since track are not normally linked to other tracks
345  break;
346  }
347  default:
348  break;
349  }
350 
351  // We could also add in the PS, HO, ..
352 
353  } // end of elementsInBlocks()
354  } // end of isCandFound
355 
356 
357 
358  return sum_ecal+sum_hcal;
359 
360 }
361 
362 void HICaloCompatibleTrackSelector::matchByDrAllowReuse(const reco::Track & trk, const edm::Handle<CaloTowerCollection> & towers, double & bestdr, double & bestpt)
363 {
364  // loop over towers
365  bestdr=1e10;
366  bestpt=0.;
367  for(unsigned int i = 0; i < towers->size(); ++i){
368  const CaloTower & tower= (*towers)[i];
369  double pt = tower.pt();
370  if (pt<towerPtMin_) continue;
371  if (fabs(tower.eta())>trkEtaMax_) continue;
372  double dr = reco::deltaR(tower,trk);
373  if (dr<bestdr) {
374  bestdr = dr;
375  bestpt = pt;
376  }
377  }
378 }
#define LogDebug(id)
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< CaloTowerCollection > srcTower_
tuple cfg
Definition: looper.py:237
int ib
Definition: cuy.py:660
std::auto_ptr< TrackingRecHitCollection > selHits_
virtual float pt() const
transverse momentum
void produce(edm::Event &evt, const edm::EventSetup &es)
process one event
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
string emptyString
Definition: archive.py:29
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:50
list elements
Definition: asciidump.py:414
edm::Ref< TrackExtraCollection > TrackExtraRef
persistent reference to a TrackExtra
Definition: TrackExtraFwd.h:17
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
std::auto_ptr< reco::TrackExtraCollection > selTrackExtras_
edm::EDGetTokenT< reco::TrackCollection > srcTracks_
source collection label
edm::EDGetTokenT< std::vector< Trajectory > > srcTrackTrajs_
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
reco::PFCandidateCollection::const_iterator CI
void add(const TrackingRecHitRef &ref)
add a reference to a RecHit
std::auto_ptr< std::vector< Trajectory > > selTrajs_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual float eta() const
momentum pseudorapidity
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:75
double pt() const
track transverse momentum
Definition: TrackBase.h:669
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:94
bool selectByTowers(TI ti, const edm::Handle< TrackCollection > hSrcTrack, const edm::Handle< CaloTowerCollection > towers, bool isTowerThere)
bool copyExtras_
copy only the tracks, not extras and rechits (for AOD)
edm::RefProd< std::vector< Trajectory > > rTrajectories_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCands_
RefProd< PROD > getRefBeforePut()
Definition: Event.h:133
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
edm::OwnVector< TrackingRecHit > TrackingRecHitCollection
collection of TrackingRecHits
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:213
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:108
double matchPFCandToTrack(const edm::Handle< PFCandidateCollection > &pfCandidates, unsigned it, double trkPt)
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:70
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:45
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:19
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:80
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:551
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
std::auto_ptr< reco::TrackCollection > selTracks_
storage
std::auto_ptr< TrajTrackAssociationCollection > selTTAss_
const TrackResiduals & residuals() const
Definition: Track.h:220
bool selectByPFCands(TI ti, const edm::Handle< TrackCollection > hSrcTrack, const edm::Handle< PFCandidateCollection > pfCandidates, bool isPFThere)
void matchByDrAllowReuse(const reco::Track &trk, const edm::Handle< CaloTowerCollection > &towers, double &bestdr, double &bestpt)
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:204
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
volatile std::atomic< bool > shutdown_flag false
virtual ParticleType particleId() const
Definition: PFCandidate.h:369
bool copyTrajectories_
copy also trajectories and trajectory-&gt;track associations
const ElementsInBlocks & elementsInBlocks() const
Definition: PFCandidate.cc:682
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
edm::EDGetTokenT< TrajTrackAssociationCollection > srcTrackTrajAssoc_
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:178