CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ChargeIsolationExtra.cc
Go to the documentation of this file.
9 
10 #include<iostream>
11 
12 namespace spr{
13 
14  double chargeIsolation(const edm::Event& iEvent, const edm::EventSetup& iSetup, CaloNavigator<DetId>& theNavigator, reco::TrackCollection::const_iterator trkItr, edm::Handle<reco::TrackCollection> trkCollection, const CaloSubdetectorGeometry* gEB, const CaloSubdetectorGeometry* gEE, TrackDetectorAssociator& associator, TrackAssociatorParameters& parameters_, int ieta, int iphi, std::string theTrackQuality, bool debug) {
15 
16  double maxNearP = -1.0;
18 
19  // const DetId anyCell,
20  reco::TrackCollection::const_iterator trkItr2;
21  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
22 
23  const reco::Track* pTrack2 = &(*trkItr2);
24 
25  bool trkQuality = pTrack2->quality(trackQuality_);
26  if ( (trkItr2 != trkItr) && trkQuality ) {
27 
28  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
29  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
30  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
31 
32  if (info2.isGoodEcal ) {
33  if (std::abs(point2.eta())<spr::etaBEEcal) {
34  const DetId anyCell = gEB->getClosestCell(point2);
35  if (debug) std::cout << "chargeIsolation:: EB cell " << (EBDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
36  if (!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
37  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
38  }
39  } else {
40  const DetId anyCell = gEE->getClosestCell(point2);
41  if (debug) std::cout << "chargeIsolation:: EE cell " << (EEDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
42  if(!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
43  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
44  }
45  }
46  } //info2.isGoodEcal
47  }
48  }
49  return maxNearP;
50  }
51 
52  //===========================================================================================================
53 
54  bool chargeIsolation(const DetId anyCell,CaloNavigator<DetId>& navigator,int ieta,int iphi) {
55 
56  bool isIsolated = false;
57 
58  DetId thisDet;
59 
60  for (int dx = -ieta; dx < ieta+1; ++dx) {
61  for (int dy = -iphi; dy < iphi+1; ++dy) {
62 
63  thisDet = navigator.offsetBy(dx, dy);
64  navigator.home();
65 
66  if (thisDet != DetId(0)) {
67  if (thisDet == anyCell) {
68  isIsolated = false;
69  return isIsolated;
70  }
71  }
72  }
73  }
74  return isIsolated;
75  }
76 
77  //===========================================================================================================
78 
79  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, std::string& theTrackQuality, bool debug) {
80 
81  const EcalBarrelGeometry *barrelGeom = (dynamic_cast< const EcalBarrelGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalBarrel)));
82  const EcalEndcapGeometry *endcapGeom = (dynamic_cast< const EcalEndcapGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalEndcap)));
83 
84  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
85  if (debug) std::cout << "chargeIsolation:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
86 
87  double maxNearP = -1.0;
89 
90  // const DetId anyCell,
91  reco::TrackCollection::const_iterator trkItr2;
92  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
93 
94  const reco::Track* pTrack2 = &(*trkItr2);
95 
96  bool trkQuality = pTrack2->quality(trackQuality_);
97  if ( (trkItr2 != trkItr) && trkQuality ) {
98 
99  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
100  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
101  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
102 
103  if (info2.isGoodEcal ) {
104  if (std::abs(point2.eta())<spr::etaBEEcal) {
105  const DetId anyCell = barrelGeom->getClosestCell(point2);
106  if (debug) std::cout << "chargeIsolation:: EB cell " << (EBDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
107  if (!spr::chargeIsolation(anyCell,vdets)) {
108  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
109  }
110  } else {
111  const DetId anyCell = endcapGeom->getClosestCell(point2);
112  if (debug) std::cout << "chargeIsolation:: EE cell " << (EEDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
113  if (!spr::chargeIsolation(anyCell,vdets)) {
114  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
115  }
116  }
117  } //info2.isGoodEcal
118  }
119  }
120  return maxNearP;
121  }
122 
123  //===========================================================================================================
124 
125  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, std::string& theTrackQuality, bool debug) {
126 
127  std::vector<DetId> dets(1,ClosestCell);
128 
129  if (debug) std::cout << (HcalDetId) ClosestCell << std::endl;
130 
131  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
132  std::vector<DetId>::iterator it;
133 
134  if (debug) {
135  for (unsigned int i=0; i<vdets.size(); i++) {
136  std::cout << "HcalDetId in " <<2*ieta+1 << "x" << 2*iphi+1 << " " << (HcalDetId) vdets[i] << std::endl;
137  }
138  }
139 
140  double maxNearP = -1.0;
141  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
142 
143  reco::TrackCollection::const_iterator trkItr2;
144  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
145 
146  const reco::Track* pTrack2 = &(*trkItr2);
147 
148  bool trkQuality = pTrack2->quality(trackQuality_);
149  if ( (trkItr2 != trkItr) && trkQuality ) {
150  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
151  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
152  const GlobalPoint point2(info2.trkGlobPosAtHcal.x(),info2.trkGlobPosAtHcal.y(),info2.trkGlobPosAtHcal.z());
153 
154  if (debug) {
155  std::cout << "Track2 (p,eta,phi) " << pTrack2->p() << " " << pTrack2->eta() << " " << pTrack2->phi() << std::endl;
156  }
157 
158  if (info2.isGoodHcal ) {
159  const DetId anyCell = gHB->getClosestCell(point2);
160  if (debug) std::cout << "chargeIsolation:: HCAL cell " << (HcalDetId)(anyCell) << " for pt " << pTrack2->p() << std::endl;
161  if (!spr::chargeIsolation(anyCell,vdets)) {
162  if(maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
163  }
164  if (debug){
165  std::cout << "maxNearP " << maxNearP << " thisCell "
166  << (HcalDetId)anyCell << " ("
167  << info2.trkGlobPosAtHcal.x() << ","
168  << info2.trkGlobPosAtHcal.y() <<","
169  << info2.trkGlobPosAtHcal.z() <<")" << std::endl;
170  }
171  }
172  }
173  }
174  return maxNearP;
175  }
176 
177 }
virtual DetId getClosestCell(const GlobalPoint &r) const
double p() const
momentum vector magnitude
Definition: TrackBase.h:610
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:45
int i
Definition: DBlmapReader.cc:9
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)
virtual DetId getClosestCell(const GlobalPoint &r) const
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:640
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:646
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:505
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 &)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
tuple cout
Definition: gather_cfg.py:145
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)