33 Double_t langaufun(Double_t *x, Double_t *par) {
47 Double_t invsq2pi = 0.3989422804014;
48 Double_t mpshift = -0.22278298;
65 mpc = par[1] - mpshift * par[0];
68 xlow = x[0] - sc * par[3];
69 xupp = x[0] + sc * par[3];
71 step = (xupp-xlow) / np;
74 for(i=1.0; i<=np/2; i++) {
75 xx = xlow + (i-.5) * step;
76 fland = TMath::Landau(xx,mpc,par[0]) / par[0];
77 sum += fland * TMath::Gaus(x[0],xx,par[3]);
79 xx = xupp - (i-.5) * step;
80 fland = TMath::Landau(xx,mpc,par[0]) / par[0];
81 sum += fland * TMath::Gaus(x[0],xx,par[3]);
84 return (par[2] * step * sum * invsq2pi / par[3]);
89 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)
111 sprintf(FunName,
"Fitfcn_%s",his->GetName());
113 TF1 *ffitold = (TF1*)gROOT->GetListOfFunctions()->FindObject(FunName);
114 if (ffitold)
delete ffitold;
116 TF1 *ffit =
new TF1(FunName,langaufun,fitrange[0],fitrange[1],4);
117 ffit->SetParameters(startvalues);
118 ffit->SetParNames(
"Width",
"MP",
"Area",
"GSigma");
120 for (i=0; i<4; i++) {
121 ffit->SetParLimits(i, parlimitslo[i], parlimitshi[i]);
124 his->Fit(FunName,
"RB0");
126 ffit->GetParameters(fitparams);
127 for (i=0; i<4; i++) {
128 fiterrors[i] = ffit->GetParError(i);
130 ChiSqr[0] = ffit->GetChisquare();
131 NDF[0] = ffit->GetNDF();
138 Int_t langaupro(Double_t *params, Double_t &maxx, Double_t &FWHM) {
145 Double_t p,x,fy,fxr,fxl;
149 Int_t MAXCALLS = 10000;
154 p = params[1] - 0.1 * params[0];
155 step = 0.05 * params[0];
160 while ( (l != lold) && (i < MAXCALLS) ) {
165 l = langaufun(&x,params);
183 p = maxx + params[0];
190 while ( (l != lold) && (i < MAXCALLS) ) {
195 l = TMath::Abs(langaufun(&x,params) - fy);
211 p = maxx - 0.5 * params[0];
217 while ( (l != lold) && (i < MAXCALLS) ) {
222 l = TMath::Abs(langaufun(&x,params) - fy);
242 Int_t data[100] = {0,0,0,0,0,0,2,6,11,18,18,55,90,141,255,323,454,563,681,
243 737,821,796,832,720,637,558,519,460,357,291,279,241,212,
244 153,164,139,106,95,91,76,80,80,59,58,51,30,49,23,35,28,23,
245 22,27,27,24,20,16,17,14,20,12,12,13,10,17,7,6,12,6,12,4,
246 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};
247 TH1F *hSNR =
new TH1F(
"snr",
"Signal-to-noise",400,0,400);
249 for (Int_t i=0; i<100; i++) hSNR->Fill(i,data[i]);
252 printf(
"Fitting...\n");
256 Double_t sv[4], pllo[4], plhi[4], fp[4], fpe[4];
257 fr[0]=0.3*hSNR->GetMean();
258 fr[1]=3.0*hSNR->GetMean();
260 pllo[0]=0.5; pllo[1]=5.0; pllo[2]=1.0; pllo[3]=0.4;
261 plhi[0]=5.0; plhi[1]=50.0; plhi[2]=1000000.0; plhi[3]=5.0;
262 sv[0]=1.8; sv[1]=20.0; sv[2]=50000.0; sv[3]=3.0;
266 TF1 *fitsnr = langaufit(hSNR,fr,sv,pllo,plhi,fp,fpe,&chisqr,&ndf);
268 Double_t SNRPeak, SNRFWHM;
269 langaupro(fp,SNRPeak,SNRFWHM);
271 printf(
"Fitting done\nPlotting results...\n");
274 gStyle->SetOptStat(1111);
275 gStyle->SetOptFit(111);
276 gStyle->SetLabelSize(0.03,
"x");
277 gStyle->SetLabelSize(0.03,
"y");
279 hSNR->GetXaxis()->SetRange(0,70);
281 fitsnr->Draw(
"lsame");