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