ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/cbrown/AnalysisFramework/Plotting/Modules/PeakFinder.C
Revision: 1.1
Committed: Wed Jun 22 11:07:37 2011 UTC (13 years, 10 months ago) by buchmann
Content type: text/plain
Branch: MAIN
Log Message:
Initial commit of Plotting tools

File Contents

# User Rev Content
1 buchmann 1.1 #include <iostream>
2     #include <sstream>
3     #include <iomanip>
4     #include <TFile.h>
5     #include <TTree.h>
6     #include <TH1.h>
7     #include <TF1.h>
8     #include <TMath.h>
9     #include <TCanvas.h>
10     #include <vector>
11     #include <TROOT.h>
12     #include <TLine.h>
13     #include <TLegend.h>
14     #include <TLatex.h>
15     #include <TRandom.h>
16     #ifndef GeneralToolBoxLoaded
17     #include "GeneralToolBox.C"
18     #endif
19     #ifndef Verbosity
20     #define Verbosity 0
21     #endif
22    
23     using namespace std;
24    
25     Double_t LogParabola(Double_t *x,Double_t *par)
26     {
27     return par[0]*TMath::Exp(-par[1]*(x[0]-par[2])*(x[0]-par[2])); // we're adding a "logarithmic parabola" :-)
28     //note: the abs() around the first parameter ensures that, when fitting, no negative values are chosen.
29     }
30    
31     bool doreject=false;
32     float low_reject=-10;
33     float hi_reject=10;
34    
35     bool dofixed=true;
36    
37    
38     bool addparabola=true;
39     float parabola_height=0;
40     float parabola_inclination=0;
41     float parabola_pointzero=0;
42    
43     float find_KM_peak(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, float &error,float &Sigma);
44     float find_Gauss_peak(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, float &error,float &Sigma,int numsig);
45    
46    
47     Double_t KrystalMallLogPar(double *x, double *par)
48     {
49     //parameters:
50     //N: the way we scale the function
51     //alpha (where the function changes)
52     //n: exponent of the power expression in the left area
53     //xbar: peak of the gaussian part (RHS)
54     //sigma: width of the gaussian part (RHS)
55     float N=par[0];
56     float alpha=par[1];
57     float n=par[2];
58     float xbar=par[3];
59     float sigma=par[4];
60     float altX=x[0];
61     float result=-999;
62     if(doreject&&x[0]>low_reject&&x[0]<hi_reject)
63     {
64     TF1::RejectPoint();
65     return 0;
66     }
67     if(((altX-xbar)/sigma>-alpha)&&((altX-xbar)/sigma<alpha)){
68     result=N*TMath::Exp(-(altX-xbar)*(altX-xbar)/(2*sigma*sigma));
69     }
70     else
71     {
72     //if we are outside the central (Gaussian) area things become more difficult ...
73     float A=TMath::Power(n/TMath::Abs(alpha),n)*TMath::Exp(-alpha*alpha/2);
74     float B=n/TMath::Abs(alpha) - TMath::Abs(alpha);
75     if((altX-xbar)/sigma<=-alpha) result=N*A*TMath::Power((B-((altX-xbar)/sigma)),-n);
76     if((altX-xbar)/sigma>=alpha) result=N*A*TMath::Power((B+((altX-xbar)/sigma)),-n);
77     }
78    
79     if(addparabola) {
80     if(dofixed) {
81     result+=parabola_height*TMath::Exp(-parabola_inclination*(x[0]-parabola_pointzero)*(x[0]-parabola_pointzero)); // we're adding a "logarithmic parabola" :-)
82     }
83     else {
84     result+=par[5]*TMath::Exp(-par[6]*(x[0]-par[7])*(x[0]-par[7])); // we're adding a "logarithmic parabola" :-)
85     if(par[5]<0) return -999; // there can be no negative ttbar contribution, so just return a value which is going to be a horrible fit.
86     if(par[6]<0) return -999; // the parabola needs to close (i.e. tend to negative values for large |jzb|, not to large positive values)
87     }
88     }
89     return result;
90     }
91    
92    
93     void do_ttbar_fit(TH1F *ttbar,TF1 *logpar, TF1 *KM)
94     {
95     logpar->SetParameters(10,2,3);
96     ttbar->Fit(logpar,"NQ");
97     ttbar->Fit(logpar,"NQ");
98     ttbar->Fit(logpar,"NQ");
99     ttbar->Fit(logpar,"NQ");
100     ttbar->SetStats(0);
101     parabola_height=logpar->GetParameter(0);
102     parabola_inclination=logpar->GetParameter(1);
103     parabola_pointzero=logpar->GetParameter(2);
104     }
105    
106     void draw_complete_fit(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, TF1 *KM)
107     {
108     TCanvas *fitsummary;
109     if(is_data) {
110     fitsummary= new TCanvas("fitsummary","Fit Summary",1000,500);
111     fitsummary->Divide(2,1);
112     }
113     else {
114     fitsummary= new TCanvas("fitsummary","Fit Summary",1200,400);
115     fitsummary->Divide(3,1);
116     }
117     TF1 *logpar = new TF1("logpar",LogParabola,minfit,maxfit,3);
118    
119     logpar->SetParameters(KM->GetParameter(5),KM->GetParameter(6),KM->GetParameter(7));
120     logpar->SetLineColor(kOrange);
121     logpar->SetLineStyle(2);
122     if(!is_data)
123     {
124     ttbar->GetXaxis()->SetTitle("JZB (GeV/c)");
125     ttbar->GetYaxis()->SetTitle("events");
126     ttbar->GetXaxis()->CenterTitle();
127     ttbar->GetYaxis()->CenterTitle();
128     ttbar->SetLineColor(kRed);
129     fitsummary->cd(1);
130     ttbar->Draw();
131     fitsummary->cd(1);
132     logpar->Draw("same");
133     TLegend *leg = new TLegend(0.3,0.25,0.65,0.4);
134     leg->AddEntry(ttbar,"t#bar{t} (mc)","l");
135     leg->AddEntry(logpar,"Fit with Log. Parabola","l");
136     leg->SetLineColor(kWhite);
137     leg->SetFillColor(kWhite);
138     leg->Draw();
139     TText *title1=write_title("t#bar{t} Distribution and Fit");
140     title1->Draw();
141     }
142     fitsummary->cd(2-int(is_data));
143     fitsummary->cd(2-int(is_data))->SetLogy(1);
144     all->GetYaxis()->SetTitle("events");
145     all->GetYaxis()->CenterTitle();
146     all->Draw();
147     ttbar->SetLineColor(kRed);
148     if(!is_data) ttbar->Draw("same");
149     KM->SetLineWidth(1);
150     KM->Draw("same");
151     logpar->SetLineWidth(1);
152     logpar->Draw("same");
153     if(!is_data)ttbar->Draw("same");
154     TLegend *leg2 = new TLegend(0.65,0.65,0.89,0.89);
155     if(is_data) leg2->AddEntry(all,"Data","l");
156     else leg2->AddEntry(all,"Stacked MC","l");
157     leg2->AddEntry(KM,"Fitted KM Function","l");
158     if(!is_data) leg2->AddEntry(ttbar,"t#bar{t} MC","l");
159     leg2->AddEntry(logpar,"t#bar{t} (Fit)","l");
160     leg2->SetFillColor(kWhite);
161     leg2->SetLineColor(kWhite);
162     leg2->Draw();
163     TText *title2=write_title("Distribution and Fits (log.)");
164     title2->Draw();
165     fitsummary->cd(3-is_data);
166     all->Draw();
167     KM->Draw("same");
168     float peaklocation=KM->GetParameter(3);
169     TLine *muline = new TLine(peaklocation,0,peaklocation,all->GetMaximum());
170     muline->SetLineColor(kBlue);
171     muline->SetLineStyle(2);
172     muline->Draw();
173     TLegend *leg = new TLegend(0.6,0.6,0.89,0.89);
174     if(is_data) leg2->AddEntry(all,"Data","l");
175     else leg->AddEntry(all,"Stacked MC","l");
176     leg->AddEntry(KM,"Fitted KM Function","l");
177     stringstream mulinelabel;
178     mulinelabel<<"Peak position at #mu="<<peaklocation;
179     leg->AddEntry(muline,mulinelabel.str().c_str(),"l");
180     leg->SetLineColor(kWhite);
181     leg->SetFillColor(kWhite);
182     leg->Draw();
183     mulinelabel<<"+/-"<<TMath::Abs(KM->GetParError(3));
184     TText *title3=write_title("Distribution and Fits");
185     title3->Draw();
186     TText *titlel=write_title_low(mulinelabel.str().c_str());
187     titlel->Draw();
188    
189     stringstream printtop;
190     printtop << "#mu="<<std::setprecision(3)<<KM->GetParameter(3)<<"+/-"<<std::setprecision(3)<<KM->GetParError(3);
191     TLatex *toptext = new TLatex(0,all->GetMaximum()*1.3,printtop.str().c_str());
192     toptext->SetTextAlign(22);
193     // toptext->Draw();
194    
195     doreject=false;
196     TF1 *wholefitfunc=(TF1*)KM->Clone();
197     doreject=true;
198     wholefitfunc->SetLineColor(kRed);
199     wholefitfunc->SetLineStyle(2);
200     wholefitfunc->Draw("same");
201    
202     fitsummary->cd(2-is_data);
203     wholefitfunc->Draw("same");
204    
205     if(is_data) CompleteSave(fitsummary, "fit/Fit_Summary_Data");
206     else CompleteSave(fitsummary,"fit/Fit_Summary_MC");
207    
208     }
209    
210     float Kostas_algorithm(TH1F *hist, float &error, float &sigma, TF1* fitFunc, float lowlimit, float highlimit,bool is_data)
211     {
212     float mean = hist->GetBinCenter( hist->GetMaximumBin());
213     float rms = hist->GetRMS();
214    
215     mean = hist->GetBinCenter( hist->GetMaximumBin());
216    
217     fitFunc->SetParameter(1,mean);
218    
219     hist->Fit(fitFunc,"QLL0","",mean-10,mean+10);
220    
221     mean=fitFunc->GetParameter(1);
222     rms=fitFunc->GetParameter(2);
223     error=fitFunc->GetParError(1);
224    
225     bool printOut = false; // print the peak estimate in the i-th iteration
226    
227     // --- perform iterations
228     int numIterations=5;
229    
230     if(printOut) std::cout << " ( ";
231     for(int i=1;i<numIterations+1;i++) //--- modify the number of iterations until peak is stable
232     {
233     hist->Fit(fitFunc,"QLLN","same",mean - lowlimit*rms, mean + highlimit*rms); // fit -2 +1 sigma from previous iteration
234     mean=fitFunc->GetParameter(1);
235     fitFunc->SetRange(mean - lowlimit*rms, mean + highlimit*rms);
236     if(printOut) std::cout << mean << ",";
237     }
238     if(printOut) std::cout << " ) ";
239     if(printOut) std::cout << endl;
240     mean=fitFunc->GetParameter(1);
241     sigma=fitFunc->GetParameter(2);
242     error=1.0*fitFunc->GetParError(1);
243    
244     // below this point we're merely doing cosmetics :-)
245     TCanvas *fitcanvas = new TCanvas("fitcanvas","Fitting Canvas");
246     fitcanvas->cd();
247     hist->SetMinimum(0);
248     if(is_data) hist->Draw("e1");
249     else hist->Draw("histo");
250     fitFunc->SetLineColor(kBlue);
251     fitFunc->SetLineWidth(1);
252     fitFunc->Draw("same");
253     hist->SetStats(0);
254     TLegend *leg;
255     if(is_data) {
256     leg= make_legend("Fit (Data)");
257     leg->AddEntry(hist,"Data","p");
258     }
259     else {
260     leg= make_legend("Fit (MC)");
261     leg->AddEntry(hist,"MC","l");
262     }
263    
264     leg->AddEntry(fitFunc,"Fit","l");
265     leg->Draw();
266    
267     TText *ftitle=write_text(0.20,0.86,"Fit results:");
268     ftitle->SetTextSize(0.03);
269     ftitle->SetTextAlign(11);
270     stringstream fitresult;
271     fitresult << "#mu=" << std::setprecision(4) << mean << "+/-" << std::setprecision(4) << error;
272     // TText *title1=write_text(0.20,0.96,fitresult.str().c_str());
273     TText *title1=write_text(0.20,0.82,fitresult.str().c_str());
274     title1->SetTextSize(0.03);
275     title1->SetTextAlign(11);
276     stringstream sigmainfo;
277     sigmainfo << "#sigma=" << std::setprecision(4) << fitFunc->GetParameter(2) << "+/-" << std::setprecision(4) << fitFunc->GetParError(2);
278     // TText *sigmatext=write_text(0.80,0.96,sigmainfo.str().c_str());
279     TText *sigmatext=write_text(0.20,0.78,sigmainfo.str().c_str());
280     sigmatext->SetTextSize(0.03);
281     sigmatext->SetTextAlign(11);
282    
283     TText* toptitle;
284     if(is_data) toptitle = write_title("Fit Result (data)");
285     else toptitle = write_title("Fit Result (MC)");
286     toptitle->Draw();
287     ftitle->Draw();
288     title1->Draw();
289     sigmatext->Draw();
290     if(!is_data) CompleteSave(fitcanvas,"fit/Fit_Summary_MC");
291     if(is_data) CompleteSave(fitcanvas,"fit/Fit_Summary_Data");
292    
293    
294    
295     // cout << "[" << fitFunc->GetParameter(1) << " , " << fitFunc->GetParError(1) << "]" << endl;
296     return mean;
297     }
298    
299    
300    
301     float find_peak(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, float &error,float &Sigma, int method)
302     {
303     float peak_position;
304     if(method==0||method>1) {
305     //looking at a gaus request
306     int numsig=1;
307     if(method>1) numsig=method;
308     peak_position=find_Gauss_peak(all,ttbar,minfit,maxfit,is_data,error,Sigma,numsig);
309     }
310     if(method==1) {
311     //looking at a KM request
312     peak_position=find_KM_peak(all,ttbar,minfit,maxfit,is_data,error,Sigma);
313     }
314     if(method==-99) { // KOSTAS!!
315     TF1 *f1 = new TF1("f1","gaus",-40,40);
316     peak_position=Kostas_algorithm(all,error,Sigma,f1,2.5,2.5,is_data);
317     }
318     return peak_position;
319     }
320    
321    
322     float get_Gaussian_peak(TH1F *hist, float &error, float &sigma, TF1* fitFunc, float lowlimit, float highlimit,bool is_data,int numsig)
323     {
324     TCanvas *fitcanvas = new TCanvas("fitcanvas","fitcanvas");
325     float mean = hist->GetBinCenter( hist->GetMaximumBin());
326     float rms = hist->GetRMS();
327    
328     mean = hist->GetBinCenter( hist->GetMaximumBin());
329    
330     fitFunc->SetParameter(1,mean);
331    
332     hist->Fit(fitFunc,"QLL0","",mean-10,mean+10);
333    
334     mean=fitFunc->GetParameter(1);
335     rms=fitFunc->GetParameter(2);
336     error=fitFunc->GetParError(1);
337    
338     bool printOut = false; // print the peak estimate in the i-th iteration
339    
340     // --- perform iterations
341     int numIterations=5;
342    
343     if(printOut) std::cout << " ( ";
344     for(int i=1;i<numIterations+1;i++) //--- modify the number of iterations until peak is stable
345     {
346     hist->Fit(fitFunc,"QLLN","same",mean - numsig*rms, mean + numsig*rms); // fit -2 +1 sigma from previous iteration
347     mean=fitFunc->GetParameter(1);
348     fitFunc->SetRange(mean - numsig*rms, mean + numsig*rms);
349     if(printOut) std::cout << mean << ",";
350     }
351     if(printOut) std::cout << " ) ";
352     if(printOut) std::cout << endl;
353     mean=fitFunc->GetParameter(1);
354     sigma=fitFunc->GetParameter(2);
355     error=1.0*fitFunc->GetParError(1);
356     fitcanvas->cd();
357     hist->SetMinimum(0);
358     if(is_data) hist->Draw("e1");
359     else hist->Draw("histo");
360     fitFunc->SetLineColor(kBlue);
361     fitFunc->SetLineWidth(1);
362     fitFunc->Draw("same");
363     hist->SetStats(0);
364     TLegend *leg;
365     if(is_data) {
366     leg= make_legend("Fit (Data)");
367     leg->AddEntry(hist,"Data","p");
368     }
369     else {
370     leg= make_legend("Fit (MC)");
371     leg->AddEntry(hist,"MC","l");
372     }
373    
374     leg->AddEntry(fitFunc,"Fit","l");
375     leg->Draw();
376    
377     TText *ftitle=write_text(0.20,0.86,"Fit results:");
378     ftitle->SetTextSize(0.03);
379     ftitle->SetTextAlign(11);
380     stringstream fitresult;
381     fitresult << "#mu=" << std::setprecision(4) << mean << "+/-" << std::setprecision(4) << error;
382     // TText *title1=write_text(0.20,0.96,fitresult.str().c_str());
383     TText *title1=write_text(0.20,0.82,fitresult.str().c_str());
384     title1->SetTextSize(0.03);
385     title1->SetTextAlign(11);
386     stringstream sigmainfo;
387     sigmainfo << "#sigma=" << std::setprecision(4) << fitFunc->GetParameter(2) << "+/-" << std::setprecision(4) << fitFunc->GetParError(2);
388     // TText *sigmatext=write_text(0.80,0.96,sigmainfo.str().c_str());
389     TText *sigmatext=write_text(0.20,0.78,sigmainfo.str().c_str());
390     sigmatext->SetTextSize(0.03);
391     sigmatext->SetTextAlign(11);
392    
393     TText* toptitle;
394     if(is_data) toptitle = write_title("Fit Result (data)");
395     else toptitle = write_title("Fit Result (MC)");
396     toptitle->Draw();
397     ftitle->Draw();
398     title1->Draw();
399     sigmatext->Draw();
400     if(!is_data) CompleteSave(fitcanvas,"fit/Fit_Summary_MC");
401     if(is_data) CompleteSave(fitcanvas,"fit/Fit_Summary_Data");
402    
403    
404     // cout << "[" << fitFunc->GetParameter(1) << " , " << fitFunc->GetParError(1) << "]" << endl;
405     return mean;
406     }
407    
408    
409     float find_Gauss_peak(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, float &error,float &Sigma,int numsig)
410     {
411     TF1 *fitfunc = new TF1("fitfunc","gaus",minfit,maxfit);
412     float peakpos = get_Gaussian_peak(all,error,Sigma,fitfunc, minfit, maxfit,is_data,numsig);
413     return peakpos;
414     }
415    
416     float find_KM_peak(TH1F *all, TH1F *ttbar, float minfit, float maxfit, bool is_data, float &error,float &Sigma)
417     {
418     all->SetLineColor(kBlue);
419     all->SetStats(0);
420     all->SetTitle("");
421     all->GetXaxis()->SetTitle("JZB (GeV/c)");
422     all->GetYaxis()->SetTitle("events");
423     all->GetXaxis()->CenterTitle();
424     all->GetYaxis()->CenterTitle();
425     TF1 *fitfunc = new TF1("fitfunc",KrystalMallLogPar,0.8*minfit,0.8*maxfit,8);
426     if(!is_data)
427     {
428     TF1 *logpar = new TF1("logpar",LogParabola,minfit,maxfit,3);
429     do_ttbar_fit(ttbar,logpar,fitfunc);
430     fitfunc->SetParameters(1000,2,2.5,-1.6,4,logpar->GetParameter(0),logpar->GetParameter(1),logpar->GetParameter(2));
431     parabola_height=logpar->GetParameter(0);
432     parabola_inclination=logpar->GetParameter(1);
433     parabola_pointzero=logpar->GetParameter(2);
434     dofixed=true;//ttbar is known so we can fix the parameters and don't need to use them for fitting!
435     }
436     else
437     {
438     fitfunc->SetParameters(1000,2,2.5,-1.6,4,5.45039,0.000324593,12.3528);
439     dofixed=false;
440     }
441    
442     vector<float> chi2values;
443     addparabola=true;
444     for (int ifit=0;ifit<100;ifit++)
445     {
446     all->Fit(fitfunc,"NQ");
447     chi2values.push_back(fitfunc->GetChisquare());
448     if(ifit>5&&chi2values[ifit-2]==chi2values[ifit]) break;
449     }
450     /*
451     The parameters represent the following quantities:
452     float N=par[0];
453     float alpha=par[1];
454     float n=par[2];
455     float xbar=par[3];
456     float sigma=par[4];
457     */
458     //we are clearing an area of two sigma to the left and to the right of the center of the function for the "real fit".
459     low_reject=-2*fitfunc->GetParameter(4)+fitfunc->GetParameter(3);
460     hi_reject=fitfunc->GetParameter(3)+2*fitfunc->GetParameter(4);
461     if(low_reject>-15) low_reject=-10;
462     if(hi_reject<15) hi_reject=10;
463     doreject=true;//activating the rejection :-)
464    
465     for (int ifit=0;ifit<100;ifit++)
466     {
467     all->Fit(fitfunc,"NQ");
468     chi2values.push_back(fitfunc->GetChisquare());
469     if(ifit>5&&chi2values[ifit-2]==chi2values[ifit]) break;
470     }
471    
472     draw_complete_fit(all,ttbar,minfit,maxfit,is_data,fitfunc);
473     doreject=true;
474     error=fitfunc->GetParError(3);
475     Sigma=fitfunc->GetParameter(4);//sigma
476    
477     return fitfunc->GetParameter(3);
478     }
479    
480     Double_t InvCrystalBall(Double_t *x,Double_t *par)
481     {
482     Double_t arg1=0,arg2=0,A=0,B=0;
483     Double_t f1=0;
484     Double_t f2=0;
485     Double_t lim=0;
486     Double_t fitval=0;
487     Double_t N=0;
488     Double_t n=par[4];
489    
490     Double_t invX = -x[0];
491    
492     if (par[2] != 0)
493     arg1 = (invX-par[1])/par[2];
494    
495     arg2 = ( -pow( TMath::Abs(par[3]) , 2 ) ) / 2 ;
496    
497     if (par[3] != 0)
498     A = pow( ( n / TMath::Abs( par[3] ) ) , n) * TMath::Exp(arg2);
499    
500     if (par[3] != 0)
501     B = n / TMath::Abs(par[3]) - TMath::Abs(par[3]);
502    
503     f1 = TMath::Exp(-0.5*arg1*arg1);
504     if (par[2] != 0)
505     f2 = A * pow( ( B - (invX - par[1])/par[2] ) , -n );
506    
507     if (par[2] != 0)
508     lim = ( par[1] - invX ) / par[2] ;
509    
510     N = par[0];
511    
512    
513    
514     if(lim < par[3])
515     fitval = N * f1;
516     if(lim >= par[3])
517     fitval = N * f2;
518    
519     return fitval;
520     }
521    
522     Double_t InvCrystalBallP(Double_t *x,Double_t *par)
523     {
524     Double_t arg1=0,arg2=0,A=0,B=0;
525     Double_t f1=0;
526     Double_t f2=0;
527     Double_t lim=0;
528     Double_t fitval=0;
529     Double_t N=0;
530     Double_t n=par[4];
531    
532     Double_t invX = -x[0];
533    
534     if (par[2] != 0)
535     arg1 = (invX-par[1])/par[2];
536    
537     arg2 = ( -pow( TMath::Abs(par[3]) , 2 ) ) / 2 ;
538    
539     if (par[3] != 0)
540     A = pow( ( n / TMath::Abs( par[3] ) ) , n) * TMath::Exp(arg2);
541    
542     if (par[3] != 0)
543     B = n / TMath::Abs(par[3]) - TMath::Abs(par[3]);
544    
545     f1 = TMath::Exp(-0.5*arg1*arg1);
546     if (par[2] != 0)
547     f2 = A * pow( ( B - (invX - par[1])/par[2] ) , -n );
548    
549     if (par[2] != 0)
550     lim = ( par[1] - invX ) / par[2] ;
551    
552     N = par[0];
553    
554    
555    
556     if(lim < par[3])
557     fitval = N * f1;
558     if(lim >= par[3])
559     fitval = N * f2;
560    
561     fitval+= statErrorP(fitval);
562     return fitval;
563     }
564    
565     Double_t InvCrystalBallN(Double_t *x,Double_t *par)
566     {
567     Double_t arg1=0,arg2=0,A=0,B=0;
568     Double_t f1=0;
569     Double_t f2=0;
570     Double_t lim=0;
571     Double_t fitval=0;
572     Double_t N=0;
573     Double_t n=par[4];
574    
575     Double_t invX = -x[0];
576    
577     if (par[2] != 0)
578     arg1 = (invX-par[1])/par[2];
579    
580     arg2 = ( -pow( TMath::Abs(par[3]) , 2 ) ) / 2 ;
581    
582     if (par[3] != 0)
583     A = pow( ( n / TMath::Abs( par[3] ) ) , n) * TMath::Exp(arg2);
584    
585     if (par[3] != 0)
586     B = n / TMath::Abs(par[3]) - TMath::Abs(par[3]);
587    
588     f1 = TMath::Exp(-0.5*arg1*arg1);
589     if (par[2] != 0)
590     f2 = A * pow( ( B - (invX - par[1])/par[2] ) , -n );
591    
592     if (par[2] != 0)
593     lim = ( par[1] - invX ) / par[2] ;
594    
595     N = par[0];
596    
597    
598    
599     if(lim < par[3])
600     fitval = N * f1;
601     if(lim >= par[3])
602     fitval = N * f2;
603    
604     fitval-= statErrorN(fitval);
605     return fitval;
606     }
607    
608