CMS 3D CMS Logo

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