CMS 3D CMS Logo

TrackingParticleRefMuonProducer.cc
Go to the documentation of this file.
10 
12 public:
14 
15  void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
16 
17 private:
20  double ptmin_;
21  double pmin_;
22 };
23 
25  : tpToken_(consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("src"))),
26  skim_(iConfig.getParameter<std::string>("skim")),
27  ptmin_(iConfig.getParameter<double>("ptmin")),
28  pmin_(iConfig.getParameter<double>("pmin")) {
29  edm::LogVerbatim("TrackingParticleRefMuonProducer")
30  << "\n constructing TrackingParticleRefMuonProducer: skim = " << skim_;
31  if (skim_ == "mu")
32  edm::LogVerbatim("TrackingParticleRefMuonProducer") << "\t ptmin = " << ptmin_ << ", pmin = " << pmin_ << "\n";
33  else if (skim_ == "track")
34  edm::LogVerbatim("TrackingParticleRefMuonProducer") << "\t ptmin = " << ptmin_ << "\n";
35  else if (skim_ == "pf")
36  edm::LogVerbatim("TrackingParticleRefMuonProducer") << "\t ptmin = " << ptmin_ << ", pmin = " << pmin_ << "\n";
37  else
38  edm::LogError("TrackingParticleRefMuonProducer") << "\t undefined skim = " << skim_ << "\n";
39 
40  produces<TrackingParticleRefVector>();
41 }
42 
45  iEvent.getByToken(tpToken_, tpH);
46 
47  auto tpskim = std::make_unique<TrackingParticleRefVector>();
48 
49  if (skim_ == "mu") {
50  for (size_t i = 0, end = tpH->size(); i < end; ++i) {
51  auto tp = TrackingParticleRef(tpH, i);
52 
53  // test if the TP is a muon with pt and p above minimum thresholds
54  bool isMu = (std::abs(tp->pdgId()) == 13);
55  bool ptpOk = (tp->pt() > ptmin_) && (tp->p() > pmin_);
56  if (isMu && ptpOk)
57  tpskim->push_back(tp);
58  else {
59  // test if the TP has muon hits
60  int n_muon_hits = tp->numberOfHits() - tp->numberOfTrackerHits();
61  if (n_muon_hits > 0)
62  tpskim->push_back(tp);
63  }
64  }
65  } else if (skim_ == "track") {
66  for (size_t i = 0, end = tpH->size(); i < end; ++i) {
67  auto tp = TrackingParticleRef(tpH, i);
68 
69  // test if the pt is above a minimum cut
70  if (tp->pt() > ptmin_)
71  tpskim->push_back(tp);
72  }
73  } else if (skim_ == "pf") {
74  for (size_t i = 0, end = tpH->size(); i < end; ++i) {
75  auto tp = TrackingParticleRef(tpH, i);
76 
77  // test if p and pt are above minimum cuts
78  if ((tp->pt() > ptmin_) && (tp->p() > pmin_))
79  tpskim->push_back(tp);
80  }
81  }
82 
83  iEvent.put(std::move(tpskim));
84 }
85 
Log< level::Info, true > LogVerbatim
TrackingParticleRefMuonProducer(const edm::ParameterSet &iConfig)
Log< level::Error, false > LogError
int iEvent
Definition: GenABIO.cc:224
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HLT enums.
std::vector< TrackingParticle > TrackingParticleCollection
edm::Ref< TrackingParticleCollection > TrackingParticleRef
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< TrackingParticleCollection > tpToken_