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.
4 
6 #include <Math/DistFunc.h>
7 #include "TMath.h"
8 
10 
11 AnalyticalTrackSelector::AnalyticalTrackSelector( const edm::ParameterSet & cfg ) : MultiTrackSelector( )
12 {
13  //Spoof the pset for each track selector!
14  //Size is always 1!!!
15  qualityToSet_.reserve(1);
16  vtxNumber_.reserve(1);
17  vertexCut_.reserve(1);
18  res_par_.reserve(1);
19  chi2n_par_.reserve(1);
20  chi2n_no1Dmod_par_.reserve(1);
21  d0_par1_.reserve(1);
22  dz_par1_.reserve(1);
23  d0_par2_.reserve(1);
24  dz_par2_.reserve(1);
25  applyAdaptedPVCuts_.reserve(1);
26  max_d0_.reserve(1);
27  max_z0_.reserve(1);
28  nSigmaZ_.reserve(1);
29  min_layers_.reserve(1);
30  min_3Dlayers_.reserve(1);
31  max_lostLayers_.reserve(1);
32  min_hits_bypass_.reserve(1);
33  applyAbsCutsIfNoPV_.reserve(1);
34  max_d0NoPV_.reserve(1);
35  max_z0NoPV_.reserve(1);
36  preFilter_.reserve(1);
37  max_relpterr_.reserve(1);
38  min_nhits_.reserve(1);
39  max_minMissHitOutOrIn_.reserve(1);
40  max_lostHitFraction_.reserve(1);
41  min_eta_.reserve(1);
42  max_eta_.reserve(1);
43 
44  produces<edm::ValueMap<float> >("MVAVals");
45  useAnyMVA_ = false;
46  forest_ = 0;
47  gbrVals_ = 0;
48 
49  src_ = consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>( "src" ));
50  beamspot_ = consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>( "beamspot" ));
51  useVertices_ = cfg.getParameter<bool>( "useVertices" );
52  useVtxError_ = cfg.getParameter<bool>( "useVtxError" );
53  if (useVertices_) vertices_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>( "vertices" ));
54  copyExtras_ = cfg.getUntrackedParameter<bool>("copyExtras", false);
55  copyTrajectories_ = cfg.getUntrackedParameter<bool>("copyTrajectories", false);
56  if (copyTrajectories_) {
57  srcTraj_ = consumes<std::vector<Trajectory> >(cfg.getParameter<edm::InputTag>( "src" ));
58  srcTass_ = consumes<TrajTrackAssociationCollection>(cfg.getParameter<edm::InputTag>( "src" ));
59  }
60 
62  // parameters for vertex selection
63  vtxNumber_.push_back( useVertices_ ? cfg.getParameter<int32_t>("vtxNumber") : 0 );
64  vertexCut_.push_back( useVertices_ ? cfg.getParameter<std::string>("vertexCut") : "");
65  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
66  res_par_.push_back(cfg.getParameter< std::vector<double> >("res_par") );
67  chi2n_par_.push_back( cfg.getParameter<double>("chi2n_par") );
68  chi2n_no1Dmod_par_.push_back( cfg.getParameter<double>("chi2n_no1Dmod_par") );
69  d0_par1_.push_back(cfg.getParameter< std::vector<double> >("d0_par1"));
70  dz_par1_.push_back(cfg.getParameter< std::vector<double> >("dz_par1"));
71  d0_par2_.push_back(cfg.getParameter< std::vector<double> >("d0_par2"));
72  dz_par2_.push_back(cfg.getParameter< std::vector<double> >("dz_par2"));
73 
74  // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
75  applyAdaptedPVCuts_.push_back(cfg.getParameter<bool>("applyAdaptedPVCuts"));
76  // Impact parameter absolute cuts.
77  max_d0_.push_back(cfg.getParameter<double>("max_d0"));
78  max_z0_.push_back(cfg.getParameter<double>("max_z0"));
79  nSigmaZ_.push_back(cfg.getParameter<double>("nSigmaZ"));
80  // Cuts on numbers of layers with hits/3D hits/lost hits.
81  min_layers_.push_back(cfg.getParameter<uint32_t>("minNumberLayers") );
82  min_3Dlayers_.push_back(cfg.getParameter<uint32_t>("minNumber3DLayers") );
83  max_lostLayers_.push_back(cfg.getParameter<uint32_t>("maxNumberLostLayers"));
84  min_hits_bypass_.push_back(cfg.getParameter<uint32_t>("minHitsToBypassChecks"));
85  max_relpterr_.push_back(cfg.getParameter<double>("max_relpterr"));
86  min_nhits_.push_back(cfg.getParameter<uint32_t>("min_nhits"));
87  max_minMissHitOutOrIn_.push_back(
88  cfg.existsAs<int32_t>("max_minMissHitOutOrIn") ?
89  cfg.getParameter<int32_t>("max_minMissHitOutOrIn") : 99);
90  max_lostHitFraction_.push_back(
91  cfg.existsAs<double>("max_lostHitFraction") ?
92  cfg.getParameter<double>("max_lostHitFraction") : 1.0);
93  min_eta_.push_back(cfg.getParameter<double>("min_eta"));
94  max_eta_.push_back(cfg.getParameter<double>("max_eta"));
95 
96  // Flag to apply absolute cuts if no PV passes the selection
97  applyAbsCutsIfNoPV_.push_back(cfg.getParameter<bool>("applyAbsCutsIfNoPV"));
98  keepAllTracks_.push_back( cfg.exists("keepAllTracks") ?
99  cfg.getParameter<bool>("keepAllTracks") :
100  false );
101 
102  setQualityBit_.push_back( false );
103  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
104 
105  if(d0_par1_[0].size()!=2 || dz_par1_[0].size()!=2 || d0_par2_[0].size()!=2 || dz_par2_[0].size()!=2)
106  {
107  edm::LogError("MisConfiguration")<<"vector of size less then 2";
108  throw;
109  }
110 
111  if (cfg.exists("qualityBit")) {
112  std::string qualityStr = cfg.getParameter<std::string>("qualityBit");
113  if (qualityStr != "") {
114  setQualityBit_[0] = true;
116  }
117  }
118 
119  if (keepAllTracks_[0] && !setQualityBit_[0]) throw cms::Exception("Configuration") <<
120  "If you set 'keepAllTracks' to true, you must specify which qualityBit to set.\n";
121  if (setQualityBit_[0] && (qualityToSet_[0] == TrackBase::undefQuality)) throw cms::Exception("Configuration") <<
122  "You can't set the quality bit " << cfg.getParameter<std::string>("qualityBit") << " as it is 'undefQuality' or unknown.\n";
123  if (applyAbsCutsIfNoPV_[0]) {
124  max_d0NoPV_.push_back(cfg.getParameter<double>("max_d0NoPV"));
125  max_z0NoPV_.push_back(cfg.getParameter<double>("max_z0NoPV"));
126  }
127  else{//dummy values
128  max_d0NoPV_.push_back(0.);
129  max_z0NoPV_.push_back(0.);
130  }
131 
132  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
133  produces<reco::TrackCollection>().setBranchAlias( alias + "Tracks");
134  if (copyExtras_) {
135  produces<reco::TrackExtraCollection>().setBranchAlias( alias + "TrackExtras");
136  produces<TrackingRecHitCollection>().setBranchAlias( alias + "RecHits");
137  }
138  if (copyTrajectories_) {
139  produces< std::vector<Trajectory> >().setBranchAlias( alias + "Trajectories");
140  produces< TrajTrackAssociationCollection >().setBranchAlias( alias + "TrajectoryTrackAssociations");
141  }
142 }
143 
145 }
146 
148 {
149  using namespace std;
150  using namespace edm;
151  using namespace reco;
152 
153  Handle<TrackCollection> hSrcTrack;
157 
158  // looking for the beam spot
160  evt.getByToken(beamspot_, hBsp);
161  reco::BeamSpot vertexBeamSpot;
162  vertexBeamSpot = *hBsp;
163 
164  // Select good primary vertices for use in subsequent track selection
166  std::vector<Point> points;
167  std::vector<float> vterr, vzerr;
168  if (useVertices_) {
169  evt.getByToken(vertices_, hVtx);
170  selectVertices(0,*hVtx, points, vterr, vzerr);
171  // Debug
172  LogDebug("SelectVertex") << points.size() << " good pixel vertices";
173  }
174 
175  // Get tracks
176  evt.getByToken( src_, hSrcTrack );
177 
178  selTracks_ = auto_ptr<TrackCollection>(new TrackCollection());
180  if (copyExtras_) {
181  selTrackExtras_ = auto_ptr<TrackExtraCollection>(new TrackExtraCollection());
182  selHits_ = auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection());
185  }
186 
187  if (copyTrajectories_) trackRefs_.resize(hSrcTrack->size());
188 
189  processMVA(evt,es);
190 
191  // Loop over tracks
192  size_t current = 0;
193  for (TrackCollection::const_iterator it = hSrcTrack->begin(), ed = hSrcTrack->end(); it != ed; ++it, ++current) {
194  const Track & trk = * it;
195  // Check if this track passes cuts
196 
197  LogTrace("TrackSelection") << "ready to check track with pt="<< trk.pt() ;
198 
199  double mvaVal = 0;
200  if(useAnyMVA_)mvaVal = mvaVals_[current];
201  bool ok = select(0,vertexBeamSpot, trk, points, vterr, vzerr,mvaVal);
202  if (!ok) {
203 
204  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " NOT selected";
205 
207  if (!keepAllTracks_[0]) continue;
208  }
209  LogTrace("TrackSelection") << "track with pt="<< trk.pt() << " selected";
210  selTracks_->push_back( Track( trk ) ); // clone and store
211  if (ok && setQualityBit_[0]) {
212  selTracks_->back().setQuality(qualityToSet_[0]);
214  selTracks_->back().setQuality(TrackBase::loose);
215  }
216  else if (qualityToSet_[0]==TrackBase::highPurity) {
217  selTracks_->back().setQuality(TrackBase::loose);
218  selTracks_->back().setQuality(TrackBase::tight);
219  }
220  if (!points.empty()) {
222  selTracks_->back().setQuality(TrackBase::looseSetWithPV);
223  }
224  else if (qualityToSet_[0]==TrackBase::highPurity) {
225  selTracks_->back().setQuality(TrackBase::looseSetWithPV);
226  selTracks_->back().setQuality(TrackBase::highPuritySetWithPV);
227  }
228  }
229  }
230  if (copyExtras_) {
231  // TrackExtras
232  selTrackExtras_->push_back( TrackExtra( trk.outerPosition(), trk.outerMomentum(), trk.outerOk(),
233  trk.innerPosition(), trk.innerMomentum(), trk.innerOk(),
234  trk.outerStateCovariance(), trk.outerDetId(),
235  trk.innerStateCovariance(), trk.innerDetId(),
236  trk.seedDirection(), trk.seedRef() ) );
237  selTracks_->back().setExtra( TrackExtraRef( rTrackExtras_, selTrackExtras_->size() - 1) );
238  TrackExtra & tx = selTrackExtras_->back();
239  tx.setResiduals(trk.residuals());
240  // TrackingRecHits
241  for( trackingRecHit_iterator hit = trk.recHitsBegin(); hit != trk.recHitsEnd(); ++ hit ) {
242  selHits_->push_back( (*hit)->clone() );
243  tx.add( TrackingRecHitRef( rHits_, selHits_->size() - 1) );
244  }
245  }
246  if (copyTrajectories_) {
247  trackRefs_[current] = TrackRef(rTracks_, selTracks_->size() - 1);
248  }
249  }
250  if ( copyTrajectories_ ) {
253  evt.getByToken(srcTass_, hTTAss);
254  evt.getByToken(srcTraj_, hTraj);
255  selTrajs_ = auto_ptr< vector<Trajectory> >(new vector<Trajectory>());
256  rTrajectories_ = evt.getRefBeforePut< vector<Trajectory> >();
257  selTTAss_ = auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
258  for (size_t i = 0, n = hTraj->size(); i < n; ++i) {
259  Ref< vector<Trajectory> > trajRef(hTraj, i);
261  if (match != hTTAss->end()) {
262  const Ref<TrackCollection> &trkRef = match->val;
263  short oldKey = static_cast<short>(trkRef.key());
264  if (trackRefs_[oldKey].isNonnull()) {
265  selTrajs_->push_back( Trajectory(*trajRef) );
266  selTTAss_->insert ( Ref< vector<Trajectory> >(rTrajectories_, selTrajs_->size() - 1), trackRefs_[oldKey] );
267  }
268  }
269  }
270  }
271 
272  static const std::string emptyString;
273  evt.put(selTracks_);
274  if (copyExtras_ ) {
275  evt.put(selTrackExtras_);
276  evt.put(selHits_);
277  }
278  if ( copyTrajectories_ ) {
279  evt.put(selTrajs_);
280  evt.put(selTTAss_);
281  }
282 }
#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
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
std::vector< std::vector< double > > dz_par1_
std::auto_ptr< reco::TrackExtraCollection > selTrackExtras_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< std::vector< double > > d0_par1_
std::vector< std::vector< double > > d0_par2_
void produce(edm::Event &evt, const edm::EventSetup &es) override
process one event
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:39
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:46
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:41
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:116
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
void processMVA(edm::Event &evt, const edm::EventSetup &es)
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:50
double pt() const
track transverse momentum
Definition: TrackBase.h:129
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:58
std::vector< reco::TrackRef > trackRefs_
std::vector< unsigned int > preFilter_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:62
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
#define LogTrace(id)
std::vector< uint32_t > min_hits_bypass_
RefProd< PROD > getRefBeforePut()
Definition: Event.h:128
void selectVertices(unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr)
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:111
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:48
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:37
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
edm::EDGetTokenT< reco::BeamSpot > beamspot_
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:52
key_type key() const
Accessor for product key.
Definition: Ref.h:266
bool select(unsigned tsNum, const reco::BeamSpot &vertexBeamSpot, const reco::Track &tk, const std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr, double mvaVal)
return class, or -1 if rejected
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?
edm::EDGetTokenT< reco::VertexCollection > vertices_
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
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:44
edm::EDGetTokenT< std::vector< Trajectory > > srcTraj_
const TrackResiduals & residuals() const
Definition: Track.h:116
edm::EDGetTokenT< TrajTrackAssociationCollection > srcTass_
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:104
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:60
std::vector< int32_t > max_minMissHitOutOrIn_
tuple size
Write out results.
std::vector< int32_t > max_lostHitFraction_
bool copyTrajectories_
copy also trajectories and trajectory-&gt;track associations
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:64
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:130