CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
FSRWeightProducer Class Reference
Inheritance diagram for FSRWeightProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FSRWeightProducer (const edm::ParameterSet &)
 
 ~FSRWeightProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

double alphaRatio (double)
 
virtual void beginJob () override
 
virtual void endJob () override
 
virtual void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

edm::EDGetTokenT
< reco::GenParticleCollection
genToken_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 20 of file FSRWeightProducer.cc.

Constructor & Destructor Documentation

FSRWeightProducer::FSRWeightProducer ( const edm::ParameterSet pset)
explicit

Definition at line 36 of file FSRWeightProducer.cc.

References genToken_, edm::ParameterSet::getUntrackedParameter(), and HLT_25ns14e33_v1_cff::InputTag.

36  {
37  genToken_ = consumes<reco::GenParticleCollection>(pset.getUntrackedParameter<edm::InputTag> ("GenTag", edm::InputTag("genParticles")));
38 
39  produces<double>();
40 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
FSRWeightProducer::~FSRWeightProducer ( )

Definition at line 43 of file FSRWeightProducer.cc.

43 {}

Member Function Documentation

double FSRWeightProducer::alphaRatio ( double  pt)
private

Definition at line 125 of file FSRWeightProducer.cc.

References HLT_25ns14e33_v1_cff::A, funct::C, create_public_lumi_plots::log, and M_PI.

Referenced by produce().

125  {
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 }
#define M_PI
void FSRWeightProducer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 46 of file FSRWeightProducer.cc.

46 {}
void FSRWeightProducer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 49 of file FSRWeightProducer.cc.

49 {}
void FSRWeightProducer::produce ( edm::Event iEvent,
const edm::EventSetup  
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 52 of file FSRWeightProducer.cc.

References funct::abs(), alphaRatio(), funct::cos(), reco::CompositeRefCandidateT< D >::daughter(), reco::CompositeRefCandidateT< D >::daughterRef(), delta, reco::Candidate::energy(), reco::LeafCandidate::energy(), genParticleCandidates2GenParticles_cfi::genParticles, genToken_, edm::Event::getByToken(), i, edm::EventBase::isRealData(), j, HLT_25ns14e33_v1_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::Candidate::phi(), reco::LeafCandidate::phi(), funct::pow(), edm::Event::put(), funct::sin(), mathSSE::sqrt(), reco::LeafCandidate::status(), reco::Candidate::theta(), reco::LeafCandidate::theta(), and histoStyle::weight.

Referenced by JSONExport.JsonExport::export(), HTMLExport.HTMLExport::export(), and HTMLExport.HTMLExportStatic::export().

52  {
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 }
dbl * delta
Definition: mlp_gen.cc:36
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
virtual double energy() const =0
energy
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
virtual int status() const
status word
virtual double mass() const =0
mass
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double alphaRatio(double)
bool isRealData() const
Definition: EventBase.h:64
virtual double mass() const
mass
virtual double energy() const
energy
virtual double theta() const =0
momentum polar angle
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
virtual size_t numberOfMothers() const
number of mothers
T sqrt(T t)
Definition: SSEVec.h:48
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 theta() const
momentum polar angle
virtual int pdgId() const =0
PDG identifier.
int weight
Definition: histoStyle.py:50
virtual double phi() const
momentum azimuthal 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

Member Data Documentation

edm::EDGetTokenT<reco::GenParticleCollection> FSRWeightProducer::genToken_
private

Definition at line 31 of file FSRWeightProducer.cc.

Referenced by FSRWeightProducer(), and produce().