CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoLocalCalo/EcalRecAlgos/interface/EcalUncalibRecHitFixedAlphaBetaAlgo.h

Go to the documentation of this file.
00001 #ifndef RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitFixedAlphaBetaAlgo_HH
00002 #define RecoLocalCalo_EcalRecAlgos_EcalUncalibRecHitFixedAlphaBetaAlgo_HH
00003 
00014 #include "CLHEP/Matrix/Vector.h"
00015 #include "CLHEP/Matrix/SymMatrix.h"
00016 //#include "CLHEP/Matrix/Matrix.h"
00017 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalUncalibRecHitRecAbsAlgo.h"
00018 
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 #include <vector>
00022 #include <string>
00023 
00024 //#include "TROOT.h"
00025 //#include "TMinuit.h"
00026 //#include "TGraph.h"
00027 //#include "TF1.h"
00028 //#include "TMatrixD.h"
00029 //#include "TVectorD.h"
00030 
00031 template<class C> class EcalUncalibRecHitFixedAlphaBetaAlgo : public EcalUncalibRecHitRecAbsAlgo<C>
00032 {
00033 
00034 
00035  private:
00036   double MinAmpl_;
00037   bool dyn_pedestal;
00038   double fAlpha_;//parameter of the shape
00039   double fBeta_;//parameter of the shape
00040   double fAmp_max_;// peak amplitude 
00041   double fTim_max_ ;// time of the peak (in 25ns units)
00042   double fPed_max_;// pedestal value
00043   double alfabeta_; 
00044 
00045   
00046   int fNb_iter_;
00047   int  fNum_samp_bef_max_ ;
00048   int fNum_samp_after_max_;
00049   
00050   float fSigma_ped; 
00051   double un_sur_sigma;
00052   //temporary solution for different alpha and beta
00053   float alpha_table_[36][1701];
00054   float beta_table_[36][1701];
00055 
00056   bool doFit_;
00057 
00058   double pulseShapeFunction(double t);
00059   float PerformAnalyticFit(double* samples, int max_sample);
00060   void InitFitParameters(double* samples, int max_sample);
00061   CLHEP::HepSymMatrix DM1_ ; CLHEP::HepVector temp_;
00062    public:
00063 
00064   EcalUncalibRecHitFixedAlphaBetaAlgo<C>():fAlpha_(0.),fBeta_(0.),fAmp_max_(-1.),fTim_max_(-1),fPed_max_(0),alfabeta_(0),fNb_iter_(4),fNum_samp_bef_max_(1),fNum_samp_after_max_(3),fSigma_ped(1.1),DM1_(3),temp_(3){
00065     un_sur_sigma = 1./double(fSigma_ped) ;
00066     for (int i=0;i<36;i++){
00067       for(int j=0;j<1701;j++){
00068         alpha_table_[i][j] = 1.2 ;
00069         beta_table_[i][j]  =  1.7 ;
00070       }
00071     }
00072     doFit_ = false;
00073     MinAmpl_ = 16;
00074     dyn_pedestal = true;
00075   }
00076   EcalUncalibRecHitFixedAlphaBetaAlgo<C>(int n_iter, int n_bef_max =1, int n_aft_max =3, float sigma_ped = 1.1):fAlpha_(0.),fBeta_(0.),fAmp_max_(-1.),fTim_max_(-1),fPed_max_(0),alfabeta_(0),DM1_(3),temp_(3){
00077     
00078    fNb_iter_ = n_iter ;
00079    fNum_samp_bef_max_   = n_bef_max ;
00080    fNum_samp_after_max_ = n_aft_max ;
00081    fSigma_ped = sigma_ped;
00082    un_sur_sigma = 1./double(fSigma_ped) ;
00083    for (int i=0;i<36;i++){
00084      for(int j=0;j<1701;j++){
00085        alpha_table_[i][j] = 1.2 ;
00086        beta_table_[i][j]  =  1.7 ;
00087      }
00088    }
00089    doFit_ = false;
00090    MinAmpl_ = 16;
00091    dyn_pedestal = true;
00092   };
00093 
00094   virtual ~EcalUncalibRecHitFixedAlphaBetaAlgo<C>() { };
00095   virtual EcalUncalibratedRecHit makeRecHit(const C& dataFrame, const double* pedestals,
00096                                             const double* gainRatios,
00097                                             const EcalWeightSet::EcalWeightMatrix** weights, 
00098                                             const EcalWeightSet::EcalChi2WeightMatrix** chi2Matrix); 
00099   void SetAlphaBeta( double alpha, double beta);
00100   void SetMinAmpl(double ampl);
00101   void SetDynamicPedestal(bool dyn_pede);
00102 };
00103 
00104 
00105 
00107 template<class C> EcalUncalibratedRecHit  EcalUncalibRecHitFixedAlphaBetaAlgo<C>::makeRecHit(const C& dataFrame, const double* pedestals,
00108                                                                                              const double* gainRatios,
00109                                                                                              const EcalWeightSet::EcalWeightMatrix** weights, 
00110                                                                                              const EcalWeightSet::EcalChi2WeightMatrix** chi2Matrix){
00111   double chi2_(-1.);
00112   
00113   //  double Gain12Equivalent[4]={0,1,2,12};
00114   double frame[C::MAXSAMPLES];// will contain the ADC values
00115   double pedestal =0;     // carries pedestal for gain12 i.e. gainId==1
00116   
00117   int gainId0 = 1;        // expected gainId at the beginning of dataFrame 
00118   int iGainSwitch = 0;    // flags whether there's any gainId other than gainId0
00119   int GainId= 0;          // stores gainId at every sample
00120   double maxsample(-1);   // ADC value of maximal ped-subtracted sample
00121   int imax(-1);           // sample number of maximal ped-subtracted sample
00122   bool external_pede = false;
00123   bool isSaturated = 0;   // flag reporting whether gain0 has been found
00124 
00125   // Get time samples checking for Gain Switch and pedestals
00126   if(pedestals){
00127     external_pede = true;
00128     if(dyn_pedestal) { pedestal = (double(dataFrame.sample(0).adc()) + double(dataFrame.sample(1).adc()))/2.;}
00129     else{ pedestal  = pedestals[0];}
00130     for(int iSample = 0; iSample < C::MAXSAMPLES; iSample++) {
00131         //create frame in adc gain 12 equivalent
00132         GainId = dataFrame.sample(iSample).gainId();
00133 
00134         // FIX-ME: warning: the vector pedestal is supposed to have in the order G12, G6 and G1
00135         // if GainId is zero treat it as 3 temporarily to protect against undefined
00136         // frame will be set to ~max of gain1
00137         if ( GainId == 0 )
00138           { 
00139             GainId = 3;
00140             isSaturated = 1;
00141           }
00142 
00143         if (GainId != gainId0) iGainSwitch = 1;
00144 
00145         if(GainId==gainId0){frame[iSample] = double(dataFrame.sample(iSample).adc())-pedestal ;}
00146         else {frame[iSample] = (double(dataFrame.sample(iSample).adc())-pedestals[GainId-1])*gainRatios[GainId-1];}
00147 
00148         if( frame[iSample]>maxsample ) {
00149           maxsample = frame[iSample];
00150           imax = iSample;
00151         }
00152       }
00153   }
00154   else {// pedestal from pre-sample
00155     external_pede = false;
00156     pedestal = (double(dataFrame.sample(0).adc()) + double(dataFrame.sample(1).adc()))/2.;
00157 
00158     for(int iSample = 0; iSample < C::MAXSAMPLES; iSample++) {
00159       //create frame in adc gain 12 equivalent
00160       GainId = dataFrame.sample(iSample).gainId();
00161       //no gain switch forseen if there is no external pedestal
00162       if ( GainId == 0 ) 
00163         {
00164           GainId = 3;
00165           isSaturated = 1;
00166         }
00167 
00168       frame[iSample] = double(dataFrame.sample(iSample).adc())-pedestal ;
00169       // if gain has switched but no pedestals are available, no much good you can do...
00170       if (GainId > gainId0) iGainSwitch = 1;
00171       if( frame[iSample]>maxsample ) {
00172         maxsample = frame[iSample];
00173         imax = iSample;
00174       }
00175     } 
00176   }
00177 
00178   if( (iGainSwitch==1 && external_pede==false) ||  // ... thus you return dummy rechit
00179       imax ==-1 ){                                 // protect against all frames being <-1
00180     return EcalUncalibratedRecHit( dataFrame.id(), -1., -100., -1. , -1.);
00181   }
00182   
00183   InitFitParameters(frame, imax);
00184   chi2_ = PerformAnalyticFit(frame,imax);
00185   uint32_t flags = 0;
00186   if (isSaturated) flags = EcalUncalibratedRecHit::kSaturated;
00187 
00188   /*    std::cout << "separate fits\nA: " << fAmp_max_  << ", ResidualPed: " <<  fPed_max_
00189               <<", pedestal: "<<pedestal << ", tPeak " << fTim_max_ << std::endl;
00190   */
00191   return EcalUncalibratedRecHit( dataFrame.id(),fAmp_max_, pedestal+fPed_max_, fTim_max_ - 5 , chi2_, flags );
00192 }
00193 
00194 template<class C> double EcalUncalibRecHitFixedAlphaBetaAlgo<C>::pulseShapeFunction(double t){
00195   if( alfabeta_ <= 0 ) return((double)0.);
00196   double dtsbeta,variable,puiss;
00197   double dt = t-fTim_max_ ;
00198   if(dt > -alfabeta_)  {
00199     dtsbeta=dt/fBeta_ ;
00200     variable=1.+dt/alfabeta_ ;
00201     puiss=pow(variable,fAlpha_);
00202     return fAmp_max_*puiss*exp(-dtsbeta)+fPed_max_ ;
00203   }
00204   return  fPed_max_ ;
00205 }
00206 
00207 template<class C> void  EcalUncalibRecHitFixedAlphaBetaAlgo<C>::InitFitParameters(double* samples, int max_sample){
00208 
00209   // in a first attempt just use the value of the maximum sample 
00210   fAmp_max_ = samples[max_sample];
00211   fTim_max_ = max_sample;
00212   fPed_max_ = 0;
00213 
00214   // amplitude too low for fit to converge
00215   // timing set correctly is assumed
00216   if(fAmp_max_ <  MinAmpl_){
00217     fAmp_max_      = samples[5];
00218     double sumA    = samples[5]+samples[4]+samples[6];
00219     if(sumA != 0) { fTim_max_ = 5+(samples[6]-samples[4])/sumA; }
00220     else{ fTim_max_ = -993; }//-999+6
00221     doFit_ = false;
00222   }
00223   // if timing very badly off, that just use max sample
00224   else if(max_sample <1 || max_sample > 7)
00225     {    doFit_ = false;}
00226   else{
00227     //y=a*(x-xM)^2+b*(x-xM)+c
00228     doFit_ = true;
00229     float a = float(samples[max_sample-1]+samples[max_sample+1]-2*samples[max_sample])/2.;
00230     if(a==0){doFit_ =false; return;}
00231     float b = float(samples[max_sample+1]-samples[max_sample-1])/2.;
00232   
00233     fTim_max_ = max_sample - b/(2*a);
00234     fAmp_max_ =  samples[max_sample] - b*b/(4*a); 
00235   }  
00236   
00237 } 
00238 
00239 template<class C> float EcalUncalibRecHitFixedAlphaBetaAlgo<C>::PerformAnalyticFit(double* samples, int max_sample){
00240   
00241   //int fValue_tim_max = max_sample;  
00246   //| the pedestal (fPed_max_)  
00247   
00248   double chi2=-1 , db[3] ;
00249   
00250 
00251   //HepSymMatrix DM1(3) ; CLHEP::HepVector temp(3) ;
00252 
00253   int num_fit_min =(int)(max_sample - fNum_samp_bef_max_ ) ;
00254   int num_fit_max =(int)(max_sample + fNum_samp_after_max_) ;
00255 
00256   if (num_fit_min<0) num_fit_min = 0 ; 
00257   //if (num_fit_max>=fNsamples-1) num_fit_max = fNsamples-2 ;
00258   if (num_fit_max>= C::MAXSAMPLES ) {num_fit_max = C::MAXSAMPLES-1 ;}
00259 
00260   if(! doFit_ ) {
00261     LogDebug("EcalUncalibRecHitFixedAlphaBetaAlgo")<<"No fit performed. The amplitude of sample 5 will be used";
00262     return -1;
00263   }
00264 
00265   double func,delta ;
00266   double variation_func_max = 0. ;double variation_tim_max = 0. ; double variation_ped_max = 0. ;
00268   for (int iter=0 ; iter < fNb_iter_ ; iter ++) {
00270     chi2 = 0. ; //PROD.Zero() ;  DM1.Zero() ;
00271 
00272      for(int i1=0 ; i1<3 ; i1++) {
00273        temp_[i1]=0;
00274         for(int j1=i1 ; j1<3 ; j1++) { 
00275           DM1_.fast(j1+1,i1+1) = 0; }
00276       }
00277 
00278     fAmp_max_ += variation_func_max ;
00279     fTim_max_ += variation_tim_max ;
00280     fPed_max_ += variation_ped_max ;
00281     
00283     for( int i = num_fit_min ; i <= num_fit_max ; i++) {
00284       //if(i>fsamp_edge_fit && i<num_fit_min) continue ; // remove front edge samples
00286       func = pulseShapeFunction( (double)i  ) ;
00288       double dt =(double)i - fTim_max_ ;
00289       if(dt > -alfabeta_)  {      
00290         double dt_sur_beta = dt/fBeta_ ;
00291         double variable = (double)1. + dt/alfabeta_ ;
00292         double expo = exp(-dt_sur_beta) ;        
00293         double puissance = pow(variable,fAlpha_) ;
00294         
00295         db[0]=un_sur_sigma*puissance*expo ;
00296         db[1]=fAmp_max_*db[0]*dt_sur_beta/(alfabeta_*variable) ;
00297       }
00298       else {
00299         db[0]=0. ; db[1]=0. ; 
00300       }
00301       db[2]=un_sur_sigma ;
00303       for(int i1=0 ; i1<3 ; i1++) {
00304         for(int j1=i1 ; j1<3 ; j1++) { 
00305           //double & fast(int row, int col);
00306           DM1_.fast(j1+1,i1+1) += db[i1]*db[j1]; }
00307       }
00309       delta = (samples[i]-func)*un_sur_sigma ;
00311       for(int ii=0 ; ii<3 ;ii++) {temp_[ii] += delta*db[ii] ;}
00312       chi2 += delta * delta ;
00313     }
00314     
00315     int fail=0 ;
00316     DM1_.invert(fail) ;
00317       if(fail != 0.) {
00318       //just a guess from the value of the parameters in the previous interaction;
00319       //printf("wH4PulseFitWithFunction =====> determinant error --> No Fit Provided !\n") ;
00320       InitFitParameters(samples,max_sample);
00321       return -101 ;
00322     }
00323 /*     for(int i1=0 ; i1<3 ; i1++) { */
00324 /*       for(int j1=0 ; j1<3 ; j1++) {  */
00325 /*      //double & fast(int row, int col); */
00326 /*      std::cout<<"inverted: "<<DM1[j1][i1]<<std::endl;;} */
00327 /*     } */
00328 /*     std::cout<<"vector temp: "<< temp[0]<<" "<<temp[1]<<" "<<temp[2]<<std::endl; */
00330     CLHEP::HepVector PROD = DM1_*temp_ ;
00331     //    std::cout<<"vector PROD: "<< PROD[0]<<" "<<PROD[1]<<" "<<PROD[2]<<std::endl;
00332 
00333     // Probably the fastest way to protect against
00334     // +-inf value in the matrix DM1_ after inversion
00335     // (which is nevertheless flagged as successfull...)
00336     if ( std::isnan( PROD[0] ) ) {
00337             InitFitParameters(samples,max_sample);
00338             return -103 ;
00339     }
00340 
00341     variation_func_max = PROD[0] ;
00342     variation_tim_max = PROD[1] ;
00343     variation_ped_max = PROD[2] ;
00344     //chi2 = chi2/((double)nsamp_used - 3.) ;
00345   }
00346   
00347 
00349   if( variation_func_max > 2000. || variation_func_max < -1000. ) {
00350     InitFitParameters(samples,max_sample);
00351     return -102;
00352   }
00353   
00354 
00356   fAmp_max_ += variation_func_max ;
00357   fTim_max_ += variation_tim_max ;
00358   fPed_max_ += variation_ped_max ;
00359 
00360 
00361   // protection against unphysical results:
00362   // ampli mismatched to MaxSample, ampli largely negative, time off preselected range
00363   if( fAmp_max_>2*samples[max_sample]  ||  fAmp_max_<-100 ||  fTim_max_<0  ||  9<fTim_max_ ) {
00364     InitFitParameters(samples,max_sample);
00365     return -104;
00366   }
00367   
00368 
00369   //std::cout <<"chi2: "<<chi2<<" ampl: "<<fAmp_max_<<" time: "<<fTim_max_<<" pede: "<<fPed_max_<<std::endl;
00370   return chi2;
00371 }
00372 
00373 template<class C> void EcalUncalibRecHitFixedAlphaBetaAlgo<C>::SetAlphaBeta( double alpha, double beta){
00374   fAlpha_ = alpha;
00375   fBeta_=  beta;
00376   alfabeta_ = fAlpha_*fBeta_;
00377 }
00378 
00379 template<class C> void EcalUncalibRecHitFixedAlphaBetaAlgo<C>::SetMinAmpl( double ampl){
00380   MinAmpl_ = ampl;
00381 }
00382 template<class C> void EcalUncalibRecHitFixedAlphaBetaAlgo<C>::SetDynamicPedestal (bool p){
00383   dyn_pedestal = p;
00384 }
00385 #endif