CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions
GenSpecificAlgo Class Reference

#include <RecoMET/METAlgorithms/interface/GenSpecificAlgo.h>

Public Member Functions

reco::GenMET addInfo (edm::Handle< edm::View< reco::Candidate > > particles, CommonMETData *met, double globalThreshold=0, bool onlyFiducial=false, bool applyFiducialThresholdForFractions=false, bool usePt=false)
 

Private Types

typedef math::XYZTLorentzVector LorentzVector
 
typedef math::XYZPoint Point
 

Private Member Functions

void fillCommonMETData (CommonMETData *met, edm::Handle< edm::View< reco::Candidate > > &particles, double globalThreshold, bool onlyFiducial, bool usePt)
 
SpecificGenMETData mkSpecificGenMETData (edm::Handle< edm::View< reco::Candidate > > &particles, double globalThreshold, bool onlyFiducial, bool applyFiducialThresholdForFractions, bool usePt)
 

Detailed Description

Description: Adds generator level HEPMC specific information to MET

Implementation: [Notes on implementation]

Definition at line 29 of file GenSpecificAlgo.h.

Member Typedef Documentation

◆ LorentzVector

Definition at line 39 of file GenSpecificAlgo.h.

◆ Point

Definition at line 40 of file GenSpecificAlgo.h.

Member Function Documentation

◆ addInfo()

reco::GenMET GenSpecificAlgo::addInfo ( edm::Handle< edm::View< reco::Candidate > >  particles,
CommonMETData met,
double  globalThreshold = 0,
bool  onlyFiducial = false,
bool  applyFiducialThresholdForFractions = false,
bool  usePt = false 
)

Definition at line 17 of file GenSpecificAlgo.cc.

References genMetCalo_cfi::applyFiducialThresholdForFractions, fillCommonMETData(), HepMCValidationHelper::genMet(), HLT_2024v13_cff::globalThreshold, BTaggingMonitor_cfi::met, mkSpecificGenMETData(), ecalTrigSettings_cff::particles, timingPdfMaker::specific, HLT_2024v13_cff::usePt, and L1BJetProducer_cff::vtx.

22  {
24 
27 
28  const LorentzVector p4(met->mex, met->mey, met->mez, met->met);
29  const Point vtx(0.0, 0.0, 0.0);
30 
31  reco::GenMET genMet(specific, met->sumet, p4, vtx);
32  return genMet;
33 }
void fillCommonMETData(CommonMETData *met, edm::Handle< edm::View< reco::Candidate > > &particles, double globalThreshold, bool onlyFiducial, bool usePt)
math::XYZTLorentzVector LorentzVector
applyFiducialThresholdForFractions
using Pt instead Et
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
SpecificGenMETData mkSpecificGenMETData(edm::Handle< edm::View< reco::Candidate > > &particles, double globalThreshold, bool onlyFiducial, bool applyFiducialThresholdForFractions, bool usePt)
Structure Point Contains parameters of Gaussian fits to DMRs.

◆ fillCommonMETData()

void GenSpecificAlgo::fillCommonMETData ( CommonMETData met,
edm::Handle< edm::View< reco::Candidate > > &  particles,
double  globalThreshold,
bool  onlyFiducial,
bool  usePt 
)
private

Definition at line 36 of file GenSpecificAlgo.cc.

References funct::cos(), MillePedeFileConverter_cfg::e, mps_fire::end, l1tnanotables_cff::et, HLT_2024v13_cff::globalThreshold, BTaggingMonitor_cfi::met, ecalTrigSettings_cff::particles, phi, funct::sin(), mathSSE::sqrt(), theta(), and HLT_2024v13_cff::usePt.

Referenced by addInfo().

40  {
41  double sum_et = 0.0;
42  double sum_ex = 0.0;
43  double sum_ey = 0.0;
44  double sum_ez = 0.0;
45 
46  for (edm::View<reco::Candidate>::const_iterator iParticle = (particles.product())->begin();
47  iParticle != (particles.product())->end();
48  ++iParticle) {
49  if ((onlyFiducial && TMath::Abs(iParticle->eta()) >= 5.0))
50  continue;
51 
52  if (iParticle->et() <= globalThreshold)
53  continue;
54 
55  if (usePt) {
56  double phi = iParticle->phi();
57  double et = iParticle->pt();
58  sum_ez += iParticle->pz();
59  sum_et += et;
60  sum_ex += et * cos(phi);
61  sum_ey += et * sin(phi);
62  } else {
63  double phi = iParticle->phi();
64  double theta = iParticle->theta();
65  double e = iParticle->energy();
66  double et = e * sin(theta);
67  sum_ez += e * cos(theta);
68  sum_et += et;
69  sum_ex += et * cos(phi);
70  sum_ey += et * sin(phi);
71  }
72  }
73 
74  met->mex = -sum_ex;
75  met->mey = -sum_ey;
76  met->mez = -sum_ez;
77  met->met = sqrt(sum_ex * sum_ex + sum_ey * sum_ey);
78  met->sumet = sum_et;
79 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:23
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
Geom::Theta< T > theta() const

◆ mkSpecificGenMETData()

SpecificGenMETData GenSpecificAlgo::mkSpecificGenMETData ( edm::Handle< edm::View< reco::Candidate > > &  particles,
double  globalThreshold,
bool  onlyFiducial,
bool  applyFiducialThresholdForFractions,
bool  usePt 
)
private

Definition at line 82 of file GenSpecificAlgo.cc.

References genMetCalo_cfi::applyFiducialThresholdForFractions, mps_fire::end, HLT_2024v13_cff::globalThreshold, ALPAKA_ACCELERATOR_NAMESPACE::caPixelDoublets::if(), ecalTrigSettings_cff::particles, EgammaValidation_cff::pdgId, DiDispStaMuonMonitor_cfi::pt, timingPdfMaker::specific, and HLT_2024v13_cff::usePt.

Referenced by addInfo().

86  {
87  const static int neutralEMpdgId[] = {22 /* photon */};
88  const static std::set<int> neutralEMpdgIdSet(neutralEMpdgId, neutralEMpdgId + sizeof(neutralEMpdgId) / sizeof(int));
89 
90  const static int chargedEMpdgId[] = {11 /* e */};
91  const static std::set<int> chargedEMpdgIdSet(chargedEMpdgId, chargedEMpdgId + sizeof(chargedEMpdgId) / sizeof(int));
92 
93  const static int muonpdgId[] = {13 /* muon */};
94  const static std::set<int> muonpdgIdSet(muonpdgId, muonpdgId + sizeof(muonpdgId) / sizeof(int));
95 
96  const static int neutralHADpdgId[] = {
97  130 /* K_long */,
98  310 /* K_short */,
99  3122 /* Lambda */,
100  2112 /* n */,
101  3222 /* Neutral Cascade */
102  };
103  const static std::set<int> neutralHADpdgIdSet(neutralHADpdgId,
104  neutralHADpdgId + sizeof(neutralHADpdgId) / sizeof(int));
105 
106  const static int chargedHADpdgId[] = {
107  211 /* pi */,
108  321 /* K+/K- */,
109  2212 /* p */,
110  3312 /* Cascade - */,
111  3112 /* Sigma - */,
112  3322 /* Sigma + */,
113  3334 /* Omega - */
114  };
115  const static std::set<int> chargedHADpdgIdSet(chargedHADpdgId,
116  chargedHADpdgId + sizeof(chargedHADpdgId) / sizeof(int));
117 
118  const static int invisiblepdgId[] = {
119  12 /* e_nu */, 14 /* mu_nu */, 16 /* tau_nu */,
120  1000022 /* Neutral ~Chi_0 */, 1000012 /* LH ~e_nu */, 1000014 /* LH ~mu_nu */,
121  1000016 /* LH ~tau_nu */, 2000012 /* RH ~e_nu */, 2000014 /* RH ~mu_nu */,
122  2000016 /* RH ~tau_nu */, 39 /* G */, 1000039 /* ~G */,
123  5100039 /* KK G */, 4000012 /* excited e_nu */, 4000014 /* excited mu_nu */,
124  4000016 /* excited tau_nu */, 9900012 /* Maj e_nu */, 9900014 /* Maj mu_nu */,
125  9900016 /* Maj tau_nu */,
126  };
127  const static std::set<int> invisiblepdgIdSet(invisiblepdgId, invisiblepdgId + sizeof(invisiblepdgId) / sizeof(int));
128 
130  double Et_unclassified = 0.0;
131 
132  for (edm::View<reco::Candidate>::const_iterator iParticle = (particles.product())->begin();
133  iParticle != (particles.product())->end();
134  ++iParticle) {
136  if (onlyFiducial && (TMath::Abs(iParticle->eta()) >= 5.0))
137  continue;
139  if (iParticle->et() <= globalThreshold)
140  continue;
141 
142  int pdgId = TMath::Abs(iParticle->pdgId());
143  double pt = (usePt) ? iParticle->pt() : iParticle->et();
144  if (neutralEMpdgIdSet.count(pdgId))
145  specific.NeutralEMEtFraction += pt;
146  else if (chargedEMpdgIdSet.count(pdgId))
147  specific.ChargedEMEtFraction += pt;
148  else if (muonpdgIdSet.count(pdgId))
149  specific.MuonEtFraction += pt;
150  else if (neutralHADpdgIdSet.count(pdgId))
151  specific.NeutralHadEtFraction += pt;
152  else if (chargedHADpdgIdSet.count(pdgId))
153  specific.ChargedHadEtFraction += pt;
154  else if (invisiblepdgIdSet.count(pdgId))
155  specific.InvisibleEtFraction += pt;
156  else
157  Et_unclassified += pt;
158  }
159 
160  double Et_Total = specific.NeutralEMEtFraction + specific.NeutralHadEtFraction + specific.ChargedEMEtFraction +
161  specific.ChargedHadEtFraction + specific.MuonEtFraction + specific.InvisibleEtFraction +
162  Et_unclassified;
163 
164  if (Et_Total) {
165  specific.NeutralEMEtFraction /= Et_Total;
166  specific.NeutralHadEtFraction /= Et_Total;
167  specific.ChargedEMEtFraction /= Et_Total;
168  specific.ChargedHadEtFraction /= Et_Total;
169  specific.MuonEtFraction /= Et_Total;
170  specific.InvisibleEtFraction /= Et_Total;
171  }
172 
173  return specific;
174 }
applyFiducialThresholdForFractions
using Pt instead Et
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:88
if(threadIdxLocalY==0 &&threadIdxLocalX==0)