CMS 3D CMS Logo

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