CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/SimG4Core/GFlash/TB/fitMaximumPoint.cc

Go to the documentation of this file.
00001 Double_t Pol2(Double_t *x, Double_t *par)
00002 {
00003   //polynoom(2)     
00004    Double_t part0  = par[0]; 
00005    Double_t part1  = par[1]*x[0];
00006    Double_t part2  = par[2]*x[0]*x[0];
00007    Double_t fitval = part0 + part1 + part2;
00008    return fitval;
00009 }
00010 
00011 Double_t Pol2_Special(Double_t *x, Double_t *par)
00012 {
00013   //polynoom(2)     
00014    Double_t part0  = par[0]; 
00015    Double_t part1  = -2.*par[2]*par[1]*x[0];
00016    Double_t part2  = par[2]*x[0]*x[0];
00017    Double_t fitval = part0 + part1 + part2;
00018    return fitval;
00019 }
00020 
00021 
00022 //====================================================
00023 //====================================================
00024 Float_t Fit_MaximumPoint( TH2F *h2, Int_t Ireturn = 0){  // xxx
00025 //====================================================
00026 //====================================================
00027 
00028   // -------------------------------------------
00029   // 1) Clean beam profile
00030   //     o all bins:  only bins >  2  entries
00031   //     o per x-bin: only bins > 75% of maximum 
00032   // -------------------------------------------
00033   // get some caracteristics
00034   Int_t   x_nbins = h2->GetXaxis()->GetNbins();  
00035   Float_t x_max   = h2->GetXaxis()->GetXmax();
00036   Float_t x_min   = h2->GetXaxis()->GetXmin();
00037   Float_t delta_x = h2->GetXaxis()->GetBinWidth(10);
00038   Int_t   y_nbins = h2->GetYaxis()->GetNbins();  
00039   Float_t y_max   = h2->GetYaxis()->GetXmax();
00040   Float_t y_min   = h2->GetYaxis()->GetXmin();
00041   Float_t delta_y = h2->GetYaxis()->GetBinWidth(10);
00042   TH2F *IHBeamProf_clean_1 = new TH2F("IHBeamProf_clean_1", "Max sam versus X (clean) ",  x_nbins, x_min, x_max, y_nbins, y_min, y_max); 
00043   TH2F *IHBeamProf_clean_2 = new TH2F("IHBeamProf_clean_2", "Max sam versus X (clean) ",  x_nbins, x_min, x_max, y_nbins, y_min, y_max); 
00044 
00045   // [1a] Only keep bins with more than 2 entries
00046   //      Also remember for each x-bin the maximum bin in y that satisfies this result
00047   Int_t y_binmax_high[500] = {1}; 
00048   for (Int_t ix=1;ix<x_nbins+1;ix++){         
00049     for (Int_t iy=1;iy<y_nbins+1;iy++){
00050        Int_t Nevts_clean = 0;
00051        if( h2->GetBinContent(ix,iy) > 1 ) {
00052            IHBeamProf_clean_1->SetBinContent(ix,iy, h2->GetBinContent(ix,iy));
00053            if( iy > y_binmax_high[ix] ){  y_binmax_high[ix] = iy; }         
00054        }                 
00055     }
00056   }
00057   // [1b] Only keep events with more than 85% of the maximum
00058   Int_t y_binmax_low[500] = {1}; 
00059   for (Int_t ix=1;ix<x_nbins+1;ix++){  
00060     y_binmax_low[ix]  = (Int_t)( 0.85*(Float_t)( y_binmax_high[ix]));
00061       for (Int_t iy=y_binmax_low[ix];iy<y_binmax_high[ix]+1;iy++){
00062         IHBeamProf_clean_2->SetBinContent(ix,iy,IHBeamProf_clean_1->GetBinContent(ix,iy));
00063       }
00064   }
00065 
00066 
00067   // -----------------------------------------------------------------------
00068   // 2) Find region to fit
00069   //     o Make profile by compuing means in every bin (error = RMS/sqrt(N))
00070   //     o Find maximum
00071   // -----------------------------------------------------------------------
00072   TH1F *IHBeamProf_clean_3 = new TH1F((h2->GetName()+(TString("_clprof"))).Data(), "Meam Energy versus X (clean) ",  x_nbins, x_min, x_max);  
00073 
00074   // [2a] Make TH1F containing the mean values
00075   TH1F *h_temp = new TH1F("h_temp"," Energy distribution for single bin in X", y_nbins, y_min, y_max);
00076   for (Int_t ix=1;ix<x_nbins+1;ix++){          
00077       double Nevt_slice = 0;
00078       //Int_t Nevt_slice = 0;
00079       for (Int_t iy=1;iy<y_nbins+1;iy++){            
00080         Int_t Nevt_tmp = (int)IHBeamProf_clean_2->GetBinContent(ix,iy); 
00081         Nevt_slice += Nevt_tmp;
00082         h_temp->SetBinContent(iy,Nevt_tmp);        
00083       }
00084       Float_t Y_mean       = h_temp->GetMean();
00085       Float_t Y_rms        = h_temp->GetRMS();
00086       Float_t Y_mean_error = (Nevt_slice>0) ? Y_rms/sqrt(Nevt_slice) : 9999.;
00087       IHBeamProf_clean_3->SetBinContent(ix,Y_mean);
00088       IHBeamProf_clean_3->SetBinError(ix,Y_mean_error);
00089       printf("%d %f %f %d %f\n",ix,Y_mean,Y_rms,Nevt_slice,Y_mean_error);
00090   }
00091 
00092   // [2b] Find maximum
00093   Float_t Y_tmp_max       = -999.;
00094   Float_t x_max_guess     = 0.;
00095   Int_t   x_max_guess_bin = 0;
00096   for (Int_t ix=1;ix<x_nbins+1;ix++){  
00097     Float_t X_tmp       = IHBeamProf_clean_3->GetBinCenter(ix);
00098     Float_t Y_tmp       = IHBeamProf_clean_3->GetBinContent(ix);
00099     Float_t Y_tmp_error = IHBeamProf_clean_3->GetBinError(ix);
00100     if( (fabs(X_tmp) < 10.) && (Y_tmp > Y_tmp_max) && (Y_tmp_error < 100.)){
00101       Y_tmp_max       = Y_tmp;
00102       x_max_guess     = X_tmp; 
00103       x_max_guess_bin = ix; 
00104       printf("Xtmp = %5.3f  Ytmp = %5.3f  Ytmp_max = %5.3f\n",X_tmp,Y_tmp,Y_tmp_max);
00105     }
00106   }
00107   printf("Fit: o Start for maximum = %5.3f\n",x_max_guess);
00108    
00109   // [2c] Define the fit range  (0.975% and error should be less than 8\%)
00110   Int_t fitbinmin =  9999;
00111   Int_t fitbinmax = -9999;
00112   Float_t x_tmp,y_tmp,e_tmp,er_tmp;
00113   for (Int_t ix = x_max_guess_bin; ix<x_nbins+1;ix++){ 
00114        fitbinmax = ix;
00115        x_tmp     = IHBeamProf_clean_3->GetBinCenter(ix);
00116        y_tmp     = IHBeamProf_clean_3->GetBinContent(ix);
00117        e_tmp     = IHBeamProf_clean_3->GetBinError(ix);
00118        er_tmp    = (y_tmp>0.) ? (e_tmp/y_tmp) : 1.00;
00119        printf("%3d %f %f %f  --   %f %f\n",ix,x_tmp,y_tmp,e_tmp,y_tmp/Y_tmp_max,er_tmp);
00120        if( y_tmp < 0.975*Y_tmp_max  && er_tmp < 0.008) break;         
00121   }  
00122   for (Int_t ix = x_max_guess_bin; ix>1        ;ix--){ 
00123        fitbinmin = ix;
00124        x_tmp     = IHBeamProf_clean_3->GetBinCenter(ix);
00125        y_tmp     = IHBeamProf_clean_3->GetBinContent(ix);
00126        e_tmp     = IHBeamProf_clean_3->GetBinError(ix);
00127        er_tmp    = (y_tmp>0.) ? (e_tmp/y_tmp) : 1.00;
00128        printf("%3d %f %f %f  --   %f %f\n",ix,x_tmp,y_tmp,e_tmp,y_tmp/Y_tmp_max,er_tmp);
00129        if( y_tmp < 0.975* Y_tmp_max  && er_tmp < 0.008) break;         
00130   }
00131 
00132   Double_t posmin_fit = x_min+fitbinmin*delta_x;
00133   Double_t posmax_fit = x_min+fitbinmax*delta_x;
00134   printf("     o Fit range = %5.3f -- %5.3f -- %5.3f\n",posmin_fit,x_max_guess,posmax_fit);
00135   if( fabs(posmax_fit - posmin_fit ) < 4. ){
00136     printf("Something is wrong with this range: returning dummy value\n");
00137     posmin_fit = x_max_guess - 6.0;
00138     posmax_fit = x_max_guess + 6.0;
00139     return -99.00;
00140   } 
00141 
00142   // -------------------------------------------------
00143   // 3) Do the actual fit
00144   //    o make clone of histogram
00145   // -------------------------------------------------
00146   TH1F *h1  = (TH1F *) IHBeamProf_clean_3->Clone();
00147   // 3a] Do the actual fit
00148   Double_t fitresults[3]          = {0.};
00149   Double_t fitresults_errors[3]   = {0.};
00150   TF1 *f1 = new TF1("f1",Pol2,-50.,50.,3);  
00151   f1->SetParameters( 1.,1.,1.);
00152   h1->Fit("f1","Q","",posmin_fit, posmax_fit);
00153   for(int i=0 ; i< 3 ; i++) { 
00154     fitresults[i]        = f1->GetParameter(i); 
00155     fitresults_errors[i] = f1->GetParError(i); 
00156   }    
00157   Float_t chi2 = f1->GetChisquare()/f1->GetNDF();  
00158   Float_t a = fitresults[2]; Float_t da = fitresults_errors[2];
00159   Float_t b = fitresults[1]; Float_t db = fitresults_errors[1];
00160   Float_t c = fitresults[0]; Float_t dc = fitresults_errors[0];
00161   Float_t x0 = (-1.*b)/(2*a);
00162   printf("  a = %7.2f   b = %7.2f   c = %7.2f\n",a,b,c);
00163   printf(" da = %7.2f  db = %7.2f  dc = %7.2f\n",da,db,dc);
00164 
00165   cout << "risultati del fit polinomiale: " << fitresults[0] << " " << fitresults[1] << " " << fitresults[2] << endl;
00166 
00167   char myProfTitle[200];
00168   sprintf(myProfTitle, h2->GetTitle());
00169   strcat(myProfTitle,"_prof");
00170   h1->Write(myProfTitle);  
00171   //h1->Write(); //write the profile
00172   // ----------------------------
00173   // 4) Compute uncertainty on x0
00174   // ----------------------------
00175   // [4a] compute dxo using the covariance matrix
00176   // covariance matrix
00177   Double_t CovMat[3][3];
00178   gMinuit->mnemat(&CovMat[0][0],3);
00179   Float_t v11 = CovMat[0][0];  Float_t v12 = CovMat[0][1]; Float_t v13 = CovMat[0][2];
00180   Float_t v21 = CovMat[1][0];  Float_t v22 = CovMat[1][1]; Float_t v23 = CovMat[1][2];
00181   Float_t v31 = CovMat[2][0];  Float_t v32 = CovMat[2][1]; Float_t v33 = CovMat[2][2];
00182   printf("Covariance Matrix:   v11 = %f     v12 = %f  v13 = %f\n",CovMat[0][0],CovMat[0][1],CovMat[0][2]);
00183   printf("                     v21 = %f     v22 = %f  v23 = %f\n",CovMat[1][0],CovMat[1][1],CovMat[1][2]);
00184   printf("                     v31 = %f     v32 = %f  v33 = %f\n",CovMat[2][0],CovMat[2][1],CovMat[2][2]);
00185   // jacobiaan
00186   Float_t j1  =  b/(2*a*a);
00187   Float_t j2  = -1./(2*a);
00188   // determinant covariance martix
00189   Float_t det =   v11*(v33*v22-v32*v23)-v21*(v33*v12-v32*v13)+v31*(v23*v12-v22*v13);
00190   printf("Determinant = %f\n",det);
00191   // inverse matrix
00192   Float_t v11_i =   v33*v22-v32*v23;   Float_t v12_i = -(v33*v12-v32*v13); Float_t v13_i =   v23*v12-v22*v13;
00193   Float_t v21_i = -(v33*v21-v31*v23);  Float_t v22_i =   v33*v11-v31*v13 ; Float_t v23_i = -(v23*v11-v21*v13) ;
00194   Float_t v31_i =   v32*v21-v31*v22;   Float_t v32_i = -(v32*v11-v31*v12); Float_t v33_i =   v22*v11-v21*v12;
00195   // variance
00196   Float_t var     = j1*(j1*v11_i+j2*v12_i)+j2*(j1*v21_i+j2*v22_i);
00197            var    /= det;
00198   Float_t dx0     = sqrt(var);    
00199   printf("Type 1 fit:  o  x0 = %f  +/- %f\n",x0,dx0);
00200 
00201 
00202   // ---------------------
00203   // 5) Second type of fit
00204   // ---------------------  
00205   // 5a] Do the actual fit
00206   TH1F *h0  = (TH1F *) IHBeamProf_clean_3->Clone();
00207   Double_t fitresults0[3]          = {0.};
00208   Double_t fitresults0_errors[3]   = {0.};
00209   TF1 *f0 = new TF1("f0",Pol2_Special,-50.,50.,3);  
00210   f0->SetParameters( -1.,x_max_guess,5000.);
00211   h0->Fit("f0","Q","",posmin_fit, posmax_fit);
00212   for(int i=0 ; i< 3 ; i++) { 
00213     fitresults0[i]        = f0->GetParameter(i); 
00214     fitresults0_errors[i] = f0->GetParError(i); 
00215   }    
00216   Float_t a0 = fitresults0[2]; Float_t da0 = fitresults0_errors[2];
00217   Float_t b0 = fitresults0[1]; Float_t db0 = fitresults0_errors[1];
00218   Float_t c0 = fitresults0[0]; Float_t dc0 = fitresults0_errors[0];
00219 
00220   Double_t CovMat0[3][3];
00221   gMinuit->mnemat(&CovMat0[0][0],3);
00222   //printf("Cov0[1][1] = %f\n",CovMat0[1][1]);
00223   Float_t db0cov = sqrt(CovMat0[1][1]);
00224   printf("Type 2 fit:  o  x0 = %f  +/- %f %f \n",b0,db0,db0cov);
00225   delete  IHBeamProf_clean_1;
00226   delete  IHBeamProf_clean_2;
00227   delete  IHBeamProf_clean_3;
00228   delete  h_temp;
00229   delete  f1;
00230   delete  f0;
00231 
00232   if(Ireturn == 10) {return x0;}
00233   if(Ireturn == 11) {return dx0;}
00234   if(Ireturn == 20) {return b0;}
00235   if(Ireturn == 21) {return db0cov;}
00236   if(Ireturn == 30) {return chi2;}
00237 
00238   if(Ireturn == 99) {return c0;}
00239 
00240 
00241 
00242 }