CMS 3D CMS Logo

ISRWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
8 
10 
15 #include <Math/VectorUtil.h>
16 
17 //
18 // class declaration
19 //
21 public:
22  explicit ISRWeightProducer(const edm::ParameterSet&);
23  ~ISRWeightProducer() override;
24 
25 private:
26  void beginJob() override;
27  void produce(edm::Event&, const edm::EventSetup&) override;
28  void endJob() override;
29 
31  std::vector<double> isrBinEdges_;
32  std::vector<double> ptWeights_;
33 };
34 
37  genToken_ = consumes<reco::GenParticleCollection>(
38  pset.getUntrackedParameter<edm::InputTag>("GenTag", edm::InputTag("genParticles")));
39 
40  // Pt bin edges
41  std::vector<double> defPtEdges;
42  defPtEdges.push_back(0.);
43  defPtEdges.push_back(999999.);
44  isrBinEdges_ = pset.getUntrackedParameter<std::vector<double> >("ISRBinEdges", defPtEdges);
45  unsigned int ninputs_expected = isrBinEdges_.size() - 1;
46 
47  // Distortions in muon momentum
48  std::vector<double> defWeights;
49  defWeights.push_back(1.);
50  ptWeights_ = pset.getUntrackedParameter<std::vector<double> >("PtWeights", defWeights);
51  if (ptWeights_.size() == 1 && ninputs_expected > 1) {
52  for (unsigned int i = 1; i < ninputs_expected; i++) {
53  ptWeights_.push_back(ptWeights_[0]);
54  }
55  }
56 
57  produces<double>();
58 }
59 
62 
65 
68 
71  if (iEvent.isRealData())
72  return;
73 
75  iEvent.getByToken(genToken_, genParticles);
76  unsigned int gensize = genParticles->size();
77 
78  std::unique_ptr<double> weight(new double);
79 
80  // Set as default weight the asymptotic value at high pt (i.e. value of last bin)
81  (*weight) = ptWeights_[ptWeights_.size() - 1];
82 
83  unsigned int nbins = isrBinEdges_.size() - 1;
84  for (unsigned int i = 0; i < gensize; ++i) {
85  const reco::GenParticle& part = (*genParticles)[i];
86  int id = part.pdgId();
87  if (id != 23 && abs(id) != 24)
88  continue;
89  int status = part.status();
90  if (status != 3)
91  continue;
92  double pt = part.pt();
93  if (pt > isrBinEdges_[0] && pt < isrBinEdges_[nbins]) {
94  for (unsigned int j = 1; j <= nbins; ++j) {
95  if (pt > isrBinEdges_[j])
96  continue;
97  (*weight) = ptWeights_[j - 1];
98  break;
99  }
100  }
101  break;
102  }
103 
104  iEvent.put(std::move(weight));
105 }
106 
int pdgId() const final
PDG identifier.
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void beginJob() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double pt() const final
transverse momentum
bool isRealData() const
Definition: EventBase.h:62
std::vector< double > isrBinEdges_
~ISRWeightProducer() override
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ISRWeightProducer(const edm::ParameterSet &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void endJob() override
part
Definition: HCALResponse.h:20
std::vector< double > ptWeights_
int status() const final
status word
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511