CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/ElectroWeakAnalysis/Utilities/src/FSRWeightProducer.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 FSRWeightProducer : public edm::EDProducer {
00021    public:
00022       explicit FSRWeightProducer(const edm::ParameterSet&);
00023       ~FSRWeightProducer();
00024 
00025    private:
00026       virtual void beginJob() ;
00027       virtual void produce(edm::Event&, const edm::EventSetup&);
00028       virtual void endJob() ;
00029       double alphaRatio(double) ;
00030 
00031       edm::InputTag genTag_;
00032 };
00033 
00034 
00035 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00037 FSRWeightProducer::FSRWeightProducer(const edm::ParameterSet& pset) {
00038       genTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles"));
00039 
00040       produces<double>();
00041 } 
00042 
00044 FSRWeightProducer::~FSRWeightProducer(){}
00045 
00047 void FSRWeightProducer::beginJob() {}
00048 
00050 void FSRWeightProducer::endJob(){}
00051 
00053 void FSRWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
00054 
00055       if (iEvent.isRealData()) return;
00056 
00057       edm::Handle<reco::GenParticleCollection> genParticles;
00058       iEvent.getByLabel(genTag_, genParticles);
00059 
00060       std::auto_ptr<double> weight (new double);
00061 
00062       // Set a default weight to start with
00063       (*weight) = 1.;
00064 
00065       unsigned int gensize = genParticles->size();
00066       for (unsigned int i = 0; i<gensize; ++i) {
00067             const reco::GenParticle& lepton = (*genParticles)[i];
00068             if (lepton.status()!=3) continue;
00069             int leptonId = lepton.pdgId();
00070             if (abs(leptonId)!=11 && abs(leptonId)!=13 && abs(leptonId)!=15) continue;
00071             if (lepton.numberOfMothers()!=1) continue;
00072             const reco::Candidate * boson = lepton.mother();
00073             int bosonId = abs(boson->pdgId());
00074             if (bosonId!=23  && bosonId!=24) continue;
00075             double bosonMass = boson->mass();
00076             double leptonMass = lepton.mass();
00077             double leptonEnergy = lepton.energy();
00078             double cosLeptonTheta = cos(lepton.theta());
00079             double sinLeptonTheta = sin(lepton.theta());
00080             double leptonPhi = lepton.phi();
00081 
00082             int trueKey = i;
00083             if (lepton.numberOfDaughters()==0) { 
00084                   continue;
00085             } else if (lepton.numberOfDaughters()==1) { 
00086                   int otherleptonKey = lepton.daughterRef(0).key();
00087                   const reco::GenParticle& otherlepton = (*genParticles)[otherleptonKey];
00088                   if (otherlepton.pdgId()!=leptonId) continue;
00089                   if (otherlepton.numberOfDaughters()<=1) continue;
00090                   trueKey = otherleptonKey;
00091             }
00092 
00093             const reco::GenParticle& trueLepton = (*genParticles)[trueKey];
00094             unsigned int nDaughters = trueLepton.numberOfDaughters();
00095 
00096             for (unsigned int j = 0; j<nDaughters; ++j) {
00097                   const reco::Candidate * photon = trueLepton.daughter(j);
00098                   if (photon->pdgId()!=22) continue;
00099                   double photonEnergy = photon->energy();
00100                   double cosPhotonTheta = cos(photon->theta());
00101                   double sinPhotonTheta = sin(photon->theta());
00102                   double photonPhi = photon->phi();
00103                   double costheta = sinLeptonTheta*sinPhotonTheta*cos(leptonPhi-photonPhi)
00104                                   + cosLeptonTheta*cosPhotonTheta;
00105                   // Missing O(alpha) terms in soft-collinear approach
00106                   // Only for W, from hep-ph/0303260
00107                   if (bosonId==24) {
00108                         double betaLepton = sqrt(1-pow(leptonMass/leptonEnergy,2));
00109                         double delta = - 8*photonEnergy *(1-betaLepton*costheta)
00110                           / pow(bosonMass,3) 
00111                           / (1-pow(leptonMass/bosonMass,2))
00112                           / (4-pow(leptonMass/bosonMass,2))
00113                           * leptonEnergy * (pow(leptonMass,2)/bosonMass+2*photonEnergy);
00114                         (*weight) *= (1 + delta);
00115                   }
00116                   // Missing NLO QED orders in QED parton shower approach
00117                   // Change coupling scale from 0 to kT to estimate this effect
00118                   (*weight) *= alphaRatio(photonEnergy*sqrt(1-pow(costheta,2)));
00119             }
00120       }
00121 
00122 
00123       iEvent.put(weight);
00124 }
00125 
00126 double FSRWeightProducer::alphaRatio(double pt) {
00127 
00128       double pigaga = 0.;
00129 
00130       // Leptonic contribution (just one loop, precise at < 0.3% level)
00131       const double alphapi = 1/137.036/M_PI;
00132       const double mass_e = 0.0005;
00133       const double mass_mu = 0.106;
00134       const double mass_tau = 1.777;
00135       const double mass_Z = 91.2;
00136       if (pt>mass_e) pigaga += alphapi * (2*log(pt/mass_e)/3.-5./9.);
00137       if (pt>mass_mu) pigaga += alphapi * (2*log(pt/mass_mu)/3.-5./9.);
00138       if (pt>mass_tau) pigaga += alphapi * (2*log(pt/mass_tau)/3.-5./9.);
00139 
00140       // Hadronic vaccum contribution
00141       // Using simple effective parametrization from Physics Letters B 513 (2001) 46.
00142       // Top contribution neglected
00143       double A = 0.; 
00144       double B = 0.; 
00145       double C = 0.; 
00146       if (pt<0.7) {
00147             A = 0.0; B = 0.0023092; C = 3.9925370;
00148       } else if (pt<2.0) {
00149             A = 0.0; B = 0.0022333; C = 4.2191779;
00150       } else if (pt<4.0) {
00151             A = 0.0; B = 0.0024402; C = 3.2496684;
00152       } else if (pt<10.0) {
00153             A = 0.0; B = 0.0027340; C = 2.0995092;
00154       } else if (pt<mass_Z) {
00155             A = 0.0010485; B = 0.0029431; C = 1.0;
00156       } else if (pt<10000.) {
00157             A = 0.0012234; B = 0.0029237; C = 1.0;
00158       } else {
00159             A = 0.0016894; B = 0.0028984; C = 1.0;
00160       }
00161       pigaga += A + B*log(1.+C*pt*pt);
00162 
00163       // Done
00164       return 1./(1.-pigaga);
00165 }
00166 
00167 DEFINE_FWK_MODULE(FSRWeightProducer);