CMS 3D CMS Logo

Classes | Functions
CommonUsefulStuff.h File Reference
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "DataFormats/DetId/interface/DetId.h"
#include "DataFormats/HcalDetId/interface/HcalDetId.h"
#include "DataFormats/EcalDetId/interface/EEDetId.h"
#include "DataFormats/EcalDetId/interface/EBDetId.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
#include "FWCore/Utilities/interface/Exception.h"
#include "DataFormats/GeometryVector/interface/GlobalPoint.h"
#include "DataFormats/GeometryVector/interface/GlobalVector.h"
#include "CommonTools/UtilAlgos/interface/DeltaR.h"

Go to the source code of this file.

Classes

struct  MaxHit_struct
 

Functions

double ecalEnergyInCone (const GlobalPoint center, double radius, const EcalRecHitCollection ecalCol, const CaloGeometry *geo)
 
double ecalEnergyInCone (const GlobalVector trackMom, const GlobalPoint center, double radius, const EcalRecHitCollection ecalCol, const CaloGeometry *geo)
 
double getDistInPlane (const GlobalVector trackDirection, const GlobalPoint caloPoint, const GlobalPoint rechitPoint, double coneHeight)
 
double getDistInPlaneSimple (const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
 
double getDistInPlaneTrackDir (const GlobalPoint caloPoint, const GlobalVector caloVector, const GlobalPoint rechitPoint)
 

Function Documentation

◆ ecalEnergyInCone() [1/2]

double ecalEnergyInCone ( const GlobalPoint  center,
double  radius,
const EcalRecHitCollection  ecalCol,
const CaloGeometry geo 
)
inline

Definition at line 156 of file CommonUsefulStuff.h.

References edm::SortedCollection< T, SORT >::begin(), EcalBarrel, EcalEndcap, edm::SortedCollection< T, SORT >::end(), getDistInPlaneSimple(), CaloGeometry::getPosition(), EBDetId::hashedIndex(), EEDetId::hashedIndex(), mps_fire::i, and CosmicsPD_Skims::radius.

159  {
160  double eECALcone = 0;
161  std::vector<int> usedHitsEcal;
162  usedHitsEcal.clear();
163 
164  for (std::vector<EcalRecHit>::const_iterator ehit = ecalCol.begin(); ehit != ecalCol.end(); ehit++) {
165  //This is a precaution for the case when hitCollection contains duplicats.
166  bool hitIsUsed = false;
167  int hitHashedIndex = -10000;
168  if (ehit->id().subdetId() == EcalBarrel) {
169  EBDetId did(ehit->id());
170  hitHashedIndex = did.hashedIndex();
171  }
172  if (ehit->id().subdetId() == EcalEndcap) {
173  EEDetId did(ehit->id());
174  hitHashedIndex = did.hashedIndex();
175  }
176  for (uint32_t i = 0; i < usedHitsEcal.size(); i++) {
177  if (usedHitsEcal[i] == hitHashedIndex)
178  hitIsUsed = true;
179  }
180  if (hitIsUsed)
181  continue;
182  usedHitsEcal.push_back(hitHashedIndex);
183  // -----------------------------------------------
184 
185  const GlobalPoint& pos = geo->getPosition((*ehit).detid());
186 
187  if (getDistInPlaneSimple(center, pos) < radius) {
188  eECALcone += ehit->energy();
189  }
190  }
191  return eECALcone;
192 }
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
const_iterator begin() const
const_iterator end() const
double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
int hashedIndex() const
Definition: EEDetId.h:183

◆ ecalEnergyInCone() [2/2]

double ecalEnergyInCone ( const GlobalVector  trackMom,
const GlobalPoint  center,
double  radius,
const EcalRecHitCollection  ecalCol,
const CaloGeometry geo 
)
inline

Definition at line 195 of file CommonUsefulStuff.h.

References edm::SortedCollection< T, SORT >::begin(), EcalBarrel, EcalEndcap, edm::SortedCollection< T, SORT >::end(), getDistInPlaneTrackDir(), CaloGeometry::getPosition(), EBDetId::hashedIndex(), EEDetId::hashedIndex(), mps_fire::i, and CosmicsPD_Skims::radius.

199  {
200  double eECALcone = 0;
201  std::vector<int> usedHitsEcal;
202  usedHitsEcal.clear();
203  for (std::vector<EcalRecHit>::const_iterator ehit = ecalCol.begin(); ehit != ecalCol.end(); ehit++) {
204  //This is a precaution for the case when hitCollection contains duplicats.
205  bool hitIsUsed = false;
206  int hitHashedIndex = -10000;
207  if (ehit->id().subdetId() == EcalBarrel) {
208  EBDetId did(ehit->id());
209  hitHashedIndex = did.hashedIndex();
210  }
211  if (ehit->id().subdetId() == EcalEndcap) {
212  EEDetId did(ehit->id());
213  hitHashedIndex = did.hashedIndex();
214  }
215  for (uint32_t i = 0; i < usedHitsEcal.size(); i++) {
216  if (usedHitsEcal[i] == hitHashedIndex)
217  hitIsUsed = true;
218  }
219  if (hitIsUsed)
220  continue;
221  usedHitsEcal.push_back(hitHashedIndex);
222  // -----------------------------------------------
223 
224  const GlobalPoint& pos = geo->getPosition((*ehit).detid());
225 
226  if (getDistInPlaneTrackDir(center, trackMom, pos) < radius) {
227  eECALcone += ehit->energy();
228  }
229  }
230  return eECALcone;
231 }
double getDistInPlaneTrackDir(const GlobalPoint caloPoint, const GlobalVector caloVector, const GlobalPoint rechitPoint)
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
const_iterator begin() const
const_iterator end() const
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
int hashedIndex() const
Definition: EEDetId.h:183

◆ getDistInPlane()

double getDistInPlane ( const GlobalVector  trackDirection,
const GlobalPoint  caloPoint,
const GlobalPoint  rechitPoint,
double  coneHeight 
)
inline

Definition at line 86 of file CommonUsefulStuff.h.

References Vector3DBase< T, FrameTag >::dot(), PV3DBase< T, PVType, FrameType >::mag(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

89  {
90  // The iso track candidate hits the Calo (Ecal or Hcal) at "caloPoint"
91  // with direction "trackDirection".
92 
93  // "rechitPoint" is the position of the rechit. We only care about
94  // the direction of the rechit.
95 
96  // Consider the rechitPoint as characterized by angles theta and phi
97  // wrt the origin which points at the calo cell of the rechit. In
98  // some sense the distance along the line theta/phi is arbitrary. A
99  // simplified choice might be to put the rechit at the surface of
100  // the Hcal. Our approach is to see whether/where this line
101  // intersects a cone of height "coneHeight" with vertex at caloPoint
102  // aligned with trackDirection.
103  // To that end, this function returns the distance between the
104  // center of the base of the cone and the intersection of the rechit
105  // line and the plane that contains the base of the cone. This
106  // distance is compared with the radius of the cone outside this
107  // function.
108 
109  // Make vector of length cone height along track direction
110  const GlobalVector heightVector = trackDirection * coneHeight;
111 
112  // Make vector from origin to point where iso track intersects the
113  // calorimeter.
114  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(), caloPoint.z());
115 
116  // Make vector from origin to point in center of base of cone
117  const GlobalVector coneBaseVector = heightVector + caloIntersectVector;
118 
119  // Make point in center of base of cone
120  const GlobalPoint coneBasePoint(coneBaseVector.x(), coneBaseVector.y(), coneBaseVector.z());
121 
122  // Make unit vector pointing at rechit.
123  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
124  const GlobalVector rechitDirection = rechitVector.unit();
125 
126  // Find distance "r" along "rechit line" (with angles theta2 and
127  // phi2) where line intersects plane defined by base of cone.
128 
129  // Definition plane of that contains base of cone:
130  // trackDirection.x() (x - coneBaseVector.x()) +
131  // trackDirection.y() (y - coneBaseVector.y()) +
132  // trackDirection.z() (z - coneBaseVector.z()) = 0
133 
134  // Point P_{rh} where rechit line intersects plane:
135  // (rechitdist sin(theta2) cos(phi2),
136  // rechitdist sin(theta2) cos(phi2),
137  // rechitdist cos(theta2))
138 
139  // Substitute P_{rh} into equation for plane and solve for rechitdist.
140  // rechitDist turns out to be the ratio of dot products:
141 
142  double rechitdist = trackDirection.dot(coneBaseVector) / trackDirection.dot(rechitDirection);
143 
144  // Now find distance between point at center of cone base and point
145  // where the "rechit line" intersects the plane defined by the base
146  // of the cone; i.e. the effectiveRecHitPoint.
147  const GlobalVector effectiveRechitVector = rechitdist * rechitDirection;
148  const GlobalPoint effectiveRechitPoint(
149  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
150 
151  GlobalVector distance_vector = effectiveRechitPoint - coneBasePoint;
152  return distance_vector.mag();
153 }
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T mag() const
Definition: PV3DBase.h:64
Vector3DBase unit() const
Definition: Vector3DBase.h:54

◆ getDistInPlaneSimple()

double getDistInPlaneSimple ( const GlobalPoint  caloPoint,
const GlobalPoint  rechitPoint 
)
inline

Definition at line 33 of file CommonUsefulStuff.h.

References Vector3DBase< T, FrameTag >::dot(), PV3DBase< T, PVType, FrameType >::mag(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by ecalEnergyInCone().

33  {
34  // Simplified version of getDistInPlane
35  // Assume track direction is origin -> point of hcal intersection
36 
37  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(), caloPoint.z());
38 
39  const GlobalVector caloIntersectUnitVector = caloIntersectVector.unit();
40 
41  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
42 
43  const GlobalVector rechitUnitVector = rechitVector.unit();
44  double dotprod = caloIntersectUnitVector.dot(rechitUnitVector);
45  double rechitdist = caloIntersectVector.mag() / dotprod;
46 
47  const GlobalVector effectiveRechitVector = rechitdist * rechitUnitVector;
48  const GlobalPoint effectiveRechitPoint(
49  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
50 
51  GlobalVector distance_vector = effectiveRechitPoint - caloPoint;
52 
53  if (dotprod > 0.) {
54  return distance_vector.mag();
55  } else {
56  return 999999.;
57  }
58 }
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T mag() const
Definition: PV3DBase.h:64
Vector3DBase unit() const
Definition: Vector3DBase.h:54

◆ getDistInPlaneTrackDir()

double getDistInPlaneTrackDir ( const GlobalPoint  caloPoint,
const GlobalVector  caloVector,
const GlobalPoint  rechitPoint 
)
inline

Definition at line 60 of file CommonUsefulStuff.h.

References Vector3DBase< T, FrameTag >::dot(), PV3DBase< T, PVType, FrameType >::mag(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by ecalEnergyInCone().

62  {
63  // Simplified version of getDistInPlane : no cone "within" Hcal, but
64  // don't assume track direction is origin -> point of hcal
65  // intersection.
66  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(),
67  caloPoint.z()); //p
68 
69  const GlobalVector caloUnitVector = caloVector.unit();
70  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
71  const GlobalVector rechitUnitVector = rechitVector.unit();
72  double dotprod_denominator = caloUnitVector.dot(rechitUnitVector);
73  double dotprod_numerator = caloUnitVector.dot(caloIntersectVector);
74  double rechitdist = dotprod_numerator / dotprod_denominator;
75  const GlobalVector effectiveRechitVector = rechitdist * rechitUnitVector;
76  const GlobalPoint effectiveRechitPoint(
77  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
78  GlobalVector distance_vector = effectiveRechitPoint - caloPoint;
79  if (dotprod_denominator > 0. && dotprod_numerator > 0.) {
80  return distance_vector.mag();
81  } else {
82  return 999999.;
83  }
84 }
T z() const
Definition: PV3DBase.h:61
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
T mag() const
Definition: PV3DBase.h:64
Vector3DBase unit() const
Definition: Vector3DBase.h:54