CMS 3D CMS Logo

TrackingFailureFilter.cc
Go to the documentation of this file.
1 
11 
12 
14 
15  public:
16 
17  explicit TrackingFailureFilter(const edm::ParameterSet & iConfig);
18  ~TrackingFailureFilter() override {}
19 
20  private:
21 
22  bool filter(edm::StreamID, edm::Event & iEvent, const edm::EventSetup & iSetup) const override;
23 
28 
29  const bool taggingMode_, debug_;
30 
31 };
32 
33 
35  : jetSrcToken_ (consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("JetSource")))
36  , trackSrcToken_ (consumes<std::vector<reco::Track> >(iConfig.getParameter<edm::InputTag>("TrackSource")))
37  , vertexSrcToken_ (consumes<std::vector<reco::Vertex> >(iConfig.getParameter<edm::InputTag>("VertexSource")))
38  , dzTrVtxMax_ (iConfig.getParameter<double>("DzTrVtxMax"))
39  , dxyTrVtxMax_ (iConfig.getParameter<double>("DxyTrVtxMax"))
40  , minSumPtOverHT_ (iConfig.getParameter<double>("MinSumPtOverHT"))
41  , taggingMode_ (iConfig.getParameter<bool>("taggingMode"))
42  , debug_ (iConfig.getParameter<bool>("debug"))
43 {
44 
45  produces<bool>();
46 }
47 
48 
50 
52  iEvent.getByToken(jetSrcToken_, jets);
54  iEvent.getByToken(trackSrcToken_, tracks);
56  iEvent.getByToken(vertexSrcToken_, vtxs);
57 
58  double ht = 0;
59  for (edm::View<reco::Jet>::const_iterator j = jets->begin(); j != jets->end(); ++j) {
60  ht += j->pt();
61  }
62  double sumpt = 0;
63  if (!vtxs->empty()) {
64 // const reco::Vertex * vtx = &((*vtxs)[0]);
65  for (std::vector<reco::Track>::const_iterator tr = tracks->begin(); tr != tracks->end(); ++tr) {
66  bool associateToPV = false;
67  for(int iv=0; iv<(int)vtxs->size(); iv++){
68  const reco::Vertex * pervtx = &((*vtxs)[iv]);
69  if( fabs(tr->dz(pervtx->position())) <= dzTrVtxMax_ && fabs(tr->dxy(pervtx->position())) <= dxyTrVtxMax_ ){
70  associateToPV = true;
71  }
72  }
73 // if (fabs(tr->dz(vtx->position())) > dzTrVtxMax_) continue;
74 // if (fabs(tr->dxy(vtx->position())) > dxyTrVtxMax_) continue;
75  if( !associateToPV ) continue;
76  sumpt += tr->pt();
77  }
78  }
79  const bool pass = (sumpt/ht) > minSumPtOverHT_;
80 
81  if( !pass && debug_ )
82  edm::LogInfo("TrackingFailureFilter")
83  << "TRACKING FAILURE: "
84  << iEvent.id().run() << " : " << iEvent.id().luminosityBlock() << " : " << iEvent.id().event()
85  << " HT=" << ht
86  << " SumPt=" << sumpt;
87 
88  iEvent.put(std::make_unique<bool>(pass));
89 
90  return taggingMode_ || pass; // return false if filtering and not enough tracks in event
91 
92 }
93 
94 
96 
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool filter(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
edm::EDGetTokenT< std::vector< reco::Track > > trackSrcToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< std::vector< reco::Vertex > > vertexSrcToken_
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
const Point & position() const
position
Definition: Vertex.h:109
TrackingFailureFilter(const edm::ParameterSet &iConfig)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: Jet.py:1
vector< PseudoJet > jets
edm::EDGetTokenT< edm::View< reco::Jet > > jetSrcToken_
edm::EventID id() const
Definition: EventBase.h:59
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86