-
Notifications
You must be signed in to change notification settings - Fork 1
/
FittingHeader.h
312 lines (254 loc) · 11.5 KB
/
FittingHeader.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
/***********************************************************************************************
*** provided by, ******
*** Friederike Bock, [email protected] ******
************************************************************************************************/
#ifndef FITTINGGENERAL
#define FITTINGGENERAL
Int_t nPeaksMultGauss = 0;
//__________________________________________________________________________________________________________
// multiple Gauss fit
//__________________________________________________________________________________________________________
Double_t multGauss(Double_t * x, Double_t * par){
Double_t result = par[0] + par[1] * x[0];
for (Int_t p = 0; p < nPeaksMultGauss; p++){
Double_t norm = par[3 * p +2]; // "height" or "area"
Double_t mean = par[3 * p + 3];
Double_t sigma = par[3 * p + 4];
#if defined(__PEAKS_C_FIT_AREAS__)
norm /= sigma * (TMath::Sqrt(TMath::TwoPi())); // "area"
#endif /* defined(__PEAKS_C_FIT_AREAS__) */
result += norm * TMath::Gaus(x[0], mean, sigma);
}
return result;
}
//__________________________________________________________________________________________________________
// Fit noise nicely
//__________________________________________________________________________________________________________
Bool_t FitNoise (TH1D* histo, TF1* &fit, Double_t &mean, Double_t &meanErr, Double_t &sigma, Double_t &sigmaErr, Int_t cb, Int_t cc, TString baseName, TString nameGain, Int_t verbosity = 0){
if (!histo) return kFALSE;
if (!(histo->GetEntries() > 0)) return kFALSE;
fit = new TF1(Form("%s_B%d_C%02d",baseName.Data(),cb,cc), "gaus", 0, 100);
// preset parames
fit->SetParameter(0, histo->GetMaximum()); // Amplitude
fit->SetParameter(1, histo->GetMean()); // Mean
fit->SetParameter(2, histo->GetStdDev()); // Sigma
histo->Fit(fit,"QRMNE0");
if (verbosity > 0) std::cout << "Channel b: " << cb << "\t c: " << cc << "\t" << nameGain.Data() << " mean: " << fit->GetParameter(1) << "\t width: " << fit->GetParameter(2) << std::endl;
mean = fit->GetParameter(1);
meanErr = fit->GetParError(1);
sigma = fit->GetParameter(2);
sigmaErr = fit->GetParError(2);
return kTRUE;
}
//__________________________________________________________________________________________________________
// Fit noise nicely
//__________________________________________________________________________________________________________
Bool_t FitNoiseWithBG (TH1D* histo, TF1* &fit, Double_t &mean, Double_t &meanErr, Double_t &sigma, Double_t &sigmaErr, Int_t cb, Int_t cc, TString baseName, TString nameGain, Int_t verbosity = 0){
if (!histo) return kFALSE;
if (!(histo->GetEntries() > 0)) return kFALSE;
fit = new TF1(Form("%s_B%d_C%02d",baseName.Data(),cb,cc), "gaus", 0, 100);
// preset parames
Double_t largestContent = 0;
Int_t minBin = histo->GetXaxis()->FindBin(0.0);
Int_t maxBin = histo->GetXaxis()->FindBin(100)+0.0001;
for (Int_t i= minBin; i < maxBin; i++){
if (largestContent < histo->GetBinContent(i)){
largestContent = histo->GetBinContent(i);
}
}
fit->SetParameter(0, 0.8* largestContent); // Amplitude
fit->SetParLimits(0, 0.1, largestContent); // Amplitude
fit->SetParameter(1, 50); // Mean
fit->SetParLimits(1, 0, 100); // Mean
fit->SetParameter(2, 20); // Sigma
fit->SetParLimits(2, 0, 80); // Sigma
histo->Fit(fit,"QRMNE0");
if (verbosity > 0) std::cout << "Channel b: " << cb << "\t c: " << cc << "\t" << nameGain.Data() << " mean: " << fit->GetParameter(1) << "\t width: " << fit->GetParameter(2) << std::endl;
histo->Fit(fit,"QRMNE0","", fit->GetParameter(1)-3*fit->GetParameter(2), fit->GetParameter(1)+0.05*fit->GetParameter(2));
if (verbosity > 0) std::cout << "2nd try Channel b: " << cb << "\t c: " << cc << "\t" << nameGain.Data() << " mean: " << fit->GetParameter(1) << "\t width: " << fit->GetParameter(2) << std::endl;
mean = fit->GetParameter(1);
meanErr = fit->GetParError(1);
sigma = fit->GetParameter(2);
sigmaErr = fit->GetParError(2);
return kTRUE;
}
//__________________________________________________________________________________________________________
// Plot & Fit gain Corr
//__________________________________________________________________________________________________________
void FitAndPlotGainCorr (TH2D* hist, TF1* &fit, TString baseNameFit, Double_t minFX, Double_t maxFX, Double_t maxPX, Double_t maxPY,
Double_t &slope, Double_t &slopeErr, Double_t &offset, Double_t &offsetErr,
Int_t cb, Int_t cc, Int_t sl, Int_t sc,
Bool_t bDetailed, TCanvas* canvas, TString baseNameOut, Double_t textSizeRel = 0.04 ){
fit = new TF1(Form("%s_B%d_C%02d",baseNameFit.Data(), cb,cc), "[0]+[1]*x", minFX, maxFX);
hist->Fit(fit,"QRMNE0");
slope = fit->GetParameter(1);
slopeErr = fit->GetParError(1);
offset = fit->GetParameter(0);
offsetErr = fit->GetParError(0);
if (!bDetailed) return;
canvas->cd();
SetStyleHistoTH2ForGraphs( hist, hist->GetXaxis()->GetTitle(), hist->GetYaxis()->GetTitle(), 0.85*textSizeRel, textSizeRel, 0.85*textSizeRel, textSizeRel,0.9, 1.25);
SetStyleFit(fit , 0, maxPX, 7, 7, kBlack);
hist->GetXaxis()->SetRangeUser(0,maxPX);
hist->GetYaxis()->SetRangeUser(0,maxPY);
hist->Draw("colz");
fit->Draw("same");
TLegend* legend = GetAndSetLegend2( 0.42, 0.15, 0.95, 0.3,textSizeRel, 1, Form("CAEN B %d, C %d, Stack L %d, C%d",cb, cc, sl, sc), 42,0.2);
legend->AddEntry(fit, "f(x) = a#upoint x + b", "l");
legend->AddEntry((TObject*)0, Form("a = %2.2f #pm %2.2f",slope, slopeErr ) , " ");
legend->AddEntry((TObject*)0, Form("b = %2.2f #pm %2.2f",offset, offsetErr ) , " ");
legend->Draw();
canvas->SaveAs(Form("%s_B%d_C%02d.pdf", baseNameOut.Data(), cb,cc));
return;
}
double langaufun(double *x, double *par) {
//Fit parameters:
//par[0]=Width (scale) parameter of Landau density
//par[1]=Most Probable (MP, location) parameter of Landau density
//par[2]=Total area (integral -inf to inf, normalization constant)
//par[3]=Width (sigma) of convoluted Gaussian function
//
//In the Landau distribution (represented by the CERNLIB approximation),
//the maximum is located at x=-0.22278298 with the location parameter=0.
//This shift is corrected within this function, so that the actual
//maximum is identical to the MP parameter.
// Numeric constants
double invsq2pi = 0.3989422804014; // (2 pi)^(-1/2)
double mpshift = -0.22278298; // Landau maximum location
// Control constants
double np = 100.0; // number of convolution steps
double sc = 5.0; // convolution extends to +-sc Gaussian sigmas
// Variables
double xx;
double mpc;
double fland;
double sum = 0.0;
double xlow,xupp;
double step;
double i;
// MP shift correction
mpc = par[1] - mpshift * par[0];
// Range of convolution integral
xlow = x[0] - sc * par[3];
xupp = x[0] + sc * par[3];
step = (xupp-xlow) / np;
// Convolution integral of Landau and Gaussian by sum
for(i=1.0; i<=np/2; i++) {
xx = xlow + (i-.5) * step;
fland = TMath::Landau(xx,mpc,par[0]) / par[0];
sum += fland * TMath::Gaus(x[0],xx,par[3]);
xx = xupp - (i-.5) * step;
fland = TMath::Landau(xx,mpc,par[0]) / par[0];
sum += fland * TMath::Gaus(x[0],xx,par[3]);
}
return (par[2] * step * sum * invsq2pi / par[3]);
}
TF1 *langaufit(TH1D *his, double *fitrange, double *startvalues, double *parlimitslo, double *parlimitshi, double *fitparams, double *fiterrors, double *ChiSqr, int *NDF, Int_t verbosity = 0)
{
// Once again, here are the Landau * Gaussian parameters:
// par[0]=Width (scale) parameter of Landau density
// par[1]=Most Probable (MP, location) parameter of Landau density
// par[2]=Total area (integral -inf to inf, normalization constant)
// par[3]=Width (sigma) of convoluted Gaussian function
//
// Variables for langaufit call:
// his histogram to fit
// fitrange[2] lo and hi boundaries of fit range
// startvalues[4] reasonable start values for the fit
// parlimitslo[4] lower parameter limits
// parlimitshi[4] upper parameter limits
// fitparams[4] returns the final fit parameters
// fiterrors[4] returns the final fit errors
// ChiSqr returns the chi square
// NDF returns ndf
int i;
char FunName[100];
sprintf(FunName,"Fitfcn_%s",his->GetName());
TF1 *ffitold = (TF1*)gROOT->GetListOfFunctions()->FindObject(FunName);
if (ffitold) delete ffitold;
TF1 *ffit = new TF1(FunName,langaufun,fitrange[0],fitrange[1],4);
ffit->SetParameters(startvalues);
ffit->SetParNames("Width","MP","Area","GSigma");
for (i=0; i<4; i++) {
ffit->SetParLimits(i, parlimitslo[i], parlimitshi[i]);
}
TString fitOption = "QRLMNE0";
if (verbosity > 1)
fitOption = "RLMNE0";
his->Fit(FunName,fitOption); // fit within specified range, use ParLimits, do not plot
ffit->GetParameters(fitparams); // obtain fit parameters
for (i=0; i<4; i++) {
fiterrors[i] = ffit->GetParError(i); // obtain fit parameter errors
}
ChiSqr[0] = ffit->GetChisquare(); // obtain chi^2
NDF[0] = ffit->GetNDF(); // obtain ndf
return (ffit); // return fit function
}
int langaupro(double *params, double &maxx, double &FWHM) {
// Searches for the location (x value) at the maximum of the
// Landau-Gaussian convolute and its full width at half-maximum.
//
// The search is probably not very efficient, but it's a first try.
double p,x,fy,fxr,fxl;
double step;
double l,lold;
int i = 0;
int MAXCALLS = 10000;
// Search for maximum
p = params[1] - 0.1 * params[0];
step = 0.05 * params[0];
lold = -2.0;
l = -1.0;
while ( (l != lold) && (i < MAXCALLS) ) {
i++;
lold = l;
x = p + step;
l = langaufun(&x,params);
if (l < lold)
step = -step/10;
p += step;
}
if (i == MAXCALLS)
return (-1);
maxx = x;
fy = l/2;
// Search for right x location of fy
p = maxx + params[0];
step = params[0];
lold = -2.0;
l = -1e300;
i = 0;
while ( (l != lold) && (i < MAXCALLS) ) {
i++;
lold = l;
x = p + step;
l = TMath::Abs(langaufun(&x,params) - fy);
if (l > lold)
step = -step/10;
p += step;
}
if (i == MAXCALLS)
return (-2);
fxr = x;
// Search for left x location of fy
p = maxx - 0.5 * params[0];
step = -params[0];
lold = -2.0;
l = -1e300;
i = 0;
while ( (l != lold) && (i < MAXCALLS) ) {
i++;
lold = l;
x = p + step;
l = TMath::Abs(langaufun(&x,params) - fy);
if (l > lold)
step = -step/10;
p += step;
}
if (i == MAXCALLS)
return (-3);
fxl = x;
FWHM = fxr - fxl;
return (0);
}
#endif