CMS 3D CMS Logo

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