CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/ElectroWeakAnalysis/Utilities/src/PdfWeightProducer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/Framework/interface/EDProducer.h"
00005 
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/MakerMacros.h"
00008 
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 
00011 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00012 #include "DataFormats/Candidate/interface/CompositeCandidate.h"
00013 #include "CommonTools/CandUtils/interface/AddFourMomenta.h"
00014 #include "CommonTools/CandUtils/interface/Booster.h"
00015 #include <Math/VectorUtil.h>
00016 
00017 //
00018 // class declaration
00019 //
00020 class PdfWeightProducer : public edm::EDProducer {
00021    public:
00022       explicit PdfWeightProducer(const edm::ParameterSet&);
00023       ~PdfWeightProducer();
00024 
00025    private:
00026       virtual void beginJob() ;
00027       virtual void produce(edm::Event&, const edm::EventSetup&);
00028       virtual void endJob() ;
00029 
00030       std::string fixPOWHEG_;
00031       edm::InputTag genTag_;
00032       edm::InputTag pdfInfoTag_;
00033       std::vector<std::string> pdfSetNames_;
00034       std::vector<std::string> pdfShortNames_;
00035 };
00036 
00037 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00038 namespace LHAPDF {
00039       void initPDFSet(int nset, const std::string& filename, int member=0);
00040       int numberPDF(int nset);
00041       void usePDFMember(int nset, int member);
00042       double xfx(int nset, double x, double Q, int fl);
00043       double getXmin(int nset, int member);
00044       double getXmax(int nset, int member);
00045       double getQ2min(int nset, int member);
00046       double getQ2max(int nset, int member);
00047       void extrapolate(bool extrapolate=true);
00048 }
00049 
00051 PdfWeightProducer::PdfWeightProducer(const edm::ParameterSet& pset) :
00052  fixPOWHEG_(pset.getUntrackedParameter<std::string> ("FixPOWHEG", "")),
00053  genTag_(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles"))),
00054  pdfInfoTag_(pset.getUntrackedParameter<edm::InputTag> ("PdfInfoTag", edm::InputTag("generator"))),
00055  pdfSetNames_(pset.getUntrackedParameter<std::vector<std::string> > ("PdfSetNames"))
00056 {
00057       if (fixPOWHEG_ != "") pdfSetNames_.insert(pdfSetNames_.begin(),fixPOWHEG_);
00058 
00059       if (pdfSetNames_.size()>3) {
00060             edm::LogWarning("") << pdfSetNames_.size() << " PDF sets requested on input. Using only the first 3 sets and ignoring the rest!!";
00061             pdfSetNames_.erase(pdfSetNames_.begin()+3,pdfSetNames_.end());
00062       }
00063 
00064       for (unsigned int k=0; k<pdfSetNames_.size(); k++) {
00065             size_t dot = pdfSetNames_[k].find_first_of('.');
00066             size_t underscore = pdfSetNames_[k].find_first_of('_');
00067             if (underscore<dot) {
00068                   pdfShortNames_.push_back(pdfSetNames_[k].substr(0,underscore));
00069             } else {
00070                   pdfShortNames_.push_back(pdfSetNames_[k].substr(0,dot));
00071             }
00072             produces<std::vector<double> >(pdfShortNames_[k].data());
00073       }
00074 } 
00075 
00077 PdfWeightProducer::~PdfWeightProducer(){}
00078 
00080 void PdfWeightProducer::beginJob() {
00081       for (unsigned int k=1; k<=pdfSetNames_.size(); k++) {
00082             LHAPDF::initPDFSet(k,pdfSetNames_[k-1]);
00083       }
00084 }
00085 
00087 void PdfWeightProducer::endJob(){}
00088 
00090 void PdfWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
00091 
00092       if (iEvent.isRealData()) return;
00093 
00094       edm::Handle<GenEventInfoProduct> pdfstuff;
00095       if (!iEvent.getByLabel(pdfInfoTag_, pdfstuff)) {
00096             edm::LogError("PDFWeightProducer") << ">>> PdfInfo not found: " << pdfInfoTag_.encode() << " !!!";
00097             return;
00098       }
00099 
00100       float Q = pdfstuff->pdf()->scalePDF;
00101 
00102       int id1 = pdfstuff->pdf()->id.first;
00103       double x1 = pdfstuff->pdf()->x.first;
00104       double pdf1 = pdfstuff->pdf()->xPDF.first;
00105 
00106       int id2 = pdfstuff->pdf()->id.second;
00107       double x2 = pdfstuff->pdf()->x.second;
00108       double pdf2 = pdfstuff->pdf()->xPDF.second; 
00109 
00110       // Ad-hoc fix for POWHEG
00111       if (fixPOWHEG_!="") {
00112             edm::Handle<reco::GenParticleCollection> genParticles;
00113             if (!iEvent.getByLabel(genTag_, genParticles)) {
00114                   edm::LogError("PDFWeightProducer") << ">>> genParticles  not found: " << genTag_.encode() << " !!!";
00115                   return;
00116             }
00117             unsigned int gensize = genParticles->size();
00118             double mboson = 0.;
00119             for(unsigned int i = 0; i<gensize; ++i) {
00120                   const reco::GenParticle& part = (*genParticles)[i];
00121                   int status = part.status();
00122                   if (status!=3) continue;
00123                   int id = part.pdgId();
00124                   if (id!=23 && abs(id)!=24) continue;
00125                   mboson = part.mass();
00126                   break;
00127             }
00128             Q = sqrt(mboson*mboson+Q*Q);
00129             LHAPDF::usePDFMember(1,0);
00130             pdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
00131             pdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
00132       }
00133 
00134       // Put PDF weights in the event
00135       for (unsigned int k=1; k<=pdfSetNames_.size(); ++k) {
00136             std::auto_ptr<std::vector<double> > weights (new std::vector<double>);
00137             unsigned int nweights = 1;
00138             if (LHAPDF::numberPDF(k)>1) nweights += LHAPDF::numberPDF(k);
00139             weights->reserve(nweights);
00140         
00141             for (unsigned int i=0; i<nweights; ++i) {
00142                   LHAPDF::usePDFMember(k,i);
00143                   double newpdf1 = LHAPDF::xfx(k, x1, Q, id1)/x1;
00144                   double newpdf2 = LHAPDF::xfx(k, x2, Q, id2)/x2;
00145                   weights->push_back(newpdf1/pdf1*newpdf2/pdf2);
00146             }
00147             iEvent.put(weights,pdfShortNames_[k-1]);
00148       }
00149 }
00150 
00151 DEFINE_FWK_MODULE(PdfWeightProducer);