CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FSRWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
8 
10 
15 #include <Math/VectorUtil.h>
16 
17 //
18 // class declaration
19 //
21  public:
22  explicit FSRWeightProducer(const edm::ParameterSet&);
24 
25  private:
26  virtual void beginJob() override ;
27  virtual void produce(edm::Event&, const edm::EventSetup&) override;
28  virtual void endJob() override ;
29  double alphaRatio(double) ;
30 
32 };
33 
34 
37  genToken_ = consumes<reco::GenParticleCollection>(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles")));
38 
39  produces<double>();
40 }
41 
44 
47 
50 
53 
54  if (iEvent.isRealData()) return;
55 
57  iEvent.getByToken(genToken_, genParticles);
58 
59  std::auto_ptr<double> weight (new double);
60 
61  // Set a default weight to start with
62  (*weight) = 1.;
63 
64  unsigned int gensize = genParticles->size();
65  for (unsigned int i = 0; i<gensize; ++i) {
66  const reco::GenParticle& lepton = (*genParticles)[i];
67  if (lepton.status()!=3) continue;
68  int leptonId = lepton.pdgId();
69  if (abs(leptonId)!=11 && abs(leptonId)!=13 && abs(leptonId)!=15) continue;
70  if (lepton.numberOfMothers()!=1) continue;
71  const reco::Candidate * boson = lepton.mother();
72  int bosonId = abs(boson->pdgId());
73  if (bosonId!=23 && bosonId!=24) continue;
74  double bosonMass = boson->mass();
75  double leptonMass = lepton.mass();
76  double leptonEnergy = lepton.energy();
77  double cosLeptonTheta = cos(lepton.theta());
78  double sinLeptonTheta = sin(lepton.theta());
79  double leptonPhi = lepton.phi();
80 
81  int trueKey = i;
82  if (lepton.numberOfDaughters()==0) {
83  continue;
84  } else if (lepton.numberOfDaughters()==1) {
85  int otherleptonKey = lepton.daughterRef(0).key();
86  const reco::GenParticle& otherlepton = (*genParticles)[otherleptonKey];
87  if (otherlepton.pdgId()!=leptonId) continue;
88  if (otherlepton.numberOfDaughters()<=1) continue;
89  trueKey = otherleptonKey;
90  }
91 
92  const reco::GenParticle& trueLepton = (*genParticles)[trueKey];
93  unsigned int nDaughters = trueLepton.numberOfDaughters();
94 
95  for (unsigned int j = 0; j<nDaughters; ++j) {
96  const reco::Candidate * photon = trueLepton.daughter(j);
97  if (photon->pdgId()!=22) continue;
98  double photonEnergy = photon->energy();
99  double cosPhotonTheta = cos(photon->theta());
100  double sinPhotonTheta = sin(photon->theta());
101  double photonPhi = photon->phi();
102  double costheta = sinLeptonTheta*sinPhotonTheta*cos(leptonPhi-photonPhi)
103  + cosLeptonTheta*cosPhotonTheta;
104  // Missing O(alpha) terms in soft-collinear approach
105  // Only for W, from hep-ph/0303260
106  if (bosonId==24) {
107  double betaLepton = sqrt(1-pow(leptonMass/leptonEnergy,2));
108  double delta = - 8*photonEnergy *(1-betaLepton*costheta)
109  / pow(bosonMass,3)
110  / (1-pow(leptonMass/bosonMass,2))
111  / (4-pow(leptonMass/bosonMass,2))
112  * leptonEnergy * (pow(leptonMass,2)/bosonMass+2*photonEnergy);
113  (*weight) *= (1 + delta);
114  }
115  // Missing NLO QED orders in QED parton shower approach
116  // Change coupling scale from 0 to kT to estimate this effect
117  (*weight) *= alphaRatio(photonEnergy*sqrt(1-pow(costheta,2)));
118  }
119  }
120 
121 
122  iEvent.put(weight);
123 }
124 
126 
127  double pigaga = 0.;
128 
129  // Leptonic contribution (just one loop, precise at < 0.3% level)
130  const double alphapi = 1/137.036/M_PI;
131  const double mass_e = 0.0005;
132  const double mass_mu = 0.106;
133  const double mass_tau = 1.777;
134  const double mass_Z = 91.2;
135  if (pt>mass_e) pigaga += alphapi * (2*log(pt/mass_e)/3.-5./9.);
136  if (pt>mass_mu) pigaga += alphapi * (2*log(pt/mass_mu)/3.-5./9.);
137  if (pt>mass_tau) pigaga += alphapi * (2*log(pt/mass_tau)/3.-5./9.);
138 
139  // Hadronic vaccum contribution
140  // Using simple effective parametrization from Physics Letters B 513 (2001) 46.
141  // Top contribution neglected
142  double A = 0.;
143  double B = 0.;
144  double C = 0.;
145  if (pt<0.7) {
146  A = 0.0; B = 0.0023092; C = 3.9925370;
147  } else if (pt<2.0) {
148  A = 0.0; B = 0.0022333; C = 4.2191779;
149  } else if (pt<4.0) {
150  A = 0.0; B = 0.0024402; C = 3.2496684;
151  } else if (pt<10.0) {
152  A = 0.0; B = 0.0027340; C = 2.0995092;
153  } else if (pt<mass_Z) {
154  A = 0.0010485; B = 0.0029431; C = 1.0;
155  } else if (pt<10000.) {
156  A = 0.0012234; B = 0.0029237; C = 1.0;
157  } else {
158  A = 0.0016894; B = 0.0028984; C = 1.0;
159  }
160  pigaga += A + B*log(1.+C*pt*pt);
161 
162  // Done
163  return 1./(1.-pigaga);
164 }
165 
dbl * delta
Definition: mlp_gen.cc:36
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual double energy() const =0
energy
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
virtual double energy() const final
energy
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual double mass() const =0
mass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
virtual double phi() const final
momentum azimuthal angle
double alphaRatio(double)
bool isRealData() const
Definition: EventBase.h:63
FSRWeightProducer(const edm::ParameterSet &)
virtual int status() const final
status word
int iEvent
Definition: GenABIO.cc:230
virtual double theta() const =0
momentum polar angle
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
virtual size_t numberOfMothers() const
number of mothers
T sqrt(T t)
Definition: SSEVec.h:18
virtual size_t numberOfDaughters() const
number of daughters
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual double mass() const final
mass
#define M_PI
virtual int pdgId() const =0
PDG identifier.
virtual void produce(edm::Event &, const edm::EventSetup &) override
virtual void beginJob() override
virtual int pdgId() const final
PDG identifier.
virtual double theta() const final
momentum polar angle
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual double phi() const =0
momentum azimuthal angle
virtual void endJob() override