CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenSimInfo.cc
Go to the documentation of this file.
6 
7 #include<iostream>
8 
9 namespace spr{
10 
11  void eGenSimInfo(const DetId& coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector<spr::propagatedGenTrackID>& trackIds, const CaloGeometry* geo, const CaloTopology* caloTopology, int ieta, int iphi, genSimInfo & info, bool debug) {
12 
13  if (debug) std::cout << "eGenSimInfo:: For track " << (*trkItr)->momentum().rho() << "/" << (*trkItr)->momentum().eta() << "/" << (*trkItr)->momentum().phi() << " with ieta:iphi " << ieta << ":" << iphi << std::endl;
14 
15  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, false);
16  if (debug) spr::debugEcalDets(0, vdets);
17  spr::cGenSimInfo(vdets, trkItr, trackIds, true, info, debug);
18  }
19 
20  void eGenSimInfo(const DetId& coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector<spr::propagatedGenTrackID>& trackIds, const CaloGeometry* geo, const CaloTopology* caloTopology, double dR, const GlobalVector& trackMom, spr::genSimInfo & info, bool debug) {
21 
22  if (debug) std::cout << "eGenSimInfo:: For track " << (*trkItr)->momentum().rho() << "/" << (*trkItr)->momentum().eta() << "/" << (*trkItr)->momentum().phi() << " with dR,tMom " << dR << " " << trackMom << std::endl;
23 
24  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, dR, trackMom, geo, caloTopology, false);
25  if (debug) spr::debugEcalDets(0, vdets);
26  spr::cGenSimInfo(vdets, trkItr, trackIds, true, info, debug);
27  }
28 
29  void eGenSimInfo(const DetId& coreDet, reco::GenParticleCollection::const_iterator trkItr, std::vector<spr::propagatedGenParticleID>& trackIds, const CaloGeometry* geo, const CaloTopology* caloTopology, int ieta, int iphi, genSimInfo & info, bool debug) {
30 
31  if (debug) std::cout << "eGenSimInfo:: For track " << trkItr->momentum().R() << "/" << trkItr->momentum().eta() << "/" << trkItr->momentum().phi() << " with ieta:iphi " << ieta << ":" << iphi << std::endl;
32 
33  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, false);
34  if (debug) spr::debugEcalDets(0, vdets);
35  spr::cGenSimInfo(vdets, trkItr, trackIds, true, info, debug);
36  }
37 
38  void eGenSimInfo(const DetId& coreDet, reco::GenParticleCollection::const_iterator trkItr, std::vector<spr::propagatedGenParticleID>& trackIds, const CaloGeometry* geo, const CaloTopology* caloTopology, double dR, const GlobalVector& trackMom, spr::genSimInfo & info, bool debug) {
39 
40  if (debug) std::cout << "eGenSimInfo:: For track " << trkItr->momentum().R() << "/" << trkItr->momentum().eta() << "/" << trkItr->momentum().phi() << " with dR,tMom " << dR << " " << trackMom << std::endl;
41 
42  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, dR, trackMom, geo, caloTopology, false);
43  if (debug) spr::debugEcalDets(0, vdets);
44  spr::cGenSimInfo(vdets, trkItr, trackIds, true, info, debug);
45  }
46 
47  void hGenSimInfo(const DetId& coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector<spr::propagatedGenTrackID>& trackIds, const HcalTopology* topology, int ieta, int iphi, genSimInfo & info, bool includeHO, bool debug) {
48 
49  if (debug) std::cout << "hGenSimInfo:: For track " << (*trkItr)->momentum().rho() << "/" << (*trkItr)->momentum().eta() << "/" << (*trkItr)->momentum().phi() << " with ieta:iphi " << ieta << ":" << iphi << std::endl;
50 
51  std::vector<DetId> dets;
52  dets.push_back(coreDet);
53  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, includeHO, false);
54  if (debug) spr::debugHcalDets(0, vdets);
55  spr::cGenSimInfo(vdets, trkItr, trackIds, false, info, debug);
56  }
57 
58  void hGenSimInfo(const DetId& coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector<spr::propagatedGenTrackID>& trackIds, const CaloGeometry* geo, const HcalTopology* topology, double dR, const GlobalVector& trackMom, spr::genSimInfo & info, bool includeHO, bool debug) {
59 
60  if (debug) std::cout << "hGenSimInfo:: For track " << (*trkItr)->momentum().rho() << "/" << (*trkItr)->momentum().eta() << "/" << (*trkItr)->momentum().phi() << " with dR,tMom " << dR << " " << trackMom << std::endl;
61 
62  std::vector<DetId> vdets = spr::matrixHCALIds(coreDet, geo, topology, dR, trackMom, includeHO, false);
63  if (debug) spr::debugHcalDets(0, vdets);
64  spr::cGenSimInfo(vdets, trkItr, trackIds, false, info, debug);
65  }
66 
67  void hGenSimInfo(const DetId& coreDet, reco::GenParticleCollection::const_iterator trkItr, std::vector<spr::propagatedGenParticleID>& trackIds, const HcalTopology* topology, int ieta, int iphi, genSimInfo & info, bool includeHO, bool debug) {
68 
69  if (debug) std::cout << "hGenSimInfo:: For track " << trkItr->momentum().R() << "/" << trkItr->momentum().eta() << "/" << trkItr->momentum().phi() << " with ieta:iphi " << ieta << ":" << iphi << std::endl;
70 
71  std::vector<DetId> dets;
72  dets.push_back(coreDet);
73  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, includeHO, false);
74  if (debug) spr::debugHcalDets(0, vdets);
75  spr::cGenSimInfo(vdets, trkItr, trackIds, false, info, debug);
76  }
77 
78  void hGenSimInfo(const DetId& coreDet, reco::GenParticleCollection::const_iterator trkItr, std::vector<spr::propagatedGenParticleID>& trackIds, const CaloGeometry* geo, const HcalTopology* topology, double dR, const GlobalVector& trackMom, spr::genSimInfo & info, bool includeHO, bool debug) {
79 
80  if (debug) std::cout << "hGenSimInfo:: For track " << trkItr->momentum().R() << "/" << trkItr->momentum().eta() << "/" << trkItr->momentum().phi() << " with dR,tMom " << dR << " " << trackMom << std::endl;
81 
82  std::vector<DetId> vdets = spr::matrixHCALIds(coreDet, geo, topology, dR, trackMom, includeHO, false);
83  if (debug) spr::debugHcalDets(0, vdets);
84  spr::cGenSimInfo(vdets, trkItr, trackIds, false, info, debug);
85  }
86 
87  void cGenSimInfo(std::vector<DetId>& vdets, HepMC::GenEvent::particle_const_iterator trkItr, std::vector<spr::propagatedGenTrackID>& trackIds, bool ifECAL, spr::genSimInfo & info, bool debug) {
88 
89  info.maxNearP=-1.0;
90  info.cHadronEne=info.nHadronEne=info.eleEne=info.muEne=info.photonEne=0.0;
91  info.isChargedIso=true;
92  for (int i=0; i<3; ++i) info.cHadronEne_[i]=0.0;
93  for (unsigned int i=0; i<trackIds.size(); ++i) {
94  HepMC::GenEvent::particle_const_iterator trkItr2 = trackIds[i].trkItr;
95  // avoid the track under consideration
96  if ( (trkItr2 != trkItr) && trackIds[i].ok) {
97  int charge = trackIds[i].charge;
98  int pdgid = trackIds[i].pdgId;
99  double p = (*trkItr2)->momentum().rho();
100  bool isolat= false;
101  if (ifECAL) {
102  const DetId anyCell = trackIds[i].detIdECAL;
103  isolat = spr::chargeIsolation(anyCell,vdets);
104  } else {
105  const DetId anyCell = trackIds[i].detIdHCAL;
106  isolat = spr::chargeIsolation(anyCell,vdets);
107  }
108  if (!isolat) spr::cGenSimInfo(charge, pdgid, p, info, debug);
109  }
110  }
111  if (debug) {
112  std::cout << "Isolation variables: isChargedIso :" << info.isChargedIso
113  << " maxNearP " << info.maxNearP << " Energy e/mu/g/ch/nh "
114  << info.eleEne << "," << info.muEne << "," << info.photonEne
115  << "," << info.cHadronEne << "," << info.nHadronEne
116  << " charge " << info.cHadronEne_[0] << ","
117  << info.cHadronEne_[1] << "," << info.cHadronEne_[2]
118  << std::endl;
119  }
120  }
121 
122  void cGenSimInfo(std::vector<DetId>& vdets, reco::GenParticleCollection::const_iterator trkItr, std::vector<spr::propagatedGenParticleID>& trackIds, bool ifECAL, spr::genSimInfo & info, bool debug) {
123 
124  info.maxNearP=-1.0;
125  info.cHadronEne=info.nHadronEne=info.eleEne=info.muEne=info.photonEne=0.0;
126  info.isChargedIso=true;
127  for (int i=0; i<3; ++i) info.cHadronEne_[i]=0.0;
128  for (unsigned int i=0; i<trackIds.size(); ++i) {
129  reco::GenParticleCollection::const_iterator trkItr2 = trackIds[i].trkItr;
130  // avoid the track under consideration
131  if ( (trkItr2 != trkItr) && trackIds[i].ok) {
132  int charge = trackIds[i].charge;
133  int pdgid = trackIds[i].pdgId;
134  double p = trkItr2->momentum().R();
135  bool isolat= false;
136  if (ifECAL) {
137  const DetId anyCell = trackIds[i].detIdECAL;
138  isolat = spr::chargeIsolation(anyCell,vdets);
139  } else {
140  const DetId anyCell = trackIds[i].detIdHCAL;
141  isolat = spr::chargeIsolation(anyCell,vdets);
142  }
143  if (!isolat) spr::cGenSimInfo(charge, pdgid, p, info, debug);
144  }
145  }
146 
147  if (debug) {
148  std::cout << "Isolation variables: isChargedIso :" << info.isChargedIso
149  << " maxNearP " << info.maxNearP << " Energy e/mu/g/ch/nh "
150  << info.eleEne << "," << info.muEne << "," << info.photonEne
151  << "," << info.cHadronEne << "," << info.nHadronEne
152  << " charge " << info.cHadronEne_[0] << ","
153  << info.cHadronEne_[1] << "," << info.cHadronEne_[2]
154  << std::endl;
155  }
156  }
157 
158  void cGenSimInfo(int charge, int pdgid, double p, spr::genSimInfo & info, bool debug) {
159 
160  if (pdgid==22 ) info.photonEne += p;
161  else if (pdgid==11) info.eleEne += p;
162  else if (pdgid==13) info.muEne += p;
163  else if (std::abs(charge)>0) {
164  info.isChargedIso = false;
165  info.cHadronEne += p;
166  if (p>1.0) info.cHadronEne_[0] += p;
167  if (p>2.0) info.cHadronEne_[1] += p;
168  if (p>3.0) info.cHadronEne_[2] += p;
169  if (info.maxNearP<p) info.maxNearP=p;
170  } else if (std::abs(charge)==0) {
171  info.nHadronEne += p;
172  }
173  }
174 }
int i
Definition: DBlmapReader.cc:9
static const TGPicture * info(bool iBackgroundIsBlack)
void debugEcalDets(unsigned int, const DetId &, bool)
Definition: DebugInfo.cc:11
double cHadronEne_[3]
Definition: GenSimInfo.h:47
double photonEne
Definition: GenSimInfo.h:45
CaloTopology const * topology(0)
bool chargeIsolation(const DetId anyCell, std::vector< DetId > &vdets)
void cGenSimInfo(std::vector< DetId > &vdets, HepMC::GenEvent::particle_const_iterator trkItr, std::vector< spr::propagatedGenTrackID > &trackIds, bool ifECAL, spr::genSimInfo &info, bool debug=false)
Definition: GenSimInfo.cc:87
void debugHcalDets(unsigned int, std::vector< DetId > &)
Definition: DebugInfo.cc:42
double charge(const std::vector< uint8_t > &Ampls)
double cHadronEne
Definition: GenSimInfo.h:45
double nHadronEne
Definition: GenSimInfo.h:45
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void eGenSimInfo(const DetId &coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector< spr::propagatedGenTrackID > &trackIds, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, spr::genSimInfo &info, bool debug=false)
Definition: GenSimInfo.cc:11
Definition: DetId.h:18
#define debug
Definition: HDRShower.cc:19
void matrixECALIds(const DetId &det, int ieta, int iphi, const CaloGeometry *geo, const CaloTopology *caloTopology, std::vector< DetId > &vdets, bool debug=false)
for(const auto &pset:thresholds)
std::vector< DetId > matrixHCALIds(std::vector< DetId > &dets, const HcalTopology *topology, int ieta, int iphi, bool includeHO=false, bool debug=false)
void hGenSimInfo(const DetId &coreDet, HepMC::GenEvent::particle_const_iterator trkItr, std::vector< spr::propagatedGenTrackID > &trackIds, const HcalTopology *topology, int ieta, int iphi, spr::genSimInfo &info, bool includeHO=false, bool debug=false)
Definition: GenSimInfo.cc:47
tuple cout
Definition: gather_cfg.py:121
double maxNearP
Definition: GenSimInfo.h:44