CMS 3D CMS Logo

WeakEffectsWeightProducer.cc
Go to the documentation of this file.
3 
5 
7 public:
10  void beginJob() override;
11  void produce(edm::Event&, const edm::EventSetup&) override;
12  void endJob() override;
13 
14 private:
16  double rhoParameter_;
17 
18  double alphaQED(double q2);
19  double sigma0_qqbarll(unsigned int quark_type, double Q, double rho);
20 };
21 
29 
32  : // genParticlesToken_(consumes<reco::GenParticleCollection>(pset.getUntrackedParameter<edm::InputTag> ("GenParticlesTag", edm::InputTag("genParticles")))),
33  genParticlesToken_(consumes<reco::GenParticleCollection>(edm::InputTag("genParticles"))),
34  rhoParameter_(pset.getUntrackedParameter<double>("RhoParameter", 1.004)) {
35  produces<double>();
36 }
37 
40 
43 
46 
49  if (iEvent.isRealData())
50  return;
51 
54  unsigned int gensize = genParticles->size();
55 
56  std::unique_ptr<double> weight(new double);
57 
58  // Set default weight
59  (*weight) = 1.;
60 
61  // Only DY implemented for the time being
62  for (unsigned int i = 0; i < gensize; ++i) {
63  const reco::GenParticle& part = (*genParticles)[i];
64  int status = part.status();
65  if (status != 3)
66  break;
67  int id = part.pdgId();
68  if (id != 23)
69  continue;
70  double Q = part.mass();
71  unsigned int nmothers = part.numberOfMothers();
72  if (nmothers <= 0)
73  continue;
74  size_t key = part.motherRef(0).key();
75  unsigned int quark_id = abs((*genParticles)[key].pdgId());
76  if (quark_id > 0 && quark_id < 6) {
77  (*weight) *= sigma0_qqbarll(quark_id, Q, rhoParameter_) / sigma0_qqbarll(quark_id, Q, 1.0);
78  }
79  break;
80  }
81 
82  //printf(" \t >>>>> WeakEffectsWeightProducer: Final weight = %f\n", (*weight));
83  iEvent.put(std::move(weight));
84 }
85 
87  double pigaga = -0.010449239475366825 - 0.0023228196282246765 * log(q2) - 0.0288 -
88  0.002980 * (log(q2 / 8464.) + 0.006307 * (q2 / 8464. - 1.));
89  return (1. / 137.0359895) / (1. + pigaga);
90 }
91 
92 double WeakEffectsWeightProducer::sigma0_qqbarll(unsigned int quark_id, double Q, double rho) {
93  double MZ = 91.188;
94  double GZ = 2.495;
95  double sin2eff = 0.232;
96 
97  double vl = -0.5 + 2. * sin2eff;
98  double al = -0.5;
99 
100  double qq = 0.;
101  double vq = 0.;
102  double aq = 0.;
103  double alphaW = 2.7e-3 * pow(log(Q * Q / 80.4 / 80.4), 2);
104  double alphaZ = 2.7e-3 * pow(log(Q * Q / MZ / MZ), 2);
105  double sudakov_factor = 1.;
106  if (quark_id % 2 == 1) {
107  qq = -1. / 3.;
108  vq = -0.5 - 2. * qq * sin2eff;
109  aq = -0.5;
110  sudakov_factor = 1 + (-2.139 + 0.864) * alphaW - 0.385 * alphaZ;
111  } else {
112  qq = 2. / 3.;
113  vq = 0.5 - 2. * qq * sin2eff;
114  aq = 0.5;
115  sudakov_factor = 1 + (-3.423 + 1.807) * alphaW - 0.557 * alphaZ;
116  }
117 
118  double alfarn = alphaQED(Q * Q);
119  double zcoupl = sqrt(2.) * 1.166389e-5 * MZ * MZ / 4. / M_PI;
120  double gll = zcoupl * MZ / 3. * (vl * vl + al * al);
121  double gdd = zcoupl * MZ / 3. * (vq * vq + aq * aq);
122  double denom = (Q * Q - MZ * MZ) * (Q * Q - MZ * MZ) + pow(Q, 4) * GZ * GZ / MZ / MZ;
123  double qed = M_PI * qq * qq * alfarn * alfarn / Q / Q;
124  double zint = rho * 2 * M_PI * zcoupl * alfarn * qq * vq * vl * (Q * Q - MZ * MZ) / denom;
125  double zonly = rho * rho * 9. * M_PI * gll * gdd / MZ / MZ * Q * Q / denom;
126 
127  return (qed + zint + zonly) * sudakov_factor;
128 }
129 
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
WeakEffectsWeightProducer
Definition: WeakEffectsWeightProducer.cc:6
genParticles2HepMC_cfi.genParticles
genParticles
Definition: genParticles2HepMC_cfi.py:4
WeakEffectsWeightProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: WeakEffectsWeightProducer.cc:48
EDProducer.h
makePileupJSON.denom
denom
Definition: makePileupJSON.py:147
ESHandle.h
reco::GenParticle
Definition: GenParticle.h:21
mps_update.status
status
Definition: mps_update.py:69
edm::EDGetTokenT< reco::GenParticleCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
mps_merge.weight
weight
Definition: mps_merge.py:88
MZ
#define MZ
Definition: ZeeCalibration.cc:63
WeakEffectsWeightProducer::genParticlesToken_
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken_
Definition: WeakEffectsWeightProducer.cc:15
ZElectronSkim_cff.rho
rho
Definition: ZElectronSkim_cff.py:38
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:85964
reco::GenParticleCollection
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
Definition: GenParticleFwd.h:13
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
WeakEffectsWeightProducer::endJob
void endJob() override
Definition: WeakEffectsWeightProducer.cc:45
edm::Handle< reco::GenParticleCollection >
class-composition.Q
Q
Definition: class-composition.py:82
GenParticle.h
MakerMacros.h
part
part
Definition: HCALResponse.h:20
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
WeakEffectsWeightProducer::WeakEffectsWeightProducer
WeakEffectsWeightProducer(const edm::ParameterSet &pset)
Definition: WeakEffectsWeightProducer.cc:31
WeakEffectsWeightProducer::beginJob
void beginJob() override
Definition: WeakEffectsWeightProducer.cc:42
q2
double q2[4]
Definition: TauolaWrapper.h:88
DDAxes::rho
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
WeakEffectsWeightProducer::rhoParameter_
double rhoParameter_
Definition: WeakEffectsWeightProducer.cc:16
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
edm::EventSetup
Definition: EventSetup.h:57
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
WeakEffectsWeightProducer::sigma0_qqbarll
double sigma0_qqbarll(unsigned int quark_type, double Q, double rho)
Definition: WeakEffectsWeightProducer.cc:92
WeakEffectsWeightProducer::alphaQED
double alphaQED(double q2)
Definition: WeakEffectsWeightProducer.cc:86
EventSetup.h
edm::EDProducer
Definition: EDProducer.h:35
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
WeakEffectsWeightProducer::~WeakEffectsWeightProducer
~WeakEffectsWeightProducer() override
Definition: WeakEffectsWeightProducer.cc:39
edm::Event
Definition: Event.h:73
crabWrapper.key
key
Definition: crabWrapper.py:19
weight
Definition: weight.py:1
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27