CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
heppy::FSRWeightAlgo Class Reference

#include <FSRWeightAlgo.h>

Public Member Functions

void addGenParticle (const reco::GenParticle &gen)
 
void clear ()
 
 FSRWeightAlgo ()
 
double weight () const
 
virtual ~FSRWeightAlgo ()
 

Private Member Functions

double alphaRatio (double) const
 

Private Attributes

std::vector< reco::GenParticlegenParticles_
 

Detailed Description

Definition at line 8 of file FSRWeightAlgo.h.

Constructor & Destructor Documentation

◆ FSRWeightAlgo()

heppy::FSRWeightAlgo::FSRWeightAlgo ( )
inline

Definition at line 10 of file FSRWeightAlgo.h.

10 {}

◆ ~FSRWeightAlgo()

virtual heppy::FSRWeightAlgo::~FSRWeightAlgo ( )
inlinevirtual

Definition at line 11 of file FSRWeightAlgo.h.

11 {}

Member Function Documentation

◆ addGenParticle()

void heppy::FSRWeightAlgo::addGenParticle ( const reco::GenParticle gen)
inline

Definition at line 12 of file FSRWeightAlgo.h.

References genParticles_.

12 { genParticles_.push_back(gen); }
std::vector< reco::GenParticle > genParticles_
Definition: FSRWeightAlgo.h:19

◆ alphaRatio()

double heppy::FSRWeightAlgo::alphaRatio ( double  pt) const
private

Definition at line 78 of file FSRWeightAlgo.cc.

References A, B, correctionTermsCaloMet_cff::C, dqm-mbProfile::log, M_PI, and DiDispStaMuonMonitor_cfi::pt.

Referenced by weight().

78  {
79  double pigaga = 0.;
80 
81  // Leptonic contribution (just one loop, precise at < 0.3% level)
82  const double alphapi = 1 / 137.036 / M_PI;
83  const double mass_e = 0.0005;
84  const double mass_mu = 0.106;
85  const double mass_tau = 1.777;
86  const double mass_Z = 91.2;
87  if (pt > mass_e)
88  pigaga += alphapi * (2 * log(pt / mass_e) / 3. - 5. / 9.);
89  if (pt > mass_mu)
90  pigaga += alphapi * (2 * log(pt / mass_mu) / 3. - 5. / 9.);
91  if (pt > mass_tau)
92  pigaga += alphapi * (2 * log(pt / mass_tau) / 3. - 5. / 9.);
93 
94  // Hadronic vaccum contribution
95  // Using simple effective parametrization from Physics Letters B 513 (2001) 46.
96  // Top contribution neglected
97  double A = 0.;
98  double B = 0.;
99  double C = 0.;
100  if (pt < 0.7) {
101  A = 0.0;
102  B = 0.0023092;
103  C = 3.9925370;
104  } else if (pt < 2.0) {
105  A = 0.0;
106  B = 0.0022333;
107  C = 4.2191779;
108  } else if (pt < 4.0) {
109  A = 0.0;
110  B = 0.0024402;
111  C = 3.2496684;
112  } else if (pt < 10.0) {
113  A = 0.0;
114  B = 0.0027340;
115  C = 2.0995092;
116  } else if (pt < mass_Z) {
117  A = 0.0010485;
118  B = 0.0029431;
119  C = 1.0;
120  } else if (pt < 10000.) {
121  A = 0.0012234;
122  B = 0.0029237;
123  C = 1.0;
124  } else {
125  A = 0.0016894;
126  B = 0.0028984;
127  C = 1.0;
128  }
129  pigaga += A + B * log(1. + C * pt * pt);
130 
131  // Done
132  return 1. / (1. - pigaga);
133  }
Definition: APVGainStruct.h:7
#define M_PI
Definition: APVGainStruct.h:7

◆ clear()

void heppy::FSRWeightAlgo::clear ( void  )
inline

Definition at line 13 of file FSRWeightAlgo.h.

References genParticles_.

13 { genParticles_.clear(); }
std::vector< reco::GenParticle > genParticles_
Definition: FSRWeightAlgo.h:19

◆ weight()

double heppy::FSRWeightAlgo::weight ( ) const

Definition at line 10 of file FSRWeightAlgo.cc.

References funct::abs(), alphaRatio(), funct::cos(), reco::CompositeRefCandidateT< D >::daughter(), reco::CompositeRefCandidateT< D >::daughterRef(), dumpMFGeometry_cfg::delta, reco::LeafCandidate::energy(), genParticles_, mps_fire::i, dqmiolumiharvest::j, HLT_2023v12_cff::leptonId, reco::Candidate::mass(), reco::LeafCandidate::mass(), reco::CompositeRefCandidateT< D >::mother(), reco::CompositeRefCandidateT< D >::numberOfDaughters(), reco::CompositeRefCandidateT< D >::numberOfMothers(), reco::Candidate::pdgId(), reco::LeafCandidate::pdgId(), reco::LeafCandidate::phi(), displacedMuons_cfi::photon, MaterialEffects_cfi::photonEnergy, HLTObjectMonitor_cfi::photonPhi, conifer::pow(), funct::sin(), mathSSE::sqrt(), reco::LeafCandidate::status(), and reco::LeafCandidate::theta().

10  {
11  double weight = 1.;
12 
13  unsigned int gensize = genParticles_.size();
14  for (unsigned int i = 0; i < gensize; ++i) {
15  const reco::GenParticle& lepton = genParticles_[i];
16  if (lepton.status() != 3)
17  continue;
18  int leptonId = lepton.pdgId();
19  if (abs(leptonId) != 11 && abs(leptonId) != 13 && abs(leptonId) != 15)
20  continue;
21  if (lepton.numberOfMothers() != 1)
22  continue;
23  const reco::Candidate* boson = lepton.mother();
24  int bosonId = abs(boson->pdgId());
25  if (bosonId != 23 && bosonId != 24)
26  continue;
27  double bosonMass = boson->mass();
28  double leptonMass = lepton.mass();
29  double leptonEnergy = lepton.energy();
30  double cosLeptonTheta = cos(lepton.theta());
31  double sinLeptonTheta = sin(lepton.theta());
32  double leptonPhi = lepton.phi();
33 
34  int trueKey = i;
35  if (lepton.numberOfDaughters() == 0) {
36  continue;
37  } else if (lepton.numberOfDaughters() == 1) {
38  int otherleptonKey = lepton.daughterRef(0).key();
39  const reco::GenParticle& otherlepton = genParticles_[otherleptonKey];
40  if (otherlepton.pdgId() != leptonId)
41  continue;
42  if (otherlepton.numberOfDaughters() <= 1)
43  continue;
44  trueKey = otherleptonKey;
45  }
46 
47  const reco::GenParticle& trueLepton = genParticles_[trueKey];
48  unsigned int nDaughters = trueLepton.numberOfDaughters();
49 
50  for (unsigned int j = 0; j < nDaughters; ++j) {
51  const reco::Candidate* photon = trueLepton.daughter(j);
52  if (photon->pdgId() != 22)
53  continue;
54  double photonEnergy = photon->energy();
55  double cosPhotonTheta = cos(photon->theta());
56  double sinPhotonTheta = sin(photon->theta());
57  double photonPhi = photon->phi();
58  double costheta =
59  sinLeptonTheta * sinPhotonTheta * cos(leptonPhi - photonPhi) + cosLeptonTheta * cosPhotonTheta;
60  // Missing O(alpha) terms in soft-collinear approach
61  // Only for W, from hep-ph/0303260
62  if (bosonId == 24) {
63  double betaLepton = sqrt(1 - pow(leptonMass / leptonEnergy, 2));
64  double delta = -8 * photonEnergy * (1 - betaLepton * costheta) / pow(bosonMass, 3) /
65  (1 - pow(leptonMass / bosonMass, 2)) / (4 - pow(leptonMass / bosonMass, 2)) * leptonEnergy *
66  (pow(leptonMass, 2) / bosonMass + 2 * photonEnergy);
67  weight *= (1 + delta);
68  }
69  // Missing NLO QED orders in QED parton shower approach
70  // Change coupling scale from 0 to kT to estimate this effect
71  weight *= alphaRatio(photonEnergy * sqrt(1 - pow(costheta, 2)));
72  }
73  }
74 
75  return weight;
76  }
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
double weight() const
std::vector< reco::GenParticle > genParticles_
Definition: FSRWeightAlgo.h:19
virtual double mass() const =0
mass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
size_t numberOfDaughters() const override
number of daughters
int status() const final
status word
Definition: weight.py:1
constexpr int pow(int x)
Definition: conifer.h:24
size_t numberOfMothers() const override
number of mothers
double alphaRatio(double) const
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
int pdgId() const final
PDG identifier.
double theta() const final
momentum polar angle
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int pdgId() const =0
PDG identifier.
double mass() const final
mass
double phi() const final
momentum azimuthal angle
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
double energy() const final
energy

Member Data Documentation

◆ genParticles_

std::vector<reco::GenParticle> heppy::FSRWeightAlgo::genParticles_
private

Definition at line 19 of file FSRWeightAlgo.h.

Referenced by addGenParticle(), clear(), and weight().