test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PdfWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
8 
10 
16 #include <Math/VectorUtil.h>
17 
18 //
19 // class declaration
20 //
22  public:
23  explicit PdfWeightProducer(const edm::ParameterSet&);
25 
26  private:
27  virtual void beginJob() override ;
28  virtual void produce(edm::Event&, const edm::EventSetup&) override;
29  virtual void endJob() override ;
30 
37  std::vector<std::string> pdfSetNames_;
38  std::vector<std::string> pdfShortNames_;
39 };
40 
41 namespace LHAPDF {
42  void initPDFSet(int nset, const std::string& filename, int member=0);
43  int numberPDF(int nset);
44  void usePDFMember(int nset, int member);
45  double xfx(int nset, double x, double Q, int fl);
46  double getXmin(int nset, int member);
47  double getXmax(int nset, int member);
48  double getQ2min(int nset, int member);
49  double getQ2max(int nset, int member);
50  void extrapolate(bool extrapolate=true);
51 }
52 
55  fixPOWHEG_(pset.getUntrackedParameter<std::string> ("FixPOWHEG", "")),
56  useFirstAsDefault_(pset.getUntrackedParameter<bool>("useFirstAsDefault",false)),
57  genTag_(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles"))),
58  genToken_(mayConsume<reco::GenParticleCollection>(genTag_)),
59  pdfInfoTag_(pset.getUntrackedParameter<edm::InputTag> ("PdfInfoTag", edm::InputTag("generator"))),
60  pdfInfoToken_(consumes<GenEventInfoProduct>(pdfInfoTag_)),
61  pdfSetNames_(pset.getUntrackedParameter<std::vector<std::string> > ("PdfSetNames"))
62 {
63  if (fixPOWHEG_ != "") pdfSetNames_.insert(pdfSetNames_.begin(),fixPOWHEG_);
64 
65  if (pdfSetNames_.size()>3) {
66  edm::LogWarning("") << pdfSetNames_.size() << " PDF sets requested on input. Using only the first 3 sets and ignoring the rest!!";
67  pdfSetNames_.erase(pdfSetNames_.begin()+3,pdfSetNames_.end());
68  }
69 
70  for (unsigned int k=0; k<pdfSetNames_.size(); k++) {
71  size_t dot = pdfSetNames_[k].find_first_of('.');
72  size_t underscore = pdfSetNames_[k].find_first_of('_');
73  if (underscore<dot) {
74  pdfShortNames_.push_back(pdfSetNames_[k].substr(0,underscore));
75  } else {
76  pdfShortNames_.push_back(pdfSetNames_[k].substr(0,dot));
77  }
78  produces<std::vector<double> >(pdfShortNames_[k].data());
79  }
80 }
81 
84 
87  for (unsigned int k=1; k<=pdfSetNames_.size(); k++) {
89  }
90 }
91 
94 
97 
98  if (iEvent.isRealData()) return;
99 
101  if (!iEvent.getByToken(pdfInfoToken_, pdfstuff)) {
102  edm::LogError("PDFWeightProducer") << ">>> PdfInfo not found: " << pdfInfoTag_.encode() << " !!!";
103  return;
104  }
105 
106  float Q = pdfstuff->pdf()->scalePDF;
107 
108  int id1 = pdfstuff->pdf()->id.first;
109  double x1 = pdfstuff->pdf()->x.first;
110  double pdf1 = pdfstuff->pdf()->xPDF.first;
111 
112  int id2 = pdfstuff->pdf()->id.second;
113  double x2 = pdfstuff->pdf()->x.second;
114  double pdf2 = pdfstuff->pdf()->xPDF.second;
115  if (useFirstAsDefault_ && pdf1 == -1. && pdf2 == -1. ) {
117  pdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
118  pdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
119  }
120 
121  // Ad-hoc fix for POWHEG
122  if (fixPOWHEG_!="") {
124  if (!iEvent.getByToken(genToken_, genParticles)) {
125  edm::LogError("PDFWeightProducer") << ">>> genParticles not found: " << genTag_.encode() << " !!!";
126  return;
127  }
128  unsigned int gensize = genParticles->size();
129  double mboson = 0.;
130  for(unsigned int i = 0; i<gensize; ++i) {
131  const reco::GenParticle& part = (*genParticles)[i];
132  int status = part.status();
133  if (status!=3) continue;
134  int id = part.pdgId();
135  if (id!=23 && abs(id)!=24) continue;
136  mboson = part.mass();
137  break;
138  }
139  Q = sqrt(mboson*mboson+Q*Q);
141  pdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
142  pdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
143  }
144 
145  // Put PDF weights in the event
146  for (unsigned int k=1; k<=pdfSetNames_.size(); ++k) {
147  std::unique_ptr<std::vector<double> > weights (new std::vector<double>);
148  unsigned int nweights = 1;
149  if (LHAPDF::numberPDF(k)>1) nweights += LHAPDF::numberPDF(k);
150  weights->reserve(nweights);
151 
152  for (unsigned int i=0; i<nweights; ++i) {
154  double newpdf1 = LHAPDF::xfx(k, x1, Q, id1)/x1;
155  double newpdf2 = LHAPDF::xfx(k, x2, Q, id2)/x2;
156  weights->push_back(newpdf1/pdf1*newpdf2/pdf2);
157  }
158  iEvent.put(std::move(weights),pdfShortNames_[k-1]);
159  }
160 }
161 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int i
Definition: DBlmapReader.cc:9
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
edm::InputTag genTag_
bool isRealData() const
Definition: EventBase.h:62
std::string encode() const
Definition: InputTag.cc:164
virtual void beginJob() override
void initPDFSet(int nset, const std::string &filename, int member=0)
std::vector< std::string > pdfSetNames_
virtual int status() const final
status word
T x() const
Cartesian x coordinate.
double getXmax(int nset, int member)
double getQ2max(int nset, int member)
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:18
void extrapolate(bool extrapolate=true)
def move
Definition: eostools.py:510
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int numberPDF(int nset)
double getXmin(int nset, int member)
virtual double mass() const final
mass
PdfWeightProducer(const edm::ParameterSet &)
edm::EDGetTokenT< GenEventInfoProduct > pdfInfoToken_
edm::InputTag pdfInfoTag_
part
Definition: HCALResponse.h:20
virtual void produce(edm::Event &, const edm::EventSetup &) override
std::vector< std::string > pdfShortNames_
virtual int pdgId() const final
PDG identifier.
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
double xfx(int nset, double x, double Q, int fl)
tuple filename
Definition: lut2db_cfg.py:20
volatile std::atomic< bool > shutdown_flag false
double getQ2min(int nset, int member)
tuple status
Definition: mps_update.py:57
void usePDFMember(int nset, int member)
virtual void endJob() override