CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AnalyticalTrackSelector.cc
Go to the documentation of this file.
3 
5 #include <Math/DistFunc.h>
6 #include "TMath.h"
7 
9 
10 AnalyticalTrackSelector::AnalyticalTrackSelector( const edm::ParameterSet & cfg ) : MultiTrackSelector( )
11 {
12  //Spoof the pset for each track selector!
13  //Size is always 1!!!
14  qualityToSet_.reserve(1);
15  vtxNumber_.reserve(1);
16  vertexCut_.reserve(1);
17  res_par_.reserve(1);
18  chi2n_par_.reserve(1);
19  chi2n_no1Dmod_par_.reserve(1);
20  d0_par1_.reserve(1);
21  dz_par1_.reserve(1);
22  d0_par2_.reserve(1);
23  dz_par2_.reserve(1);
24  applyAdaptedPVCuts_.reserve(1);
25  max_d0_.reserve(1);
26  max_z0_.reserve(1);
27  nSigmaZ_.reserve(1);
28  min_layers_.reserve(1);
29  min_3Dlayers_.reserve(1);
30  max_lostLayers_.reserve(1);
31  min_hits_bypass_.reserve(1);
32  applyAbsCutsIfNoPV_.reserve(1);
33  max_d0NoPV_.reserve(1);
34  max_z0NoPV_.reserve(1);
35  preFilter_.reserve(1);
36  max_relpterr_.reserve(1);
37  min_nhits_.reserve(1);
38 
39  src_ = cfg.getParameter<edm::InputTag>( "src" );
40  beamspot_ = cfg.getParameter<edm::InputTag>( "beamspot" );
41  useVertices_ = cfg.getParameter<bool>( "useVertices" );
42  useVtxError_ = cfg.getParameter<bool>( "useVtxError" );
43  vertices_ = useVertices_ ? cfg.getParameter<edm::InputTag>( "vertices" ) : edm::InputTag("NONE");
44  copyExtras_ = cfg.getUntrackedParameter<bool>("copyExtras", false);
45  copyTrajectories_ = cfg.getUntrackedParameter<bool>("copyTrajectories", false);
46  minEta_ = cfg.getParameter<double>("min_eta");
47  maxEta_ = cfg.getParameter<double>("max_eta");
48 
50  // parameters for vertex selection
51  vtxNumber_.push_back( useVertices_ ? cfg.getParameter<int32_t>("vtxNumber") : 0 );
52  vertexCut_.push_back( useVertices_ ? cfg.getParameter<std::string>("vertexCut") : "");
53  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
54  res_par_.push_back(cfg.getParameter< std::vector<double> >("res_par") );
55  chi2n_par_.push_back( cfg.getParameter<double>("chi2n_par") );
56  chi2n_no1Dmod_par_.push_back( cfg.getParameter<double>("chi2n_no1Dmod_par") );
57  d0_par1_.push_back(cfg.getParameter< std::vector<double> >("d0_par1"));
58  dz_par1_.push_back(cfg.getParameter< std::vector<double> >("dz_par1"));
59  d0_par2_.push_back(cfg.getParameter< std::vector<double> >("d0_par2"));
60  dz_par2_.push_back(cfg.getParameter< std::vector<double> >("dz_par2"));
61 
62  // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
63  applyAdaptedPVCuts_.push_back(cfg.getParameter<bool>("applyAdaptedPVCuts"));
64  // Impact parameter absolute cuts.
65  max_d0_.push_back(cfg.getParameter<double>("max_d0"));
66  max_z0_.push_back(cfg.getParameter<double>("max_z0"));
67  nSigmaZ_.push_back(cfg.getParameter<double>("nSigmaZ"));
68  // Cuts on numbers of layers with hits/3D hits/lost hits.
69  min_layers_.push_back(cfg.getParameter<uint32_t>("minNumberLayers") );
70  min_3Dlayers_.push_back(cfg.getParameter<uint32_t>("minNumber3DLayers") );
71  max_lostLayers_.push_back(cfg.getParameter<uint32_t>("maxNumberLostLayers"));
72  min_hits_bypass_.push_back(cfg.getParameter<uint32_t>("minHitsToBypassChecks"));
73  max_relpterr_.push_back(cfg.getParameter<double>("max_relpterr"));
74  min_nhits_.push_back(cfg.getParameter<uint32_t>("min_nhits"));
75 
76  // Flag to apply absolute cuts if no PV passes the selection
77  applyAbsCutsIfNoPV_.push_back(cfg.getParameter<bool>("applyAbsCutsIfNoPV"));
78  keepAllTracks_.push_back( cfg.exists("keepAllTracks") ?
79  cfg.getParameter<bool>("keepAllTracks") :
80  false );
81 
82  setQualityBit_.push_back( false );
83  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
84 
85  if(d0_par1_[0].size()!=2 || dz_par1_[0].size()!=2 || d0_par2_[0].size()!=2 || dz_par2_[0].size()!=2)
86  {
87  edm::LogError("MisConfiguration")<<"vector of size less then 2";
88  throw;
89  }
90 
91  if (cfg.exists("qualityBit")) {
92  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
93  if (qualityStr != "") {
94  setQualityBit_[0] = true;
95  qualityToSet_ [0] = TrackBase::qualityByName(cfg.getParameter<std::string>("qualityBit"));
96  }
97  }
98 
99  if (keepAllTracks_[0] && !setQualityBit_[0]) throw cms::Exception("Configuration") <<
100  "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
101  if (setQualityBit_[0] && (qualityToSet_[0] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
102  "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
103  if (applyAbsCutsIfNoPV_[0]) {
104  max_d0NoPV_.push_back(cfg.getParameter<double>("max_d0NoPV"));
105  max_z0NoPV_.push_back(cfg.getParameter<double>("max_z0NoPV"));
106  }
107  else{//dummy values
108  max_d0NoPV_.push_back(0.);
109  max_z0NoPV_.push_back(0.);
110  }
111 
112  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
113  produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
114  if (copyExtras_) {
115  produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
116  produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
117  }
118  if (copyTrajectories_) {
119  produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
120  produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
121  }
122 }
123 
125 }
126 
128 {
129  using namespace std;
130  using namespace edm;
131  using namespace reco;
132 
133  Handle<TrackCollection> hSrcTrack;
137 
138  // looking for the beam spot
140  evt.getByLabel(beamspot_, hBsp);
141  reco::BeamSpot vertexBeamSpot;
142  vertexBeamSpot = *hBsp;
143 
144  // Select good primary vertices for use in subsequent track selection
146  std::vector<Point> points;
147  std::vector<double> vterr, vzerr;
148  if (useVertices_) {
149  evt.getByLabel(vertices_, hVtx);
150  selectVertices(0,*hVtx, points, vterr, vzerr);
151  // Debug
152  LogDebug("SelectVertex") << points.size() << " good pixel vertices";
153  }
154 
155  // Get tracks
156  evt.getByLabel( src_, hSrcTrack );
157 
158  selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
160  if (copyExtras_) {
161  selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
162  selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
165  }
166 
167  if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
168 
169  // Loop over tracks
170  size_t current = 0;
171  for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
172  const Track & trk = * it;
173  // Check if this track passes cuts
174 
175  LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
176 
177  bool ok = trk.eta()>minEta_ && trk.eta()<maxEta_ && select(0,vertexBeamSpot, trk, points, vterr, vzerr);
178  if (!ok) {
179 
180  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
181 
183  if (!keepAllTracks_[0]) continue;
184  }
185  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
186  selTracks_->push_back( Track( trk ) ); // clone and store
187  if (ok && setQualityBit_[0]) {
188  selTracks_->back().setQuality(qualityToSet_[0]);
189  if (!points.empty()) {
192  }
193  }
194  if (copyExtras_) {
195  // TrackExtras
196  selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
197  trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
198  trk.outerStateCovariance(), trk.outerDetId(),
199  trk.innerStateCovariance(), trk.innerDetId(),
200  trk.seedDirection(), trk.seedRef() ) );
201  selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
202  TrackExtra & tx = selTrackExtras_->back();
203  tx.setResiduals(trk.residuals());
204  // TrackingRecHits
205  for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
206  selHits_->push_back( (*hit)->clone() );
207  tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
208  }
209  }
210  if (copyTrajectories_) {
211  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
212  }
213  }
214  if ( copyTrajectories_ ) {
217  evt.getByLabel(src_, hTTAss);
218  evt.getByLabel(src_, hTraj);
219  selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>());
220  rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
221  selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
222  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
223  Ref< vector<Trajectory> > trajRef(hTraj, i);
225  if (match != hTTAss->end()) {
226  const Ref<TrackCollection> &trkRef = match->val;
227  short oldKey = static_cast<short>(trkRef.key());
228  if (trackRefs_[oldKey].isNonnull()) {
229  selTrajs_->push_back( Trajectory(*trajRef) );
230  selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
231  }
232  }
233  }
234  }
235 
236  static const std::string emptyString;
237  evt.put(selTracks_);
238  if (copyExtras_ ) {
239  evt.put(selTrackExtras_);
240  evt.put(selHits_);
241  }
242  if ( copyTrajectories_ ) {
243  evt.put(selTrajs_);
244  evt.put(selTTAss_);
245  }
246 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > max_d0_
Impact parameter absolute cuts.
int i
Definition: DBlmapReader.cc:9
void selectVertices(unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< double > &vterr, std::vector< double > &vzerr)
std::vector< std::vector< double > > dz_par1_
std::auto_ptr< reco::TrackExtraCollection > selTrackExtras_
std::vector< std::vector< double > > d0_par1_
std::vector< std::vector< double > > d0_par2_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::auto_ptr< std::vector< Trajectory > > selTrajs_
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:40
std::vector< uint32_t > min_3Dlayers_
edm::Ref< TrackExtraCollection > TrackExtraRef
persistent reference to a TrackExtra
Definition: TrackExtraFwd.h:13
std::vector< int32_t > vtxNumber_
vertex cuts
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
bool copyExtras_
copy only the tracks, not extras and rechits (for AOD)
edm::RefProd< std::vector< Trajectory > > rTrajectories_
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
std::auto_ptr< TrajTrackAssociationCollection > selTTAss_
std::vector< std::vector< double > > res_par_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
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:131
std::vector< uint32_t > min_nhits_
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
std::vector< reco::TrackRef > trackRefs_
std::vector< unsigned int > preFilter_
bool select(unsigned tsNum, const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk, const std::vector< Point > &points, std::vector< double > &vterr, std::vector< double > &vzerr)
return class, or -1 if rejected
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
#define LogTrace(id)
std::vector< uint32_t > min_hits_bypass_
RefProd< PROD > getRefBeforePut()
Definition: Event.h:97
edm::InputTag src_
source collection label
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::vector< double > chi2n_no1Dmod_par_
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:49
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
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
std::vector< TrackBase::TrackQuality > qualityToSet_
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
std::vector< bool > setQualityBit_
do I have to set a quality bit?
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
static std::string const emptyString("")
std::auto_ptr< reco::TrackCollection > selTracks_
storage
void produce(edm::Event &evt, const edm::EventSetup &es)
process one event
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
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
std::vector< uint32_t > max_lostLayers_
std::vector< std::vector< double > > dz_par2_
std::auto_ptr< TrackingRecHitCollection > selHits_
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
tuple size
Write out results.
bool copyTrajectories_
copy also trajectories and trajectory-&gt;track associations
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