CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalUncalibRecHitRatioMethodAlgo.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
2 #define RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
3 
11 #include "Math/SVector.h"
12 #include "Math/SMatrix.h"
15 #include <vector>
16 
17 template < class C > class EcalUncalibRecHitRatioMethodAlgo {
18  public:
19  struct Ratio {
20  unsigned int index;
21  unsigned int step;
22  double value;
23  double error;
24  };
25  struct Tmax {
26  unsigned int index;
27  unsigned int step;
28  double value;
29  double error;
30  double amplitude;
31  double chi2;
32  };
34  double amplitudeMax;
35  double timeMax;
36  double timeError;
37  double chi2;
38  };
39 
41  virtual EcalUncalibratedRecHit makeRecHit(const C & dataFrame,
42  const EcalSampleMask & sampleMask,
43  const double *pedestals,
44  const double* pedestalRMSes,
45  const double *gainRatios,
46  std::vector < double >&timeFitParameters,
47  std::vector < double >&amplitudeFitParameters,
48  std::pair < double, double >&timeFitLimits);
49 
50  // more function to be able to compute
51  // amplitude and time separately
52  void init( const C &dataFrame, const EcalSampleMask &sampleMask, const double * pedestals, const double * pedestalRMSes, const double * gainRatios );
53  void computeTime(std::vector < double >&timeFitParameters, std::pair < double, double >&timeFitLimits, std::vector< double > &amplitudeFitParameters);
54  void computeAmplitude( std::vector< double > &amplitudeFitParameters );
55  CalculatedRecHit getCalculatedRecHit() { return calculatedRechit_; };
56  bool fixMGPAslew( const C &dataFrame );
57 
58  protected:
59 
62  std::vector < double > amplitudes_;
63  std::vector < double > amplitudeErrors_;
64  std::vector < Ratio > ratios_;
65  std::vector < Tmax > times_;
66  std::vector < Tmax > timesAB_;
67 
68  double pedestal_;
69  int num_;
70  double ampMaxError_;
71 
72  CalculatedRecHit calculatedRechit_;
73 };
74 
75 template <class C>
76 void EcalUncalibRecHitRatioMethodAlgo<C>::init( const C &dataFrame, const EcalSampleMask &sampleMask,
77  const double * pedestals, const double * pedestalRMSes, const double * gainRatios )
78 {
79  sampleMask_ = sampleMask;
80  theDetId_ = DetId(dataFrame.id().rawId());
81 
82  calculatedRechit_.timeMax = 5;
83  calculatedRechit_.amplitudeMax = 0;
84  calculatedRechit_.timeError = -999;
85  amplitudes_.clear();
86  amplitudes_.reserve(C::MAXSAMPLES);
87  amplitudeErrors_.clear();
88  amplitudeErrors_.reserve(C::MAXSAMPLES);
89  ratios_.clear();
90  ratios_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
91  times_.clear();
92  times_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
93  timesAB_.clear();
94  timesAB_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
95 
96  // to obtain gain 12 pedestal:
97  // -> if it's in gain 12, use first sample
98  // --> average it with second sample if in gain 12 and 3-sigma-noise compatible (better LF noise cancellation)
99  // -> else use pedestal from database
100  pedestal_ = 0;
101  num_ = 0;
102  if (dataFrame.sample(0).gainId() == 1 &&
103  sampleMask_.useSample(0, theDetId_ )
104  ) {
105  pedestal_ += double (dataFrame.sample(0).adc());
106  num_++;
107  }
108  if (num_!=0 &&
109  dataFrame.sample(1).gainId() == 1 &&
110  sampleMask_.useSample(1, theDetId_) &&
111  fabs(dataFrame.sample(1).adc()-dataFrame.sample(0).adc())<3*pedestalRMSes[0]) {
112  pedestal_ += double (dataFrame.sample(1).adc());
113  num_++;
114  }
115  if (num_ != 0)
116  pedestal_ /= num_;
117  else
118  pedestal_ = pedestals[0];
119 
120  // fill vector of amplitudes, pedestal subtracted and vector
121  // of amplitude uncertainties Also, find the uncertainty of a
122  // sample with max amplitude. We will use it later.
123 
124  ampMaxError_ = 0;
125  double ampMaxValue = -1000;
126 
127  // ped-subtracted and gain-renormalized samples. It is VERY
128  // IMPORTANT to have samples one clock apart which means to
129  // have vector size equal to MAXSAMPLES
130  double sample;
131  double sampleError;
132  int GainId;
133  for (int iSample = 0; iSample < C::MAXSAMPLES; iSample++) {
134 
135 
136 
137  GainId = dataFrame.sample(iSample).gainId();
138 
139 
140  // only use normally samples which are desired; if sample not to be used
141  // inflate error so won't generate ratio considered for the measurement
142  if (!sampleMask_.useSample(iSample, theDetId_ ) ) {
143  sample = 1e-9;
144  sampleError = 1e+9;
145  }
146  else if (GainId == 1) {
147  sample = double (dataFrame.sample(iSample).adc() - pedestal_);
148  sampleError = pedestalRMSes[0];
149  }
150  else if (GainId == 2 || GainId == 3){
151  sample = (double (dataFrame.sample(iSample).adc() - pedestals[GainId - 1])) *gainRatios[GainId - 1];
152  sampleError = pedestalRMSes[GainId-1]*gainRatios[GainId-1];
153  }
154  else {
155  sample = 1e-9; // GainId=0 case falls here, from saturation
156  sampleError = 1e+9; // inflate error so won't generate ratio considered for the measurement
157  }
158 
159 
160  if(sampleError>0){
161  amplitudes_.push_back(sample);
162  amplitudeErrors_.push_back(sampleError);
163  if(ampMaxValue < sample){
164  ampMaxValue = sample;
165  ampMaxError_ = sampleError;
166  }
167  }else{
168  // inflate error for useless samples
169  amplitudes_.push_back(sample);
170  amplitudeErrors_.push_back(1e+9);
171  }
172  }
173 }
174 
175 template <class C>
177 {
178 
179  // This fuction finds sample(s) preceeding gain switching and
180  // inflates errors on this sample, therefore, making this sample
181  // invisible for Ratio Method. Only gain switching DOWN is
182  // considered Only gainID=1,2,3 are considered. In case of the
183  // saturation (gainID=0), we keep "distorted" sample because it is
184  // the only chance to make time measurement; the qualilty of it will
185  // be bad anyway.
186 
187  bool result = false;
188 
189  int GainIdPrev;
190  int GainIdNext;
191  for (int iSample = 1; iSample < C::MAXSAMPLES; iSample++) {
192 
193  // only use samples which are desired
194  if (!sampleMask_.useSample(iSample, theDetId_) ) continue;
195 
196  GainIdPrev = dataFrame.sample(iSample-1).gainId();
197  GainIdNext = dataFrame.sample(iSample).gainId();
198  if( GainIdPrev>=1 && GainIdPrev<=3 &&
199  GainIdNext>=1 && GainIdNext<=3 &&
200  GainIdPrev<GainIdNext ){
201  amplitudes_[iSample-1]=1e-9;
202  amplitudeErrors_[iSample-1]=1e+9;
203  result = true;
204  }
205  }
206  return result;
207 
208 }
209 
210 template<class C>
211 void EcalUncalibRecHitRatioMethodAlgo<C>::computeTime(std::vector < double >&timeFitParameters,
212  std::pair < double, double >&timeFitLimits, std::vector < double >&amplitudeFitParameters)
213 {
215  // //
216  // RATIO METHOD FOR TIME STARTS HERE //
217  // //
219  double ampMaxAlphaBeta = 0;
220  double tMaxAlphaBeta = 5;
221  double tMaxErrorAlphaBeta = 999;
222  double tMaxRatio = 5;
223  double tMaxErrorRatio = 999;
224 
225  double sumAA = 0;
226  double sumA = 0;
227  double sum1 = 0;
228  double sum0 = 0;
229  double sumAf = 0;
230  double sumff = 0;
231  double NullChi2 = 0;
232 
233  // null hypothesis = no pulse, pedestal only
234  for(unsigned int i = 0; i < amplitudes_.size(); i++){
235  double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
236  sum0 += 1;
237  sum1 += 1/err2;
238  sumA += amplitudes_[i]/err2;
239  sumAA += amplitudes_[i]*amplitudes_[i]/err2;
240  }
241  if(sum0>0){
242  NullChi2 = (sumAA - sumA*sumA/sum1)/sum0;
243  }else{
244  // not enough samples to reconstruct the pulse
245  return;
246  }
247 
248  // Make all possible Ratio's based on any pair of samples i and j
249  // (j>i) with positive amplitudes_
250  //
251  // Ratio[k] = Amp[i]/Amp[j]
252  // where Amp[i] is pedestal subtracted ADC value in a time sample [i]
253  //
254  double alphabeta = amplitudeFitParameters[0]*amplitudeFitParameters[1];
255  double alpha = amplitudeFitParameters[0];
256  double beta = amplitudeFitParameters[1];
257 
258  for(unsigned int i = 0; i < amplitudes_.size()-1; i++){
259  for(unsigned int j = i+1; j < amplitudes_.size(); j++){
260 
261  if(amplitudes_[i]>1 && amplitudes_[j]>1){
262 
263  // ratio
264  double Rtmp = amplitudes_[i]/amplitudes_[j];
265 
266  // error^2 due to stat fluctuations of time samples
267  // (uncorrelated for both samples)
268 
269  double err1 = Rtmp*Rtmp*( (amplitudeErrors_[i]*amplitudeErrors_[i]/(amplitudes_[i]*amplitudes_[i])) + (amplitudeErrors_[j]*amplitudeErrors_[j]/(amplitudes_[j]*amplitudes_[j])) );
270 
271  // error due to fluctuations of pedestal (common to both samples)
272  double stat;
273  if(num_>0) stat = num_; // num presampeles used to compute pedestal
274  else stat = 1; // pedestal from db
275  double err2 = amplitudeErrors_[j]*(amplitudes_[i]-amplitudes_[j])/(amplitudes_[j]*amplitudes_[j])/sqrt(stat);
276 
277  //error due to integer round-down. It is relevant to low
278  //amplitudes_ in gainID=1 and negligible otherwise.
279  double err3 = 0.289/amplitudes_[j];
280 
281  double totalError = sqrt(err1 + err2*err2 +err3*err3);
282 
283 
284  // don't include useless ratios
285  if(totalError < 1.0
286  && Rtmp>0.001
287  && Rtmp<exp(double(j-i)/beta)-0.001
288  ){
289  Ratio currentRatio = { i, (j-i), Rtmp, totalError };
290  ratios_.push_back(currentRatio);
291  }
292 
293  }
294 
295  }
296  }
297 
298  // No useful ratios, return zero amplitude and no time measurement
299  if(!ratios_.size() >0)
300  return;
301 
302  // make a vector of Tmax measurements that correspond to each ratio
303  // and based on Alpha-Beta parameterization of the pulse shape
304 
305  for(unsigned int i = 0; i < ratios_.size(); i++){
306 
307  double stepOverBeta = double(ratios_[i].step)/beta;
308  double offset = double(ratios_[i].index) + alphabeta;
309 
310  double Rmin = ratios_[i].value - ratios_[i].error;
311  if(Rmin<0.001) Rmin=0.001;
312 
313  double Rmax = ratios_[i].value + ratios_[i].error;
314  double RLimit = exp(stepOverBeta)-0.001;
315  if( Rmax > RLimit ) Rmax = RLimit;
316 
317  double time1 = offset - ratios_[i].step/(exp((stepOverBeta-log(Rmin))/alpha)-1.0);
318  double time2 = offset - ratios_[i].step/(exp((stepOverBeta-log(Rmax))/alpha)-1.0);
319 
320  // this is the time measurement based on the ratios[i]
321  double tmax = 0.5 * (time1 + time2);
322  double tmaxerr = 0.5 * sqrt( (time1 - time2)*(time1 - time2) );
323 
324  // calculate chi2
325  sumAf = 0;
326  sumff = 0;
327  for(unsigned int it = 0; it < amplitudes_.size(); it++){
328  double err2 = amplitudeErrors_[it]*amplitudeErrors_[it];
329  double offset = (double(it) - tmax)/alphabeta;
330  double term1 = 1.0 + offset;
331  if(term1>1e-6){
332  double f = exp( alpha*(log(1.0+offset) - offset) );
333  sumAf += amplitudes_[it]*f/err2;
334  sumff += f*f/err2;
335  }
336  }
337 
338  double chi2 = sumAA;
339  double amp = 0;
340  if( sumff > 0 ){
341  chi2 = sumAA - sumAf*sumAf/sumff;
342  amp = sumAf/sumff;
343  }
344  chi2 /= sum0;
345 
346  // choose reasonable measurements. One might argue what is
347  // reasonable and what is not.
348  if(chi2 > 0 && tmaxerr > 0 && tmax > 0){
349  Tmax currentTmax={ ratios_[i].index, ratios_[i].step, tmax, tmaxerr, amp, chi2 };
350  timesAB_.push_back(currentTmax);
351  }
352  }
353 
354  // no reasonable time measurements!
355  if( !(timesAB_.size()> 0))
356  return;
357 
358  // find minimum chi2
359  double chi2min = 1.0e+9;
360  //double timeMinimum = 5;
361  //double errorMinimum = 999;
362  for(unsigned int i = 0; i < timesAB_.size(); i++){
363  if( timesAB_[i].chi2 <= chi2min ){
364  chi2min = timesAB_[i].chi2;
365  //timeMinimum = timesAB_[i].value;
366  //errorMinimum = timesAB_[i].error;
367  }
368  }
369 
370  // make a weighted average of tmax measurements with "small" chi2
371  // (within 1 sigma of statistical uncertainty :-)
372  double chi2Limit = chi2min + 1.0;
373  double time_max = 0;
374  double time_wgt = 0;
375  for(unsigned int i = 0; i < timesAB_.size(); i++){
376  if( timesAB_[i].chi2 < chi2Limit ){
377  double inverseSigmaSquared = 1.0/(timesAB_[i].error*timesAB_[i].error);
378  time_wgt += inverseSigmaSquared;
379  time_max += timesAB_[i].value*inverseSigmaSquared;
380  }
381  }
382 
383  tMaxAlphaBeta = time_max/time_wgt;
384  tMaxErrorAlphaBeta = 1.0/sqrt(time_wgt);
385 
386  // find amplitude and chi2
387  sumAf = 0;
388  sumff = 0;
389  for(unsigned int i = 0; i < amplitudes_.size(); i++){
390  double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
391  double offset = (double(i) - tMaxAlphaBeta)/alphabeta;
392  double term1 = 1.0 + offset;
393  if(term1>1e-6){
394  double f = exp( alpha*(log(1.0+offset) - offset) );
395  sumAf += amplitudes_[i]*f/err2;
396  sumff += f*f/err2;
397  }
398  }
399 
400  if( sumff > 0 ){
401  ampMaxAlphaBeta = sumAf/sumff;
402  double chi2AlphaBeta = (sumAA - sumAf*sumAf/sumff)/sum0;
403  if(chi2AlphaBeta > NullChi2){
404  // null hypothesis is better
405  return;
406  }
407 
408  }else{
409  // no visible pulse here
410  return;
411  }
412 
413  // if we got to this point, we have a reconstructied Tmax
414  // using RatioAlphaBeta Method. To summarize:
415  //
416  // tMaxAlphaBeta - Tmax value
417  // tMaxErrorAlphaBeta - error on Tmax, but I would not trust it
418  // ampMaxAlphaBeta - amplitude of the pulse
419  // ampMaxError_ - uncertainty of the time sample with max amplitude
420  //
421 
422 
423 
424  // Do Ratio's Method with "large" pulses
425  if( ampMaxAlphaBeta/ampMaxError_ > 5.0 ){
426 
427  // make a vector of Tmax measurements that correspond to each
428  // ratio. Each measurement have it's value and the error
429 
430  double time_max = 0;
431  double time_wgt = 0;
432 
433 
434  for (unsigned int i = 0; i < ratios_.size(); i++) {
435 
436  if(ratios_[i].step == 1
437  && ratios_[i].value >= timeFitLimits.first
438  && ratios_[i].value <= timeFitLimits.second
439  ){
440 
441  double time_max_i = ratios_[i].index;
442 
443  // calculate polynomial for Tmax
444 
445  double u = timeFitParameters[timeFitParameters.size() - 1];
446  for (int k = timeFitParameters.size() - 2; k >= 0; k--) {
447  u = u * ratios_[i].value + timeFitParameters[k];
448  }
449 
450  // calculate derivative for Tmax error
451  double du =
452  (timeFitParameters.size() -
453  1) * timeFitParameters[timeFitParameters.size() - 1];
454  for (int k = timeFitParameters.size() - 2; k >= 1; k--) {
455  du = du * ratios_[i].value + k * timeFitParameters[k];
456  }
457 
458 
459  // running sums for weighted average
460  double errorsquared =
461  ratios_[i].error * ratios_[i].error * du * du;
462  if (errorsquared > 0) {
463 
464  time_max += (time_max_i - u) / errorsquared;
465  time_wgt += 1.0 / errorsquared;
466  Tmax currentTmax =
467  { ratios_[i].index, 1, (time_max_i - u),
468  sqrt(errorsquared),0,1 };
469  times_.push_back(currentTmax);
470 
471  }
472  }
473  }
474 
475 
476  // calculate weighted average of all Tmax measurements
477  if (time_wgt > 0) {
478  tMaxRatio = time_max/time_wgt;
479  tMaxErrorRatio = 1.0/sqrt(time_wgt);
480 
481  // combine RatioAlphaBeta and Ratio Methods
482 
483  if( ampMaxAlphaBeta/ampMaxError_ > 10.0 ){
484 
485  // use pure Ratio Method
486  calculatedRechit_.timeMax = tMaxRatio;
487  calculatedRechit_.timeError = tMaxErrorRatio;
488 
489  }else{
490 
491  // combine two methods
492  calculatedRechit_.timeMax = ( tMaxAlphaBeta*(10.0-ampMaxAlphaBeta/ampMaxError_) + tMaxRatio*(ampMaxAlphaBeta/ampMaxError_ - 5.0) )/5.0;
493  calculatedRechit_.timeError = ( tMaxErrorAlphaBeta*(10.0-ampMaxAlphaBeta/ampMaxError_) + tMaxErrorRatio*(ampMaxAlphaBeta/ampMaxError_ - 5.0) )/5.0;
494 
495  }
496 
497  }else{
498 
499  // use RatioAlphaBeta Method
500  calculatedRechit_.timeMax = tMaxAlphaBeta;
501  calculatedRechit_.timeError = tMaxErrorAlphaBeta;
502 
503  }
504 
505  }else{
506 
507  // use RatioAlphaBeta Method
508  calculatedRechit_.timeMax = tMaxAlphaBeta;
509  calculatedRechit_.timeError = tMaxErrorAlphaBeta;
510 
511  }
512 }
513 
514 template<class C>
515 void EcalUncalibRecHitRatioMethodAlgo<C>::computeAmplitude( std::vector< double > &amplitudeFitParameters )
516 {
518  // //
519  // CALCULATE AMPLITUDE //
520  // //
522 
523 
524  double alpha = amplitudeFitParameters[0];
525  double beta = amplitudeFitParameters[1];
526 
527  // calculate pedestal, again
528 
529  double pedestalLimit = calculatedRechit_.timeMax - (alpha * beta) - 1.0;
530  double sumA = 0;
531  double sumF = 0;
532  double sumAF = 0;
533  double sumFF = 0;
534  double sum1 = 0;
535  for (unsigned int i = 0; i < amplitudes_.size(); i++) {
536  double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
537  double f = 0;
538  double termOne = 1 + (i - calculatedRechit_.timeMax) / (alpha * beta);
539  if (termOne > 1.e-5) f = exp(alpha * log(termOne)) * exp(-(i - calculatedRechit_.timeMax) / beta);
540 
541  // apply range of interesting samples
542 
543  if ( (i < pedestalLimit)
544  || (f > 0.6 && i <= calculatedRechit_.timeMax)
545  || (f > 0.4 && i >= calculatedRechit_.timeMax)) {
546  sum1 += 1/err2;
547  sumA += amplitudes_[i]/err2;
548  sumF += f/err2;
549  sumAF += f*amplitudes_[i]/err2;
550  sumFF += f*f/err2;
551  }
552  }
553 
554  calculatedRechit_.amplitudeMax = 0;
555  if(sum1 > 0){
556  double denom = sumFF*sum1 - sumF*sumF;
557  if(fabs(denom)>1.0e-20){
558  calculatedRechit_.amplitudeMax = (sumAF*sum1 - sumA*sumF)/denom;
559  }
560  }
561 }
562 
563 
564 
565 template < class C > EcalUncalibratedRecHit
567  const EcalSampleMask & sampleMask,
568  const double *pedestals,
569  const double *pedestalRMSes,
570  const double *gainRatios,
571  std::vector < double >&timeFitParameters,
572  std::vector < double >&amplitudeFitParameters,
573  std::pair < double, double >&timeFitLimits)
574 {
575 
576  init( dataFrame, sampleMask, pedestals, pedestalRMSes, gainRatios );
577  computeTime( timeFitParameters, timeFitLimits, amplitudeFitParameters );
578  computeAmplitude( amplitudeFitParameters );
579 
580  // 1st parameters is id
581  //
582  // 2nd parameters is amplitude. It is calculated by this method.
583  //
584  // 3rd parameter is pedestal. It is not calculated. This method
585  // relies on input parameters for pedestals and gain ratio. Return
586  // zero.
587  //
588  // 4th parameter is jitter which is a bad choice to call Tmax. It is
589  // calculated by this method (in 25 nsec clock units)
590  //
591  // GF subtract 5 so that jitter==0 corresponds to synchronous hit
592  //
593  //
594  // 5th parameter is chi2. It is possible to calculate chi2 for
595  // Tmax. It is possible to calculate chi2 for Amax. However, these
596  // values are not very useful without TmaxErr and AmaxErr. This
597  // method can return one value for chi2 but there are 4 different
598  // parameters that have useful information about the quality of Amax
599  // ans Tmax. For now we can return TmaxErr. The quality of Tmax and
600  // Amax can be judged from the magnitude of TmaxErr
601 
602  return EcalUncalibratedRecHit(dataFrame.id(),
603  calculatedRechit_.amplitudeMax,
604  pedestal_,
605  calculatedRechit_.timeMax - 5,
606  calculatedRechit_.timeError);
607 }
608 #endif
const double beta
int i
Definition: DBlmapReader.cc:9
float alpha
Definition: AMPTWrapper.h:95
void computeAmplitude(std::vector< double > &amplitudeFitParameters)
list step
Definition: launcher.py:15
void computeTime(std::vector< double > &timeFitParameters, std::pair< double, double > &timeFitLimits, std::vector< double > &amplitudeFitParameters)
int init
Definition: HydjetWrapper.h:63
T sqrt(T t)
Definition: SSEVec.h:46
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
double f[11][100]
virtual EcalUncalibratedRecHit makeRecHit(const C &dataFrame, const EcalSampleMask &sampleMask, const double *pedestals, const double *pedestalRMSes, const double *gainRatios, std::vector< double > &timeFitParameters, std::vector< double > &amplitudeFitParameters, std::pair< double, double > &timeFitLimits)
unsigned int offset(bool)
int k[5][pyjets_maxn]
static const double tmax[3]
Definition: DetId.h:20
const double Rmax[kNumberCalorimeter]
void init(const C &dataFrame, const EcalSampleMask &sampleMask, const double *pedestals, const double *pedestalRMSes, const double *gainRatios)
const double Rmin[kNumberCalorimeter]