CMS 3D CMS Logo

CosmicTrackSelector.cc
Go to the documentation of this file.
1 
11 #include <utility>
12 #include <vector>
13 #include <memory>
14 #include <algorithm>
15 #include <map>
21 
29 
30 using namespace reco;
31 
33 private:
34 public:
35  // constructor
36  explicit CosmicTrackSelector(const edm::ParameterSet &cfg);
37  // destructor
38  ~CosmicTrackSelector() override;
39 
40 private:
42  // process one event
43  void produce(edm::Event &evt, const edm::EventSetup &es) override;
44  // return class, or -1 if rejected
45  bool select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk);
46  // source collection label
49  // copy only the tracks, not extras and rechits (for AOD)
51  // copy also trajectories and trajectory->track associations
55 
56  // save all the tracks
58  // do I have to set a quality bit?
61 
62  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
63  std::vector<double> res_par_;
64  double chi2n_par_;
65 
66  // Impact parameter absolute cuts
67  double max_d0_;
68  double max_z0_;
69  // Trackk parameter cuts
70  double min_pt_;
71  double max_eta_;
72  // Cut on number of valid hits
73  uint32_t min_nHit_;
74  // Cut on number of valid Pixel hits
75  uint32_t min_nPixelHit_;
76  // Cuts on numbers of layers with hits/3D hits/lost hits.
77  uint32_t min_layers_;
78  uint32_t min_3Dlayers_;
79  uint32_t max_lostLayers_;
80 
81  // storage
82  std::unique_ptr<reco::TrackCollection> selTracks_;
83  std::unique_ptr<reco::TrackExtraCollection> selTrackExtras_;
84  std::unique_ptr<TrackingRecHitCollection> selHits_;
85  std::unique_ptr<std::vector<Trajectory>> selTrajs_;
86  std::unique_ptr<std::vector<const Trajectory *>> selTrajPtrs_;
87  std::unique_ptr<TrajTrackAssociationCollection> selTTAss_;
92  std::vector<reco::TrackRef> trackRefs_;
93 };
94 
95 #include <Math/DistFunc.h>
96 #include "TMath.h"
97 
99  : src_(consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>("src"))),
100  beamspot_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamspot"))),
101  copyExtras_(cfg.getUntrackedParameter<bool>("copyExtras", false)),
102  copyTrajectories_(cfg.getUntrackedParameter<bool>("copyTrajectories", false)),
103  keepAllTracks_(cfg.exists("keepAllTracks") ? cfg.getParameter<bool>("keepAllTracks")
104  : false), // as this is what you expect from a well behaved selector
105  setQualityBit_(false),
106  qualityToSet_(TrackBase::undefQuality),
107  chi2n_par_(cfg.getParameter<double>("chi2n_par")),
108  // Impact parameter absolute cuts.
109  max_d0_(cfg.getParameter<double>("max_d0")),
110  max_z0_(cfg.getParameter<double>("max_z0")),
111  // Track parameter cuts.
112  min_pt_(cfg.getParameter<double>("min_pt")),
113  max_eta_(cfg.getParameter<double>("max_eta")),
114  // Cut on number of valid hits
115  min_nHit_(cfg.getParameter<uint32_t>("min_nHit")),
116  // Cut on number of valid hits
117  min_nPixelHit_(cfg.getParameter<uint32_t>("min_nPixelHit")),
118  // Cuts on numbers of layers with hits/3D hits/lost hits.
119  min_layers_(cfg.getParameter<uint32_t>("minNumberLayers")),
120  min_3Dlayers_(cfg.getParameter<uint32_t>("minNumber3DLayers")),
121  max_lostLayers_(cfg.getParameter<uint32_t>("maxNumberLostLayers")) {
122  if (cfg.exists("qualityBit")) {
123  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
124  if (!qualityStr.empty()) {
125  setQualityBit_ = true;
126  qualityToSet_ = TrackBase::qualityByName(cfg.getParameter<std::string>("qualityBit"));
127  }
128  }
130  throw cms::Exception("Configuration")
131  << "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
132  if (setQualityBit_ && (qualityToSet_ == TrackBase::undefQuality))
133  throw cms::Exception("Configuration")
134  << "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit")
135  << " as it is 'undefQuality' or unknown.\n";
136 
137  std::string alias(cfg.getParameter<std::string>("@module_label"));
138  produces<reco::TrackCollection>().setBranchAlias(alias + "Tracks");
139  if (copyExtras_) {
140  produces<reco::TrackExtraCollection>().setBranchAlias(alias + "TrackExtras");
141  produces<TrackingRecHitCollection>().setBranchAlias(alias + "RecHits");
142  }
143  if (copyTrajectories_) {
144  srcTraj_ = consumes<std::vector<Trajectory>>(cfg.getParameter<edm::InputTag>("src"));
145  srcTass_ = consumes<TrajTrackAssociationCollection>(cfg.getParameter<edm::InputTag>("src"));
146  produces<std::vector<Trajectory>>().setBranchAlias(alias + "Trajectories");
147  produces<TrajTrackAssociationCollection>().setBranchAlias(alias + "TrajectoryTrackAssociations");
148  }
149 }
150 
152 
154  using namespace std;
155  using namespace edm;
156  using namespace reco;
157 
158  Handle<TrackCollection> hSrcTrack;
162 
163  // looking for the beam spot
165  evt.getByToken(beamspot_, hBsp);
166  reco::BeamSpot vertexBeamSpot;
167  vertexBeamSpot = *hBsp;
168 
169  // Get tracks
170  evt.getByToken(src_, hSrcTrack);
171 
172  selTracks_ = std::make_unique<TrackCollection>();
174  if (copyExtras_) {
175  selTrackExtras_ = std::make_unique<TrackExtraCollection>();
176  selHits_ = std::make_unique<TrackingRecHitCollection>();
179  }
180 
181  if (copyTrajectories_)
182  trackRefs_.resize(hSrcTrack->size());
183 
184  // Loop over tracks
185  size_t current = 0;
186  for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
187  const Track &trk = *it;
188  // Check if this track passes cuts
189  bool ok = select(vertexBeamSpot, trk);
190  if (!ok) {
191  if (copyTrajectories_)
192  trackRefs_[current] = reco::TrackRef();
193  if (!keepAllTracks_)
194  continue;
195  }
196  selTracks_->push_back(Track(trk)); // clone and store
197  if (ok && setQualityBit_)
198  selTracks_->back().setQuality(qualityToSet_);
199  if (copyExtras_) {
200  // TrackExtras
201  selTrackExtras_->push_back(TrackExtra(trk.outerPosition(),
202  trk.outerMomentum(),
203  trk.outerOk(),
204  trk.innerPosition(),
205  trk.innerMomentum(),
206  trk.innerOk(),
207  trk.outerStateCovariance(),
208  trk.outerDetId(),
209  trk.innerStateCovariance(),
210  trk.innerDetId(),
211  trk.seedDirection(),
212  trk.seedRef()));
213  selTracks_->back().setExtra(TrackExtraRef(rTrackExtras_, selTrackExtras_->size() - 1));
214  TrackExtra &tx = selTrackExtras_->back();
215  tx.setResiduals(trk.residuals());
216  // TrackingRecHits
217  auto const firstHitIndex = selHits_->size();
218  for (trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++hit) {
219  selHits_->push_back((*hit)->clone());
220  }
221  tx.setHits(rHits_, firstHitIndex, selHits_->size() - firstHitIndex);
222  tx.setTrajParams(trk.extra()->trajParams(), trk.extra()->chi2sX5());
223  }
224  if (copyTrajectories_) {
225  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
226  }
227  }
228  if (copyTrajectories_) {
231  evt.getByToken(srcTass_, hTTAss);
232  evt.getByToken(srcTraj_, hTraj);
233  selTrajs_ = std::make_unique<std::vector<Trajectory>>();
234  rTrajectories_ = evt.getRefBeforePut<vector<Trajectory>>();
235  selTTAss_ = std::make_unique<TrajTrackAssociationCollection>(&evt.productGetter());
236  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
237  Ref<vector<Trajectory>> trajRef(hTraj, i);
239  if (match != hTTAss->end()) {
240  const Ref<TrackCollection> &trkRef = match->val;
241  short oldKey = static_cast<short>(trkRef.key());
242  if (trackRefs_[oldKey].isNonnull()) {
243  selTrajs_->push_back(Trajectory(*trajRef));
244  selTTAss_->insert(Ref<vector<Trajectory>>(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey]);
245  }
246  }
247  }
248  }
249 
250  static const std::string emptyString;
251  evt.put(std::move(selTracks_));
252  if (copyExtras_) {
254  evt.put(std::move(selHits_));
255  }
256  if (copyTrajectories_) {
257  evt.put(std::move(selTrajs_));
258  evt.put(std::move(selTTAss_));
259  }
260 }
261 
262 bool CosmicTrackSelector::select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk) {
263  // Decide if the given track passes selection cuts.
264 
265  using namespace std;
266 
267  // Cuts on numbers of layers with hits/3D hits/lost hits.
268  uint32_t nlayers = tk.hitPattern().trackerLayersWithMeasurement();
269  uint32_t nlayers3D =
272 
273  // Get the number of valid hits and PixelHits
274  uint32_t nHit = 0;
275  uint32_t nPixelHit = 0;
277  if (!((*recHit)->isValid()))
278  continue;
279  ++nHit;
280  DetId id((*recHit)->geographicalId());
281  if ((unsigned int)id.subdetId() == PixelSubdetector::PixelBarrel ||
282  (unsigned int)id.subdetId() == PixelSubdetector::PixelEndcap)
283  ++nPixelHit;
284  }
285 
286  // Cut on the number of valid hits
287  if (nHit < min_nHit_)
288  return false;
289  // Cut on the number of valid Pixel hits
290  if (nPixelHit < min_nPixelHit_)
291  return false;
292  if (nlayers < min_layers_)
293  return false;
294  if (nlayers3D < min_3Dlayers_)
295  return false;
296  if (nlayersLost > max_lostLayers_)
297  return false;
298 
299  // Get track parameters
300  double pt = tk.pt(), eta = tk.eta(), chi2n = tk.normalizedChi2();
301  double d0 = -tk.dxy(vertexBeamSpot.position()), dz = tk.dz();
302 
303  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
304  if (abs(d0) > max_d0_)
305  return false;
306  if (abs(dz) > max_z0_)
307  return false;
308 
309  // optimized cuts adapted to the track eta, pt and chiquare/ndof
310  if (abs(eta) > max_eta_)
311  return false;
312  if (pt < min_pt_)
313  return false;
314  if (chi2n > chi2n_par_ * nlayers)
315  return false;
316 
317  else
318  return true;
319 }
320 
323 
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 Point & position() const
position
Definition: BeamSpot.h:59
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:65
Quality qualityByName(std::string const &name)
void produce(edm::Event &evt, const edm::EventSetup &es) override
TrackQuality
track quality
Definition: TrackBase.h:150
TrackBase::TrackQuality qualityToSet_
edm::RefProd< std::vector< Trajectory > > rTrajectories_
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
reco::TrackExtraRefProd rTrackExtras_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::BeamSpot > beamspot_
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< std::vector< Trajectory > > srcTraj_
key_type key() const
Accessor for product key.
Definition: Ref.h:250
const_iterator find(const key_type &k) const
find element with specified reference key
std::unique_ptr< TrackingRecHitCollection > selHits_
const_iterator end() const
last iterator over the map (read only)
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:71
EDProductGetter const & productGetter() const
Definition: Event.cc:106
double pt() const
track transverse momentum
Definition: TrackBase.h:637
std::vector< double > res_par_
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:369
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
std::unique_ptr< std::vector< const Trajectory * > > selTrajPtrs_
std::unique_ptr< reco::TrackExtraCollection > selTrackExtras_
TrackingRecHitRefProd rHits_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:91
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
RefProd< PROD > getRefBeforePut()
Definition: Event.h:158
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:53
#define dso_hidden
Definition: Visibility.h:12
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
Definition: DetId.h:17
std::unique_ptr< TrajTrackAssociationCollection > selTTAss_
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 constexpr float d0
std::unique_ptr< std::vector< Trajectory > > selTrajs_
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:553
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
bool select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk)
CosmicTrackSelector(const edm::ParameterSet &cfg)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
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
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
static std::string const emptyString("")
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::TrackCollection > src_
edm::EDGetTokenT< TrajTrackAssociationCollection > srcTass_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:56
reco::TrackRefProd rTracks_
std::vector< reco::TrackRef > trackRefs_
def move(src, dest)
Definition: eostools.py:511
std::unique_ptr< reco::TrackCollection > selTracks_
void setTrajParams(TrajParams tmps, Chi2sFive chi2s)
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:50
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:132