1 #ifndef RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
2 #define RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
11 #include "Math/SVector.h"
12 #include "Math/SMatrix.h"
27 inline float fast_expf(
float x) {
return unsafe_expf<6>(x); }
28 inline float fast_logf(
float x) {
return unsafe_logf<7>(x); }
57 const double *pedestals,
58 const double *pedestalRMSes,
59 const double *gainRatios,
60 std::vector<double> &timeFitParameters,
61 std::vector<double> &litudeFitParameters,
62 std::pair<double, double> &timeFitLimits);
66 void init(
const C &dataFrame,
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> &litudeFitParameters);
98 const double *pedestals,
99 const double *pedestalRMSes,
100 const double *gainRatios) {
101 sampleMask_ = sampleMask;
102 theDetId_ =
DetId(dataFrame.id().rawId());
104 calculatedRechit_.timeMax = 5;
105 calculatedRechit_.amplitudeMax = 0;
106 calculatedRechit_.timeError = -999;
115 if (dataFrame.sample(0).gainId() == 1 && sampleMask_.useSample(0, theDetId_)) {
116 pedestal_ += double(dataFrame.sample(0).adc());
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());
127 pedestal_ = pedestals[0];
134 double ampMaxValue = -1000;
142 for (
int iSample = 0; iSample < C::MAXSAMPLES; iSample++) {
143 GainId = dataFrame.sample(iSample).gainId();
148 if (!sampleMask_.useSample(iSample, theDetId_)) {
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];
165 useless_[iSample] = (sampleError <= 0) | bad;
166 amplitudes_[iSample] =
sample;
168 amplitudeErrors_[iSample] = useless_[iSample] ? 1
e+9 : sampleError;
169 amplitudeIE2_[iSample] = useless_[iSample] ? 0 : 1 / (amplitudeErrors_[iSample] * amplitudeErrors_[iSample]);
170 if (sampleError > 0) {
171 if (ampMaxValue < sample) {
173 ampMaxError_ = sampleError;
192 for (
int iSample = 1; iSample < C::MAXSAMPLES; iSample++) {
194 if (!sampleMask_.useSample(iSample, theDetId_))
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] = 1
e+9;
202 amplitudeIE2_[iSample - 1] = 0;
203 useless_[iSample - 1] =
true;
212 std::pair<double, double> &timeFitLimits,
213 std::vector<double> &litudeFitParameters) {
219 double ampMaxAlphaBeta = 0;
220 double tMaxAlphaBeta = 5;
221 double tMaxErrorAlphaBeta = 999;
222 double tMaxRatio = 5;
223 double tMaxErrorRatio = 999;
234 for (
unsigned int i = 0;
i < amplitudes_.size();
i++) {
237 double inverr2 = amplitudeIE2_[
i];
240 sumA += amplitudes_[
i] * inverr2;
241 sumAA += amplitudes_[
i] * (amplitudes_[
i] * inverr2);
244 NullChi2 = (sumAA - sumA * sumA / sum1) / sum0;
256 double alphabeta = amplitudeFitParameters[0] * amplitudeFitParameters[1];
257 double invalphabeta = 1. / alphabeta;
258 double alpha = amplitudeFitParameters[0];
259 double beta = amplitudeFitParameters[1];
261 Ratio ratios_[C::MAXSAMPLES * (C::MAXSAMPLES - 1) / 2];
262 unsigned int ratios_size = 0;
264 double Rlim[amplitudes_.size()];
265 for (
unsigned int k = 1;
k != amplitudes_.size(); ++
k)
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]);
275 for (
unsigned int i = 0;
i < amplitudes_.size() - 1;
i++) {
278 for (
unsigned int j = i + 1;
j < amplitudes_.size();
j++) {
282 if (amplitudes_[i] > 1 && amplitudes_[j] > 1) {
284 double Rtmp = amplitudes_[
i] / amplitudes_[
j];
289 double err1 = Rtmp * Rtmp * (relErr2[
i] + relErr2[
j]);
298 amplitudeErrors_[
j] * (amplitudes_[
i] - amplitudes_[
j]) * (invampl[j] * invampl[j]);
303 double err3 = (0.289 * 0.289) * (invampl[
j] * invampl[
j]);
305 double totalError =
std::sqrt(err1 + err2 + err3);
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;
317 if (0 == ratios_size)
321 Tmax timesAB_[C::MAXSAMPLES * (C::MAXSAMPLES - 1) / 2];
322 unsigned int timesAB_size = 0;
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;
331 double Rmin = ratios_[
i].value - ratios_[
i].error;
335 double Rmax = ratios_[
i].value + ratios_[
i].error;
336 double RLimit = Rlim[ratios_[
i].step];
346 double tmax = 0.5 * (time1 + time2);
347 double tmaxerr = 0.5 *
std::sqrt((time1 - time2) * (time1 - time2));
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;
358 double inverr2 = amplitudeIE2_[it];
360 double term1 = 1.0 + loffset;
363 sumAf += amplitudes_[it] * (f * inverr2);
364 sumff += f * (f * inverr2);
370 chi2 = sumAA - sumAf * (sumAf / sumff);
371 amp = (sumAf / sumff);
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;
384 if (0 == timesAB_size)
388 double chi2min = 1.0e+9;
391 for (
unsigned int i = 0; i < timesAB_size; i++) {
392 if (timesAB_[i].
chi2 < chi2min) {
393 chi2min = timesAB_[
i].
chi2;
401 double chi2Limit = chi2min + 1.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;
412 tMaxAlphaBeta = time_max / time_wgt;
413 tMaxErrorAlphaBeta = 1.0 /
sqrt(time_wgt);
418 for (
unsigned int i = 0; i < amplitudes_.size(); i++) {
421 double inverr2 = amplitudeIE2_[
i];
422 double offset = (double(i) - tMaxAlphaBeta) * invalphabeta;
423 double term1 = 1.0 +
offset;
426 sumAf += amplitudes_[
i] * (f * inverr2);
427 sumff += f * (f * inverr2);
432 ampMaxAlphaBeta = sumAf / sumff;
433 double chi2AlphaBeta = (sumAA - sumAf * sumAf / sumff) / sum0;
434 if (chi2AlphaBeta > NullChi2) {
454 if (ampMaxAlphaBeta / ampMaxError_ > 5.0) {
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;
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];
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];
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;
493 tMaxRatio = time_max / time_wgt;
494 tMaxErrorRatio = 1.0 /
sqrt(time_wgt);
498 if (ampMaxAlphaBeta / ampMaxError_ > 10.0) {
500 calculatedRechit_.timeMax = tMaxRatio;
501 calculatedRechit_.timeError = tMaxErrorRatio;
505 calculatedRechit_.timeMax = (tMaxAlphaBeta * (10.0 - ampMaxAlphaBeta / ampMaxError_) +
506 tMaxRatio * (ampMaxAlphaBeta / ampMaxError_ - 5.0)) /
508 calculatedRechit_.timeError = (tMaxErrorAlphaBeta * (10.0 - ampMaxAlphaBeta / ampMaxError_) +
509 tMaxErrorRatio * (ampMaxAlphaBeta / ampMaxError_ - 5.0)) /
515 calculatedRechit_.timeMax = tMaxAlphaBeta;
516 calculatedRechit_.timeError = tMaxErrorAlphaBeta;
521 calculatedRechit_.timeMax = tMaxAlphaBeta;
522 calculatedRechit_.timeError = tMaxErrorAlphaBeta;
528 calculatedRechit_.amplitudeMax = computeAmplitudeImpl(amplitudeFitParameters, 1., 1.);
541 double alpha = amplitudeFitParameters[0];
542 double beta = amplitudeFitParameters[1];
546 double pedestalLimit = calculatedRechit_.timeMax - (alpha *
beta) - 1.0;
552 for (
unsigned int i = 0;
i < amplitudes_.size();
i++) {
555 double inverr2 = amplitudeIE2_[
i];
557 double termOne = 1 + (i - calculatedRechit_.timeMax) / (alpha * beta);
563 if ((i < pedestalLimit) || (f > 0.6 * corr6 && i <= calculatedRechit_.timeMax) ||
564 (f > 0.4 * corr4 && i >= calculatedRechit_.timeMax)) {
566 sumA += amplitudes_[
i] * inverr2;
568 sumAF += (f * inverr2) * amplitudes_[i];
569 sumFF += f * (f * inverr2);
573 double amplitudeMax = 0;
575 double denom = sumFF * sum1 - sumF * sumF;
577 amplitudeMax = (sumAF * sum1 - sumA * sumF) / denom;
586 const double *pedestals,
587 const double *pedestalRMSes,
588 const double *gainRatios,
589 std::vector<double> &timeFitParameters,
590 std::vector<double> &litudeFitParameters,
591 std::pair<double, double> &timeFitLimits) {
592 init(dataFrame, sampleMask, pedestals, pedestalRMSes, gainRatios);
593 computeTime(timeFitParameters, timeFitLimits, amplitudeFitParameters);
594 computeAmplitude(amplitudeFitParameters);
619 calculatedRechit_.amplitudeMax,
621 calculatedRechit_.timeMax - 5,
622 calculatedRechit_.timeError);
void computeAmplitude(std::vector< double > &litudeFitParameters)
std::array< double, C::MAXSAMPLES > amplitudeErrors_
CalculatedRecHit getCalculatedRecHit()
void computeTime(std::vector< double > &timeFitParameters, std::pair< double, double > &timeFitLimits, std::vector< double > &litudeFitParameters)
CalculatedRecHit calculatedRechit_
EcalSampleMask sampleMask_
Abs< T >::type abs(const T &t)
EcalUncalibratedRecHit makeRecHit(const C &dataFrame, const EcalSampleMask &sampleMask, const double *pedestals, const double *pedestalRMSes, const double *gainRatios, std::vector< double > &timeFitParameters, std::vector< double > &litudeFitParameters, std::pair< double, double > &timeFitLimits)
std::array< double, C::MAXSAMPLES > amplitudeIE2_
std::array< double, C::MAXSAMPLES > amplitudes_
double computeAmplitudeImpl(std::vector< double > &litudeFitParameters, double, double)
static const double tmax[3]
bool fixMGPAslew(const C &dataFrame)
std::array< bool, C::MAXSAMPLES > useless_
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]