CMS 3D CMS Logo

CommonUsefulStuff.h
Go to the documentation of this file.
1 #ifndef Calibration_HcalCalibALgos_CommonUsefulStuff_h
2 #define Calibration_HcalCalibALgos_CommonUsefulStuff_h
3 
10 
19 
22 
24 
27 
31 
32 #include "TROOT.h"
33 #include "TFile.h"
34 #include "TTree.h"
35 #include "TH1F.h"
36 #include "TString.h"
37 
38 /* getDist functions by Jim:
39 /uscms/home/jhirsch/IsoTracks_314/CMSSW_3_1_4/src/JetMETCorrections/IsolatedParticles/interface/FindCaloHit.icc
40 */
41 
42 struct MaxHit_struct {
43  int iphihitm;
44  int ietahitm;
45  int depthhit;
46  float hitenergy;
47  float dr;
49  MaxHit_struct() : iphihitm(0), ietahitm(0), depthhit(0), hitenergy(-100), dr(0) {}
50 };
51 
52 inline double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint) {
53  // Simplified version of getDistInPlane
54  // Assume track direction is origin -> point of hcal intersection
55 
56  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(), caloPoint.z());
57 
58  const GlobalVector caloIntersectUnitVector = caloIntersectVector.unit();
59 
60  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
61 
62  const GlobalVector rechitUnitVector = rechitVector.unit();
63  double dotprod = caloIntersectUnitVector.dot(rechitUnitVector);
64  double rechitdist = caloIntersectVector.mag() / dotprod;
65 
66  const GlobalVector effectiveRechitVector = rechitdist * rechitUnitVector;
67  const GlobalPoint effectiveRechitPoint(
68  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
69 
70  GlobalVector distance_vector = effectiveRechitPoint - caloPoint;
71 
72  if (dotprod > 0.) {
73  return distance_vector.mag();
74  } else {
75  return 999999.;
76  }
77 }
78 
79 inline double getDistInPlaneTrackDir(const GlobalPoint caloPoint,
80  const GlobalVector caloVector,
81  const GlobalPoint rechitPoint) {
82  // Simplified version of getDistInPlane : no cone "within" Hcal, but
83  // don't assume track direction is origin -> point of hcal
84  // intersection.
85  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(),
86  caloPoint.z()); //p
87 
88  const GlobalVector caloUnitVector = caloVector.unit();
89  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
90  const GlobalVector rechitUnitVector = rechitVector.unit();
91  double dotprod_denominator = caloUnitVector.dot(rechitUnitVector);
92  double dotprod_numerator = caloUnitVector.dot(caloIntersectVector);
93  double rechitdist = dotprod_numerator / dotprod_denominator;
94  const GlobalVector effectiveRechitVector = rechitdist * rechitUnitVector;
95  const GlobalPoint effectiveRechitPoint(
96  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
97  GlobalVector distance_vector = effectiveRechitPoint - caloPoint;
98  if (dotprod_denominator > 0. && dotprod_numerator > 0.) {
99  return distance_vector.mag();
100  } else {
101  return 999999.;
102  }
103 }
104 
105 inline double getDistInPlane(const GlobalVector trackDirection,
106  const GlobalPoint caloPoint,
107  const GlobalPoint rechitPoint,
108  double coneHeight) {
109  // The iso track candidate hits the Calo (Ecal or Hcal) at "caloPoint"
110  // with direction "trackDirection".
111 
112  // "rechitPoint" is the position of the rechit. We only care about
113  // the direction of the rechit.
114 
115  // Consider the rechitPoint as characterized by angles theta and phi
116  // wrt the origin which points at the calo cell of the rechit. In
117  // some sense the distance along the line theta/phi is arbitrary. A
118  // simplified choice might be to put the rechit at the surface of
119  // the Hcal. Our approach is to see whether/where this line
120  // intersects a cone of height "coneHeight" with vertex at caloPoint
121  // aligned with trackDirection.
122  // To that end, this function returns the distance between the
123  // center of the base of the cone and the intersection of the rechit
124  // line and the plane that contains the base of the cone. This
125  // distance is compared with the radius of the cone outside this
126  // function.
127 
128  // Make vector of length cone height along track direction
129  const GlobalVector heightVector = trackDirection * coneHeight;
130 
131  // Make vector from origin to point where iso track intersects the
132  // calorimeter.
133  const GlobalVector caloIntersectVector(caloPoint.x(), caloPoint.y(), caloPoint.z());
134 
135  // Make vector from origin to point in center of base of cone
136  const GlobalVector coneBaseVector = heightVector + caloIntersectVector;
137 
138  // Make point in center of base of cone
139  const GlobalPoint coneBasePoint(coneBaseVector.x(), coneBaseVector.y(), coneBaseVector.z());
140 
141  // Make unit vector pointing at rechit.
142  const GlobalVector rechitVector(rechitPoint.x(), rechitPoint.y(), rechitPoint.z());
143  const GlobalVector rechitDirection = rechitVector.unit();
144 
145  // Find distance "r" along "rechit line" (with angles theta2 and
146  // phi2) where line intersects plane defined by base of cone.
147 
148  // Definition plane of that contains base of cone:
149  // trackDirection.x() (x - coneBaseVector.x()) +
150  // trackDirection.y() (y - coneBaseVector.y()) +
151  // trackDirection.z() (z - coneBaseVector.z()) = 0
152 
153  // Point P_{rh} where rechit line intersects plane:
154  // (rechitdist sin(theta2) cos(phi2),
155  // rechitdist sin(theta2) cos(phi2),
156  // rechitdist cos(theta2))
157 
158  // Substitute P_{rh} into equation for plane and solve for rechitdist.
159  // rechitDist turns out to be the ratio of dot products:
160 
161  double rechitdist = trackDirection.dot(coneBaseVector) / trackDirection.dot(rechitDirection);
162 
163  // Now find distance between point at center of cone base and point
164  // where the "rechit line" intersects the plane defined by the base
165  // of the cone; i.e. the effectiveRecHitPoint.
166  const GlobalVector effectiveRechitVector = rechitdist * rechitDirection;
167  const GlobalPoint effectiveRechitPoint(
168  effectiveRechitVector.x(), effectiveRechitVector.y(), effectiveRechitVector.z());
169 
170  GlobalVector distance_vector = effectiveRechitPoint - coneBasePoint;
171  return distance_vector.mag();
172 }
173 
174 /* Function to calculate Ecal Energy in Cone (given in cm) */
175 inline double ecalEnergyInCone(const GlobalPoint center,
176  double radius,
177  const EcalRecHitCollection ecalCol,
178  const CaloGeometry* geo) {
179  double eECALcone = 0;
180  std::vector<int> usedHitsEcal;
181  usedHitsEcal.clear();
182 
183  for (std::vector<EcalRecHit>::const_iterator ehit = ecalCol.begin(); ehit != ecalCol.end(); ehit++) {
184  //This is a precaution for the case when hitCollection contains duplicats.
185  bool hitIsUsed = false;
186  int hitHashedIndex = -10000;
187  if (ehit->id().subdetId() == EcalBarrel) {
188  EBDetId did(ehit->id());
189  hitHashedIndex = did.hashedIndex();
190  }
191  if (ehit->id().subdetId() == EcalEndcap) {
192  EEDetId did(ehit->id());
193  hitHashedIndex = did.hashedIndex();
194  }
195  for (uint32_t i = 0; i < usedHitsEcal.size(); i++) {
196  if (usedHitsEcal[i] == hitHashedIndex)
197  hitIsUsed = true;
198  }
199  if (hitIsUsed)
200  continue;
201  usedHitsEcal.push_back(hitHashedIndex);
202  // -----------------------------------------------
203 
204  const GlobalPoint& pos = geo->getPosition((*ehit).detid());
205 
206  if (getDistInPlaneSimple(center, pos) < radius) {
207  eECALcone += ehit->energy();
208  }
209  }
210  return eECALcone;
211 }
212 
213 /* This is another version of ecalEnergy calculation using the getDistInPlaneTrackDir() */
214 inline double ecalEnergyInCone(const GlobalVector trackMom,
215  const GlobalPoint center,
216  double radius,
217  const EcalRecHitCollection ecalCol,
218  const CaloGeometry* geo) {
219  double eECALcone = 0;
220  std::vector<int> usedHitsEcal;
221  usedHitsEcal.clear();
222  for (std::vector<EcalRecHit>::const_iterator ehit = ecalCol.begin(); ehit != ecalCol.end(); ehit++) {
223  //This is a precaution for the case when hitCollection contains duplicats.
224  bool hitIsUsed = false;
225  int hitHashedIndex = -10000;
226  if (ehit->id().subdetId() == EcalBarrel) {
227  EBDetId did(ehit->id());
228  hitHashedIndex = did.hashedIndex();
229  }
230  if (ehit->id().subdetId() == EcalEndcap) {
231  EEDetId did(ehit->id());
232  hitHashedIndex = did.hashedIndex();
233  }
234  for (uint32_t i = 0; i < usedHitsEcal.size(); i++) {
235  if (usedHitsEcal[i] == hitHashedIndex)
236  hitIsUsed = true;
237  }
238  if (hitIsUsed)
239  continue;
240  usedHitsEcal.push_back(hitHashedIndex);
241  // -----------------------------------------------
242 
243  const GlobalPoint& pos = geo->getPosition((*ehit).detid());
244 
245  if (getDistInPlaneTrackDir(center, trackMom, pos) < radius) {
246  eECALcone += ehit->energy();
247  }
248  }
249  return eECALcone;
250 }
251 
252 #endif
Vector3DBase
Definition: Vector3DBase.h:8
EBDetId::hashedIndex
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
mps_fire.i
i
Definition: mps_fire.py:428
MaxHit_struct
Definition: CommonUsefulStuff.h:42
getDistInPlane
double getDistInPlane(const GlobalVector trackDirection, const GlobalPoint caloPoint, const GlobalPoint rechitPoint, double coneHeight)
Definition: CommonUsefulStuff.h:105
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
MaxHit_struct::hitenergy
float hitenergy
Definition: CommonUsefulStuff.h:46
CaloGeometry::getPosition
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
EBDetId
Definition: EBDetId.h:17
pos
Definition: PixelAliasList.h:18
EBDetId.h
EEDetId.h
MaxHit_struct::dr
float dr
Definition: CommonUsefulStuff.h:47
edm::SortedCollection< EcalRecHit >
EDAnalyzer.h
EcalBarrel
Definition: EcalSubdetector.h:10
EcalRecHitCollections.h
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
MakerMacros.h
CaloGeometry
Definition: CaloGeometry.h:21
Vector3DBase::unit
Vector3DBase unit() const
Definition: Vector3DBase.h:54
MaxHit_struct::iphihitm
int iphihitm
Definition: CommonUsefulStuff.h:43
edm::SortedCollection::begin
const_iterator begin() const
Definition: SortedCollection.h:262
Service.h
MaxHit_struct::ietahitm
int ietahitm
Definition: CommonUsefulStuff.h:44
Point3DBase< float, GlobalTag >
EEDetId
Definition: EEDetId.h:14
CaloGeometryRecord.h
EcalEndcap
Definition: EcalSubdetector.h:10
DeltaR.h
TFileService.h
CaloSubdetectorGeometry.h
getDistInPlaneSimple
double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
Definition: CommonUsefulStuff.h:52
Event.h
edm::SortedCollection::end
const_iterator end() const
Definition: SortedCollection.h:267
HcalDetId.h
EEDetId::hashedIndex
int hashedIndex() const
Definition: EEDetId.h:183
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
Vector3DBase::dot
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
CaloCellGeometry.h
DetId.h
Frameworkfwd.h
CaloGeometry.h
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
GlobalVector.h
EventSetup.h
MaxHit_struct::MaxHit_struct
MaxHit_struct()
Definition: CommonUsefulStuff.h:49
HcalRecHitCollections.h
Exception.h
ecalEnergyInCone
double ecalEnergyInCone(const GlobalPoint center, double radius, const EcalRecHitCollection ecalCol, const CaloGeometry *geo)
Definition: CommonUsefulStuff.h:175
ParameterSet.h
GlobalPoint.h
MaxHit_struct::posMax
GlobalPoint posMax
Definition: CommonUsefulStuff.h:48
MaxHit_struct::depthhit
int depthhit
Definition: CommonUsefulStuff.h:45
getDistInPlaneTrackDir
double getDistInPlaneTrackDir(const GlobalPoint caloPoint, const GlobalVector caloVector, const GlobalPoint rechitPoint)
Definition: CommonUsefulStuff.h:79