CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/ElectroWeakAnalysis/Utilities/src/ISRWeightProducer.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 ISRWeightProducer : public edm::EDProducer {
00021    public:
00022       explicit ISRWeightProducer(const edm::ParameterSet&);
00023       ~ISRWeightProducer();
00024 
00025    private:
00026       virtual void beginJob() ;
00027       virtual void produce(edm::Event&, const edm::EventSetup&);
00028       virtual void endJob() ;
00029 
00030       edm::InputTag genTag_;
00031       std::vector<double> isrBinEdges_;
00032       std::vector<double> ptWeights_;
00033 };
00034 
00035 
00036 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00038 ISRWeightProducer::ISRWeightProducer(const edm::ParameterSet& pset) {
00039       genTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genPArticles"));
00040 
00041   // Pt bin edges
00042       std::vector<double> defPtEdges;
00043       defPtEdges.push_back(0.);
00044       defPtEdges.push_back(999999.);
00045       isrBinEdges_ = pset.getUntrackedParameter<std::vector<double> > ("ISRBinEdges",defPtEdges);
00046       unsigned int ninputs_expected = isrBinEdges_.size()-1;
00047 
00048   // Distortions in muon momentum
00049       std::vector<double> defWeights;
00050       defWeights.push_back(1.);
00051       ptWeights_ = pset.getUntrackedParameter<std::vector<double> > ("PtWeights",defWeights);
00052       if (ptWeights_.size()==1 && ninputs_expected>1) {
00053             for (unsigned int i=1; i<ninputs_expected; i++){ ptWeights_.push_back(ptWeights_[0]);}
00054       }
00055 
00056       produces<double>();
00057 } 
00058 
00060 ISRWeightProducer::~ISRWeightProducer(){}
00061 
00063 void ISRWeightProducer::beginJob() {}
00064 
00066 void ISRWeightProducer::endJob(){}
00067 
00069 void ISRWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
00070 
00071       if (iEvent.isRealData()) return;
00072 
00073       edm::Handle<reco::GenParticleCollection> genParticles;
00074       iEvent.getByLabel(genTag_, genParticles);
00075       unsigned int gensize = genParticles->size();
00076 
00077       std::auto_ptr<double> weight (new double);
00078 
00079       // Set as default weight the asymptotic value at high pt (i.e. value of last bin)
00080       (*weight) = ptWeights_[ptWeights_.size()-1];
00081 
00082       unsigned int nbins = isrBinEdges_.size()-1;
00083       for(unsigned int i = 0; i<gensize; ++i) {
00084             const reco::GenParticle& part = (*genParticles)[i];
00085             int id = part.pdgId();
00086             if (id!=23 && abs(id)!=24) continue;
00087             int status = part.status();
00088             if (status!=3) continue;
00089             double pt = part.pt();
00090             if (pt>isrBinEdges_[0] && pt<isrBinEdges_[nbins]) {
00091                   for (unsigned int j=1; j<=nbins; ++j) {
00092                         if (pt>isrBinEdges_[j]) continue;
00093                         (*weight) = ptWeights_[j-1];
00094                         break;
00095                   }
00096             }
00097             break;
00098       }
00099 
00100       iEvent.put(weight);
00101 }
00102 
00103 DEFINE_FWK_MODULE(ISRWeightProducer);