00001 #ifndef RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
00002 #define RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitRatioMethodAlgo_HH
00003
00014 #include "Math/SVector.h"
00015 #include "Math/SMatrix.h"
00016 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalUncalibRecHitRecAbsAlgo.h"
00017 #include "CondFormats/EcalObjects/interface/EcalWeightSet.h"
00018 #include <vector>
00019
00020 template < class C > class EcalUncalibRecHitRatioMethodAlgo {
00021 public:
00022 struct Ratio {
00023 int index;
00024 int step;
00025 double value;
00026 double error;
00027 };
00028 struct Tmax {
00029 int index;
00030 int step;
00031 double value;
00032 double error;
00033 double amplitude;
00034 double chi2;
00035 };
00036 struct CalculatedRecHit {
00037 double amplitudeMax;
00038 double timeMax;
00039 double timeError;
00040 double chi2;
00041 };
00042
00043 virtual ~ EcalUncalibRecHitRatioMethodAlgo < C > () { };
00044 virtual EcalUncalibratedRecHit makeRecHit(const C & dataFrame,
00045 const double *pedestals,
00046 const double* pedestalRMSes,
00047 const double *gainRatios,
00048 std::vector < double >&timeFitParameters,
00049 std::vector < double >&litudeFitParameters,
00050 std::pair < double, double >&timeFitLimits);
00051
00052
00053
00054 void init( const C &dataFrame, const double * pedestals, const double * pedestalRMSes, const double * gainRatios );
00055 void computeTime(std::vector < double >&timeFitParameters, std::pair < double, double >&timeFitLimits, std::vector< double > &litudeFitParameters);
00056 void computeAmplitude( std::vector< double > &litudeFitParameters );
00057 CalculatedRecHit getCalculatedRecHit() { return calculatedRechit_; };
00058 bool fixMGPAslew( const C &dataFrame );
00059
00060 protected:
00061 std::vector < double > amplitudes_;
00062 std::vector < double > amplitudeErrors_;
00063 std::vector < Ratio > ratios_;
00064 std::vector < Tmax > times_;
00065 std::vector < Tmax > timesAB_;
00066
00067 double pedestal_;
00068 int num_;
00069 double ampMaxError_;
00070
00071 CalculatedRecHit calculatedRechit_;
00072 };
00073
00074 template <class C>
00075 void EcalUncalibRecHitRatioMethodAlgo<C>::init( const C &dataFrame, const double * pedestals, const double * pedestalRMSes, const double * gainRatios )
00076 {
00077 calculatedRechit_.timeMax = 5;
00078 calculatedRechit_.amplitudeMax = 0;
00079 calculatedRechit_.timeError = -999;
00080 amplitudes_.clear();
00081 amplitudes_.reserve(C::MAXSAMPLES);
00082 amplitudeErrors_.clear();
00083 amplitudeErrors_.reserve(C::MAXSAMPLES);
00084 ratios_.clear();
00085 ratios_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
00086 times_.clear();
00087 times_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
00088 timesAB_.clear();
00089 timesAB_.reserve(C::MAXSAMPLES*(C::MAXSAMPLES-1)/2);
00090
00091
00092
00093
00094
00095 pedestal_ = 0;
00096 num_ = 0;
00097 if (dataFrame.sample(0).gainId() == 1) {
00098 pedestal_ += double (dataFrame.sample(0).adc());
00099 num_++;
00100 }
00101 if (num_!=0 &&
00102 dataFrame.sample(1).gainId() == 1 &&
00103 fabs(dataFrame.sample(1).adc()-dataFrame.sample(0).adc())<3*pedestalRMSes[0]) {
00104 pedestal_ += double (dataFrame.sample(1).adc());
00105 num_++;
00106 }
00107 if (num_ != 0)
00108 pedestal_ /= num_;
00109 else
00110 pedestal_ = pedestals[0];
00111
00112
00113
00114
00115
00116 ampMaxError_ = 0;
00117 double ampMaxValue = -1000;
00118
00119
00120
00121
00122 double sample;
00123 double sampleError;
00124 int GainId;
00125 for (int iSample = 0; iSample < C::MAXSAMPLES; iSample++) {
00126 GainId = dataFrame.sample(iSample).gainId();
00127
00128 if (GainId == 1) {
00129 sample = double (dataFrame.sample(iSample).adc() - pedestal_);
00130 sampleError = pedestalRMSes[0];
00131 } else if (GainId == 2 || GainId == 3){
00132 sample = (double (dataFrame.sample(iSample).adc() - pedestals[GainId - 1])) *gainRatios[GainId - 1];
00133 sampleError = pedestalRMSes[GainId-1]*gainRatios[GainId-1];
00134 } else {
00135 sample = 1e-9;
00136 sampleError = 1e+9;
00137 }
00138
00139
00140 if(sampleError>0){
00141 amplitudes_.push_back(sample);
00142 amplitudeErrors_.push_back(sampleError);
00143 if(ampMaxValue < sample){
00144 ampMaxValue = sample;
00145 ampMaxError_ = sampleError;
00146 }
00147 }else{
00148
00149 amplitudes_.push_back(sample);
00150 amplitudeErrors_.push_back(1e+9);
00151 }
00152 }
00153 }
00154
00155 template <class C>
00156 bool EcalUncalibRecHitRatioMethodAlgo<C>::fixMGPAslew( const C &dataFrame )
00157 {
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167 bool result = false;
00168
00169 int GainIdPrev;
00170 int GainIdNext;
00171 for (int iSample = 1; iSample < C::MAXSAMPLES; iSample++) {
00172 GainIdPrev = dataFrame.sample(iSample-1).gainId();
00173 GainIdNext = dataFrame.sample(iSample).gainId();
00174 if( GainIdPrev>=1 && GainIdPrev<=3 &&
00175 GainIdNext>=1 && GainIdNext<=3 &&
00176 GainIdPrev<GainIdNext ){
00177 amplitudes_[iSample-1]=1e-9;
00178 amplitudeErrors_[iSample-1]=1e+9;
00179 result = true;
00180 }
00181 }
00182 return result;
00183
00184 }
00185
00186 template<class C>
00187 void EcalUncalibRecHitRatioMethodAlgo<C>::computeTime(std::vector < double >&timeFitParameters,
00188 std::pair < double, double >&timeFitLimits, std::vector < double >&litudeFitParameters)
00189 {
00191
00192
00193
00195 double ampMaxAlphaBeta = 0;
00196 double tMaxAlphaBeta = 5;
00197 double tMaxErrorAlphaBeta = 999;
00198 double tMaxRatio = 5;
00199 double tMaxErrorRatio = 999;
00200
00201 double sumAA = 0;
00202 double sumA = 0;
00203 double sum1 = 0;
00204 double sum0 = 0;
00205 double sumAf = 0;
00206 double sumff = 0;
00207 double NullChi2 = 0;
00208
00209
00210 for(unsigned int i = 0; i < amplitudes_.size(); i++){
00211 double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
00212 sum0 += 1;
00213 sum1 += 1/err2;
00214 sumA += amplitudes_[i]/err2;
00215 sumAA += amplitudes_[i]*amplitudes_[i]/err2;
00216 }
00217 if(sum0>0){
00218 NullChi2 = (sumAA - sumA*sumA/sum1)/sum0;
00219 }else{
00220
00221 return;
00222 }
00223
00224
00225
00226
00227
00228
00229
00230 double alphabeta = amplitudeFitParameters[0]*amplitudeFitParameters[1];
00231 double alpha = amplitudeFitParameters[0];
00232 double beta = amplitudeFitParameters[1];
00233
00234 for(unsigned int i = 0; i < amplitudes_.size()-1; i++){
00235 for(unsigned int j = i+1; j < amplitudes_.size(); j++){
00236
00237 if(amplitudes_[i]>1 && amplitudes_[j]>1){
00238
00239
00240 double Rtmp = amplitudes_[i]/amplitudes_[j];
00241
00242
00243
00244
00245 double err1 = Rtmp*Rtmp*( (amplitudeErrors_[i]*amplitudeErrors_[i]/(amplitudes_[i]*amplitudes_[i])) + (amplitudeErrors_[j]*amplitudeErrors_[j]/(amplitudes_[j]*amplitudes_[j])) );
00246
00247
00248 double stat;
00249 if(num_>0) stat = num_;
00250 else stat = 1;
00251 double err2 = amplitudeErrors_[j]*(amplitudes_[i]-amplitudes_[j])/(amplitudes_[j]*amplitudes_[j])/sqrt(stat);
00252
00253
00254
00255 double err3 = 0.289/amplitudes_[j];
00256
00257 double totalError = sqrt(err1 + err2*err2 +err3*err3);
00258
00259
00260
00261 if(totalError < 1.0
00262 && Rtmp>0.001
00263 && Rtmp<exp(double(j-i)/beta)-0.001
00264 ){
00265 Ratio currentRatio = { i, (j-i), Rtmp, totalError };
00266 ratios_.push_back(currentRatio);
00267 }
00268
00269 }
00270
00271 }
00272 }
00273
00274
00275 if(!ratios_.size() >0)
00276 return;
00277
00278
00279
00280
00281 for(unsigned int i = 0; i < ratios_.size(); i++){
00282
00283 double stepOverBeta = double(ratios_[i].step)/beta;
00284 double offset = double(ratios_[i].index) + alphabeta;
00285
00286 double Rmin = ratios_[i].value - ratios_[i].error;
00287 if(Rmin<0.001) Rmin=0.001;
00288
00289 double Rmax = ratios_[i].value + ratios_[i].error;
00290 double RLimit = exp(stepOverBeta)-0.001;
00291 if( Rmax > RLimit ) Rmax = RLimit;
00292
00293 double time1 = offset - ratios_[i].step/(exp((stepOverBeta-log(Rmin))/alpha)-1.0);
00294 double time2 = offset - ratios_[i].step/(exp((stepOverBeta-log(Rmax))/alpha)-1.0);
00295
00296
00297 double tmax = 0.5 * (time1 + time2);
00298 double tmaxerr = 0.5 * sqrt( (time1 - time2)*(time1 - time2) );
00299
00300
00301 sumAf = 0;
00302 sumff = 0;
00303 for(unsigned int it = 0; it < amplitudes_.size(); it++){
00304 double err2 = amplitudeErrors_[it]*amplitudeErrors_[it];
00305 double offset = (double(it) - tmax)/alphabeta;
00306 double term1 = 1.0 + offset;
00307 if(term1>1e-6){
00308 double f = exp( alpha*(log(1.0+offset) - offset) );
00309 sumAf += amplitudes_[it]*f/err2;
00310 sumff += f*f/err2;
00311 }
00312 }
00313
00314 double chi2 = sumAA;
00315 double amp = 0;
00316 if( sumff > 0 ){
00317 chi2 = sumAA - sumAf*sumAf/sumff;
00318 amp = sumAf/sumff;
00319 }
00320 chi2 /= sum0;
00321
00322
00323
00324 if(chi2 > 0 && tmaxerr > 0 && tmax > 0){
00325 Tmax currentTmax={ ratios_[i].index, ratios_[i].step, tmax, tmaxerr, amp, chi2 };
00326 timesAB_.push_back(currentTmax);
00327 }
00328 }
00329
00330
00331 if( !(timesAB_.size()> 0))
00332 return;
00333
00334
00335 double chi2min = 1.0e+9;
00336
00337
00338 for(unsigned int i = 0; i < timesAB_.size(); i++){
00339 if( timesAB_[i].chi2 <= chi2min ){
00340 chi2min = timesAB_[i].chi2;
00341
00342
00343 }
00344 }
00345
00346
00347
00348 double chi2Limit = chi2min + 1.0;
00349 double time_max = 0;
00350 double time_wgt = 0;
00351 for(unsigned int i = 0; i < timesAB_.size(); i++){
00352 if( timesAB_[i].chi2 < chi2Limit ){
00353 double inverseSigmaSquared = 1.0/(timesAB_[i].error*timesAB_[i].error);
00354 time_wgt += inverseSigmaSquared;
00355 time_max += timesAB_[i].value*inverseSigmaSquared;
00356 }
00357 }
00358
00359 tMaxAlphaBeta = time_max/time_wgt;
00360 tMaxErrorAlphaBeta = 1.0/sqrt(time_wgt);
00361
00362
00363 sumAf = 0;
00364 sumff = 0;
00365 for(unsigned int i = 0; i < amplitudes_.size(); i++){
00366 double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
00367 double offset = (double(i) - tMaxAlphaBeta)/alphabeta;
00368 double term1 = 1.0 + offset;
00369 if(term1>1e-6){
00370 double f = exp( alpha*(log(1.0+offset) - offset) );
00371 sumAf += amplitudes_[i]*f/err2;
00372 sumff += f*f/err2;
00373 }
00374 }
00375
00376 if( sumff > 0 ){
00377 ampMaxAlphaBeta = sumAf/sumff;
00378 double chi2AlphaBeta = (sumAA - sumAf*sumAf/sumff)/sum0;
00379 if(chi2AlphaBeta > NullChi2){
00380
00381 return;
00382 }
00383
00384 }else{
00385
00386 return;
00387 }
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401 if( ampMaxAlphaBeta/ampMaxError_ > 5.0 ){
00402
00403
00404
00405
00406 double time_max = 0;
00407 double time_wgt = 0;
00408
00409
00410 for (unsigned int i = 0; i < ratios_.size(); i++) {
00411
00412 if(ratios_[i].step == 1
00413 && ratios_[i].value >= timeFitLimits.first
00414 && ratios_[i].value <= timeFitLimits.second
00415 ){
00416
00417 double time_max_i = ratios_[i].index;
00418
00419
00420
00421 double u = timeFitParameters[timeFitParameters.size() - 1];
00422 for (int k = timeFitParameters.size() - 2; k >= 0; k--) {
00423 u = u * ratios_[i].value + timeFitParameters[k];
00424 }
00425
00426
00427 double du =
00428 (timeFitParameters.size() -
00429 1) * timeFitParameters[timeFitParameters.size() - 1];
00430 for (int k = timeFitParameters.size() - 2; k >= 1; k--) {
00431 du = du * ratios_[i].value + k * timeFitParameters[k];
00432 }
00433
00434
00435
00436 double errorsquared =
00437 ratios_[i].error * ratios_[i].error * du * du;
00438 if (errorsquared > 0) {
00439
00440 time_max += (time_max_i - u) / errorsquared;
00441 time_wgt += 1.0 / errorsquared;
00442 Tmax currentTmax =
00443 { ratios_[i].index, 1, (time_max_i - u),
00444 sqrt(errorsquared),0,1 };
00445 times_.push_back(currentTmax);
00446
00447 }
00448 }
00449 }
00450
00451
00452
00453 if (time_wgt > 0) {
00454 tMaxRatio = time_max/time_wgt;
00455 tMaxErrorRatio = 1.0/sqrt(time_wgt);
00456
00457
00458
00459 if( ampMaxAlphaBeta/ampMaxError_ > 10.0 ){
00460
00461
00462 calculatedRechit_.timeMax = tMaxRatio;
00463 calculatedRechit_.timeError = tMaxErrorRatio;
00464
00465 }else{
00466
00467
00468 calculatedRechit_.timeMax = ( tMaxAlphaBeta*(10.0-ampMaxAlphaBeta/ampMaxError_) + tMaxRatio*(ampMaxAlphaBeta/ampMaxError_ - 5.0) )/5.0;
00469 calculatedRechit_.timeError = ( tMaxErrorAlphaBeta*(10.0-ampMaxAlphaBeta/ampMaxError_) + tMaxErrorRatio*(ampMaxAlphaBeta/ampMaxError_ - 5.0) )/5.0;
00470
00471 }
00472
00473 }else{
00474
00475
00476 calculatedRechit_.timeMax = tMaxAlphaBeta;
00477 calculatedRechit_.timeError = tMaxErrorAlphaBeta;
00478
00479 }
00480
00481 }else{
00482
00483
00484 calculatedRechit_.timeMax = tMaxAlphaBeta;
00485 calculatedRechit_.timeError = tMaxErrorAlphaBeta;
00486
00487 }
00488 }
00489
00490 template<class C>
00491 void EcalUncalibRecHitRatioMethodAlgo<C>::computeAmplitude( std::vector< double > &litudeFitParameters )
00492 {
00494
00495
00496
00498
00499
00500 double alpha = amplitudeFitParameters[0];
00501 double beta = amplitudeFitParameters[1];
00502
00503
00504
00505 double pedestalLimit = calculatedRechit_.timeMax - (alpha * beta) - 1.0;
00506 double sumA = 0;
00507 double sumF = 0;
00508 double sumAF = 0;
00509 double sumFF = 0;
00510 double sum1 = 0;
00511 for (unsigned int i = 0; i < amplitudes_.size(); i++) {
00512 double err2 = amplitudeErrors_[i]*amplitudeErrors_[i];
00513 double f = 0;
00514 double termOne = 1 + (i - calculatedRechit_.timeMax) / (alpha * beta);
00515 if (termOne > 1.e-5) f = exp(alpha * log(termOne)) * exp(-(i - calculatedRechit_.timeMax) / beta);
00516
00517
00518
00519 if ( (i < pedestalLimit)
00520 || (f > 0.6 && i <= calculatedRechit_.timeMax)
00521 || (f > 0.4 && i >= calculatedRechit_.timeMax)) {
00522 sum1 += 1/err2;
00523 sumA += amplitudes_[i]/err2;
00524 sumF += f/err2;
00525 sumAF += f*amplitudes_[i]/err2;
00526 sumFF += f*f/err2;
00527 }
00528 }
00529
00530 calculatedRechit_.amplitudeMax = 0;
00531 if(sum1 > 0){
00532 double denom = sumFF*sum1 - sumF*sumF;
00533 if(fabs(denom)>1.0e-20){
00534 calculatedRechit_.amplitudeMax = (sumAF*sum1 - sumA*sumF)/denom;
00535 }
00536 }
00537 }
00538
00539
00540
00541 template < class C > EcalUncalibratedRecHit
00542 EcalUncalibRecHitRatioMethodAlgo < C >::makeRecHit(const C & dataFrame,
00543 const double *pedestals,
00544 const double *pedestalRMSes,
00545 const double *gainRatios,
00546 std::vector < double >&timeFitParameters,
00547 std::vector < double >&litudeFitParameters,
00548 std::pair < double, double >&timeFitLimits)
00549 {
00550
00551 init( dataFrame, pedestals, pedestalRMSes, gainRatios );
00552 computeTime( timeFitParameters, timeFitLimits, amplitudeFitParameters );
00553 computeAmplitude( amplitudeFitParameters );
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577 return EcalUncalibratedRecHit(dataFrame.id(),
00578 calculatedRechit_.amplitudeMax,
00579 pedestal_,
00580 calculatedRechit_.timeMax - 5,
00581 calculatedRechit_.timeError);
00582 }
00583 #endif