CMS 3D CMS Logo

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