CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/Calibration/HcalCalibAlgos/interface/CommonUsefulStuff.h

Go to the documentation of this file.
00001 #ifndef CommonUsefulStuff_h
00002 #define CommonUsefulStuff_h
00003 
00004 // $Id: CommonUsefulStuff.h,v 1.3 2010/03/16 23:25:48 andrey Exp $
00005 
00006 #include "FWCore/Framework/interface/Frameworkfwd.h"
00007 #include "FWCore/Framework/interface/EDAnalyzer.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00015 #include "DataFormats/DetId/interface/DetId.h"
00016 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00017 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00018 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00019 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00020 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00021 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00022 
00023 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00024 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00025 
00026 #include "FWCore/Utilities/interface/Exception.h"
00027 
00028 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00029 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00030 
00031 #include "FWCore/ServiceRegistry/interface/Service.h"
00032 #include "CommonTools/UtilAlgos/interface/DeltaR.h"
00033 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00034 
00035 #include "TROOT.h"
00036 #include "TFile.h"
00037 #include "TTree.h"
00038 #include "TH1F.h"
00039 #include "TString.h"
00040 
00041 /* getDist functions by Jim:
00042 /uscms/home/jhirsch/IsoTracks_314/CMSSW_3_1_4/src/JetMETCorrections/IsolatedParticles/interface/FindCaloHit.icc
00043 */
00044 
00045 
00046 struct MaxHit_struct
00047 {
00048   int iphihitm;
00049   int ietahitm;
00050   int depthhit;
00051   float hitenergy;
00052   float dr;
00053   GlobalPoint posMax;
00054   MaxHit_struct():iphihitm(0),ietahitm(0),
00055                   depthhit(0),hitenergy(-100),dr(0){}
00056 };
00057 
00058 
00059 inline double getDistInPlaneSimple(const GlobalPoint caloPoint,
00060                             const GlobalPoint rechitPoint)
00061 {
00062 
00063   // Simplified version of getDistInPlane
00064   // Assume track direction is origin -> point of hcal intersection
00065 
00066   const GlobalVector caloIntersectVector(caloPoint.x(),
00067                                          caloPoint.y(),
00068                                          caloPoint.z());
00069 
00070   const GlobalVector caloIntersectUnitVector = caloIntersectVector.unit();
00071 
00072   const GlobalVector rechitVector(rechitPoint.x(),
00073                                   rechitPoint.y(),
00074                                   rechitPoint.z());
00075 
00076   const GlobalVector rechitUnitVector = rechitVector.unit();
00077   double dotprod = caloIntersectUnitVector.dot(rechitUnitVector);
00078   double rechitdist = caloIntersectVector.mag()/dotprod;
00079 
00080 
00081   const GlobalVector effectiveRechitVector = rechitdist*rechitUnitVector;
00082   const GlobalPoint effectiveRechitPoint(effectiveRechitVector.x(),
00083                                          effectiveRechitVector.y(),
00084                                          effectiveRechitVector.z());
00085 
00086 
00087   GlobalVector distance_vector = effectiveRechitPoint-caloPoint;
00088 
00089   if (dotprod > 0.)
00090   {
00091     return distance_vector.mag();
00092   }
00093   else
00094   {
00095     return 999999.;
00096 
00097   }
00098 
00099 }
00100 
00101 inline double getDistInPlaneTrackDir(const GlobalPoint  caloPoint,
00102                               const GlobalVector caloVector,
00103                               const GlobalPoint  rechitPoint)
00104 {
00105 
00106   // Simplified version of getDistInPlane : no cone "within" Hcal, but
00107   // don't assume track direction is origin -> point of hcal
00108   // intersection.
00109   const GlobalVector caloIntersectVector(caloPoint.x(),
00110                                          caloPoint.y(),
00111                                          caloPoint.z()); //p
00112 
00113   const GlobalVector caloUnitVector = caloVector.unit();
00114   const GlobalVector rechitVector(rechitPoint.x(),
00115                                   rechitPoint.y(),
00116                                   rechitPoint.z());
00117   const GlobalVector rechitUnitVector = rechitVector.unit();
00118   double dotprod_denominator = caloUnitVector.dot(rechitUnitVector);
00119   double dotprod_numerator   = caloUnitVector.dot(caloIntersectVector);
00120   double rechitdist = dotprod_numerator/dotprod_denominator;
00121   const GlobalVector effectiveRechitVector = rechitdist*rechitUnitVector;
00122   const GlobalPoint effectiveRechitPoint(effectiveRechitVector.x(),
00123                                          effectiveRechitVector.y(),
00124                                          effectiveRechitVector.z());
00125   GlobalVector distance_vector = effectiveRechitPoint-caloPoint;
00126   if (dotprod_denominator > 0. && dotprod_numerator > 0.)
00127   {
00128 
00129     return distance_vector.mag();
00130   }
00131   else
00132   {
00133     return 999999.;
00134   }
00135 }
00136 
00137 
00138 
00139 inline double getDistInPlane(const GlobalVector trackDirection,
00140                       const GlobalPoint caloPoint,
00141                       const GlobalPoint rechitPoint,
00142                       double coneHeight)
00143 {
00144 
00145   // The iso track candidate hits the Calo (Ecal or Hcal) at "caloPoint"
00146   // with direction "trackDirection".
00147 
00148   // "rechitPoint" is the position of the rechit.  We only care about
00149   // the direction of the rechit.
00150 
00151   // Consider the rechitPoint as characterized by angles theta and phi
00152   // wrt the origin which points at the calo cell of the rechit.  In
00153   // some sense the distance along the line theta/phi is arbitrary. A
00154   // simplified choice might be to put the rechit at the surface of
00155   // the Hcal.  Our approach is to see whether/where this line
00156   // intersects a cone of height "coneHeight" with vertex at caloPoint
00157   // aligned with trackDirection.
00158  // To that end, this function returns the distance between the
00159   // center of the base of the cone and the intersection of the rechit
00160   // line and the plane that contains the base of the cone.  This
00161   // distance is compared with the radius of the cone outside this
00162   // function.
00163 
00164 
00165   // Make vector of length cone height along track direction
00166   const GlobalVector heightVector = trackDirection*coneHeight;
00167 
00168   // Make vector from origin to point where iso track intersects the
00169   // calorimeter.
00170   const GlobalVector caloIntersectVector(caloPoint.x(),
00171                                          caloPoint.y(),
00172                                          caloPoint.z());
00173 
00174   // Make vector from origin to point in center of base of cone
00175   const GlobalVector coneBaseVector = heightVector+caloIntersectVector;
00176 
00177 // Make point in center of base of cone
00178   const GlobalPoint coneBasePoint(coneBaseVector.x(),
00179                                   coneBaseVector.y(),
00180                                   coneBaseVector.z());
00181 
00182   // Make unit vector pointing at rechit.
00183   const GlobalVector rechitVector(rechitPoint.x(),
00184                                   rechitPoint.y(),
00185                                   rechitPoint.z());
00186   const GlobalVector rechitDirection = rechitVector.unit();
00187 
00188   // Find distance "r" along "rechit line" (with angles theta2 and
00189   // phi2) where line intersects plane defined by base of cone.
00190 
00191   // Definition plane of that contains base of cone:
00192   // trackDirection.x() (x - coneBaseVector.x()) +
00193   // trackDirection.y() (y - coneBaseVector.y()) +
00194   // trackDirection.z() (z - coneBaseVector.z()) = 0
00195 
00196   // Point P_{rh} where rechit line intersects plane:
00197   // (rechitdist sin(theta2) cos(phi2),
00198   //  rechitdist sin(theta2) cos(phi2),
00199   //  rechitdist cos(theta2))
00200 
00201   // Substitute P_{rh} into equation for plane and solve for rechitdist.
00202   // rechitDist turns out to be the ratio of dot products:
00203 
00204   double rechitdist = trackDirection.dot(coneBaseVector)/trackDirection.dot(rechitDirection);
00205 
00206   // Now find distance between point at center of cone base and point
00207   // where the "rechit line" intersects the plane defined by the base
00208   // of the cone; i.e. the effectiveRecHitPoint.
00209   const GlobalVector effectiveRechitVector = rechitdist*rechitDirection;
00210   const GlobalPoint effectiveRechitPoint(effectiveRechitVector.x(),
00211                                          effectiveRechitVector.y(),
00212                                          effectiveRechitVector.z());
00213 
00214 
00215   GlobalVector distance_vector = effectiveRechitPoint-coneBasePoint;
00216   return distance_vector.mag();
00217 }
00218 
00219 
00220 /*  Function to calculate Ecal Energy in Cone (given in cm) */
00221 inline double ecalEnergyInCone(const GlobalPoint center, double radius, const EcalRecHitCollection ecalCol, const CaloGeometry *geo)
00222 {
00223   double eECALcone = 0;
00224   std::vector<int> usedHitsEcal; 
00225   usedHitsEcal.clear();
00226   
00227   for (std::vector<EcalRecHit>::const_iterator ehit=ecalCol.begin(); ehit!=ecalCol.end(); ehit++)
00228     {
00229       //This is a precaution for the case when hitCollection contains duplicats.
00230       bool hitIsUsed=false;
00231       int hitHashedIndex=-10000;
00232       if (ehit->id().subdetId()==EcalBarrel)
00233         {
00234           EBDetId did(ehit->id());
00235           hitHashedIndex=did.hashedIndex();
00236         }
00237       if (ehit->id().subdetId()==EcalEndcap)
00238         {
00239           EEDetId did(ehit->id());
00240           hitHashedIndex=did.hashedIndex();
00241         }
00242       for (uint32_t i=0; i<usedHitsEcal.size(); i++)
00243         {
00244           if (usedHitsEcal[i]==hitHashedIndex) hitIsUsed=true;
00245         }
00246       if (hitIsUsed) continue;
00247       usedHitsEcal.push_back(hitHashedIndex);
00248       // -----------------------------------------------
00249       
00250       GlobalPoint pos = geo->getPosition((*ehit).detid());
00251   
00252       if (getDistInPlaneSimple(center,pos) < radius)
00253         {
00254           eECALcone += ehit->energy();
00255         }
00256     }
00257   return eECALcone;
00258 }
00259 
00260 /*  This is another version of ecalEnergy calculation using the getDistInPlaneTrackDir()  */
00261 inline double ecalEnergyInCone(const GlobalVector trackMom, const GlobalPoint center, double radius, const EcalRecHitCollection ecalCol, const CaloGeometry *geo)
00262 {
00263   
00264   double eECALcone = 0;
00265   std::vector<int> usedHitsEcal; 
00266   usedHitsEcal.clear();
00267   for (std::vector<EcalRecHit>::const_iterator ehit=ecalCol.begin(); ehit!=ecalCol.end(); ehit++)
00268     {
00269       //This is a precaution for the case when hitCollection contains duplicats.
00270       bool hitIsUsed=false;
00271       int hitHashedIndex=-10000;
00272       if (ehit->id().subdetId()==EcalBarrel)
00273         {
00274           EBDetId did(ehit->id());
00275           hitHashedIndex=did.hashedIndex();
00276         }
00277       if (ehit->id().subdetId()==EcalEndcap)
00278         {
00279           EEDetId did(ehit->id());
00280           hitHashedIndex=did.hashedIndex();
00281         }
00282       for (uint32_t i=0; i<usedHitsEcal.size(); i++)
00283         {
00284           if (usedHitsEcal[i]==hitHashedIndex) hitIsUsed=true;
00285         }
00286       if (hitIsUsed) continue;
00287       usedHitsEcal.push_back(hitHashedIndex);
00288       // -----------------------------------------------
00289 
00290       GlobalPoint pos = geo->getPosition((*ehit).detid());
00291   
00292       if (getDistInPlaneTrackDir(center, trackMom ,pos) < radius)
00293         {
00294           eECALcone += ehit->energy();
00295         }
00296       }
00297   return eECALcone;
00298 }
00299 
00300 
00301 
00302 #endif
00303