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 ISRGammaWeightProducer : public edm::EDProducer { 00021 public: 00022 explicit ISRGammaWeightProducer(const edm::ParameterSet&); 00023 ~ISRGammaWeightProducer(); 00024 00025 virtual void beginJob() ; 00026 virtual void produce(edm::Event&, const edm::EventSetup&) override; 00027 virtual void endJob() ; 00028 00029 private: 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 ISRGammaWeightProducer::ISRGammaWeightProducer(const edm::ParameterSet& pset) { 00039 genTag_ = pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles")); 00040 00041 produces<double>(); 00042 } 00043 00045 ISRGammaWeightProducer::~ISRGammaWeightProducer(){} 00046 00048 void ISRGammaWeightProducer::beginJob() {} 00049 00051 void ISRGammaWeightProducer::endJob(){} 00052 00054 void ISRGammaWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) { 00055 00056 if (iEvent.isRealData()) return; 00057 00058 edm::Handle<reco::GenParticleCollection> genParticles; 00059 iEvent.getByLabel(genTag_, genParticles); 00060 unsigned int gensize = genParticles->size(); 00061 00062 std::auto_ptr<double> weight (new double); 00063 00064 // Set a default weight to start with 00065 (*weight) = 1.; 00066 00067 // Find the boson at the hard scattering level 00068 const reco::GenParticle* boson = 0; 00069 int parton1Key = -1; 00070 int parton2Key = -1; 00071 for (unsigned int i = 0; i<gensize; ++i) { 00072 const reco::GenParticle& part = (*genParticles)[i]; 00073 int status = abs(part.status()); 00074 if (status!=3) continue; 00075 if (part.numberOfMothers()!=2) continue; 00076 int partId = abs(part.pdgId()); 00077 if (status==3 && (partId==23||abs(partId)==24)) { 00078 boson = &(*genParticles)[i]; 00079 parton1Key = part.motherRef(0).key(); 00080 parton2Key = part.motherRef(1).key(); 00081 break; 00082 } 00083 } 00084 00085 // Consider only photons near the hard-scattering process 00086 const reco::GenParticle* photon = 0; 00087 if (boson) { 00088 for (unsigned int i = 0; i<gensize; ++i) { 00089 photon = 0; 00090 const reco::GenParticle& part = (*genParticles)[i]; 00091 int status = abs(part.status()); 00092 if (status!=1) continue; 00093 int partId = abs(part.pdgId()); 00094 if (partId!=22) continue; 00095 if (part.numberOfMothers()!=1) continue; 00096 int keyM = part.motherRef(0).key(); 00097 const reco::GenParticle* mother = &(*genParticles)[keyM]; 00098 if (mother->status()!=3) continue; 00099 int mId = mother->pdgId(); 00100 if (abs(mId)>6 && mId!=2212) continue; 00101 for (unsigned int j=0; j<mother->numberOfDaughters(); ++j){ 00102 int keyD = mother->daughterRef(j).key(); 00103 if (keyD==parton1Key || keyD==parton2Key) { 00104 photon = ∂ 00105 break; 00106 } 00107 } 00108 if (photon) break; 00109 } 00110 } 00111 00112 if (boson && photon) { 00113 math::XYZTLorentzVector smom = boson->p4() + photon->p4(); 00114 double s = smom.M2(); 00115 double sqrts = smom.M(); 00116 00117 // Go to CM using the boost direction of the boson+photon system 00118 ROOT::Math::Boost cmboost(smom.BoostToCM()); 00119 math::XYZTLorentzVector photonCM(cmboost(photon->p4())); 00120 double pcostheta = ( smom.x()*photonCM.x() 00121 + smom.y()*photonCM.y() 00122 + smom.z()*photonCM.z() ) / smom.P(); 00123 00124 // Determine kinematic invariants 00125 double t = - sqrts * (photonCM.t()-pcostheta); 00126 double MV = boson->mass(); 00127 double u = MV*MV - s - t; 00128 (*weight) = 1. - 2*t*u/(s*s+MV*MV*MV*MV); 00129 //printf(">>>>>>>>> s %f t %f u %f, MV %f, weight = %f\n", s, t, u, MV, (*weight)); 00130 } 00131 00132 iEvent.put(weight); 00133 } 00134 00135 DEFINE_FWK_MODULE(ISRGammaWeightProducer);