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.
8 
9 #include<iostream>
10 
11 namespace spr{
12 
13  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) {
14 
15  double maxNearP = -1.0;
17 
18  // const DetId anyCell,
19  reco::TrackCollection::const_iterator trkItr2;
20  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
21 
22  const reco::Track* pTrack2 = &(*trkItr2);
23 
24  bool trkQuality = pTrack2->quality(trackQuality_);
25  if ( (trkItr2 != trkItr) && trkQuality ) {
26 
27  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
28  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
29  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
30 
31  if (info2.isGoodEcal ) {
32  if (std::abs(point2.eta())<1.479) {
33  const DetId anyCell = gEB->getClosestCell(point2);
34  if (!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
35  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
36  }
37  } else {
38  const DetId anyCell = gEE->getClosestCell(point2);
39  if(!spr::chargeIsolation(anyCell,theNavigator,ieta, iphi)) {
40  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
41  }
42  }
43  } //info2.isGoodEcal
44  }
45  }
46  return maxNearP;
47  }
48 
49  //===========================================================================================================
50 
51  bool chargeIsolation(const DetId anyCell,CaloNavigator<DetId>& navigator,int ieta,int iphi) {
52 
53  bool isIsolated = false;
54 
55  DetId thisDet;
56 
57  for (int dx = -ieta; dx < ieta+1; ++dx) {
58  for (int dy = -iphi; dy < iphi+1; ++dy) {
59 
60  thisDet = navigator.offsetBy(dx, dy);
61  navigator.home();
62 
63  if (thisDet != DetId(0)) {
64  if (thisDet == anyCell) {
65  isIsolated = false;
66  return isIsolated;
67  }
68  }
69  }
70  }
71  return isIsolated;
72  }
73 
74  //===========================================================================================================
75 
76  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) {
77 
78  const EcalBarrelGeometry *barrelGeom = (dynamic_cast< const EcalBarrelGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalBarrel)));
79  const EcalEndcapGeometry *endcapGeom = (dynamic_cast< const EcalEndcapGeometry *> (geo->getSubdetectorGeometry(DetId::Ecal,EcalEndcap)));
80 
81  std::vector<DetId> vdets = spr::matrixECALIds(coreDet, ieta, iphi, geo, caloTopology, debug);
82  if (debug) std::cout << "chargeIsolation:: eta/phi/dets " << ieta << " " << iphi << " " << vdets.size() << std::endl;
83 
84  double maxNearP = -1.0;
86 
87  // const DetId anyCell,
88  reco::TrackCollection::const_iterator trkItr2;
89  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
90 
91  const reco::Track* pTrack2 = &(*trkItr2);
92 
93  bool trkQuality = pTrack2->quality(trackQuality_);
94  if ( (trkItr2 != trkItr) && trkQuality ) {
95 
96  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
97  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
98  const GlobalPoint point2(info2.trkGlobPosAtEcal.x(),info2.trkGlobPosAtEcal.y(),info2.trkGlobPosAtEcal.z());
99 
100  if (info2.isGoodEcal ) {
101  if (std::abs(point2.eta())<1.479) {
102  const DetId anyCell = barrelGeom->getClosestCell(point2);
103  if (!spr::chargeIsolation(anyCell,vdets)) {
104  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
105  }
106  } else {
107  const DetId anyCell = endcapGeom->getClosestCell(point2);
108  if (!spr::chargeIsolation(anyCell,vdets)) {
109  if (maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
110  }
111  }
112  } //info2.isGoodEcal
113  }
114  }
115  return maxNearP;
116  }
117 
118  //===========================================================================================================
119 
120  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) {
121 
122  std::vector<DetId> dets(1,ClosestCell);
123 
124  if (debug) std::cout << (HcalDetId) ClosestCell << std::endl;
125 
126  std::vector<DetId> vdets = spr::matrixHCALIds(dets, topology, ieta, iphi, false, debug);
127  std::vector<DetId>::iterator it;
128 
129  if (debug) {
130  for (unsigned int i=0; i<vdets.size(); i++) {
131  std::cout << "HcalDetId in " <<2*ieta+1 << "x" << 2*iphi+1 << " " << (HcalDetId) vdets[i] << std::endl;
132  }
133  }
134 
135  double maxNearP = -1.0;
136  reco::TrackBase::TrackQuality trackQuality_=reco::TrackBase::qualityByName(theTrackQuality);
137 
138  reco::TrackCollection::const_iterator trkItr2;
139  for (trkItr2 = trkCollection->begin(); trkItr2 != trkCollection->end(); ++trkItr2) {
140 
141  const reco::Track* pTrack2 = &(*trkItr2);
142 
143  bool trkQuality = pTrack2->quality(trackQuality_);
144  if ( (trkItr2 != trkItr) && trkQuality ) {
145  const FreeTrajectoryState fts2 = associator.getFreeTrajectoryState(iSetup, *pTrack2);
146  TrackDetMatchInfo info2 = associator.associate(iEvent, iSetup, fts2, parameters_);
147  const GlobalPoint point2(info2.trkGlobPosAtHcal.x(),info2.trkGlobPosAtHcal.y(),info2.trkGlobPosAtHcal.z());
148 
149  if (debug) {
150  std::cout << "Track2 (p,eta,phi) " << pTrack2->p() << " " << pTrack2->eta() << " " << pTrack2->phi() << std::endl;
151  }
152 
153  if (info2.isGoodHcal ) {
154  const DetId anyCell = gHB->getClosestCell(point2);
155  if (!spr::chargeIsolation(anyCell,vdets)) {
156  if(maxNearP<pTrack2->p()) maxNearP=pTrack2->p();
157  }
158  if (debug){
159  std::cout << "maxNearP " << maxNearP << " thisCell "
160  << (HcalDetId)anyCell << " ("
161  << info2.trkGlobPosAtHcal.x() << ","
162  << info2.trkGlobPosAtHcal.y() <<","
163  << info2.trkGlobPosAtHcal.z() <<")" << std::endl;
164  }
165  }
166  }
167  }
168  return maxNearP;
169  }
170 
171 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:591
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
CaloTopology const * topology(0)
bool chargeIsolation(const DetId anyCell, std::vector< DetId > &vdets)
TrackQuality
track quality
Definition: TrackBase.h:113
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:621
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:627
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual DetId getClosestCell(const GlobalPoint &r) const
Definition: DetId.h:18
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:102
#define debug
Definition: HDRShower.cc:19
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:479
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:121
double chargeIsolationHcal(unsigned int trkIndex, std::vector< spr::propagatedTrackID > &vdetIds, const HcalTopology *topology, int ieta, int iphi, bool debug=false)