CMS 3D CMS Logo

ChargeIsolation.cc
Go to the documentation of this file.
8 
9 #include<iostream>
10 
11 //#define EDM_ML_DEBUG
12 
13 namespace spr{
14 
15 
16  double chargeIsolationEcal(unsigned int trkIndex, std::vector<spr::propagatedTrackID>& vdetIds, const CaloGeometry* geo, const CaloTopology* caloTopology, int ieta, int iphi, bool debug) {
17 
18  const DetId coreDet = vdetIds[trkIndex].detIdECAL;
19 #ifdef EDM_ML_DEBUG
20  if (debug) {
21  if (coreDet.subdetId() == EcalBarrel)
22  std::cout << "DetId " << (EBDetId)(coreDet) << " Flag " << vdetIds[trkIndex].okECAL << std::endl;
23  else
24  std::cout << "DetId " << (EEDetId)(coreDet) << " Flag " << vdetIds[trkIndex].okECAL << std::endl;
25  }
26 #endif
27  double maxNearP = -1.0;
28  if (vdetIds[trkIndex].okECAL) {
29  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
30 #ifdef EDM_ML_DEBUG
31  if (debug) std::cout << "chargeIsolationEcal:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
32 #endif
33 
34  for (unsigned int indx=0; indx<vdetIds.size(); ++indx) {
35  if (indx != trkIndex && vdetIds[indx].ok && vdetIds[indx].okECAL) {
36  const DetId anyCell = vdetIds[indx].detIdECAL;
37  if (!spr::chargeIsolation(anyCell,vdets)) {
38  const reco::Track* pTrack = &(*(vdetIds[indx].trkItr));
39  if (maxNearP < pTrack->p()) maxNearP = pTrack->p();
40  }
41  }
42  }
43  }
44  return maxNearP;
45  }
46 
47  //===========================================================================================================
48 
49  double chargeIsolationEcal(const DetId& coreDet, reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, const CaloGeometry* geo, const CaloTopology* caloTopology, const MagneticField* bField, int ieta, int iphi, std::string& theTrackQuality, bool debug) {
50 
53 
54  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
55 #ifdef EDM_ML_DEBUG
56  if (debug) std::cout << "chargeIsolation:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
57 #endif
58  double maxNearP = -1.0;
60 
61  // const DetId anyCell,
62  reco::TrackCollection::const_iterator trkItr2;
63  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
64 
65  const reco::Track* pTrack2 = &(*trkItr2);
66 
67  bool trkQuality = pTrack2->quality(trackQuality_);
68  if ( (trkItr2 != trkItr) && trkQuality ) {
69 
70  std::pair<math::XYZPoint,bool> info = spr::propagateECAL(pTrack2,bField);
71  const GlobalPoint point2(info.first.x(),info.first.y(),info.first.z());
72 
73  if (info.second) {
74  if (std::abs(point2.eta())<spr::etaBEEcal) {
75  const DetId anyCell = barrelGeom->getClosestCell(point2);
76  if (!spr::chargeIsolation(anyCell,vdets)) {
77 #ifdef EDM_ML_DEBUG
78  if (debug) std::cout << "chargeIsolationEcal Cell " << (EBDetId)(anyCell) << " pt " << pTrack2->p() << std::endl;
79 #endif
80  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
81  }
82  } else {
83  if (endcapGeom) {
84  const DetId anyCell = endcapGeom->getClosestCell(point2);
85  if (!spr::chargeIsolation(anyCell,vdets)) {
86 #ifdef EDM_ML_DEBUG
87  if (debug) std::cout << "chargeIsolationEcal Cell " << (EEDetId)(anyCell) << " pt " << pTrack2->p() << std::endl;
88 #endif
89  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
90  }
91  }
92  }
93  } //info.second
94  }
95  }
96  return maxNearP;
97  }
98 
99  //===========================================================================================================
100 
101  double chargeIsolationHcal(unsigned int trkIndex, std::vector<spr::propagatedTrackID> & vdetIds, const HcalTopology* topology, int ieta, int iphi, bool debug) {
102 
103  std::vector<DetId> dets(1,vdetIds[trkIndex].detIdHCAL);
104 #ifdef EDM_ML_DEBUG
105  if (debug) {
106  std::cout << "DetId " << (HcalDetId)(dets[0]) << " Flag " << vdetIds[trkIndex].okHCAL << std::endl;
107  }
108 #endif
109  double maxNearP = -1.0;
110  if (vdetIds[trkIndex].okHCAL) {
111  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
112 #ifdef EDM_ML_DEBUG
113  if (debug) std::cout << "chargeIsolationHcal:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
114 #endif
115  for (unsigned indx = 0; indx<vdetIds.size(); ++indx) {
116  if (indx != trkIndex && vdetIds[indx].ok && vdetIds[indx].okHCAL) {
117  const DetId anyCell = vdetIds[indx].detIdHCAL;
118  if (!spr::chargeIsolation(anyCell,vdets)) {
119  const reco::Track* pTrack = &(*(vdetIds[indx].trkItr));
120 #ifdef EDM_ML_DEBUG
121  if (debug) std::cout << "chargeIsolationHcal Cell " << (HcalDetId)(anyCell) << " pt " << pTrack->p() << std::endl;
122 #endif
123  if (maxNearP < pTrack->p()) maxNearP = pTrack->p();
124  }
125  }
126  }
127  }
128  return maxNearP;
129  }
130 
131  //===========================================================================================================
132 
133  double chargeIsolationHcal(reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, const DetId ClosestCell, const HcalTopology* topology, const CaloSubdetectorGeometry* gHB, const MagneticField* bField, int ieta, int iphi, std::string& theTrackQuality, bool debug) {
134 
135  std::vector<DetId> dets(1,ClosestCell);
136 #ifdef EDM_ML_DEBUG
137  if (debug) std::cout << (HcalDetId) ClosestCell << std::endl;
138 #endif
139  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
140  std::vector<DetId>::iterator it;
141 
142 #ifdef EDM_ML_DEBUG
143  if (debug) {
144  for (unsigned int i=0; i<vdets.size(); i++) {
145  std::cout << "HcalDetId in " <<2*ieta+1 << "x" << 2*iphi+1 << " " << (HcalDetId) vdets[i] << std::endl;
146  }
147  }
148 #endif
149  double maxNearP = -1.0;
150  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
151 
152  reco::TrackCollection::const_iterator trkItr2;
153  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
154 
155  const reco::Track* pTrack2 = &(*trkItr2);
156 
157  bool trkQuality = pTrack2->quality(trackQuality_);
158  if ( (trkItr2 != trkItr) && trkQuality ) {
159  std::pair<math::XYZPoint,bool> info = spr::propagateHCAL(pTrack2,bField);
160  const GlobalPoint point2(info.first.x(),info.first.y(),info.first.z());
161 
162 #ifdef EDM_ML_DEBUG
163  if (debug) {
164  std::cout << "Track2 (p,eta,phi) " << pTrack2->p() << " " << pTrack2->eta() << " " << pTrack2->phi() << std::endl;
165  }
166 #endif
167  if (info.second) {
168  const DetId anyCell = gHB->getClosestCell(point2);
169  if (!spr::chargeIsolation(anyCell,vdets)) {
170  if(maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
171  }
172 #ifdef EDM_ML_DEBUG
173  if (debug){
174  std::cout << "maxNearP " << maxNearP << " thisCell "
175  << (HcalDetId)anyCell << " ("
176  << info.first.x() << "," << info.first.y() <<","
177  << info.first.z() << ")" << std::endl;
178  }
179 #endif
180  }
181  }
182  }
183  return maxNearP;
184  }
185 
186  //===========================================================================================================
187 
188  bool chargeIsolation(const DetId anyCell, std::vector<DetId>& vdets) {
189  bool isIsolated = true;
190  for (unsigned int i=0; i<vdets.size(); i++){
191  if (anyCell == vdets[i] ) {
192  isIsolated = false;
193  break;
194  }
195  }
196  return isIsolated;
197  }
198 
199  //===========================================================================================================
200 
201  double coneChargeIsolation(const edm::Event& iEvent, const edm::EventSetup& iSetup, reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, TrackDetectorAssociator& associator, TrackAssociatorParameters& parameters_, std::string theTrackQuality, int &nNearTRKs, int &nLayers_maxNearP, int &trkQual_maxNearP, double &maxNearP_goodTrk, const GlobalPoint& hpoint1, const GlobalVector& trackMom, double dR) {
202 
203  nNearTRKs=0;
204  nLayers_maxNearP=0;
205  trkQual_maxNearP=-1;
206  maxNearP_goodTrk = -999.0;
207  double maxNearP = -999.0;
208  reco::TrackBase::TrackQuality trackQuality_= reco::TrackBase::qualityByName(theTrackQuality);
209 
210  // Iterate over tracks
211  reco::TrackCollection::const_iterator trkItr2;
212  for( trkItr2 = trkCollection->begin();
213  trkItr2 != trkCollection->end(); ++trkItr2){
214 
215  // Get track
216  const reco::Track* pTrack2 = &(*trkItr2);
217 
218  // Get track qual, nlayers, and hit pattern
219  if (pTrack2->quality(trackQuality_)) trkQual_maxNearP = 1;
220  const reco::HitPattern& hitp = pTrack2->hitPattern();
221  nLayers_maxNearP = hitp.trackerLayersWithMeasurement() ;
222 
223  // Skip if the neighboring track candidate is the iso-track
224  // candidate
225  if (trkItr2 != trkItr) {
226 
227  // Get propagator
228  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
229  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
230 
231  // Make sure it reaches Hcal
232  if (info2.isGoodHcal ) {
233 
234  const GlobalPoint point2(info2.trkGlobPosAtHcal.x(),
235  info2.trkGlobPosAtHcal.y(),
236  info2.trkGlobPosAtHcal.z());
237 
238  int isConeChargedIso = spr::coneChargeIsolation(hpoint1, point2, trackMom, dR);
239 
240  if (isConeChargedIso==0) {
241  nNearTRKs++;
242  if(maxNearP<pTrack2->p()) {
243  maxNearP=pTrack2->p();
244  if (trkQual_maxNearP>0 && nLayers_maxNearP>7 && maxNearP_goodTrk<pTrack2->p()) {
245  maxNearP_goodTrk=pTrack2->p();
246  }
247  }
248  }
249  }
250  }
251  } // Iterate over track loop
252 
253  return maxNearP;
254  }
255 
256  double chargeIsolationCone(unsigned int trkIndex, std::vector<spr::propagatedTrackDirection> & trkDirs, double dR, int & nNearTRKs, bool debug) {
257 
258  double maxNearP = -1.0;
259  nNearTRKs = 0;
260  if (trkDirs[trkIndex].okHCAL) {
261 #ifdef EDM_ML_DEBUG
262  if (debug) std::cout << "chargeIsolationCone with " << trkDirs.size() << " tracks " << std::endl;
263 #endif
264  for (unsigned int indx=0; indx<trkDirs.size(); ++indx) {
265  if (indx != trkIndex && trkDirs[indx].ok && trkDirs[indx].okHCAL) {
266  int isConeChargedIso = spr::coneChargeIsolation(trkDirs[trkIndex].pointHCAL, trkDirs[indx].pointHCAL, trkDirs[trkIndex].directionHCAL, dR);
267  if (isConeChargedIso==0) {
268  nNearTRKs++;
269  const reco::Track* pTrack = &(*(trkDirs[indx].trkItr));
270  if (maxNearP < pTrack->p()) maxNearP = pTrack->p();
271  }
272  }
273  }
274  }
275 #ifdef EDM_ML_DEBUG
276  if (debug) std::cout << "chargeIsolationCone Track " << trkDirs[trkIndex].okHCAL << " maxNearP " << maxNearP << std::endl;
277 #endif
278  return maxNearP;
279  }
280 
281 
282  int coneChargeIsolation(const GlobalPoint& hpoint1, const GlobalPoint& point2, const GlobalVector& trackMom, double dR) {
283 
284  int isIsolated = 1;
285  if (spr::getDistInPlaneTrackDir(hpoint1, trackMom, point2) > dR) isIsolated = 1;
286  else isIsolated = 0;
287  return isIsolated;
288  }
289 
290 }
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:44
static const double etaBEEcal
Definition: CaloConstants.h:12
static const TGPicture * info(bool iBackgroundIsBlack)
double getDistInPlaneTrackDir(const GlobalPoint &caloPoint, const GlobalVector &caloVector, const GlobalPoint &rechitPoint, bool debug=false)
Definition: FindDistCone.cc: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
std::pair< math::XYZPoint, bool > propagateHCAL(const reco::Track *, const MagneticField *, bool debug=false)
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
double chargeIsolationCone(unsigned int trkIndex, std::vector< spr::propagatedTrackDirection > &trkDirs, double dR, int &nNearTRKs, bool debug=false)
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:518
math::XYZPoint trkGlobPosAtHcal
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
std::pair< math::XYZPoint, bool > propagateECAL(const reco::Track *, const MagneticField *, bool debug=false)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
virtual DetId getClosestCell(const GlobalPoint &r) const
Definition: DetId.h:18
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
#define debug
Definition: HDRShower.cc:19
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:446
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)
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)
double coneChargeIsolation(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::TrackCollection::const_iterator trkItr, edm::Handle< reco::TrackCollection > trkCollection, TrackDetectorAssociator &associator, TrackAssociatorParameters &parameters_, std::string theTrackQuality, int &nNearTRKs, int &nLayers_maxNearP, int &trkQual_maxNearP, double &maxNearP_goodTrk, const GlobalPoint &hpoint1, const GlobalVector &trackMom, double dR)