CbmRoot
Loading...
Searching...
No Matches
PairAnalysisSignalFit.cxx
Go to the documentation of this file.
1/* Copyright (C) 2015 Justus-Liebig-Universitaet Giessen, Giessen
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Julian Book [committer] */
4
6// Dielectron SignalFit //
7// //
8// //
9/*
10
11 Class used for extracting the signal from an invariant mass spectrum.
12 It implements the PairAnalysisSignalBase and -Ext classes and it uses user provided
13 functions to fit the spectrum with a combined signa+background fit.
14 Used invariant mass spectra are provided via an array of histograms. There are serveral method
15 to estimate the background and to extract the raw yield from the background subtracted spectra.
16
17 Example usage:
18
19 PairAnalysisSignalFit *sig = new PairAnalysisSignalFit();
20
21
22 1) invariant mass input spectra
23
24 1.1) Assuming a PairAnalysisCF container as data format (check class for more details)
25 PairAnalysisCFdraw *cf = new PairAnalysisCFdraw("path/to/the/output/file.root");
26 TObjArray *arrHists = cf->CollectMinvProj(cf->FindStep("Config"));
27
28 1.2) Assuming a PairAnalysisHF grid as data format (check class for more details)
29 PairAnalysisHFhelper *hf = new PairAnalysisHFhelper("path/to/the/output/file.root", "ConfigName");
30 TObjArray *arrHists = hf->CollectHistos(PairAnalysisVarManager::kM);
31
32 1.3) Assuming a single histograms
33 TObjArray *histoArray = new TObjArray();
34 arrHists->Add(signalPP); // add the spectrum histograms to the array
35 arrHists->Add(signalPM); // the order is important !!!
36 arrHists->Add(signalMM);
37
38
39 2) background estimation
40
41 2.1) set the method for the background estimation (methods can be found in PairAnalysisSignalBase)
42 sig->SetMethod(PairAnalysisSignalBase::kFitted);
43 2.2) rebin the spectras if needed
44 // sig->SetRebin(2);
45 2.3) add any background function you like
46 TF1 *fB = new TF1("fitBgrd","pol3",minFit,maxFit);
47
48
49 3) configure the signal extraction
50
51 3.1) chose one of the signal functions (MCshape, CrystalBall, Gauss)
52 TF1 *fS = new TF1("fitSign",PairAnalysisFunction::PeakFunCB,minFit,maxFit,5); // has 5 parameters
53 // TF1 *fS = new TF1("fitSign",PairAnalysisFunction::PeakFunGaus,minFit,maxFit,3); // has 3 parameters
54 // sig->SetMCSignalShape(hMCsign);
55 // TF1 *fS = new TF1("fitSign",PairAnalysisFunction::PeakFunMC,minFit,maxFit,1); // requires a MC shape
56 3.2) set the method for the signal extraction (methods can be found in PairAnalysisSignalBase)
57 depending on the method serveral inputs are needed (e.g. MC shape, PDG code of the particle of interest)
58 // sig->SetParticleOfInterest(443); //default is jpsi
59 // sig->SetMCSignalShape(signalMC);
60 // sig->SetIntegralRange(minInt, maxInt);
61 sig->SetExtractionMethod(PairAnalysisSignal::BinCounting); // this is the default
62
63
64 4) combined fit of bgrd+signal
65
66 4.1) combine the two functions
67 sig->CombineFunc(fS,fB);
68 4.2) apply fitting ranges and the fit options
69 sig->SetFitRange(minFit, maxFit);
70 sig->SetFitOption("NR");
71
72
73 5) start the processing
74
75 sig->Process(arrHists);
76 sig->Print(""); // print values and errors extracted
77
78
79 6) access the spectra and values created
80
81 6.1) standard spectra as provided filled in PairAnalysisSignalExt
82 TH1F *hsign = (TH1F*) sig->GetUnlikeSignHistogram(); // same as the input (rebinned)
83 TH1F *hbgrd = (TH1F*) sig->GetBackgroundHistogram(); // filled histogram with fitBgrd
84 TH1F *hextr = (TH1F*) sig->GetSignalHistogram(); // after backgound extraction (rebinned)
85 TObject *oPeak = (TObject*) sig->GetPeakShape(); // can be a TF1 or TH1 depending on the method
86 6.2) flow spectra
87 TF1 *fFitSign = sig->GetCombinedFunction(); // combined fit function
88 TF1 *fFitExtr = sig->GetSignalFunction(); // signal function
89 TF1 *fFitBgrd = sig->GetBackgroundFunction(); // background function
90 6.3) access the extracted values and errors
91 sig->GetValues(); or GetErrors(); // yield extraction
92
93*/
94// //
96
97#include <TF1.h>
98#include <TFitResult.h>
99#include <TGraph.h>
100#include <TH1.h>
101#include <TList.h>
102#include <TMath.h>
103#include <TPaveText.h>
104#include <TString.h>
105//#include <../hist/hist/src/TF1Helper.h>
106
108
110
111 //TH1F* PairAnalysisSignalFunc::fgHistSimPM=0x0;
112 // // TF1* PairAnalysisSignalFunc::fFuncSignal=0x0;
113 // // TF1* PairAnalysisSignalFunc::fFuncBackground=0x0;
114 // // Int_t PairAnalysisSignalFunc::fNparPeak=0;
115 // // Int_t PairAnalysisSignalFunc::fNparBgnd=0;
116
117
120// ,PairAnalysisFunction()
121{
122 //
123 // Default Constructor
124 //
125}
126
127//______________________________________________
128PairAnalysisSignalFit::PairAnalysisSignalFit(const char* name, const char* title)
129 : PairAnalysisSignalExt(name, title) //,
130// PairAnalysisFunction()
131{
132 //
133 // Named Constructor
134 //
135}
136
137//______________________________________________
139{
140 //
141 // Default Destructor
142 //
143}
144
145
146//______________________________________________
147void PairAnalysisSignalFit::Process(TObjArray* const arrhist)
148{
149 //
150 // Fit the invariant mass histograms and retrieve the signal and background
151 //
152 switch (fMethod) {
153 case kFitted: ProcessFit(arrhist); break;
154
155 case kLikeSignFit: ProcessFitLS(arrhist); break;
156
157 case kEventMixingFit: ProcessFitEM(arrhist); break;
158
159 case kLikeSign:
160 case kLikeSignArithm:
161 case kLikeSignRcorr:
163 case kEventMixing:
164 case kRotation: PairAnalysisSignalExt::Process(arrhist); break;
165
166 default: Error("Process", "Background substraction method not known!");
167 }
168}
169
170//______________________________________________
171void PairAnalysisSignalFit::ProcessFit(TObjArray* const arrhist)
172{
173 //
174 // Fit the +- invariant mass distribution only
175 // Here we assume that the combined fit function is a sum of the signal and background functions
176 // and that the signal function is always the first term of this sum
177 //
178
179 fHistDataPM = (TH1F*) (arrhist->At(1))->Clone("histPM"); // +- SE
180 fHistDataPM->Sumw2();
181 if (fRebin > 1) fHistDataPM->Rebin(fRebin);
182
183 fHistSignal = new TH1F("HistSignal", "Fit substracted signal", fHistDataPM->GetXaxis()->GetNbins(),
184 fHistDataPM->GetXaxis()->GetXmin(), fHistDataPM->GetXaxis()->GetXmax());
185 fHistBackground = new TH1F("HistBackground", "Fit contribution", fHistDataPM->GetXaxis()->GetNbins(),
186 fHistDataPM->GetXaxis()->GetXmin(), fHistDataPM->GetXaxis()->GetXmax());
187
188 // the starting parameters of the fit function and their limits can be tuned
189 // by the user in its macro
190 fHistDataPM->Fit(fFuncSigBack, fFitOpt.Data(), "", fFitMin, fFitMax);
191 TFitResultPtr pmFitPtr = fHistDataPM->Fit(fFuncSigBack, fFitOpt.Data(), "", fFitMin, fFitMax);
192 //TFitResult *pmFitResult = pmFitPtr.Get(); // used only with TF1Helper
193 //fFuncBackground->SetParameters(fFuncSigBack->GetParameters());
194 fFuncSignal->SetParameters(fFuncSigBack->GetParameters());
195 fFuncBackground->SetParameters(fFuncSigBack->GetParameters() + fFuncSignal->GetNpar());
196
197 // fill the background spectrum
199 // set the error for the background histogram
201 Double_t inte = fFuncBackground->IntegralError(fIntMin, fIntMax) / fHistDataPM->GetBinWidth(1);
202 Double_t binte = inte / TMath::Sqrt((fHistDataPM->FindBin(fIntMax) - fHistDataPM->FindBin(fIntMin)) + 1);
203 for (Int_t iBin = fHistDataPM->FindBin(fIntMin); iBin <= fHistDataPM->FindBin(fIntMax); iBin++) {
204 fHistBackground->SetBinError(iBin, binte);
205 }
206
207 for (Int_t iBin = 1; iBin <= fHistDataPM->GetXaxis()->GetNbins(); iBin++) {
208 // Double_t m = fHistDataPM->GetBinCenter(iBin);
209 Double_t pm = fHistDataPM->GetBinContent(iBin);
210 Double_t epm = fHistDataPM->GetBinError(iBin);
211 Double_t bknd = fHistBackground->GetBinContent(iBin);
212 Double_t ebknd = fHistBackground->GetBinError(iBin);
213 Double_t signal = pm - bknd;
214 Double_t error = TMath::Sqrt(epm * epm + ebknd);
215 // theres no signal extraction outside the fit region
216 if (fHistDataPM->GetBinLowEdge(iBin) > fFitMax
217 || fHistDataPM->GetBinLowEdge(iBin) + fHistDataPM->GetBinWidth(iBin) < fFitMin) {
218 signal = 0.0;
219 error = 0.0;
220 }
221 fHistSignal->SetBinContent(iBin, signal);
222 fHistSignal->SetBinError(iBin, error);
223 }
224
225 if (fUseIntegral) {
226 // signal
227 fValues(0) = fFuncSignal->Integral(fIntMin, fIntMax) / fHistDataPM->GetBinWidth(1);
228 fErrors(0) = fFuncSignal->IntegralError(fIntMin, fIntMax) / fHistDataPM->GetBinWidth(1);
229 // fErrors(0) = 0;
230 // background
231 fValues(1) = fFuncBackground->Integral(fIntMin, fIntMax) / fHistDataPM->GetBinWidth(1);
232 fErrors(1) = fFuncBackground->IntegralError(fIntMin, fIntMax) / fHistDataPM->GetBinWidth(1);
233 }
234 else {
235 // signal
236 fValues(0) =
237 fHistSignal->IntegralAndError(fHistSignal->FindBin(fIntMin), fHistSignal->FindBin(fIntMax), fErrors(0));
238 // background
239 fValues(1) = fHistBackground->IntegralAndError(fHistBackground->FindBin(fIntMin), fHistBackground->FindBin(fIntMax),
240 fErrors(1));
241 }
242 // S/B and significance
244 fValues(4) = fFuncSigBack->GetParameter(fParMass);
245 fErrors(4) = fFuncSigBack->GetParError(fParMass);
246 fValues(5) = fFuncSigBack->GetParameter(fParMassWidth);
247 fErrors(5) = fFuncSigBack->GetParError(fParMassWidth);
248
249 // flag
250 fProcessed = kTRUE;
251
252 fHistBackground->GetListOfFunctions()->Add(fFuncBackground);
253}
254
255//______________________________________________
256void PairAnalysisSignalFit::ProcessFitLS(TObjArray* const arrhist)
257{
258 //
259 // Substract background using the like-sign spectrum
260 //
261 fHistDataPP = (TH1F*) (arrhist->At(0))->Clone("histPP"); // ++ SE
262 fHistDataPM = (TH1F*) (arrhist->At(1))->Clone("histPM"); // +- SE
263 fHistDataMM = (TH1F*) (arrhist->At(2))->Clone("histMM"); // -- SE
264 if (fRebin > 1) {
265 fHistDataPP->Rebin(fRebin);
266 fHistDataPM->Rebin(fRebin);
267 fHistDataMM->Rebin(fRebin);
268 }
269 fHistDataPP->Sumw2();
270 fHistDataPM->Sumw2();
271 fHistDataMM->Sumw2();
272
273 fHistSignal = new TH1F("HistSignal", "Like-Sign substracted signal", fHistDataPM->GetXaxis()->GetNbins(),
274 fHistDataPM->GetXaxis()->GetXmin(), fHistDataPM->GetXaxis()->GetXmax());
275 fHistBackground = new TH1F("HistBackground", "Like-sign contribution", fHistDataPM->GetXaxis()->GetNbins(),
276 fHistDataPM->GetXaxis()->GetXmin(), fHistDataPM->GetXaxis()->GetXmax());
277
278 // fit the +- mass distribution
279 fHistDataPM->Fit(fFuncSigBack, fFitOpt.Data(), "", fFitMin, fFitMax);
280 fHistDataPM->Fit(fFuncSigBack, fFitOpt.Data(), "", fFitMin, fFitMax);
281 // declare the variables where the like-sign fit results will be stored
282 // TFitResult *ppFitResult = 0x0;
283 // TFitResult *mmFitResult = 0x0;
284 // fit the like sign background
285 TF1* funcClonePP = (TF1*) fFuncBackground->Clone("funcClonePP");
286 TF1* funcCloneMM = (TF1*) fFuncBackground->Clone("funcCloneMM");
287 fHistDataPP->Fit(funcClonePP, fFitOpt.Data(), "", fFitMin, fFitMax);
288 // TFitResultPtr ppFitPtr = fHistDataPP->Fit(funcClonePP, fFitOpt.Data(), "", fFitMin, fFitMax);
289 // ppFitResult = ppFitPtr.Get();
290 fHistDataMM->Fit(funcCloneMM, fFitOpt.Data(), "", fFitMin, fFitMax);
291 // TFitResultPtr mmFitPtr = fHistDataMM->Fit(funcCloneMM, fFitOpt.Data(), "", fFitMin, fFitMax);
292 // mmFitResult = mmFitPtr.Get();
293
294 for (Int_t iBin = 1; iBin <= fHistDataPM->GetXaxis()->GetNbins(); iBin++) {
295 Double_t m = fHistDataPM->GetBinCenter(iBin);
296 Double_t pm = fHistDataPM->GetBinContent(iBin);
297 Double_t pp = funcClonePP->Eval(m);
298 Double_t mm = funcCloneMM->Eval(m);
299 Double_t epm = fHistDataPM->GetBinError(iBin);
300 // TODO: use TFitResults for errors?
301 Double_t epp = 0;
302 Double_t emm = 0;
303
304 Double_t signal = pm - 2.0 * TMath::Sqrt(pp * mm);
305 Double_t background = 2.0 * TMath::Sqrt(pp * mm);
306
307 // error propagation on the signal calculation above
308 Double_t esignal = TMath::Sqrt(epm * epm + (mm / pp) * epp + (pp / mm) * emm);
309 Double_t ebackground = TMath::Sqrt((mm / pp) * epp + (pp / mm) * emm);
310 fHistSignal->SetBinContent(iBin, signal);
311 fHistSignal->SetBinError(iBin, esignal);
312 fHistBackground->SetBinContent(iBin, background);
313 fHistBackground->SetBinError(iBin, ebackground);
314 }
315
316 // signal
317 fValues(0) = fHistSignal->IntegralAndError(fHistSignal->FindBin(fIntMin), fHistSignal->FindBin(fIntMax), fErrors(0));
318 // background
319 fValues(1) =
320 fHistBackground->IntegralAndError(fHistBackground->FindBin(fIntMin), fHistBackground->FindBin(fIntMax), fErrors(1));
321 // S/B and significance
323 fValues(4) = fFuncSigBack->GetParameter(fParMass);
324 fErrors(4) = fFuncSigBack->GetParError(fParMass);
325 fValues(5) = fFuncSigBack->GetParameter(fParMassWidth);
326 fErrors(5) = fFuncSigBack->GetParError(fParMassWidth);
327
328 // flag
329 fProcessed = kTRUE;
330}
331
332//______________________________________________
333void PairAnalysisSignalFit::ProcessFitEM(TObjArray* const arrhist)
334{
335 //
336 // Substract background with the event mixing technique
337 //
338 arrhist->GetEntries(); // just to avoid the unused parameter warning
339 Error("ProcessFitEM", "Event mixing for background substraction method not yet implemented!");
340}
341
342//______________________________________________
343void PairAnalysisSignalFit::Draw(const Option_t* option)
344{
345 //
346 // Draw the fitted function
347 //
348 // TODO: update
349 TString drawOpt(option);
350 drawOpt.ToLower();
351
352 Bool_t optStat = drawOpt.Contains("stat");
353
354 fFuncSigBack->SetNpx(200);
355 fFuncSigBack->SetRange(fIntMin, fIntMax);
356 fFuncBackground->SetNpx(200);
357 fFuncBackground->SetRange(fIntMin, fIntMax);
358
359 TGraph* grSig = new TGraph(fFuncSigBack);
360 grSig->SetFillColor(kGreen);
361 grSig->SetFillStyle(3001);
362
363 TGraph* grBack = new TGraph(fFuncBackground);
364 grBack->SetFillColor(kRed);
365 grBack->SetFillStyle(3001);
366
367 grSig->SetPoint(0, grBack->GetX()[0], grBack->GetY()[0]);
368 grSig->SetPoint(grSig->GetN() - 1, grBack->GetX()[grBack->GetN() - 1], grBack->GetY()[grBack->GetN() - 1]);
369
370 grBack->SetPoint(0, grBack->GetX()[0], 0.);
371 grBack->SetPoint(grBack->GetN() - 1, grBack->GetX()[grBack->GetN() - 1], 0.);
372
373 fFuncSigBack->SetRange(fFitMin, fFitMax);
374 fFuncBackground->SetRange(fFitMin, fFitMax);
375
376 if (!drawOpt.Contains("same")) {
377 if (fHistDataPM) {
378 fHistDataPM->Draw();
379 grSig->Draw("f");
380 }
381 else {
382 grSig->Draw("af");
383 }
384 }
385 else {
386 grSig->Draw("f");
387 }
388 if (fMethod == kFitted || fMethod == kFittedMC) grBack->Draw("f");
389 fFuncSigBack->Draw("same");
390 fFuncSigBack->SetLineWidth(2);
391 if (fMethod == kLikeSign) {
392 fHistDataPP->SetLineWidth(2);
393 fHistDataPP->SetLineColor(6);
394 fHistDataPP->Draw("same");
395 fHistDataMM->SetLineWidth(2);
396 fHistDataMM->SetLineColor(8);
397 fHistDataMM->Draw("same");
398 }
399
400 if (fMethod == kFitted || fMethod == kFittedMC) fFuncBackground->Draw("same");
401
402 if (optStat) DrawStats();
403}
404
405//______________________________________________
406void PairAnalysisSignalFit::Print(Option_t* /*option*/) const
407{
408 //
409 // Print the statistics
410 //
411 PairAnalysisSignalBase::Print();
412 printf("Fit int. : %.5g - %.5g \n", fFitMin, fFitMax);
413 printf("Fit chi/%d: %.5g \n", fDof, fChi2Dof);
414}
ClassImp(PairAnalysisSignalFit) PairAnalysisSignalFit
void Process(TObjArray *const arrhist)
TPaveText * DrawStats(Double_t x1=0., Double_t y1=0., Double_t x2=0., Double_t y2=0., TString opt="pRnbsSmrc")
virtual void Print(Option_t *option="") const
virtual void Process(TObjArray *const arrhist)
void ProcessFitLS(TObjArray *const arrhist)
void ProcessFitEM(TObjArray *const arrhist)
virtual void Draw(const Option_t *option="")
void ProcessFit(TObjArray *const arrhist)