CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/QCDAnalysis/ChargedHadronSpectra/src/EnergyLossPlain.cc

Go to the documentation of this file.
00001 #include "QCDAnalysis/ChargedHadronSpectra/interface/EnergyLossPlain.h"
00002 
00003 #include "DataFormats/GeometryVector/interface/GlobalTag.h"
00004 #include "DataFormats/GeometryVector/interface/Vector2DBase.h"
00005 typedef Vector2DBase<float,GlobalTag> Global2DVector;
00006 
00007 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00008 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00009 
00010 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00011 #include "Geometry/TrackerGeometryBuilder/interface/RectangularPixelTopology.h"
00012 
00013 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00014 
00015 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
00016 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00017 
00018 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00019 
00020 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00021 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00022 
00023 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00024 
00025 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00026 
00027 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/ClusterShape.h"
00028 #include "RecoPixelVertexing/PixelLowPtUtilities/interface/ClusterData.h"
00029 
00030 #include "Geometry/CommonDetUnit/interface/GeomDetEnumerators.h"
00031 
00032 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00033 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00034 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00035 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00036 
00037 #include "FWCore/ParameterSet/interface/FileInPath.h"
00038 #include <fstream>
00039 
00040 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00041 
00042 // m_e c^2
00043 #define Tmax 1.
00044 
00045 using namespace std;
00046 
00047 bool  EnergyLossPlain::isFirst = true;
00048 float EnergyLossPlain::optimalWeight[61][61];
00049 
00050 /*****************************************************************************/
00051 EnergyLossPlain::EnergyLossPlain
00052   (const TrackerGeometry* theTracker_,
00053    double pixelToStripMultiplier_,
00054    double pixelToStripExponent_  ) : theTracker(theTracker_),
00055       pixelToStripMultiplier(pixelToStripMultiplier_),
00056       pixelToStripExponent  (pixelToStripExponent_  )
00057 {
00058   // Load data
00059   if(isFirst == true)
00060   {
00061     loadOptimalWeights();
00062     isFirst = false;
00063   }
00064 }
00065 
00066 /*****************************************************************************/
00067 EnergyLossPlain::~EnergyLossPlain()
00068 {
00069 }
00070 
00071 
00072 /*****************************************************************************/
00073 void EnergyLossPlain::loadOptimalWeights()
00074 {
00075   edm::FileInPath
00076     fileInPath("QCDAnalysis/ChargedHadronSpectra/data/energyWeights.dat");
00077   ifstream inFile(fileInPath.fullPath().c_str());
00078 
00079   while(inFile.eof() == false)
00080   {
00081     int i; float w; int n;
00082     inFile >> i;
00083     inFile >> w;
00084     inFile >> n;
00085 
00086     EnergyLossPlain::optimalWeight[n][i] = w;
00087   }
00088 
00089   inFile.close();
00090 
00091   LogTrace("MinBiasTracking")
00092     << " [EnergyLossEstimator] optimal weights loaded";
00093 }
00094 
00095 /*****************************************************************************/
00096 double EnergyLossPlain::average(std::vector<pair<double,double> >& values)
00097 {
00098   int num = values.size();
00099   double sum[2] = {0.,0.};
00100 
00101   for(int i = 0; i < num; i++)
00102   {
00103     sum[1] += values[i].first;
00104     sum[0] += 1;
00105   }
00106 
00107   return sum[1] / sum[0];
00108 }
00109 
00110 /*****************************************************************************/
00111 double EnergyLossPlain::logTruncate(std::vector<pair<double,double> >& values_)
00112 {
00113   std::vector<double> values;
00114   for(std::vector<pair<double,double> >::iterator
00115       v = values_.begin(); v!= values_.end(); v++)
00116     values.push_back((*v).first);
00117 
00118   sort(values.begin(), values.end());
00119 
00120   int num = values.size();
00121   double sum[2] = {0.,1.};
00122 
00123   for(int i = 0; i < (num+1)/2; i++)
00124   { 
00125     double weight = 1.;
00126 
00127     if(num%2 == 1 && i == (num-1)/2) weight = 1./2;
00128   
00129     sum[1] *= pow(values[i], weight); 
00130     sum[0] += weight;
00131 
00132   }
00133   
00134   return pow(sum[1], 1./sum[0]);
00135 } 
00136 
00137 /*****************************************************************************/
00138 double EnergyLossPlain::truncate(std::vector<pair<double,double> >& values_)
00139 {
00140   std::vector<double> values;
00141   for(std::vector<pair<double,double> >::iterator 
00142       v = values_.begin(); v!= values_.end(); v++)
00143     values.push_back((*v).first);
00144 
00145   sort(values.begin(), values.end());
00146 
00147   int num = values.size();
00148   double sum[2] = {0.,0.};
00149 
00150   for(int i = 0; i < (num+1)/2; i++)
00151   {
00152     double weight = 1.;
00153 
00154     if(num%2 == 1 && i == (num-1)/2) weight = 1./2;
00155 
00156     sum[1] += weight * values[i];
00157     sum[0] += weight;
00158   }
00159 
00160   return sum[1] / sum[0];
00161 }
00162 
00163 /*****************************************************************************/
00164 double EnergyLossPlain::optimal(std::vector<pair<double,double> >& values_)
00165 {
00166   std::vector<double> values;
00167   for(std::vector<pair<double,double> >::iterator
00168       v = values_.begin(); v!= values_.end(); v++)
00169     values.push_back((*v).first);
00170 
00171   int num = values.size();
00172   sort(values.begin(), values.end());
00173 
00174   // First guess
00175   double sum = 0.;
00176   for(int i = 0; i < num; i++)
00177   {
00178     double weight = optimalWeight[num][i];
00179     sum += weight * values[i];
00180   }
00181 
00182   // Correct every deposit with log(path)
00183   for(int i = 0; i < num; i++)
00184     values_[i].first -= 0.178*log(values_[i].second/0.03) * 0.38 * sum;
00185 
00186   // Sort again 
00187   values.clear();
00188   for(std::vector<pair<double,double> >::iterator
00189       v = values_.begin(); v!= values_.end(); v++)
00190     values.push_back((*v).first);
00191   sort(values.begin(), values.end()); 
00192 
00193   // Corrected estimate
00194   sum = 0.;
00195   for(int i = 0; i < num; i++)
00196   {
00197     double weight = optimalWeight[num][i];
00198     sum += weight * values[i];
00199   }
00200 
00201   return sum;
00202 }
00203 
00204 /*****************************************************************************/
00205 double EnergyLossPlain::expected(double Delta1, double Delta2)
00206 {
00207   return log(Delta2/Delta1) / (1/Delta1 - 1/Delta2);
00208 //    return 1 + (Delta2*log(Delta1) - Delta1*log(Delta2)) / (Delta2 - Delta1);
00209 }
00210 
00211 /*****************************************************************************/
00212 void EnergyLossPlain::process
00213   (LocalVector ldir, const SiPixelRecHit* recHit, std::vector<pair<double,double> >& values)
00214 {
00215   DetId id = recHit->geographicalId();
00216   const PixelGeomDetUnit* pixelDet =
00217     dynamic_cast<const PixelGeomDetUnit*> (theTracker->idToDet(id));
00218 
00219   // Check if cluster does not touch the edge
00220   if(recHit->cluster()->minPixelRow() == 0 ||
00221      recHit->cluster()->maxPixelRow() ==
00222        pixelDet->specificTopology().nrows() - 1 ||
00223      recHit->cluster()->minPixelCol() == 0 ||
00224      recHit->cluster()->maxPixelCol() == 
00225        pixelDet->specificTopology().ncolumns() - 1)
00226     return;
00227 
00228   double MeVperElec = 3.61e-6;
00229 
00230   // Collect adc
00231   double elec = recHit->cluster()->charge();
00232   double Delta = elec * MeVperElec;
00233 
00234   // Length
00235   double x = ldir.mag()/fabsf(ldir.z()) *
00236              pixelDet->surface().bounds().thickness();
00237 
00238   double pix = Delta/x;
00239 
00240   // MeV/cm, only if not low deposit
00241   if(pix > 1.5 * 0.795)
00242     values.push_back(std::pair<double,double>(pix, x)); 
00243 }
00244 
00245 /*****************************************************************************/
00246 void EnergyLossPlain::process
00247   (LocalVector ldir, const SiStripRecHit2D* recHit, std::vector<pair<double,double> >& values)
00248 {
00249   DetId id = recHit->geographicalId();
00250   const StripGeomDetUnit* stripDet =
00251     dynamic_cast<const StripGeomDetUnit*> (theTracker->idToDet(id)); 
00252 
00253   // Check if cluster does not touch the edge
00254   if(recHit->cluster()->firstStrip() == 0 ||
00255      int(recHit->cluster()->firstStrip() +
00256          recHit->cluster()->amplitudes().size()) ==
00257        stripDet->specificTopology().nstrips())
00258     return;
00259 
00260   double MeVperAdc = 313 * 3.61e-6;
00261 
00262   // Collect adc
00263   double Delta = 0;
00264 //  for(std::vector<uint16_t>::const_iterator
00265   for(std::vector<uint8_t>::const_iterator
00266     i = (recHit->cluster()->amplitudes()).begin();
00267     i!= (recHit->cluster()->amplitudes()).end(); i++)
00268   {
00269     // MeV
00270     double delta = (*i + 0.5) * MeVperAdc;
00271 
00272     if(*i >= 254)
00273     {
00274       double delta1, delta2; 
00275       if(*i == 254) { delta1 = 254 * MeVperAdc; delta2 = 512 * MeVperAdc; }
00276                else { delta1 = 512 * MeVperAdc; delta2 = Tmax; }
00277 
00278       delta = expected(delta1,delta2);
00279     }
00280 
00281     Delta += delta;
00282   }
00283 
00284   // Length
00285   double x = ldir.mag()/fabsf(ldir.z()) *
00286              stripDet->surface().bounds().thickness();
00287 
00288   double str = Delta/x;
00289 
00290   // MeV/cm, only if not low deposit
00291   if(str > 1.5)
00292     values.push_back(std::pair<double,double>(str, x));
00293 }
00294 
00295 /*****************************************************************************/
00296 int EnergyLossPlain::estimate
00297   (const Trajectory* trajectory, std::vector<pair<int,double> >& arithmeticMean,
00298                                  std::vector<pair<int,double> >& truncatedMean)
00299 {
00300   // (dE/dx, dx)
00301   std::vector<pair<double,double> > vpix,vstr; 
00302 
00303   for(std::vector<TrajectoryMeasurement>::const_iterator
00304         meas = trajectory->measurements().begin();
00305         meas!= trajectory->measurements().end(); meas++)
00306   {
00307     const TrackingRecHit* recHit = meas->recHit()->hit();
00308     DetId id = recHit->geographicalId();
00309 
00310     if(recHit->isValid())
00311     {
00312       LocalVector ldir = meas->updatedState().localDirection();
00313   
00314       if(theTracker->idToDet(id)->subDetector() ==
00315            GeomDetEnumerators::PixelBarrel ||
00316          theTracker->idToDet(id)->subDetector() ==
00317            GeomDetEnumerators::PixelEndcap)
00318       { 
00319         // Pixel
00320         const SiPixelRecHit* pixelRecHit =
00321           dynamic_cast<const SiPixelRecHit *>(recHit);
00322   
00323         if(pixelRecHit != 0)
00324           process(ldir,pixelRecHit , vpix);
00325       }
00326       else
00327       {
00328         // Strip
00329         const SiStripMatchedRecHit2D* stripMatchedRecHit =
00330           dynamic_cast<const SiStripMatchedRecHit2D *>(recHit);
00331         const ProjectedSiStripRecHit2D* stripProjectedRecHit =
00332           dynamic_cast<const ProjectedSiStripRecHit2D *>(recHit);
00333         const SiStripRecHit2D* stripRecHit =
00334           dynamic_cast<const SiStripRecHit2D *>(recHit);
00335 
00336         std::pair<double,double> v;
00337   
00338         if(stripMatchedRecHit != 0)
00339         {
00340           auto m = stripMatchedRecHit->monoHit();
00341           auto s = stripMatchedRecHit->stereoHit();
00342           process(ldir,&m, vstr);
00343           process(ldir,&s, vstr);
00344         } 
00345 
00346         if(stripProjectedRecHit != 0)
00347           process(ldir,&(stripProjectedRecHit->originalHit()), vstr);
00348   
00349         if(stripRecHit != 0)
00350           process(ldir,stripRecHit, vstr);
00351       }
00352     }
00353   }
00354 
00355   // Transform
00356   std::vector<pair<double,double> > vall;
00357 
00358   for(unsigned int i = 0; i < vpix.size(); i++)
00359   {
00360     float a = 0.795;
00361     float s = 10.1;
00362 
00363     std::pair<double,double> str(vpix[i]);
00364 
00365     double y = str.first / a / s;
00366     if(y > 0.9999) y =  0.9999;
00367     if(y <-0.9999) y = -0.9999;
00368 
00369     str.first = s * atanh(y);
00370 
00371     vall.push_back(str);
00372   }
00373 
00374   for(unsigned int i = 0; i < vstr.size(); i++) vall.push_back(vstr[i]);
00375 
00376   // Arithmetic mean
00377   arithmeticMean.push_back(std::pair<int,double>(vpix.size(), average(vpix)));
00378   arithmeticMean.push_back(std::pair<int,double>(vstr.size(), average(vstr)));
00379   arithmeticMean.push_back(std::pair<int,double>(vall.size(), average(vall)));
00380 
00381   // Wighted mean
00382   truncatedMean.push_back(std::pair<int,double>(vpix.size(), optimal(vpix)));
00383   truncatedMean.push_back(std::pair<int,double>(vstr.size(), optimal(vstr)));
00384   truncatedMean.push_back(std::pair<int,double>(vall.size(), optimal(vall)));
00385 
00386   return vall.size();
00387 }
00388