CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/DQM/SiStripHistoricInfoClient/interface/fitUtilities.h

Go to the documentation of this file.
00001 #ifndef DQM_SiStripHistoricInfoClient_fitUtilities_H
00002 #define DQM_SiStripHistoricInfoClient_fitUtilities_H
00003 
00004 #include <vector>
00005 #include <cstring>
00006 #include <iostream>
00007 #include <sstream>
00008 #include <string>
00009 #include "TH1.h"
00010 #include "TF1.h"
00011 #include "TMath.h"
00012 #include "DQMServices/Core/interface/MonitorElement.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 
00024 double langaufun(double *x, double *par);
00025 int32_t langaupro(double *params, double &maxx, double &FWHM);
00026 double Gauss(double *x, double *par);
00027 
00028 class fitUtilities{
00029 
00030  public: 
00031 
00032   fitUtilities();
00033   ~fitUtilities();
00034  
00035   void   init();
00036   double doLanGaussFit(MonitorElement* ME){return doLanGaussFit(ME->getTH1F());}
00037   double doLanGaussFit(TH1F*);
00038 
00039   double doGaussFit(MonitorElement* ME){return doGaussFit(ME->getTH1F());}
00040   double doGaussFit(TH1F*);
00041   
00042   double getLanGaussPar(std::string s)   ;
00043   double getLanGaussParErr(std::string s);
00044   double getLanGaussConv(std::string s)  ;
00045   
00046   double getGaussPar(std::string s)   ;
00047   double getGaussParErr(std::string s);
00048 
00049   double getFitChi()  {return chi2GausS;}
00050   int    getFitnDof() {return nDofGausS;}
00051  
00052  private:
00053   
00054   double pLanGausS[4], epLanGausS[4];
00055   double pGausS[3], epGausS[3];
00056   double pLanConv[2];
00057   double chi2GausS;
00058   int32_t nDofGausS;
00059   TF1 *langausFit; 
00060   TF1 *gausFit;  
00061 };
00062 
00063 #endif // DQM_SiStripHistoricInfoClient_fitUtilities_H
00064 
00065 //-----------------------------------------------------------------------------------------------
00066 double langaufun(double *x, double *par){
00067   //Fit parameters:
00068   //par[0]=Width (scale) parameter of Landau density
00069   //par[1]=Most Probable (MP, location) parameter of Landau density
00070   //par[2]=Total area (integral -inf to inf, normalization constant)
00071   //par[3]=Width (sigma) of convoluted Gaussian function
00072   //
00073   //In the Landau distribution (represented by the CERNLIB approximation), 
00074   //the maximum is located at x=-0.22278298 with the location parameter=0.
00075   //This shift is corrected within this function, so that the actual
00076   //maximum is identical to the MP parameter.
00077 
00078   // Numeric constants
00079   double invsq2pi = 0.3989422804014;   // (2 pi)^(-1/2)
00080   double mpshift  = -0.22278298;       // Landau maximum location
00081 
00082   // Control constants
00083   double np = 100.0;      // number of convolution steps
00084   double sc =   5.0;      // convolution extends to +-sc Gaussian sigmas
00085 
00086   // Variables
00087   double xx;
00088   double mpc;
00089   double fland;
00090   double sum = 0.0;
00091   double xlow,xupp;
00092   double step;
00093   double i;
00094 
00095 
00096   // MP shift correction
00097   mpc = par[1] - mpshift * par[0]; 
00098 
00099   // Range of convolution integral
00100   xlow = x[0] - sc * par[3];
00101   xupp = x[0] + sc * par[3];
00102 
00103   step = (xupp-xlow) / np;
00104 
00105   // Landau Distribution Production
00106   for(i=1.0; i<=np/2; i++) {
00107     xx = xlow + (i-.5) * step;
00108     fland = TMath::Landau(xx,mpc,par[0]) / par[0];
00109     sum += fland * TMath::Gaus(x[0],xx,par[3]);
00110 
00111     xx = xupp - (i-.5) * step;
00112     fland = TMath::Landau(xx,mpc,par[0]) / par[0];
00113     sum += fland * TMath::Gaus(x[0],xx,par[3]);
00114   }
00115 
00116   return (par[2] * step * sum * invsq2pi / par[3]);
00117 }
00118 
00119 //-----------------------------------------------------------------------------------------------
00120 int32_t langaupro(double *params, double &maxx, double &FWHM) {
00121    edm::LogInfo("fitUtility") << "inside langaupro " << std::endl;
00122   // Seaches for the location (x value) at the maximum of the 
00123   // Landau and its full width at half-maximum.
00124   //
00125   // The search is probably not very efficient, but it's a first try.
00126 
00127   double p,x,fy,fxr,fxl;
00128   double step;
00129   double l,lold,dl;
00130   int32_t i = 0;
00131   const int32_t MAXCALLS = 10000;
00132   const double dlStop = 1e-3; // relative change < .001
00133 
00134   // Search for maximum
00135   p = params[1] - 0.1 * params[0];
00136   step = 0.05 * params[0];
00137   lold = -2.0;
00138   l    = -1.0;
00139 
00140   dl = (l-lold)/lold;    // FIXME catch divide by zero
00141   while ( (TMath::Abs(dl)>dlStop ) && (i < MAXCALLS) ) {
00142     i++;
00143 
00144     lold = l;
00145     x = p + step;
00146     l = langaufun(&x,params);
00147     dl = (l-lold)/lold; // FIXME catch divide by zero
00148         
00149     if (l < lold)
00150       step = -step/10;
00151  
00152     p += step;
00153   }
00154 
00155   if (i == MAXCALLS)
00156     return (-1);
00157 
00158   maxx = x;
00159 
00160   fy = l/2;
00161 
00162 
00163   // Search for right x location of fy
00164   p = maxx + params[0];
00165   step = params[0];
00166   lold = -2.0;
00167   l    = -1e300;
00168   i    = 0;
00169 
00170   dl = (l-lold)/lold;   // FIXME catch divide by zero
00171   while ( ( TMath::Abs(dl)>dlStop ) && (i < MAXCALLS) ) {
00172     i++;
00173 
00174     lold = l;
00175     x = p + step;
00176     l = TMath::Abs(langaufun(&x,params) - fy);
00177     dl = (l-lold)/lold; // FIXME catch divide by zero
00178  
00179     if (l > lold)
00180       step = -step/10;
00181  
00182     p += step;
00183   }
00184 
00185   if (i == MAXCALLS)
00186     return (-2);
00187 
00188   fxr = x;
00189 
00190 
00191   // Search for left x location of fy
00192   p = maxx - 0.5 * params[0];
00193   step = -params[0];
00194   lold = -2.0;
00195   l    = -1e300;
00196   i    = 0;
00197 
00198   dl = (l-lold)/lold;    // FIXME catch divide by zero
00199   while ( ( TMath::Abs(dl)>dlStop ) && (i < MAXCALLS) ) {
00200     i++;
00201 
00202     lold = l;
00203     x = p + step;
00204     l = TMath::Abs(langaufun(&x,params) - fy);
00205     dl = (l-lold)/lold; // FIXME catch divide by zero
00206  
00207     if (l > lold)
00208       step = -step/10;
00209  
00210     p += step;
00211   }
00212 
00213   if (i == MAXCALLS)
00214     return (-3);
00215 
00216 
00217   fxl = x;
00218 
00219   FWHM = fxr - fxl;
00220   return (0);
00221 }
00222 
00223 
00224 //-----------------------------------------------------------------------------------------------
00225 double Gauss(double *x, double *par)
00226 {
00227   // The noise function: a gaussian
00228 
00229   double arg = 0;
00230   if (par[2]) arg = (x[0] - par[1])/par[2];
00231 
00232   double noise = par[0]*TMath::Exp(-0.5*arg*arg);
00233   return noise;
00234 }
00235 
00236 
00237 
00238 //-----------------------------------------------------------------------------------------------
00239 fitUtilities::fitUtilities():langausFit(0),gausFit(0){
00240   init();
00241 }
00242 
00243 //-----------------------------------------------------------------------------------------------
00244 void fitUtilities::init(){
00245   pLanGausS[0]=0 ;   pLanGausS[1]=0;   pLanGausS[2]=0;   pLanGausS[3]=0;
00246   epLanGausS[0]=0;  epLanGausS[1]=0;  epLanGausS[2]=0;  epLanGausS[3]=0;
00247   pGausS[0]=0 ;   pGausS[1]=0;   pGausS[2]=0;
00248   epGausS[0]=0;  epGausS[1]=0;  epGausS[2]=0;
00249   pLanConv[0]=0;  pLanConv[1]=0;
00250 
00251   //if ( langausFit!=0 ) delete langausFit;
00252   //if ( gausFit!=0 ) delete gausFit;
00253 
00254   chi2GausS = -9999.;
00255   nDofGausS = -9999;
00256 }
00257 
00258 //-----------------------------------------------------------------------------------------------
00259 fitUtilities::~fitUtilities(){
00260   if ( langausFit!=0 ) delete langausFit;
00261   if ( gausFit!=0 ) delete gausFit;
00262 }
00263 
00264   
00265  
00266 //-----------------------------------------------------------------------------------------------
00267 double fitUtilities::doLanGaussFit(TH1F* htoFit){
00268   init();
00269  
00270   if (htoFit->GetEntries()!=0) {
00271      edm::LogInfo("fitUtility")<<"Fitting "<< htoFit->GetTitle() <<std::endl;
00272     // Setting fit range and start values
00273     double fr[2];
00274     double sv[4], pllo[4], plhi[4];
00275     fr[0]=0.5*htoFit->GetMean();
00276     fr[1]=3.0*htoFit->GetMean();
00277               
00278     // (EM) parameters setting good for signal only 
00279     int32_t imax=htoFit->GetMaximumBin();
00280     double xmax=htoFit->GetBinCenter(imax);
00281     double ymax=htoFit->GetBinContent(imax);
00282     int32_t i[2];
00283     int32_t iArea[2];
00284               
00285     i[0]=htoFit->GetXaxis()->FindBin(fr[0]);
00286     i[1]=htoFit->GetXaxis()->FindBin(fr[1]);
00287                   
00288     iArea[0]=htoFit->GetXaxis()->FindBin(fr[0]);
00289     iArea[1]=htoFit->GetXaxis()->FindBin(fr[1]);
00290     double AreaFWHM=htoFit->Integral(iArea[0],iArea[1],"width");
00291                   
00292     sv[1]=xmax;
00293     sv[2]=htoFit->Integral(i[0],i[1],"width");
00294     sv[3]=AreaFWHM/(4*ymax);
00295     sv[0]=sv[3];
00296                   
00297     plhi[0]=25.0; plhi[1]=200.0; plhi[2]=1000000.0; plhi[3]=50.0;
00298     pllo[0]=1.5 ; pllo[1]=10.0 ; pllo[2]=1.0      ; pllo[3]= 1.0;
00299                   
00300     Char_t FunName[100];
00301     sprintf(FunName,"FitfcnLG_%s",htoFit->GetName());  
00302                   
00303     langausFit = new TF1(FunName,langaufun,fr[0],fr[1],4);
00304     langausFit->SetParameters(sv);
00305     langausFit->SetParNames("Width","MP","Area","GSigma");
00306                   
00307     for (int32_t i=0; i<4; i++) {
00308       langausFit->SetParLimits(i,pllo[i],plhi[i]);
00309     }  
00310                   
00311     try{
00312       htoFit->Fit(langausFit,"R0");  // "R" fit in a range,"0" quiet fit
00313       
00314       langausFit->SetRange(fr[0],fr[1]);
00315       langausFit->GetParameters(pLanGausS);
00316       std::memcpy((void*) epLanGausS, (void*) langausFit->GetParErrors(), 4*sizeof(double));
00317       
00318       chi2GausS =langausFit->GetChisquare();  // obtain chi^2
00319       nDofGausS = langausFit->GetNDF();           // obtain ndf
00320       
00321       double sPeak, sFWHM;
00322       langaupro(pLanGausS,sPeak,sFWHM);
00323       pLanConv[0]=sPeak;
00324       pLanConv[1]=sFWHM;
00325       edm::LogInfo("fitUtility") << "langaupro:  max  " << sPeak << std::endl;
00326       edm::LogInfo("fitUtility") << "langaupro:  FWHM " << sFWHM << std::endl;
00327     }
00328     catch(...){
00329       edm::LogError("fitUtility") << "problem in fitting " << htoFit->GetTitle() << " \n\tDefault values of the parameters will be used";
00330       pLanGausS[0]=-9999; pLanGausS[1]=-9999; pLanGausS[2]=-9999; pLanGausS[3]=-9999;
00331       epLanGausS[0]=-9999; epLanGausS[1]=-9999; epLanGausS[2]=-9999; epLanGausS[3]=-9999;
00332       pLanConv[0]=-9999;   pLanConv[1]=-9999;   
00333       chi2GausS=-9999;  nDofGausS=-9999;    
00334     }
00335   }
00336   else {  
00337     pLanGausS[0]=-9999; pLanGausS[1]=-9999; pLanGausS[2]=-9999; pLanGausS[3]=-9999;
00338     epLanGausS[0]=-9999; epLanGausS[1]=-9999; epLanGausS[2]=-9999; epLanGausS[3]=-9999;
00339     pLanConv[0]=-9999;   pLanConv[1]=-9999;   
00340     chi2GausS=-9999;  nDofGausS=-9999;    
00341   }
00342 
00343   return htoFit->GetEntries();
00344 
00345 }
00346 
00347   
00348 //-----------------------------------------------------------------------------------------------
00349 double fitUtilities::doGaussFit(TH1F* htoFit){
00350   init();
00351   if (htoFit->GetEntries()!=0) {
00352     
00353     // Setting fit range and start values
00354     double fr[2];
00355     double sv[3], pllo[3], plhi[3];
00356     fr[0]=htoFit->GetMean()-5*htoFit->GetRMS();
00357     fr[1]=htoFit->GetMean()+5*htoFit->GetRMS();
00358                   
00359     int32_t imax=htoFit->GetMaximumBin();
00360     double xmax=htoFit->GetBinCenter(imax);
00361     double ymax=htoFit->GetBinContent(imax);
00362     int32_t i[2];
00363     int32_t iArea[2];
00364                   
00365     i[0]=htoFit->GetXaxis()->FindBin(fr[0]);
00366     i[1]=htoFit->GetXaxis()->FindBin(fr[1]);
00367                   
00368     iArea[0]=htoFit->GetXaxis()->FindBin(fr[0]);
00369     iArea[1]=htoFit->GetXaxis()->FindBin(fr[1]);
00370     double AreaFWHM=htoFit->Integral(iArea[0],iArea[1],"width");
00371                   
00372     sv[2]=AreaFWHM/(4*ymax);
00373     sv[1]=xmax;
00374     sv[0]=htoFit->Integral(i[0],i[1],"width");
00375                   
00376     plhi[0]=1000000.0; plhi[1]=10.0; plhi[2]=10.;
00377     pllo[0]=1.5 ; pllo[1]=0.1; pllo[2]=0.3;
00378     Char_t FunName[100];
00379     sprintf(FunName,"FitfcnLG_%s",htoFit->GetName());
00380     gausFit = new TF1(FunName,Gauss,fr[0],fr[1],3);
00381     gausFit->SetParameters(sv);
00382     gausFit->SetParNames("Constant","GaussPeak","Sigma");
00383                   
00384     for (int32_t i=0; i<3; i++) {
00385       gausFit->SetParLimits(i,pllo[i],plhi[i]);
00386     }
00387  
00388     try{
00389       htoFit->Fit(gausFit,"R0");
00390       
00391       gausFit->SetRange(fr[0],fr[1]);
00392       gausFit->GetParameters(pGausS);
00393       std::memcpy((void*) epGausS, (void*) gausFit->GetParErrors(), 3*sizeof(double));
00394       
00395       chi2GausS =langausFit->GetChisquare(); // obtain chi^2
00396       nDofGausS = langausFit->GetNDF();// obtain ndf
00397     }
00398     catch(...){
00399       edm::LogError("fitUtility") << "problem in fitting " << htoFit->GetTitle() << " \n\tDefault values of the parameters will be used";
00400       pGausS[0]=-9999; pGausS[1]=-9999; pGausS[2]=-9999;
00401       epGausS[0]=-9999; epGausS[1]=-9999; epGausS[2]=-9999;
00402       chi2GausS=-9999;  nDofGausS=-9999;    
00403     }
00404                  
00405   }
00406   else {
00407     pGausS[0]=-9999; pGausS[1]=-9999; pGausS[2]=-9999;
00408     epGausS[0]=-9999; epGausS[1]=-9999; epGausS[2]=-9999;
00409     chi2GausS=-9999;  nDofGausS=-9999;    
00410   }
00411   
00412   return htoFit->GetEntries();
00413 }
00414 
00415 //-----------------------------------------------------------------------------------------------
00416 double fitUtilities::getLanGaussPar(std::string s){
00417   if(s=="landau_width")
00418     return pLanGausS[0];
00419   else if(s=="mpv")
00420     return pLanGausS[1];
00421   else if(s=="area")
00422     return pLanGausS[2];
00423   else if(s=="gauss_sigma")
00424     return pLanGausS[3];
00425   else
00426     return -99999;
00427 }
00428 
00429 double fitUtilities::getLanGaussParErr(std::string s){
00430   if(s=="landau_width")
00431     return epLanGausS[0];
00432   else if(s=="mpv")
00433     return epLanGausS[1];
00434   else if(s=="area")
00435     return epLanGausS[2];
00436   else if(s=="gauss_sigma")
00437     return epLanGausS[3];
00438   else
00439     return -99999;
00440 }
00441 
00442 double fitUtilities::getLanGaussConv(std::string s) {
00443   if(s=="mpv")
00444     return pLanConv[0];
00445   else if(s=="fwhm")
00446     return pLanConv[1];
00447   else
00448     return -99999;
00449 }
00450 
00451 double fitUtilities::getGaussPar(std::string s) {
00452   if(s=="area")
00453     return pGausS[0];
00454   else if(s=="mean")
00455     return pGausS[1];
00456   else if(s=="sigma")
00457     return pGausS[2];
00458   else
00459     return -99999;
00460 }
00461 
00462 double fitUtilities::getGaussParErr(std::string s) {
00463   if(s=="area")
00464     return epGausS[0];
00465   else if(s=="mean")
00466     return epGausS[1];
00467   else if(s=="sigma")
00468     return epGausS[2];
00469   else
00470     return -99999;
00471 }