test
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.
2 
3 #include <Math/DistFunc.h>
4 #include "TMath.h"
5 
7 
8 CosmicTrackSelector::CosmicTrackSelector( const edm::ParameterSet & cfg ) :
9  src_( consumes<reco::TrackCollection>( cfg.getParameter<edm::InputTag>( "src" ) ) ),
10  beamspot_( consumes<reco::BeamSpot>( cfg.getParameter<edm::InputTag>( "beamspot" ) ) ),
11  copyExtras_(cfg.getUntrackedParameter<bool>("copyExtras", false)),
12  copyTrajectories_(cfg.getUntrackedParameter<bool>("copyTrajectories", false)),
13  keepAllTracks_( cfg.exists("keepAllTracks") ?
14  cfg.getParameter<bool>("keepAllTracks") :
15  false ), // as this is what you expect from a well behaved selector
16  setQualityBit_( false ),
17  qualityToSet_( TrackBase::undefQuality ),
18  chi2n_par_( cfg.getParameter<double>("chi2n_par") ),
19  // Impact parameter absolute cuts.
20  max_d0_(cfg.getParameter<double>("max_d0")),
21  max_z0_(cfg.getParameter<double>("max_z0")),
22  // Track parameter cuts.
23  min_pt_(cfg.getParameter<double>("min_pt")),
24  max_eta_(cfg.getParameter<double>("max_eta")),
25  // Cut on number of valid hits
26  min_nHit_(cfg.getParameter<uint32_t>("min_nHit")),
27  // Cut on number of valid hits
28  min_nPixelHit_(cfg.getParameter<uint32_t>("min_nPixelHit")),
29  // Cuts on numbers of layers with hits/3D hits/lost hits.
30  min_layers_(cfg.getParameter<uint32_t>("minNumberLayers") ),
31  min_3Dlayers_(cfg.getParameter<uint32_t>("minNumber3DLayers") ),
32  max_lostLayers_(cfg.getParameter<uint32_t>("maxNumberLostLayers") )
33 {
34  if (cfg.exists("qualityBit")) {
35  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
36  if (qualityStr != "") {
37  setQualityBit_ = true;
39  }
40  }
41  if (keepAllTracks_ && !setQualityBit_) throw cms::Exception("Configuration") <<
42  "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
43  if (setQualityBit_ && (qualityToSet_ == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
44  "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
45 
46  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
47  produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
48  if (copyExtras_) {
49  produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
50  produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
51  }
52  if (copyTrajectories_) {
53  srcTraj_ = consumes<std::vector<Trajectory> >(cfg.getParameter<edm::InputTag>( "src" ));
54  srcTass_ = consumes<TrajTrackAssociationCollection>(cfg.getParameter<edm::InputTag>( "src" ));
55  produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
56  produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
57  }
58 
59 }
60 
62 }
63 
65 {
66  using namespace std;
67  using namespace edm;
68  using namespace reco;
69 
70  Handle<TrackCollection> hSrcTrack;
74 
75  // looking for the beam spot
77  evt.getByToken(beamspot_, hBsp);
78  reco::BeamSpot vertexBeamSpot;
79  vertexBeamSpot = *hBsp;
80 
81  // Get tracks
82  evt.getByToken( src_, hSrcTrack );
83 
84  selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
86  if (copyExtras_) {
87  selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
88  selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
91  }
92 
93  if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
94 
95  // Loop over tracks
96  size_t current = 0;
97  for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
98  const Track & trk = * it;
99  // Check if this track passes cuts
100  bool ok = select(vertexBeamSpot, trk);
101  if (!ok) {
103  if (!keepAllTracks_) continue;
104  }
105  selTracks_->push_back( Track( trk ) ); // clone and store
106  if (ok && setQualityBit_) selTracks_->back().setQuality(qualityToSet_);
107  if (copyExtras_) {
108  // TrackExtras
109  selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
110  trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
111  trk.outerStateCovariance(), trk.outerDetId(),
112  trk.innerStateCovariance(), trk.innerDetId(),
113  trk.seedDirection(), trk.seedRef() ) );
114  selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
115  TrackExtra & tx = selTrackExtras_->back();
116  tx.setResiduals(trk.residuals());
117  // TrackingRecHits
118  for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
119  selHits_->push_back( (*hit)->clone() );
120  tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
121  }
122  }
123  if (copyTrajectories_) {
124  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
125  }
126  }
127  if ( copyTrajectories_ ) {
130  evt.getByToken(srcTass_, hTTAss);
131  evt.getByToken(srcTraj_, hTraj);
132  selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>());
133  rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
134  selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
135  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
136  Ref< vector<Trajectory> > trajRef(hTraj, i);
138  if (match != hTTAss->end()) {
139  const Ref<TrackCollection> &trkRef = match->val;
140  short oldKey = static_cast<short>(trkRef.key());
141  if (trackRefs_[oldKey].isNonnull()) {
142  selTrajs_->push_back( Trajectory(*trajRef) );
143  selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
144  }
145  }
146  }
147  }
148 
149  static const std::string emptyString;
150  evt.put(selTracks_);
151  if (copyExtras_ ) {
152  evt.put(selTrackExtras_);
153  evt.put(selHits_);
154  }
155  if ( copyTrajectories_ ) {
156  evt.put(selTrajs_);
157  evt.put(selTTAss_);
158  }
159 }
160 
161 
162 bool CosmicTrackSelector::select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk) {
163  // Decide if the given track passes selection cuts.
164 
165  using namespace std;
166 
167  // Cuts on numbers of layers with hits/3D hits/lost hits.
168  uint32_t nlayers = tk.hitPattern().trackerLayersWithMeasurement();
169  uint32_t nlayers3D = tk.hitPattern().pixelLayersWithMeasurement() +
171  uint32_t nlayersLost = tk.hitPattern().trackerLayersWithoutMeasurement();
172 
173  // Get the number of valid hits and PixelHits
174  uint32_t nHit = 0;
175  uint32_t nPixelHit = 0;
176  for ( trackingRecHit_iterator recHit = tk.recHitsBegin(); recHit != tk.recHitsEnd(); ++recHit ) {
177  if ( !((*recHit)->isValid()) ) continue;
178  ++nHit;
179  DetId id((*recHit)->geographicalId());
180  if ( (unsigned int)id.subdetId() == PixelSubdetector::PixelBarrel
181  || (unsigned int)id.subdetId() == PixelSubdetector::PixelEndcap )
182  ++nPixelHit;
183  }
184 
185  // Cut on the number of valid hits
186  if (nHit < min_nHit_) return false;
187  // Cut on the number of valid Pixel hits
188  if (nPixelHit < min_nPixelHit_) return false;
189  if (nlayers < min_layers_) return false;
190  if (nlayers3D < min_3Dlayers_) return false;
191  if (nlayersLost > max_lostLayers_) return false;
192 
193  // Get track parameters
194  double pt = tk.pt(),eta = tk.eta(), chi2n = tk.normalizedChi2();
195  double d0 = -tk.dxy(vertexBeamSpot.position()), dz = tk.dz();
196 
197  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
198  if (abs(d0) > max_d0_) return false;
199  if (abs(dz) > max_z0_) return false;
200 
201  // optimized cuts adapted to the track eta, pt and chiquare/ndof
202  if (abs(eta) > max_eta_) return false;
203  if (pt < min_pt_) return false;
204  if (chi2n > chi2n_par_*nlayers) return false;
205 
206 
207  else
208  return true;
209 
210 }
211 
212 
T getParameter(std::string const &) const
std::auto_ptr< reco::TrackExtraCollection > selTrackExtras_
int i
Definition: DBlmapReader.cc:9
edm::EDGetTokenT< std::vector< Trajectory > > srcTraj_
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:109
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::RefProd< std::vector< Trajectory > > rTrajectories_
string emptyString
Definition: archive.py:29
std::auto_ptr< std::vector< Trajectory > > selTrajs_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::TrackExtraRefProd rTrackExtras_
int pixelLayersWithMeasurement() const
Definition: HitPattern.h:745
std::vector< reco::TrackRef > trackRefs_
int trackerLayersWithoutMeasurement() const
Definition: HitPattern.h:758
edm::EDGetTokenT< reco::TrackCollection > src_
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:40
edm::Ref< TrackExtraCollection > TrackExtraRef
persistent reference to a TrackExtra
Definition: TrackExtraFwd.h:13
T eta() const
void produce(edm::Event &evt, const edm::EventSetup &es) override
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
edm::EDGetTokenT< reco::BeamSpot > beamspot_
TrackBase::TrackQuality qualityToSet_
std::auto_ptr< TrajTrackAssociationCollection > selTTAss_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
int trackerLayersWithMeasurement() const
Definition: HitPattern.h:740
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:51
double pt() const
track transverse momentum
Definition: TrackBase.h:129
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
std::auto_ptr< TrackingRecHitCollection > selHits_
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
RefProd< PROD > getRefBeforePut()
Definition: Event.h:128
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:125
Definition: DetId.h:18
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
edm::OwnVector< TrackingRecHit > TrackingRecHitCollection
collection of TrackingRecHits
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:112
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
std::auto_ptr< reco::TrackCollection > selTracks_
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:49
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:38
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:53
key_type key() const
Accessor for product key.
Definition: Ref.h:266
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
bool select(const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk)
edm::EDGetTokenT< TrajTrackAssociationCollection > srcTass_
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
const TrackResiduals & residuals() const
Definition: Track.h:117
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:105
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:6
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:119
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
int numberOfValidStripLayersWithMonoAndStereo(uint32_t stripdet, uint32_t layer) const
Definition: HitPattern.cc:221
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:131