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 
35 
38  genToken_ = consumes<reco::GenParticleCollection>(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++){ ptWeights_.push_back(ptWeights_[0]);}
53  }
54 
55  produces<double>();
56 }
57 
60 
63 
66 
69 
70  if (iEvent.isRealData()) return;
71 
73  iEvent.getByToken(genToken_, genParticles);
74  unsigned int gensize = genParticles->size();
75 
76  std::unique_ptr<double> weight (new double);
77 
78  // Set as default weight the asymptotic value at high pt (i.e. value of last bin)
79  (*weight) = ptWeights_[ptWeights_.size()-1];
80 
81  unsigned int nbins = isrBinEdges_.size()-1;
82  for(unsigned int i = 0; i<gensize; ++i) {
83  const reco::GenParticle& part = (*genParticles)[i];
84  int id = part.pdgId();
85  if (id!=23 && abs(id)!=24) continue;
86  int status = part.status();
87  if (status!=3) continue;
88  double pt = part.pt();
89  if (pt>isrBinEdges_[0] && pt<isrBinEdges_[nbins]) {
90  for (unsigned int j=1; j<=nbins; ++j) {
91  if (pt>isrBinEdges_[j]) continue;
92  (*weight) = ptWeights_[j-1];
93  break;
94  }
95  }
96  break;
97  }
98 
99  iEvent.put(std::move(weight));
100 }
101 
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:125
void beginJob() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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