CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  virtual ~CosmicTrackSelector() ;
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::auto_ptr<reco::TrackCollection> selTracks_;
83  std::auto_ptr<reco::TrackExtraCollection> selTrackExtras_;
84  std::auto_ptr< TrackingRecHitCollection> selHits_;
85  std::auto_ptr< std::vector<Trajectory> > selTrajs_;
86  std::auto_ptr< std::vector<const Trajectory *> > selTrajPtrs_;
87  std::auto_ptr< TrajTrackAssociationCollection > selTTAss_;
92  std::vector<reco::TrackRef> trackRefs_;
93 
94  };
95 
96 #include <Math/DistFunc.h>
97 #include "TMath.h"
98 
100  src_( consumes<reco::TrackCollection>( cfg.getParameter<edm::InputTag>( "src" ) ) ),
101  beamspot_( consumes<reco::BeamSpot>( cfg.getParameter<edm::InputTag>( "beamspot" ) ) ),
102  copyExtras_(cfg.getUntrackedParameter<bool>("copyExtras", false)),
103  copyTrajectories_(cfg.getUntrackedParameter<bool>("copyTrajectories", false)),
104  keepAllTracks_( cfg.exists("keepAllTracks") ?
105  cfg.getParameter<bool>("keepAllTracks") :
106  false ), // as this is what you expect from a well behaved selector
107  setQualityBit_( false ),
108  qualityToSet_( TrackBase::undefQuality ),
109  chi2n_par_( cfg.getParameter<double>("chi2n_par") ),
110  // Impact parameter absolute cuts.
111  max_d0_(cfg.getParameter<double>("max_d0")),
112  max_z0_(cfg.getParameter<double>("max_z0")),
113  // Track parameter cuts.
114  min_pt_(cfg.getParameter<double>("min_pt")),
115  max_eta_(cfg.getParameter<double>("max_eta")),
116  // Cut on number of valid hits
117  min_nHit_(cfg.getParameter<uint32_t>("min_nHit")),
118  // Cut on number of valid hits
119  min_nPixelHit_(cfg.getParameter<uint32_t>("min_nPixelHit")),
120  // Cuts on numbers of layers with hits/3D hits/lost hits.
121  min_layers_(cfg.getParameter<uint32_t>("minNumberLayers") ),
122  min_3Dlayers_(cfg.getParameter<uint32_t>("minNumber3DLayers") ),
123  max_lostLayers_(cfg.getParameter<uint32_t>("maxNumberLostLayers") )
124 {
125  if (cfg.exists("qualityBit")) {
126  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
127  if (qualityStr != "") {
128  setQualityBit_ = true;
129  qualityToSet_ = TrackBase::qualityByName(cfg.getParameter<std::string>("qualityBit"));
130  }
131  }
132  if (keepAllTracks_ && !setQualityBit_) throw cms::Exception("Configuration") <<
133  "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
134  if (setQualityBit_ && (qualityToSet_ == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
135  "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " 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 }
151 
153 }
154 
156 {
157  using namespace std;
158  using namespace edm;
159  using namespace reco;
160 
161  Handle<TrackCollection> hSrcTrack;
165 
166  // looking for the beam spot
168  evt.getByToken(beamspot_, hBsp);
169  reco::BeamSpot vertexBeamSpot;
170  vertexBeamSpot = *hBsp;
171 
172  // Get tracks
173  evt.getByToken( src_, hSrcTrack );
174 
175  selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
177  if (copyExtras_) {
178  selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
179  selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
182  }
183 
184  if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
185 
186  // Loop over tracks
187  size_t current = 0;
188  for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
189  const Track & trk = * it;
190  // Check if this track passes cuts
191  bool ok = select(vertexBeamSpot, trk);
192  if (!ok) {
193  if (copyTrajectories_) trackRefs_[current] = reco::TrackRef();
194  if (!keepAllTracks_) continue;
195  }
196  selTracks_->push_back( Track( trk ) ); // clone and store
197  if (ok && setQualityBit_) selTracks_->back().setQuality(qualityToSet_);
198  if (copyExtras_) {
199  // TrackExtras
200  selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
201  trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
202  trk.outerStateCovariance(), trk.outerDetId(),
203  trk.innerStateCovariance(), trk.innerDetId(),
204  trk.seedDirection(), trk.seedRef() ) );
205  selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
206  TrackExtra & tx = selTrackExtras_->back();
207  tx.setResiduals(trk.residuals());
208  // TrackingRecHits
209  auto const firstHitIndex = selHits_->size();
210  for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
211  selHits_->push_back( (*hit)->clone() );
212  }
213  tx.setHits( rHits_, firstHitIndex, selHits_->size() - firstHitIndex );
214  }
215  if (copyTrajectories_) {
216  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
217  }
218  }
219  if ( copyTrajectories_ ) {
222  evt.getByToken(srcTass_, hTTAss);
223  evt.getByToken(srcTraj_, hTraj);
224  selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>());
225  rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
226  selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection(&evt.productGetter()));
227  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
228  Ref< vector<Trajectory> > trajRef(hTraj, i);
230  if (match != hTTAss->end()) {
231  const Ref<TrackCollection> &trkRef = match->val;
232  short oldKey = static_cast<short>(trkRef.key());
233  if (trackRefs_[oldKey].isNonnull()) {
234  selTrajs_->push_back( Trajectory(*trajRef) );
235  selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
236  }
237  }
238  }
239  }
240 
241  static const std::string emptyString;
242  evt.put(selTracks_);
243  if (copyExtras_ ) {
244  evt.put(selTrackExtras_);
245  evt.put(selHits_);
246  }
247  if ( copyTrajectories_ ) {
248  evt.put(selTrajs_);
249  evt.put(selTTAss_);
250  }
251 }
252 
253 
254 bool CosmicTrackSelector::select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk) {
255  // Decide if the given track passes selection cuts.
256 
257  using namespace std;
258 
259  // Cuts on numbers of layers with hits/3D hits/lost hits.
261  uint32_t nlayers3D = tk.hitPattern().pixelLayersWithMeasurement() +
264 
265  // Get the number of valid hits and PixelHits
266  uint32_t nHit = 0;
267  uint32_t nPixelHit = 0;
268  for ( trackingRecHit_iterator recHit = tk.recHitsBegin(); recHit != tk.recHitsEnd(); ++recHit ) {
269  if ( !((*recHit)->isValid()) ) continue;
270  ++nHit;
271  DetId id((*recHit)->geographicalId());
272  if ( (unsigned int)id.subdetId() == PixelSubdetector::PixelBarrel
273  || (unsigned int)id.subdetId() == PixelSubdetector::PixelEndcap )
274  ++nPixelHit;
275  }
276 
277  // Cut on the number of valid hits
278  if (nHit < min_nHit_) return false;
279  // Cut on the number of valid Pixel hits
280  if (nPixelHit < min_nPixelHit_) return false;
281  if (nlayers < min_layers_) return false;
282  if (nlayers3D < min_3Dlayers_) return false;
283  if (nlayersLost > max_lostLayers_) return false;
284 
285  // Get track parameters
286  double pt = tk.pt(),eta = tk.eta(), chi2n = tk.normalizedChi2();
287  double d0 = -tk.dxy(vertexBeamSpot.position()), dz = tk.dz();
288 
289  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
290  if (abs(d0) > max_d0_) return false;
291  if (abs(dz) > max_z0_) return false;
292 
293  // optimized cuts adapted to the track eta, pt and chiquare/ndof
294  if (abs(eta) > max_eta_) return false;
295  if (pt < min_pt_) return false;
296  if (chi2n > chi2n_par_*nlayers) return false;
297 
298 
299  else
300  return true;
301 
302 }
303 
304 
307 
309 
T getParameter(std::string const &) const
#define dso_hidden
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:293
std::auto_ptr< TrajTrackAssociationCollection > selTTAss_
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:556
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
void produce(edm::Event &evt, const edm::EventSetup &es) override
TrackQuality
track quality
Definition: TrackBase.h:151
TrackBase::TrackQuality qualityToSet_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
EDProductGetter const & productGetter() const
Definition: Event.cc:55
std::auto_ptr< reco::TrackCollection > selTracks_
reco::TrackExtraRefProd rTrackExtras_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::BeamSpot > beamspot_
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:50
edm::Ref< TrackExtraCollection > TrackExtraRef
persistent reference to a TrackExtra
Definition: TrackExtraFwd.h:17
key_type key() const
Accessor for product key.
Definition: Ref.h:264
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:508
edm::EDGetTokenT< std::vector< Trajectory > > srcTraj_
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:527
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
edm::RefProd< std::vector< Trajectory > > rTrajectories_
std::vector< double > res_par_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:646
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:350
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:75
double pt() const
track transverse momentum
Definition: TrackBase.h:616
TrackingRecHitRefProd rHits_
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:94
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< TrajTrackAssociationCollection > srcTass_
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
std::auto_ptr< reco::TrackExtraCollection > selTrackExtras_
RefProd< PROD > getRefBeforePut()
Definition: Event.h:141
std::auto_ptr< std::vector< const Trajectory * > > selTrajPtrs_
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:604
Definition: DetId.h:18
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
std::auto_ptr< std::vector< Trajectory > > selTrajs_
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
bool select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk)
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:445
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
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:80
std::auto_ptr< TrackingRecHitCollection > selHits_
static std::string const emptyString("")
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:544
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
const TrackResiduals & residuals() const
Definition: Track.h:220
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:204
edm::EDGetTokenT< reco::TrackCollection > src_
const Point & position() const
position
Definition: BeamSpot.h:62
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
volatile std::atomic< bool > shutdown_flag false
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:586
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
reco::TrackRefProd rTracks_
std::vector< reco::TrackRef > trackRefs_
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
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:179