CMS 3D CMS Logo

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