CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingFailureFilter.cc
Go to the documentation of this file.
1 
10 
11 
13 
14  public:
15 
16  explicit TrackingFailureFilter(const edm::ParameterSet & iConfig);
18 
19  private:
20 
21  virtual bool filter(edm::Event & iEvent, const edm::EventSetup & iSetup) override;
22 
27 
28  const bool taggingMode_, debug_;
29 
30 };
31 
32 
34  : jetSrcToken_ (consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("JetSource")))
35  , trackSrcToken_ (consumes<std::vector<reco::Track> >(iConfig.getParameter<edm::InputTag>("TrackSource")))
36  , vertexSrcToken_ (consumes<std::vector<reco::Vertex> >(iConfig.getParameter<edm::InputTag>("VertexSource")))
37  , dzTrVtxMax_ (iConfig.getParameter<double>("DzTrVtxMax"))
38  , dxyTrVtxMax_ (iConfig.getParameter<double>("DxyTrVtxMax"))
39  , minSumPtOverHT_ (iConfig.getParameter<double>("MinSumPtOverHT"))
40  , taggingMode_ (iConfig.getParameter<bool>("taggingMode"))
41  , debug_ (iConfig.getParameter<bool>("debug"))
42 {
43 
44  produces<bool>();
45 }
46 
47 
49 
51  iEvent.getByToken(jetSrcToken_, jets);
53  iEvent.getByToken(trackSrcToken_, tracks);
55  iEvent.getByToken(vertexSrcToken_, vtxs);
56 
57  double ht = 0;
58  for (edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j) {
59  ht += j->pt();
60  }
61  double sumpt = 0;
62  if (vtxs->size() > 0) {
63 // const reco::Vertex * vtx = &((*vtxs)[0]);
64  for (std::vector<reco::Track>::const_iterator tr = tracks->begin(); tr != tracks->end(); ++tr) {
65  bool associateToPV = false;
66  for(int iv=0; iv<(int)vtxs->size(); iv++){
67  const reco::Vertex * pervtx = &((*vtxs)[iv]);
68  if( fabs(tr->dz(pervtx->position())) <= dzTrVtxMax_ && fabs(tr->dxy(pervtx->position())) <= dxyTrVtxMax_ ){
69  associateToPV = true;
70  }
71  }
72 // if (fabs(tr->dz(vtx->position())) > dzTrVtxMax_) continue;
73 // if (fabs(tr->dxy(vtx->position())) > dxyTrVtxMax_) continue;
74  if( !associateToPV ) continue;
75  sumpt += tr->pt();
76  }
77  }
78  const bool pass = (sumpt/ht) > minSumPtOverHT_;
79 
80  if( !pass && debug_ )
81  std::cout << "TRACKING FAILURE: "
82  << iEvent.id().run() << " : " << iEvent.id().luminosityBlock() << " : " << iEvent.id().event()
83  << " HT=" << ht
84  << " SumPt=" << sumpt
85  << std::endl;
86 
87  iEvent.put( std::auto_ptr<bool>(new bool(pass)) );
88 
89  return taggingMode_ || pass; // return false if filtering and not enough tracks in event
90 
91 }
92 
93 
95 
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
virtual bool filter(edm::Event &iEvent, const edm::EventSetup &iSetup) override
edm::EDGetTokenT< std::vector< reco::Track > > trackSrcToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexSrcToken_
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
const Point & position() const
position
Definition: Vertex.h:106
TrackingFailureFilter(const edm::ParameterSet &iConfig)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
int j
Definition: DBlmapReader.cc:9
tuple tracks
Definition: testEve_cfg.py:39
edm::EDGetTokenT< edm::View< reco::Jet > > jetSrcToken_
edm::EventID id() const
Definition: EventBase.h:60
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
tuple cout
Definition: gather_cfg.py:121