00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #include "TH1.h"
00018 #include "TF1.h"
00019 #include "TROOT.h"
00020 #include "TStyle.h"
00021 #include "TMath.h"
00022
00023 Double_t langaufun(Double_t *x, Double_t *par) {
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037 Double_t invsq2pi = 0.3989422804014;
00038 Double_t mpshift = -0.22278298;
00039
00040
00041 Double_t np = 100.0;
00042 Double_t sc = 5.0;
00043
00044
00045 Double_t xx;
00046 Double_t mpc;
00047 Double_t fland;
00048 Double_t sum = 0.0;
00049 Double_t xlow,xupp;
00050 Double_t step;
00051 Double_t i;
00052
00053
00054
00055 mpc = par[1] - mpshift * par[0];
00056
00057
00058 xlow = x[0] - sc * par[3];
00059 xupp = x[0] + sc * par[3];
00060
00061 step = (xupp-xlow) / np;
00062
00063
00064 for(i=1.0; i<=np/2; i++) {
00065 xx = xlow + (i-.5) * step;
00066 fland = TMath::Landau(xx,mpc,par[0]) / par[0];
00067 sum += fland * TMath::Gaus(x[0],xx,par[3]);
00068
00069 xx = xupp - (i-.5) * step;
00070 fland = TMath::Landau(xx,mpc,par[0]) / par[0];
00071 sum += fland * TMath::Gaus(x[0],xx,par[3]);
00072 }
00073
00074 return (par[2] * step * sum * invsq2pi / par[3]);
00075 }
00076
00077
00078
00079 TF1 *langaufit(TH1F *his, Double_t *fitrange, Double_t *startvalues, Double_t *parlimitslo, Double_t *parlimitshi, Double_t *fitparams, Double_t *fiterrors, Double_t *ChiSqr, Int_t *NDF)
00080 {
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098 Int_t i;
00099 Char_t FunName[100];
00100
00101 sprintf(FunName,"Fitfcn_%s",his->GetName());
00102
00103 TF1 *ffitold = (TF1*)gROOT->GetListOfFunctions()->FindObject(FunName);
00104 if (ffitold) delete ffitold;
00105
00106 TF1 *ffit = new TF1(FunName,langaufun,fitrange[0],fitrange[1],4);
00107 ffit->SetParameters(startvalues);
00108 ffit->SetParNames("Width","MP","Area","GSigma");
00109
00110 for (i=0; i<4; i++) {
00111 ffit->SetParLimits(i, parlimitslo[i], parlimitshi[i]);
00112 }
00113
00114 his->Fit(FunName,"RB0");
00115
00116 ffit->GetParameters(fitparams);
00117 for (i=0; i<4; i++) {
00118 fiterrors[i] = ffit->GetParError(i);
00119 }
00120 ChiSqr[0] = ffit->GetChisquare();
00121 NDF[0] = ffit->GetNDF();
00122
00123 return (ffit);
00124
00125 }
00126
00127
00128 Int_t langaupro(Double_t *params, Double_t &maxx, Double_t &FWHM) {
00129
00130
00131
00132
00133
00134
00135 Double_t p,x,fy,fxr,fxl;
00136 Double_t step;
00137 Double_t l,lold;
00138 Int_t i = 0;
00139 Int_t MAXCALLS = 10000;
00140
00141
00142
00143
00144 p = params[1] - 0.1 * params[0];
00145 step = 0.05 * params[0];
00146 lold = -2.0;
00147 l = -1.0;
00148
00149
00150 while ( (l != lold) && (i < MAXCALLS) ) {
00151 i++;
00152
00153 lold = l;
00154 x = p + step;
00155 l = langaufun(&x,params);
00156
00157 if (l < lold)
00158 step = -step/10;
00159
00160 p += step;
00161 }
00162
00163 if (i == MAXCALLS)
00164 return (-1);
00165
00166 maxx = x;
00167
00168 fy = l/2;
00169
00170
00171
00172
00173 p = maxx + params[0];
00174 step = params[0];
00175 lold = -2.0;
00176 l = -1e300;
00177 i = 0;
00178
00179
00180 while ( (l != lold) && (i < MAXCALLS) ) {
00181 i++;
00182
00183 lold = l;
00184 x = p + step;
00185 l = TMath::Abs(langaufun(&x,params) - fy);
00186
00187 if (l > lold)
00188 step = -step/10;
00189
00190 p += step;
00191 }
00192
00193 if (i == MAXCALLS)
00194 return (-2);
00195
00196 fxr = x;
00197
00198
00199
00200
00201 p = maxx - 0.5 * params[0];
00202 step = -params[0];
00203 lold = -2.0;
00204 l = -1e300;
00205 i = 0;
00206
00207 while ( (l != lold) && (i < MAXCALLS) ) {
00208 i++;
00209
00210 lold = l;
00211 x = p + step;
00212 l = TMath::Abs(langaufun(&x,params) - fy);
00213
00214 if (l > lold)
00215 step = -step/10;
00216
00217 p += step;
00218 }
00219
00220 if (i == MAXCALLS)
00221 return (-3);
00222
00223
00224 fxl = x;
00225
00226 FWHM = fxr - fxl;
00227 return (0);
00228 }
00229
00230 void langaus() {
00231
00232 Int_t data[100] = {0,0,0,0,0,0,2,6,11,18,18,55,90,141,255,323,454,563,681,
00233 737,821,796,832,720,637,558,519,460,357,291,279,241,212,
00234 153,164,139,106,95,91,76,80,80,59,58,51,30,49,23,35,28,23,
00235 22,27,27,24,20,16,17,14,20,12,12,13,10,17,7,6,12,6,12,4,
00236 9,9,10,3,4,5,2,4,1,5,5,1,7,1,6,3,3,3,4,5,4,4,2,2,7,2,4};
00237 TH1F *hSNR = new TH1F("snr","Signal-to-noise",400,0,400);
00238
00239 for (Int_t i=0; i<100; i++) hSNR->Fill(i,data[i]);
00240
00241
00242 printf("Fitting...\n");
00243
00244
00245 Double_t fr[2];
00246 Double_t sv[4], pllo[4], plhi[4], fp[4], fpe[4];
00247 fr[0]=0.3*hSNR->GetMean();
00248 fr[1]=3.0*hSNR->GetMean();
00249
00250 pllo[0]=0.5; pllo[1]=5.0; pllo[2]=1.0; pllo[3]=0.4;
00251 plhi[0]=5.0; plhi[1]=50.0; plhi[2]=1000000.0; plhi[3]=5.0;
00252 sv[0]=1.8; sv[1]=20.0; sv[2]=50000.0; sv[3]=3.0;
00253
00254 Double_t chisqr;
00255 Int_t ndf;
00256 TF1 *fitsnr = langaufit(hSNR,fr,sv,pllo,plhi,fp,fpe,&chisqr,&ndf);
00257
00258 Double_t SNRPeak, SNRFWHM;
00259 langaupro(fp,SNRPeak,SNRFWHM);
00260
00261 printf("Fitting done\nPlotting results...\n");
00262
00263
00264 gStyle->SetOptStat(1111);
00265 gStyle->SetOptFit(111);
00266 gStyle->SetLabelSize(0.03,"x");
00267 gStyle->SetLabelSize(0.03,"y");
00268
00269 hSNR->GetXaxis()->SetRange(0,70);
00270 hSNR->Draw();
00271 fitsnr->Draw("lsame");
00272 }
00273