CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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&);
24 
25  private:
26  virtual void beginJob() override ;
27  virtual void produce(edm::Event&, const edm::EventSetup&) override;
28  virtual void endJob() override ;
29 
31  std::vector<double> isrBinEdges_;
32  std::vector<double> ptWeights_;
33 };
34 
35 
39  genTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genPArticles"));
40 
41  // Pt bin edges
42  std::vector<double> defPtEdges;
43  defPtEdges.push_back(0.);
44  defPtEdges.push_back(999999.);
45  isrBinEdges_ = pset.getUntrackedParameter<std::vector<double> > ("ISRBinEdges",defPtEdges);
46  unsigned int ninputs_expected = isrBinEdges_.size()-1;
47 
48  // Distortions in muon momentum
49  std::vector<double> defWeights;
50  defWeights.push_back(1.);
51  ptWeights_ = pset.getUntrackedParameter<std::vector<double> > ("PtWeights",defWeights);
52  if (ptWeights_.size()==1 && ninputs_expected>1) {
53  for (unsigned int i=1; i<ninputs_expected; i++){ ptWeights_.push_back(ptWeights_[0]);}
54  }
55 
56  produces<double>();
57 }
58 
61 
64 
67 
70 
71  if (iEvent.isRealData()) return;
72 
74  iEvent.getByLabel(genTag_, genParticles);
75  unsigned int gensize = genParticles->size();
76 
77  std::auto_ptr<double> weight (new double);
78 
79  // Set as default weight the asymptotic value at high pt (i.e. value of last bin)
80  (*weight) = ptWeights_[ptWeights_.size()-1];
81 
82  unsigned int nbins = isrBinEdges_.size()-1;
83  for(unsigned int i = 0; i<gensize; ++i) {
84  const reco::GenParticle& part = (*genParticles)[i];
85  int id = part.pdgId();
86  if (id!=23 && abs(id)!=24) continue;
87  int status = part.status();
88  if (status!=3) continue;
89  double pt = part.pt();
90  if (pt>isrBinEdges_[0] && pt<isrBinEdges_[nbins]) {
91  for (unsigned int j=1; j<=nbins; ++j) {
92  if (pt>isrBinEdges_[j]) continue;
93  (*weight) = ptWeights_[j-1];
94  break;
95  }
96  }
97  break;
98  }
99 
100  iEvent.put(weight);
101 }
102 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void beginJob() override
virtual int pdgId() const GCC11_FINAL
PDG identifier.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool isRealData() const
Definition: EventBase.h:60
std::vector< double > isrBinEdges_
int iEvent
Definition: GenABIO.cc:243
virtual int status() const GCC11_FINAL
status word
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
ISRWeightProducer(const edm::ParameterSet &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual void endJob() override
part
Definition: HCALResponse.h:20
std::vector< double > ptWeights_
virtual void produce(edm::Event &, const edm::EventSetup &) override
edm::InputTag genTag_
int weight
Definition: histoStyle.py:50
tuple status
Definition: ntuplemaker.py:245
virtual float pt() const GCC11_FINAL
transverse momentum