CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/ElectroWeakAnalysis/Utilities/src/WeakEffectsWeightProducer.cc

Go to the documentation of this file.
00001 
00002 #include "FWCore/Framework/interface/EDProducer.h"
00003 #include "FWCore/Utilities/interface/InputTag.h"
00004 
00005 class WeakEffectsWeightProducer: public edm::EDProducer {
00006 public:
00007       WeakEffectsWeightProducer(const edm::ParameterSet& pset);
00008       virtual ~WeakEffectsWeightProducer();
00009       virtual void beginJob() ;
00010       virtual void produce(edm::Event &, const edm::EventSetup&);
00011       virtual void endJob() ;
00012 private:
00013       edm::InputTag genParticlesTag_;
00014       double rhoParameter_;
00015 
00016       double alphaQED(double q2);
00017       double sigma0_qqbarll(unsigned int quark_type, double Q, double rho);
00018 
00019 };
00020 
00022 #include "FWCore/Framework/interface/MakerMacros.h"
00023 #include "FWCore/Framework/interface/Frameworkfwd.h"
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/Framework/interface/EventSetup.h"
00026 #include "FWCore/Framework/interface/ESHandle.h"
00027 #include "DataFormats/Common/interface/Handle.h"
00028 
00029 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00030 
00032 WeakEffectsWeightProducer::WeakEffectsWeightProducer(const edm::ParameterSet& pset) :
00033       genParticlesTag_(pset.getUntrackedParameter<edm::InputTag> ("GenParticlesTag", edm::InputTag("genParticles"))), 
00034       rhoParameter_(pset.getUntrackedParameter<double> ("RhoParameter", 1.004)) 
00035 {
00036       produces<double>();
00037 }
00038 
00040 WeakEffectsWeightProducer::~WeakEffectsWeightProducer(){}
00041 
00043 void WeakEffectsWeightProducer::beginJob(){
00044 }
00045 
00047 void WeakEffectsWeightProducer::endJob(){
00048 }
00049 
00051 void WeakEffectsWeightProducer::produce(edm::Event & iEvent, const edm::EventSetup&){
00052       if (iEvent.isRealData()) return;
00053 
00054       edm::Handle<reco::GenParticleCollection> genParticles;
00055       iEvent.getByLabel("genParticles", genParticles);
00056       unsigned int gensize = genParticles->size();
00057 
00058       std::auto_ptr<double> weight (new double);
00059 
00060       // Set default weight
00061       (*weight) = 1.;
00062 
00063       // Only DY implemented for the time being
00064       for(unsigned int i = 0; i<gensize; ++i) {
00065             const reco::GenParticle& part = (*genParticles)[i];
00066             int status = part.status();
00067             if (status!=3) break;
00068             int id = part.pdgId();
00069             if (id!=23) continue;
00070             double Q = part.mass();
00071             unsigned int nmothers = part.numberOfMothers();
00072             if (nmothers<=0) continue;
00073             size_t key = part.motherRef(0).key();
00074             unsigned int quark_id = abs((*genParticles)[key].pdgId());
00075             if (quark_id>0 && quark_id<6) {
00076                   (*weight) *= sigma0_qqbarll(quark_id, Q, rhoParameter_)
00077                              / sigma0_qqbarll(quark_id, Q, 1.0);
00078             }
00079             break;
00080       }
00081 
00082       //printf(" \t >>>>> WeakEffectsWeightProducer: Final weight = %f\n", (*weight)); 
00083       iEvent.put(weight);
00084 }
00085 
00086 double WeakEffectsWeightProducer::alphaQED(double q2) {
00087       double pigaga = -0.010449239475366825 - 0.0023228196282246765*log(q2)- 0.0288 - 0.002980*(log(q2/8464.)+0.006307*(q2/8464.-1.));
00088       return (1./137.0359895) / (1.+pigaga);
00089 }
00090 
00091 double WeakEffectsWeightProducer::sigma0_qqbarll(unsigned int quark_id, double Q, double rho) {
00092       double MZ = 91.188;
00093       double GZ = 2.495;
00094       double sin2eff = 0.232;
00095 
00096       double vl = -0.5 + 2.*sin2eff;
00097       double al = -0.5;
00098 
00099       double qq = 0.;
00100       double vq = 0.;
00101       double aq = 0.;
00102       double alphaW = 2.7e-3 * pow(log(Q*Q/80.4/80.4),2);
00103       double alphaZ = 2.7e-3 * pow(log(Q*Q/MZ/MZ),2);
00104       double sudakov_factor = 1.;
00105       if (abs(quark_id)%2==1) {
00106             qq = -1./3.;
00107             vq = -0.5 - 2.*qq*sin2eff;
00108             aq = -0.5;
00109             sudakov_factor = 1 + (-2.139 + 0.864)*alphaW - 0.385*alphaZ;
00110       } else {
00111             qq = 2./3.;
00112             vq = 0.5 - 2.*qq*sin2eff;
00113             aq = 0.5;
00114             sudakov_factor = 1 + (-3.423 + 1.807)*alphaW - 0.557*alphaZ;
00115       }
00116 
00117       double alfarn = alphaQED(Q*Q);
00118       double zcoupl = sqrt(2.) * 1.166389e-5 * MZ*MZ / 4. / M_PI;
00119       double gll = zcoupl * MZ/3. * (vl*vl + al*al);
00120       double gdd = zcoupl * MZ/3. * (vq*vq + aq*aq);
00121       double denom =  (Q*Q-MZ*MZ)*(Q*Q-MZ*MZ)+ pow(Q,4)*GZ*GZ/MZ/MZ;
00122       double qed = M_PI * qq*qq * alfarn*alfarn / Q/Q;
00123       double zint = rho * 2*M_PI * zcoupl * alfarn * qq * vq*vl * (Q*Q-MZ*MZ) / denom;
00124       double zonly = rho * rho * 9.*M_PI * gll * gdd / MZ/MZ * Q*Q / denom;
00125 
00126       return (qed + zint + zonly) * sudakov_factor;
00127 }
00128 
00129 DEFINE_FWK_MODULE(WeakEffectsWeightProducer);