CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PdfWeightProducer Class Reference

Inheritance diagram for PdfWeightProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 PdfWeightProducer (const edm::ParameterSet &)
 ~PdfWeightProducer ()

Private Member Functions

virtual void beginJob ()
virtual void endJob ()
virtual void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

std::string fixPOWHEG_
edm::InputTag genTag_
edm::InputTag pdfInfoTag_
std::vector< std::string > pdfSetNames_
std::vector< std::string > pdfShortNames_

Detailed Description

Definition at line 20 of file PdfWeightProducer.cc.


Constructor & Destructor Documentation

PdfWeightProducer::PdfWeightProducer ( const edm::ParameterSet pset) [explicit]

Definition at line 51 of file PdfWeightProducer.cc.

References dot(), fixPOWHEG_, gen::k, pdfSetNames_, and pdfShortNames_.

                                                                :
 fixPOWHEG_(pset.getUntrackedParameter<std::string> ("FixPOWHEG", "")),
 genTag_(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles"))),
 pdfInfoTag_(pset.getUntrackedParameter<edm::InputTag> ("PdfInfoTag", edm::InputTag("generator"))),
 pdfSetNames_(pset.getUntrackedParameter<std::vector<std::string> > ("PdfSetNames"))
{
      if (fixPOWHEG_ != "") pdfSetNames_.insert(pdfSetNames_.begin(),fixPOWHEG_);

      if (pdfSetNames_.size()>3) {
            edm::LogWarning("") << pdfSetNames_.size() << " PDF sets requested on input. Using only the first 3 sets and ignoring the rest!!";
            pdfSetNames_.erase(pdfSetNames_.begin()+3,pdfSetNames_.end());
      }

      for (unsigned int k=0; k<pdfSetNames_.size(); k++) {
            size_t dot = pdfSetNames_[k].find_first_of('.');
            size_t underscore = pdfSetNames_[k].find_first_of('_');
            if (underscore<dot) {
                  pdfShortNames_.push_back(pdfSetNames_[k].substr(0,underscore));
            } else {
                  pdfShortNames_.push_back(pdfSetNames_[k].substr(0,dot));
            }
            produces<std::vector<double> >(pdfShortNames_[k].data());
      }
} 
PdfWeightProducer::~PdfWeightProducer ( )

Definition at line 77 of file PdfWeightProducer.cc.

{}

Member Function Documentation

void PdfWeightProducer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 80 of file PdfWeightProducer.cc.

References LHAPDF::initPDFSet(), gen::k, and pdfSetNames_.

                                 {
      for (unsigned int k=1; k<=pdfSetNames_.size(); k++) {
            LHAPDF::initPDFSet(k,pdfSetNames_[k-1]);
      }
}
void PdfWeightProducer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDProducer.

Definition at line 87 of file PdfWeightProducer.cc.

{}
void PdfWeightProducer::produce ( edm::Event iEvent,
const edm::EventSetup  
) [private, virtual]

Implements edm::EDProducer.

Definition at line 90 of file PdfWeightProducer.cc.

References abs, edm::InputTag::encode(), fixPOWHEG_, genParticleCandidates2GenParticles_cfi::genParticles, genTag_, edm::Event::getByLabel(), i, edm::HandleBase::id(), edm::EventBase::isRealData(), gen::k, reco::LeafCandidate::mass(), LHAPDF::numberPDF(), RooFit::pdf1(), pdfInfoTag_, pdfSetNames_, pdfShortNames_, reco::LeafCandidate::pdgId(), edm::Event::put(), mathSSE::sqrt(), ntuplemaker::status, reco::LeafCandidate::status(), LHAPDF::usePDFMember(), create_public_pileup_plots::weights, and LHAPDF::xfx().

                                                                      {

      if (iEvent.isRealData()) return;

      edm::Handle<GenEventInfoProduct> pdfstuff;
      if (!iEvent.getByLabel(pdfInfoTag_, pdfstuff)) {
            edm::LogError("PDFWeightProducer") << ">>> PdfInfo not found: " << pdfInfoTag_.encode() << " !!!";
            return;
      }

      float Q = pdfstuff->pdf()->scalePDF;

      int id1 = pdfstuff->pdf()->id.first;
      double x1 = pdfstuff->pdf()->x.first;
      double pdf1 = pdfstuff->pdf()->xPDF.first;

      int id2 = pdfstuff->pdf()->id.second;
      double x2 = pdfstuff->pdf()->x.second;
      double pdf2 = pdfstuff->pdf()->xPDF.second; 

      // Ad-hoc fix for POWHEG
      if (fixPOWHEG_!="") {
            edm::Handle<reco::GenParticleCollection> genParticles;
            if (!iEvent.getByLabel(genTag_, genParticles)) {
                  edm::LogError("PDFWeightProducer") << ">>> genParticles  not found: " << genTag_.encode() << " !!!";
                  return;
            }
            unsigned int gensize = genParticles->size();
            double mboson = 0.;
            for(unsigned int i = 0; i<gensize; ++i) {
                  const reco::GenParticle& part = (*genParticles)[i];
                  int status = part.status();
                  if (status!=3) continue;
                  int id = part.pdgId();
                  if (id!=23 && abs(id)!=24) continue;
                  mboson = part.mass();
                  break;
            }
            Q = sqrt(mboson*mboson+Q*Q);
            LHAPDF::usePDFMember(1,0);
            pdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
            pdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
      }

      // Put PDF weights in the event
      for (unsigned int k=1; k<=pdfSetNames_.size(); ++k) {
            std::auto_ptr<std::vector<double> > weights (new std::vector<double>);
            unsigned int nweights = 1;
            if (LHAPDF::numberPDF(k)>1) nweights += LHAPDF::numberPDF(k);
            weights->reserve(nweights);
        
            for (unsigned int i=0; i<nweights; ++i) {
                  LHAPDF::usePDFMember(k,i);
                  double newpdf1 = LHAPDF::xfx(k, x1, Q, id1)/x1;
                  double newpdf2 = LHAPDF::xfx(k, x2, Q, id2)/x2;
                  weights->push_back(newpdf1/pdf1*newpdf2/pdf2);
            }
            iEvent.put(weights,pdfShortNames_[k-1]);
      }
}

Member Data Documentation

std::string PdfWeightProducer::fixPOWHEG_ [private]

Definition at line 30 of file PdfWeightProducer.cc.

Referenced by PdfWeightProducer(), and produce().

Definition at line 31 of file PdfWeightProducer.cc.

Referenced by produce().

Definition at line 32 of file PdfWeightProducer.cc.

Referenced by produce().

std::vector<std::string> PdfWeightProducer::pdfSetNames_ [private]

Definition at line 33 of file PdfWeightProducer.cc.

Referenced by beginJob(), PdfWeightProducer(), and produce().

std::vector<std::string> PdfWeightProducer::pdfShortNames_ [private]

Definition at line 34 of file PdfWeightProducer.cc.

Referenced by PdfWeightProducer(), and produce().