CMS 3D CMS Logo

FindCaloHitCone.cc
Go to the documentation of this file.
7 #include <iostream>
8 
9 namespace spr {
10 
11  //Ecal Endcap OR Barrel RecHits
12  std::vector<EcalRecHitCollection::const_iterator> findCone(const CaloGeometry* geo,
14  const GlobalPoint& hpoint1,
15  const GlobalPoint& point1,
16  double dR,
17  const GlobalVector& trackMom,
18  bool debug) {
19  std::vector<EcalRecHitCollection::const_iterator> hit;
20 
21  for (EcalRecHitCollection::const_iterator j = hits->begin(); j != hits->end(); j++) {
22  bool keepHit = false;
23 
24  if (j->id().subdetId() == EcalEndcap) {
25  EEDetId EEid = EEDetId(j->id());
26  const GlobalPoint& rechitPoint = geo->getPosition(EEid);
27  if (spr::getDistInPlaneTrackDir(point1, trackMom, rechitPoint, debug) < dR)
28  keepHit = true;
29  } else if (j->id().subdetId() == EcalBarrel) {
30  EBDetId EBid = EBDetId(j->id());
31  const GlobalPoint& rechitPoint = geo->getPosition(EBid);
32  if (spr::getDistInPlaneTrackDir(point1, trackMom, rechitPoint, debug) < dR)
33  keepHit = true;
34  }
35 
36  if (keepHit)
37  hit.push_back(j);
38  }
39  return hit;
40  }
41 
42  // Ecal Endcap AND Barrel RecHits
43  std::vector<EcalRecHitCollection::const_iterator> findCone(const CaloGeometry* geo,
46  const GlobalPoint& hpoint1,
47  const GlobalPoint& point1,
48  double dR,
49  const GlobalVector& trackMom,
50  bool debug) {
51  std::vector<EcalRecHitCollection::const_iterator> hit;
52 
53  // Only check both barrel and endcap when track is near transition
54  // region: 1.479-2*0.087 < trkEta < 1.479+2*0.087
55 
56  bool doBarrel = false, doEndcap = false;
57  if (std::abs(point1.eta()) < (spr::etaBEEcal + 2 * spr::deltaEta))
58  doBarrel = true; // 1.479+2*0.087
59  if (std::abs(point1.eta()) > (spr::etaBEEcal - 2 * spr::deltaEta))
60  doEndcap = true; // 1.479-2*0.087
61 
62  if (doBarrel) {
63  for (EcalRecHitCollection::const_iterator j = barrelhits->begin(); j != barrelhits->end(); j++) {
64  bool keepHit = false;
65  if (j->id().subdetId() == EcalBarrel) {
66  EBDetId EBid = EBDetId(j->id());
67  const GlobalPoint& rechitPoint = geo->getPosition(EBid);
68  if (spr::getDistInPlaneTrackDir(point1, trackMom, rechitPoint, debug) < dR)
69  keepHit = true;
70  } else {
71  edm::LogWarning("IsoTrack") << "PROBLEM : Endcap RecHits in Barrel Collection!?";
72  }
73  if (keepHit)
74  hit.push_back(j);
75  }
76  } // doBarrel
77 
78  if (doEndcap) {
79  for (EcalRecHitCollection::const_iterator j = endcaphits->begin(); j != endcaphits->end(); j++) {
80  bool keepHit = false;
81 
82  if (j->id().subdetId() == EcalEndcap) {
83  EEDetId EEid = EEDetId(j->id());
84  const GlobalPoint& rechitPoint = geo->getPosition(EEid);
85  if (spr::getDistInPlaneTrackDir(point1, trackMom, rechitPoint, debug) < dR)
86  keepHit = true;
87  } else {
88  edm::LogWarning("IsoTrack") << "PROBLEM : Barrel RecHits in Endcap Collection!?";
89  }
90  if (keepHit)
91  hit.push_back(j);
92  }
93  } // doEndcap
94 
95  return hit;
96  }
97 
98  //HBHE RecHits
99  std::vector<HBHERecHitCollection::const_iterator> findCone(const CaloGeometry* geo,
101  const GlobalPoint& hpoint1,
102  const GlobalPoint& point1,
103  double dR,
104  const GlobalVector& trackMom,
105  bool debug) {
106  std::vector<HBHERecHitCollection::const_iterator> hit;
107  // Loop over Hcal RecHits
108  for (HBHERecHitCollection::const_iterator j = hits->begin(); j != hits->end(); j++) {
109  DetId detId(j->id());
110  const GlobalPoint rechitPoint =
111  (static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(detId)))->getPosition(detId);
112  if (spr::getDistInPlaneTrackDir(hpoint1, trackMom, rechitPoint, debug) < dR)
113  hit.push_back(j);
114  }
115  return hit;
116  }
117 
118  // PCalo SimHits
119  std::vector<edm::PCaloHitContainer::const_iterator> findCone(const CaloGeometry* geo,
121  const GlobalPoint& hpoint1,
122  const GlobalPoint& point1,
123  double dR,
124  const GlobalVector& trackMom,
125  bool debug) {
126  std::vector<edm::PCaloHitContainer::const_iterator> hit;
127  edm::PCaloHitContainer::const_iterator ihit;
128  for (ihit = hits->begin(); ihit != hits->end(); ihit++) {
129  DetId detId(ihit->id());
130  const GlobalPoint rechitPoint =
131  (detId.det() == DetId::Hcal)
132  ? (static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(detId)))->getPosition(detId)
133  : geo->getPosition(detId);
134  if (spr::getDistInPlaneTrackDir(hpoint1, trackMom, rechitPoint, debug) < dR) {
135  hit.push_back(ihit);
136  }
137  }
138  return hit;
139  }
140 
141 } // namespace spr
static const double etaBEEcal
Definition: CaloConstants.h:12
double getDistInPlaneTrackDir(const GlobalPoint &caloPoint, const GlobalVector &caloVector, const GlobalPoint &rechitPoint, bool debug=false)
Definition: FindDistCone.cc:9
T eta() const
Definition: PV3DBase.h:73
std::vector< EcalRecHit >::const_iterator const_iterator
static const double deltaEta
Definition: CaloConstants.h:8
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const_iterator begin() const
const_iterator end() const
Definition: DetId.h:17
#define debug
Definition: HDRShower.cc:19
Log< level::Warning, false > LogWarning
std::vector< EcalRecHitCollection::const_iterator > findCone(const CaloGeometry *geo, edm::Handle< EcalRecHitCollection > &hits, const GlobalPoint &hpoint1, const GlobalPoint &point1, double dR, bool debug=false)
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34