CMS 3D CMS Logo

ChargeIsolationExtra.cc
Go to the documentation of this file.
9 
10 #include<iostream>
11 
12 //#define EDM_ML_DEBUG
13 
14 namespace spr{
15 
17  const edm::EventSetup& iSetup,
18  CaloNavigator<DetId>& theNavigator,
19  reco::TrackCollection::const_iterator trkItr,
21  const CaloSubdetectorGeometry* gEB,
22  const CaloSubdetectorGeometry* gEE,
24  TrackAssociatorParameters& parameters_, int ieta,
25  int iphi, const std::string & theTrackQuality, bool
26 #ifdef EDM_ML_DEBUG
27  debug
28 #endif
29  ) {
30 
31  double maxNearP = -1.0;
33 
34  // const DetId anyCell,
35  reco::TrackCollection::const_iterator trkItr2;
36  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
37 
38  const reco::Track* pTrack2 = &(*trkItr2);
39 
40  bool trkQuality = pTrack2->quality(trackQuality_);
41  if ( (trkItr2 != trkItr) && trkQuality ) {
42 
43  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
44  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
45  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
46 
47  if (info2.isGoodEcal ) {
48  if (std::abs(point2.eta())<spr::etaBEEcal) {
49  const DetId anyCell = gEB->getClosestCell(point2);
50 #ifdef EDM_ML_DEBUG
51  if (debug) std::cout << "chargeIsolation:: EB cell " << (EBDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
52 #endif
53  if (!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
54  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
55  }
56  } else {
57  const DetId anyCell = gEE->getClosestCell(point2);
58 #ifdef EDM_ML_DEBUG
59  if (debug) std::cout << "chargeIsolation:: EE cell " << (EEDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
60 #endif
61  if(!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
62  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
63  }
64  }
65  } //info2.isGoodEcal
66  }
67  }
68  return maxNearP;
69  }
70 
71  //===========================================================================================================
72 
73  bool chargeIsolation(const DetId anyCell,CaloNavigator<DetId>& navigator,int ieta,int iphi) {
74 
75  bool isIsolated = false;
76 
77  DetId thisDet;
78 
79  for (int dx = -ieta; dx < ieta+1; ++dx) {
80  for (int dy = -iphi; dy < iphi+1; ++dy) {
81 
82  thisDet = navigator.offsetBy(dx, dy);
83  navigator.home();
84 
85  if (thisDet != DetId(0)) {
86  if (thisDet == anyCell) {
87  isIsolated = false;
88  return isIsolated;
89  }
90  }
91  }
92  }
93  return isIsolated;
94  }
95 
96  //===========================================================================================================
97 
98  double chargeIsolationEcal(const edm::Event& iEvent, const edm::EventSetup& iSetup, const DetId& coreDet, reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, const CaloGeometry* geo, const CaloTopology* caloTopology, TrackDetectorAssociator& associator, TrackAssociatorParameters& parameters_, int ieta, int iphi, const std::string& theTrackQuality, bool debug) {
99 
102 
103  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
104 #ifdef EDM_ML_DEBUG
105  if (debug) std::cout << "chargeIsolation:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
106 #endif
107  double maxNearP = -1.0;
108  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
109 
110  // const DetId anyCell,
111  reco::TrackCollection::const_iterator trkItr2;
112  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
113 
114  const reco::Track* pTrack2 = &(*trkItr2);
115 
116  bool trkQuality = pTrack2->quality(trackQuality_);
117  if ( (trkItr2 != trkItr) && trkQuality ) {
118 
119  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
120  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
121  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
122 
123  if (info2.isGoodEcal ) {
124  if (std::abs(point2.eta())<spr::etaBEEcal) {
125  const DetId anyCell = barrelGeom->getClosestCell(point2);
126 #ifdef EDM_ML_DEBUG
127  if (debug) std::cout << "chargeIsolation:: EB cell " << (EBDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
128 #endif
129  if (!spr::chargeIsolation(anyCell,vdets)) {
130  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
131  }
132  } else {
133  const DetId anyCell = endcapGeom->getClosestCell(point2);
134 #ifdef EDM_ML_DEBUG
135  if (debug) std::cout << "chargeIsolation:: EE cell " << (EEDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
136 #endif
137  if (!spr::chargeIsolation(anyCell,vdets)) {
138  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
139  }
140  }
141  } //info2.isGoodEcal
142  }
143  }
144  return maxNearP;
145  }
146 
147  //===========================================================================================================
148 
149  double chargeIsolationHcal(const edm::Event& iEvent, const edm::EventSetup& iSetup, reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, const DetId ClosestCell, const HcalTopology* topology, const CaloSubdetectorGeometry* gHB, TrackDetectorAssociator& associator, TrackAssociatorParameters& parameters_, int ieta, int iphi, const std::string& theTrackQuality, bool debug) {
150 
151  std::vector<DetId> dets(1,ClosestCell);
152 
153 #ifdef EDM_ML_DEBUG
154  if (debug) std::cout << (HcalDetId) ClosestCell << std::endl;
155 #endif
156  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
157  std::vector<DetId>::iterator it;
158 
159 #ifdef EDM_ML_DEBUG
160  if (debug) {
161  for (unsigned int i=0; i<vdets.size(); i++) {
162  std::cout << "HcalDetId in " <<2*ieta+1 << "x" << 2*iphi+1 << " " << (HcalDetId) vdets[i] << std::endl;
163  }
164  }
165 #endif
166  double maxNearP = -1.0;
167  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
168 
169  reco::TrackCollection::const_iterator trkItr2;
170  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
171 
172  const reco::Track* pTrack2 = &(*trkItr2);
173 
174  bool trkQuality = pTrack2->quality(trackQuality_);
175  if ( (trkItr2 != trkItr) && trkQuality ) {
176  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
177  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
178  const GlobalPoint point2(info2.trkGlobPosAtHcal.x(),info2.trkGlobPosAtHcal.y(),info2.trkGlobPosAtHcal.z());
179 
180 #ifdef EDM_ML_DEBUG
181  if (debug) {
182  std::cout << "Track2 (p,eta,phi) " << pTrack2->p() << " " << pTrack2->eta() << " " << pTrack2->phi() << std::endl;
183  }
184 #endif
185  if (info2.isGoodHcal ) {
186  const DetId anyCell = gHB->getClosestCell(point2);
187 #ifdef EDM_ML_DEBUG
188  if (debug) std::cout << "chargeIsolation:: HCAL cell " << (HcalDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
189 #endif
190  if (!spr::chargeIsolation(anyCell,vdets)) {
191  if(maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
192  }
193 #ifdef EDM_ML_DEBUG
194  if (debug){
195  std::cout << "maxNearP " << maxNearP << " thisCell "
196  << (HcalDetId)anyCell << " ("
197  << info2.trkGlobPosAtHcal.x() << ","
198  << info2.trkGlobPosAtHcal.y() <<","
199  << info2.trkGlobPosAtHcal.z() <<")" << std::endl;
200  }
201 #endif
202  }
203  }
204  }
205  return maxNearP;
206  }
207 
208 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:615
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
static const double etaBEEcal
Definition: CaloConstants.h:12
CaloTopology const * topology(0)
void matrixECALIds(const DetId &det, int ieta, int iphi, const CaloGeometry *geo, const CaloTopology *caloTopology, std::vector< DetId > &vdets, bool debug=false, bool igNoreTransition=true)
bool chargeIsolation(const DetId anyCell, std::vector< DetId > &vdets)
TrackQuality
track quality
Definition: TrackBase.h:151
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
math::XYZPoint trkGlobPosAtHcal
T offsetBy(int deltaX, int deltaY) const
Free movement of arbitray steps.
Definition: CaloNavigator.h:80
double chargeIsolationEcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, bool debug=false)
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual DetId getClosestCell(const GlobalPoint &r) const
void home() const
move the navigator back to the starting point
Definition: DetId.h:18
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
#define debug
Definition: HDRShower.cc:19
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:510
std::vector< DetId > matrixHCALIds(std::vector< DetId > &dets, const HcalTopology *topology, int ieta, int iphi, bool includeHO=false, bool debug=false)
#define EDM_ML_DEBUG
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)