CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EnergyLossPlain.cc
Go to the documentation of this file.
2 
6 
9 
12 
14 
17 
19 
22 
24 
26 
29 
31 
36 
38 #include <fstream>
39 
41 
42 // m_e c^2
43 #define Tmax 1.
44 
45 using namespace std;
46 
47 bool EnergyLossPlain::isFirst = true;
48 float EnergyLossPlain::optimalWeight[61][61];
49 
50 /*****************************************************************************/
52  (const TrackerGeometry* theTracker_,
53  double pixelToStripMultiplier_,
54  double pixelToStripExponent_ ) : theTracker(theTracker_),
55  pixelToStripMultiplier(pixelToStripMultiplier_),
56  pixelToStripExponent (pixelToStripExponent_ )
57 {
58  // Load data
59  if(isFirst == true)
60  {
61  loadOptimalWeights();
62  isFirst = false;
63  }
64 }
65 
66 /*****************************************************************************/
68 {
69 }
70 
71 
72 /*****************************************************************************/
74 {
76  fileInPath("QCDAnalysis/ChargedHadronSpectra/data/energyWeights.dat");
77  ifstream inFile(fileInPath.fullPath().c_str());
78 
79  while(inFile.eof() == false)
80  {
81  int i; float w; int n;
82  inFile >> i;
83  inFile >> w;
84  inFile >> n;
85 
87  }
88 
89  inFile.close();
90 
91  LogTrace("MinBiasTracking")
92  << " [EnergyLossEstimator] optimal weights loaded";
93 }
94 
95 /*****************************************************************************/
96 double EnergyLossPlain::average(std::vector<pair<double,double> >& values)
97 {
98  int num = values.size();
99  double sum[2] = {0.,0.};
100 
101  for(int i = 0; i < num; i++)
102  {
103  sum[1] += values[i].first;
104  sum[0] += 1;
105  }
106 
107  return sum[1] / sum[0];
108 }
109 
110 /*****************************************************************************/
111 double EnergyLossPlain::logTruncate(std::vector<pair<double,double> >& values_)
112 {
113  std::vector<double> values;
114  for(std::vector<pair<double,double> >::iterator
115  v = values_.begin(); v!= values_.end(); v++)
116  values.push_back((*v).first);
117 
118  sort(values.begin(), values.end());
119 
120  int num = values.size();
121  double sum[2] = {0.,1.};
122 
123  for(int i = 0; i < (num+1)/2; i++)
124  {
125  double weight = 1.;
126 
127  if(num%2 == 1 && i == (num-1)/2) weight = 1./2;
128 
129  sum[1] *= pow(values[i], weight);
130  sum[0] += weight;
131 
132  }
133 
134  return pow(sum[1], 1./sum[0]);
135 }
136 
137 /*****************************************************************************/
138 double EnergyLossPlain::truncate(std::vector<pair<double,double> >& values_)
139 {
140  std::vector<double> values;
141  for(std::vector<pair<double,double> >::iterator
142  v = values_.begin(); v!= values_.end(); v++)
143  values.push_back((*v).first);
144 
145  sort(values.begin(), values.end());
146 
147  int num = values.size();
148  double sum[2] = {0.,0.};
149 
150  for(int i = 0; i < (num+1)/2; i++)
151  {
152  double weight = 1.;
153 
154  if(num%2 == 1 && i == (num-1)/2) weight = 1./2;
155 
156  sum[1] += weight * values[i];
157  sum[0] += weight;
158  }
159 
160  return sum[1] / sum[0];
161 }
162 
163 /*****************************************************************************/
164 double EnergyLossPlain::optimal(std::vector<pair<double,double> >& values_)
165 {
166  std::vector<double> values;
167  for(std::vector<pair<double,double> >::iterator
168  v = values_.begin(); v!= values_.end(); v++)
169  values.push_back((*v).first);
170 
171  int num = values.size();
172  sort(values.begin(), values.end());
173 
174  // First guess
175  double sum = 0.;
176  for(int i = 0; i < num; i++)
177  {
178  double weight = optimalWeight[num][i];
179  sum += weight * values[i];
180  }
181 
182  // Correct every deposit with log(path)
183  for(int i = 0; i < num; i++)
184  values_[i].first -= 0.178*log(values_[i].second/0.03) * 0.38 * sum;
185 
186  // Sort again
187  values.clear();
188  for(std::vector<pair<double,double> >::iterator
189  v = values_.begin(); v!= values_.end(); v++)
190  values.push_back((*v).first);
191  sort(values.begin(), values.end());
192 
193  // Corrected estimate
194  sum = 0.;
195  for(int i = 0; i < num; i++)
196  {
197  double weight = optimalWeight[num][i];
198  sum += weight * values[i];
199  }
200 
201  return sum;
202 }
203 
204 /*****************************************************************************/
205 double EnergyLossPlain::expected(double Delta1, double Delta2)
206 {
207  return log(Delta2/Delta1) / (1/Delta1 - 1/Delta2);
208 // return 1 + (Delta2*log(Delta1) - Delta1*log(Delta2)) / (Delta2 - Delta1);
209 }
210 
211 /*****************************************************************************/
213  (LocalVector ldir, const SiPixelRecHit* recHit, std::vector<pair<double,double> >& values)
214 {
215  DetId id = recHit->geographicalId();
216  const PixelGeomDetUnit* pixelDet =
217  dynamic_cast<const PixelGeomDetUnit*> (theTracker->idToDet(id));
218 
219  // Check if cluster does not touch the edge
220  if(recHit->cluster()->minPixelRow() == 0 ||
221  recHit->cluster()->maxPixelRow() ==
222  pixelDet->specificTopology().nrows() - 1 ||
223  recHit->cluster()->minPixelCol() == 0 ||
224  recHit->cluster()->maxPixelCol() ==
225  pixelDet->specificTopology().ncolumns() - 1)
226  return;
227 
228  double MeVperElec = 3.61e-6;
229 
230  // Collect adc
231  double elec = recHit->cluster()->charge();
232  double Delta = elec * MeVperElec;
233 
234  // Length
235  double x = ldir.mag()/fabsf(ldir.z()) *
236  pixelDet->surface().bounds().thickness();
237 
238  double pix = Delta/x;
239 
240  // MeV/cm, only if not low deposit
241  if(pix > 1.5 * 0.795)
242  values.push_back(std::pair<double,double>(pix, x));
243 }
244 
245 /*****************************************************************************/
247  (LocalVector ldir, const SiStripRecHit2D* recHit, std::vector<pair<double,double> >& values)
248 {
249  DetId id = recHit->geographicalId();
250  const StripGeomDetUnit* stripDet =
251  dynamic_cast<const StripGeomDetUnit*> (theTracker->idToDet(id));
252 
253  // Check if cluster does not touch the edge
254  if(recHit->cluster()->firstStrip() == 0 ||
255  int(recHit->cluster()->firstStrip() +
256  recHit->cluster()->amplitudes().size()) ==
257  stripDet->specificTopology().nstrips())
258  return;
259 
260  double MeVperAdc = 313 * 3.61e-6;
261 
262  // Collect adc
263  double Delta = 0;
264 // for(std::vector<uint16_t>::const_iterator
265  for(std::vector<uint8_t>::const_iterator
266  i = (recHit->cluster()->amplitudes()).begin();
267  i!= (recHit->cluster()->amplitudes()).end(); i++)
268  {
269  // MeV
270  double delta = (*i + 0.5) * MeVperAdc;
271 
272  if(*i >= 254)
273  {
274  double delta1, delta2;
275  if(*i == 254) { delta1 = 254 * MeVperAdc; delta2 = 512 * MeVperAdc; }
276  else { delta1 = 512 * MeVperAdc; delta2 = Tmax; }
277 
278  delta = expected(delta1,delta2);
279  }
280 
281  Delta += delta;
282  }
283 
284  // Length
285  double x = ldir.mag()/fabsf(ldir.z()) *
286  stripDet->surface().bounds().thickness();
287 
288  double str = Delta/x;
289 
290  // MeV/cm, only if not low deposit
291  if(str > 1.5)
292  values.push_back(std::pair<double,double>(str, x));
293 }
294 
295 /*****************************************************************************/
297  (const Trajectory* trajectory, std::vector<pair<int,double> >& arithmeticMean,
298  std::vector<pair<int,double> >& truncatedMean)
299 {
300  // (dE/dx, dx)
301  std::vector<pair<double,double> > vpix,vstr;
302 
303  for(std::vector<TrajectoryMeasurement>::const_iterator
304  meas = trajectory->measurements().begin();
305  meas!= trajectory->measurements().end(); meas++)
306  {
307  const TrackingRecHit* recHit = meas->recHit()->hit();
308  DetId id = recHit->geographicalId();
309 
310  if(recHit->isValid())
311  {
312  LocalVector ldir = meas->updatedState().localDirection();
313 
314  if(theTracker->idToDet(id)->subDetector() ==
316  theTracker->idToDet(id)->subDetector() ==
318  {
319  // Pixel
320  const SiPixelRecHit* pixelRecHit =
321  dynamic_cast<const SiPixelRecHit *>(recHit);
322 
323  if(pixelRecHit != 0)
324  process(ldir,pixelRecHit , vpix);
325  }
326  else
327  {
328  // Strip
329  const SiStripMatchedRecHit2D* stripMatchedRecHit =
330  dynamic_cast<const SiStripMatchedRecHit2D *>(recHit);
331  const ProjectedSiStripRecHit2D* stripProjectedRecHit =
332  dynamic_cast<const ProjectedSiStripRecHit2D *>(recHit);
333  const SiStripRecHit2D* stripRecHit =
334  dynamic_cast<const SiStripRecHit2D *>(recHit);
335 
336  std::pair<double,double> v;
337 
338  if(stripMatchedRecHit != 0)
339  {
340  auto m = stripMatchedRecHit->monoHit();
341  auto s = stripMatchedRecHit->stereoHit();
342  process(ldir,&m, vstr);
343  process(ldir,&s, vstr);
344  }
345 
346  if(stripProjectedRecHit != 0)
347  process(ldir,&(stripProjectedRecHit->originalHit()), vstr);
348 
349  if(stripRecHit != 0)
350  process(ldir,stripRecHit, vstr);
351  }
352  }
353  }
354 
355  // Transform
356  std::vector<pair<double,double> > vall;
357 
358  for(unsigned int i = 0; i < vpix.size(); i++)
359  {
360  float a = 0.795;
361  float s = 10.1;
362 
363  std::pair<double,double> str(vpix[i]);
364 
365  double y = str.first / a / s;
366  if(y > 0.9999) y = 0.9999;
367  if(y <-0.9999) y = -0.9999;
368 
369  str.first = s * atanh(y);
370 
371  vall.push_back(str);
372  }
373 
374  for(unsigned int i = 0; i < vstr.size(); i++) vall.push_back(vstr[i]);
375 
376  // Arithmetic mean
377  arithmeticMean.push_back(std::pair<int,double>(vpix.size(), average(vpix)));
378  arithmeticMean.push_back(std::pair<int,double>(vstr.size(), average(vstr)));
379  arithmeticMean.push_back(std::pair<int,double>(vall.size(), average(vall)));
380 
381  // Wighted mean
382  truncatedMean.push_back(std::pair<int,double>(vpix.size(), optimal(vpix)));
383  truncatedMean.push_back(std::pair<int,double>(vstr.size(), optimal(vstr)));
384  truncatedMean.push_back(std::pair<int,double>(vall.size(), optimal(vall)));
385 
386  return vall.size();
387 }
388 
dbl * delta
Definition: mlp_gen.cc:36
virtual int nstrips() const =0
int i
Definition: DBlmapReader.cc:9
bool isFirst(HepMC::GenParticle *x)
double truncate(std::vector< std::pair< double, double > > &values)
#define Tmax
virtual int ncolumns() const =0
double optimal(std::vector< std::pair< double, double > > &values)
virtual int nrows() const =0
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
void process(LocalVector ldir, const SiPixelRecHit *recHit, std::vector< std::pair< double, double > > &values)
int estimate(const Trajectory *trajectory, std::vector< std::pair< int, double > > &arithmeticMean, std::vector< std::pair< int, double > > &truncatedMean)
Vector2DBase< float, GlobalTag > Global2DVector
virtual float thickness() const =0
U second(std::pair< T, U > const &p)
EnergyLossPlain(const TrackerGeometry *theTracker_, double pixelToStripMultiplier, double pixelToStripExponent)
DataContainer const & measurements() const
Definition: Trajectory.h:203
T mag() const
Definition: PV3DBase.h:66
T z() const
Definition: PV3DBase.h:63
#define end
Definition: vmac.h:38
bool first
Definition: L1TdeRCT.cc:94
static bool isFirst
#define LogTrace(id)
Definition: DetId.h:20
const Bounds & bounds() const
Definition: BoundSurface.h:89
double average(std::vector< std::pair< double, double > > &values)
double expected(double Delta1, double Delta2)
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
long long int num
Definition: procUtils.cc:71
bool isValid() const
static float optimalWeight[61][61]
int average
Definition: PDRates.py:137
#define begin
Definition: vmac.h:31
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double a
Definition: hdecay.h:121
double logTruncate(std::vector< std::pair< double, double > > &values)
DetId geographicalId() const
std::string fullPath() const
Definition: FileInPath.cc:171
Definition: DDAxes.h:10
tuple process
Definition: LaserDQM_cfg.py:3
const SiStripRecHit2D & originalHit() const
mathSSE::Vec4< T > v
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
T w() const
Pixel Reconstructed Hit.