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 
51  const EcalBarrelGeometry *barrelGeom = (dynamic_cast< const EcalBarrelGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalBarrel)));
52  const EcalEndcapGeometry *endcapGeom = (dynamic_cast< const EcalEndcapGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalEndcap)));
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  const DetId anyCell = endcapGeom->getClosestCell(point2);
84  if (!spr::chargeIsolation(anyCell,vdets)) {
85 #ifdef EDM_ML_DEBUG
86  if (debug) std::cout << "chargeIsolationEcal Cell " << (EEDetId)(anyCell) << " pt " << pTrack2->p() << std::endl;
87 #endif
88  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
89  }
90  }
91  } //info.second
92  }
93  }
94  return maxNearP;
95  }
96 
97  //===========================================================================================================
98 
99  double chargeIsolationHcal(unsigned int trkIndex, std::vector<spr::propagatedTrackID> & vdetIds, const HcalTopology* topology, int ieta, int iphi, bool debug) {
100 
101  std::vector<DetId> dets(1,vdetIds[trkIndex].detIdHCAL);
102 #ifdef EDM_ML_DEBUG
103  if (debug) {
104  std::cout << "DetId " << (HcalDetId)(dets[0]) << " Flag " << vdetIds[trkIndex].okHCAL << std::endl;
105  }
106 #endif
107  double maxNearP = -1.0;
108  if (vdetIds[trkIndex].okHCAL) {
109  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
110 #ifdef EDM_ML_DEBUG
111  if (debug) std::cout << "chargeIsolationHcal:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
112 #endif
113  for (unsigned indx = 0; indx<vdetIds.size(); ++indx) {
114  if (indx != trkIndex && vdetIds[indx].ok && vdetIds[indx].okHCAL) {
115  const DetId anyCell = vdetIds[indx].detIdHCAL;
116  if (!spr::chargeIsolation(anyCell,vdets)) {
117  const reco::Track* pTrack = &(*(vdetIds[indx].trkItr));
118 #ifdef EDM_ML_DEBUG
119  if (debug) std::cout << "chargeIsolationHcal Cell " << (HcalDetId)(anyCell) << " pt " << pTrack->p() << std::endl;
120 #endif
121  if (maxNearP < pTrack->p()) maxNearP = pTrack->p();
122  }
123  }
124  }
125  }
126  return maxNearP;
127  }
128 
129  //===========================================================================================================
130 
131  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) {
132 
133  std::vector<DetId> dets(1,ClosestCell);
134 #ifdef EDM_ML_DEBUG
135  if (debug) std::cout << (HcalDetId) ClosestCell << std::endl;
136 #endif
137  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
138  std::vector<DetId>::iterator it;
139 
140 #ifdef EDM_ML_DEBUG
141  if (debug) {
142  for (unsigned int i=0; i<vdets.size(); i++) {
143  std::cout << "HcalDetId in " <<2*ieta+1 << "x" << 2*iphi+1 << " " << (HcalDetId) vdets[i] << std::endl;
144  }
145  }
146 #endif
147  double maxNearP = -1.0;
148  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
149 
150  reco::TrackCollection::const_iterator trkItr2;
151  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
152 
153  const reco::Track* pTrack2 = &(*trkItr2);
154 
155  bool trkQuality = pTrack2->quality(trackQuality_);
156  if ( (trkItr2 != trkItr) && trkQuality ) {
157  std::pair<math::XYZPoint,bool> info = spr::propagateHCAL(pTrack2,bField);
158  const GlobalPoint point2(info.first.x(),info.first.y(),info.first.z());
159 
160 #ifdef EDM_ML_DEBUG
161  if (debug) {
162  std::cout << "Track2 (p,eta,phi) " << pTrack2->p() << " " << pTrack2->eta() << " " << pTrack2->phi() << std::endl;
163  }
164 #endif
165  if (info.second) {
166  const DetId anyCell = gHB->getClosestCell(point2);
167  if (!spr::chargeIsolation(anyCell,vdets)) {
168  if(maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
169  }
170 #ifdef EDM_ML_DEBUG
171  if (debug){
172  std::cout << "maxNearP " << maxNearP << " thisCell "
173  << (HcalDetId)anyCell << " ("
174  << info.first.x() << "," << info.first.y() <<","
175  << info.first.z() << ")" << std::endl;
176  }
177 #endif
178  }
179  }
180  }
181  return maxNearP;
182  }
183 
184  //===========================================================================================================
185 
186  bool chargeIsolation(const DetId anyCell, std::vector<DetId>& vdets) {
187  bool isIsolated = true;
188  for (unsigned int i=0; i<vdets.size(); i++){
189  if (anyCell == vdets[i] ) {
190  isIsolated = false;
191  break;
192  }
193  }
194  return isIsolated;
195  }
196 
197  //===========================================================================================================
198 
199  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) {
200 
201  nNearTRKs=0;
202  nLayers_maxNearP=0;
203  trkQual_maxNearP=-1;
204  maxNearP_goodTrk = -999.0;
205  double maxNearP = -999.0;
206  reco::TrackBase::TrackQuality trackQuality_= reco::TrackBase::qualityByName(theTrackQuality);
207 
208  // Iterate over tracks
209  reco::TrackCollection::const_iterator trkItr2;
210  for( trkItr2 = trkCollection->begin();
211  trkItr2 != trkCollection->end(); ++trkItr2){
212 
213  // Get track
214  const reco::Track* pTrack2 = &(*trkItr2);
215 
216  // Get track qual, nlayers, and hit pattern
217  if (pTrack2->quality(trackQuality_)) trkQual_maxNearP = 1;
218  const reco::HitPattern& hitp = pTrack2->hitPattern();
219  nLayers_maxNearP = hitp.trackerLayersWithMeasurement() ;
220 
221  // Skip if the neighboring track candidate is the iso-track
222  // candidate
223  if (trkItr2 != trkItr) {
224 
225  // Get propagator
226  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
227  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
228 
229  // Make sure it reaches Hcal
230  if (info2.isGoodHcal ) {
231 
232  const GlobalPoint point2(info2.trkGlobPosAtHcal.x(),
233  info2.trkGlobPosAtHcal.y(),
234  info2.trkGlobPosAtHcal.z());
235 
236  int isConeChargedIso = spr::coneChargeIsolation(hpoint1, point2, trackMom, dR);
237 
238  if (isConeChargedIso==0) {
239  nNearTRKs++;
240  if(maxNearP<pTrack2->p()) {
241  maxNearP=pTrack2->p();
242  if (trkQual_maxNearP>0 && nLayers_maxNearP>7 && maxNearP_goodTrk<pTrack2->p()) {
243  maxNearP_goodTrk=pTrack2->p();
244  }
245  }
246  }
247  }
248  }
249  } // Iterate over track loop
250 
251  return maxNearP;
252  }
253 
254  double chargeIsolationCone(unsigned int trkIndex, std::vector<spr::propagatedTrackDirection> & trkDirs, double dR, int & nNearTRKs, bool debug) {
255 
256  double maxNearP = -1.0;
257  nNearTRKs = 0;
258  if (trkDirs[trkIndex].okHCAL) {
259 #ifdef EDM_ML_DEBUG
260  if (debug) std::cout << "chargeIsolationCone with " << trkDirs.size() << " tracks " << std::endl;
261 #endif
262  for (unsigned int indx=0; indx<trkDirs.size(); ++indx) {
263  if (indx != trkIndex && trkDirs[indx].ok && trkDirs[indx].okHCAL) {
264  int isConeChargedIso = spr::coneChargeIsolation(trkDirs[trkIndex].pointHCAL, trkDirs[indx].pointHCAL, trkDirs[trkIndex].directionHCAL, dR);
265  if (isConeChargedIso==0) {
266  nNearTRKs++;
267  const reco::Track* pTrack = &(*(trkDirs[indx].trkItr));
268  if (maxNearP < pTrack->p()) maxNearP = pTrack->p();
269  }
270  }
271  }
272  }
273 #ifdef EDM_ML_DEBUG
274  if (debug) std::cout << "chargeIsolationCone Track " << trkDirs[trkIndex].okHCAL << " maxNearP " << maxNearP << std::endl;
275 #endif
276  return maxNearP;
277  }
278 
279 
280  int coneChargeIsolation(const GlobalPoint& hpoint1, const GlobalPoint& point2, const GlobalVector& trackMom, double dR) {
281 
282  int isIsolated = 1;
283  if (spr::getDistInPlaneTrackDir(hpoint1, trackMom, point2) > dR) isIsolated = 1;
284  else isIsolated = 0;
285  return isIsolated;
286  }
287 
288 }
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:45
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
DetId getClosestCell(const GlobalPoint &r) const override
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:512
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 &)
DetId getClosestCell(const GlobalPoint &r) const override
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)