source: trunk/MagicSoft/Mars/mhist/MHFindSignificance.cc@ 2336

Last change on this file since 2336 was 2318, checked in by wittek, 22 years ago
*** empty log message ***
File size: 53.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Wolfgang Wittek, July 2003 <mailto:wittek@mppmu.mpg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHFindSignificance
28//
29// determines the significance of a gamma signal in an |alpha| plot
30//
31//
32// Input : TH1 histogram of |alpha| : with 0 < |alpha| < 90 degrees
33// alphamin, alphamax : defining the background region
34// alphasig : defining the signal region for which
35// the significance is calculated
36// degree : the degree of the polynomial to be fitted to the background
37// ( a0 + a1*x + a2*x**2 + a3*x**3 + ...)
38//
39// Output :
40//
41// - polynomial which describes the background in the background region
42// - the number of events in the signal region (Non)
43// the number of background events in the signal region (Nbg)
44// - the number of excess events in the signal region (Nex = Non - Nbg)
45// - thew effective number of background events (Noff), and gamma :
46// Nbg = gamma * Noff
47// - the significance of the gamma signal according to Li & Ma
48//
49//
50// call member function 'FindSigma'
51// to fit the background and to determine the significance
52//
53// call the member function 'SigmaVsAlpha'
54// to determine the significance as a function of alphasig
55//
56/////////////////////////////////////////////////////////////////////////////
57#include "MHFindSignificance.h"
58
59#include <fstream>
60#include <math.h>
61
62#include <TArrayD.h>
63#include <TArrayI.h>
64#include <TH1.h>
65#include <TF1.h>
66#include <TCanvas.h>
67#include <TFitter.h>
68#include <TMinuit.h>
69#include <TPaveText.h>
70#include <TStyle.h>
71
72#include "MLog.h"
73#include "MLogManip.h"
74#include "MMinuitInterface.h"
75
76
77ClassImp(MHFindSignificance);
78
79using namespace std;
80
81const TString MHFindSignificance::gsDefName = "MHFindSignificance";
82const TString MHFindSignificance::gsDefTitle = "Find Significance in alpha plot";
83
84
85
86// --------------------------------------------------------------------------
87//
88// fcnpoly
89//
90// calculates the chi2 for the fit of the polynomial function 'poly'
91// to the histogram 'fhist'
92//
93// it is called by CallMinuit() (which is called in FitPolynomial())
94//
95// bins of fhist with huge errors are ignored in the calculation of the chi2
96// (the huge errors were set in 'FitPolynomial()')
97//
98
99static void fcnpoly(Int_t &npar, Double_t *gin, Double_t &f,
100 Double_t *par, Int_t iflag)
101{
102 TH1 *fhist = (TH1*)gMinuit->GetObjectFit();
103 TF1 *fpoly = fhist->GetFunction("Poly");
104
105
106 //-------------------------------------------
107
108 Double_t chi2 = 0.0;
109
110 Int_t nbins = fhist->GetNbinsX();
111 Int_t mbins = 0;
112 for (Int_t i=1; i<=nbins; i++)
113 {
114 Double_t content = fhist->GetBinContent(i);
115 Double_t error = fhist->GetBinError(i);
116 Double_t center = fhist->GetBinCenter(i);
117
118 //-----------------------------
119 // ignore unwanted points
120 if (error > 1.e19)
121 continue;
122
123 if (content <= 0.0)
124 {
125 gLog << "fcnpoly : bin with zero content; i, content, error = "
126 << i << ", " << content << ", " << error << endl;
127 continue;
128 }
129
130 if (error <= 0.0)
131 {
132 gLog << "fcnpoly : bin with zero error; i, content, error = "
133 << i << ", " << content << ", " << error << endl;
134 continue;
135 }
136
137 //-----------------------------
138 mbins++;
139
140 Double_t fu;
141 fu = fpoly->EvalPar(&center, par);
142
143 // the fitted function must not be negative
144 if (fu <= 0.0)
145 {
146 chi2 = 1.e10;
147 break;
148 }
149
150 Double_t temp = (content - fu) / error;
151 chi2 += temp*temp;
152 }
153
154 //-------------------------------------------
155
156 f = chi2;
157
158 //-------------------------------------------
159 // final calculations
160 //if (iflag == 3)
161 //{
162 //}
163
164 //-------------------------------------------------------------
165}
166
167
168// --------------------------------------------------------------------------
169//
170// fcnpolygauss
171//
172// calculates the chi2 for the fit of the (polynomial+Gauss) function
173// 'PolyGauss' to the histogram 'fhist'
174//
175// it is called by CallMinuit() (which is called in FitGaussPoly())
176//
177// bins of fhist with huge errors are ignored in the calculation of the chi2
178// (the huge errors were set in 'FitGaussPoly()')
179//
180
181static void fcnpolygauss(Int_t &npar, Double_t *gin, Double_t &f,
182 Double_t *par, Int_t iflag)
183{
184 TH1 *fhist = (TH1*)gMinuit->GetObjectFit();
185 TF1 *fpolygauss = fhist->GetFunction("PolyGauss");
186
187
188 //-------------------------------------------
189
190 Double_t chi2 = 0.0;
191
192 Int_t nbins = fhist->GetNbinsX();
193 Int_t mbins = 0;
194 for (Int_t i=1; i<=nbins; i++)
195 {
196 Double_t content = fhist->GetBinContent(i);
197 Double_t error = fhist->GetBinError(i);
198 Double_t center = fhist->GetBinCenter(i);
199
200 //-----------------------------
201 // ignore unwanted points
202 if (error > 1.e19)
203 continue;
204
205 if (content <= 0.0)
206 {
207 gLog << "fcnpolygauss : bin with zero content; i, content, error = "
208 << i << ", " << content << ", " << error << endl;
209 continue;
210 }
211
212 if (error <= 0.0)
213 {
214 gLog << "fcnpolygauss : bin with zero error; i, content, error = "
215 << i << ", " << content << ", " << error << endl;
216 continue;
217 }
218
219 //-----------------------------
220 mbins++;
221
222 Double_t fu;
223 fu = fpolygauss->EvalPar(&center, par);
224
225 // the fitted function must not be negative
226 if (fu <= 0.0)
227 {
228 chi2 = 1.e10;
229 break;
230 }
231
232 Double_t temp = (content - fu) / error;
233 chi2 += temp*temp;
234 }
235
236 //-------------------------------------------
237
238 f = chi2;
239
240 //-------------------------------------------
241 // final calculations
242 //if (iflag == 3)
243 //{
244 //}
245
246 //-------------------------------------------------------------
247}
248
249
250
251// --------------------------------------------------------------------------
252//
253// Constructor
254//
255MHFindSignificance::MHFindSignificance(const char *name, const char *title)
256{
257 fName = name ? name : gsDefName.Data();
258 fTitle = title ? title : gsDefTitle.Data();
259
260 fSigVsAlpha = NULL;
261
262 fPoly = NULL;
263 fGPoly = NULL;
264 fGBackg = NULL;
265
266 fHist = NULL;
267 fHistOrig = NULL;
268
269 // allow rebinning of the alpha plot
270 fRebin = kTRUE;
271
272 // allow reducing the degree of the polynomial
273 fReduceDegree = kTRUE;
274
275 fCanvas = NULL;
276}
277
278// --------------------------------------------------------------------------
279//
280// Destructor.
281//
282// =====> it is not clear why one obtains sometimes a segmentation violation
283// when the destructor is active <=======================
284//
285// therefore the 'return'statement
286//
287
288MHFindSignificance::~MHFindSignificance()
289{
290 return;
291
292 *fLog << "destructor of MHFindSignificance is called" << endl;
293
294 //delete fHist;
295
296 delete fSigVsAlpha;
297 delete fPoly;
298 delete fGPoly;
299 delete fGBackg;
300 //delete fCanvas;
301}
302
303// --------------------------------------------------------------------------
304//
305// Set flag fRebin
306//
307// if flag is kTRUE rebinning of the alpha plot is allowed
308//
309//
310void MHFindSignificance::SetRebin(Bool_t b)
311{
312 fRebin = b;
313
314 *fLog << "MHFindSignificance::SetRebin; flag fRebin set to "
315 << (b? "kTRUE" : "kFALSE") << endl;
316}
317
318// --------------------------------------------------------------------------
319//
320// Set flag fReduceDegree
321//
322// if flag is kTRUE reducing of the degree of the polynomial is allowed
323//
324//
325void MHFindSignificance::SetReduceDegree(Bool_t b)
326{
327 fReduceDegree = b;
328
329 *fLog << "MHFindSignificance::SetReduceDegree; flag fReduceDegree set to "
330 << (b? "kTRUE" : "kFALSE") << endl;
331}
332
333// --------------------------------------------------------------------------
334//
335// FindSigma
336//
337// calls FitPolynomial to fit the background in the background region
338// calls DetExcess to determine the number of excess events
339// using an extrapolation of the polynomial
340// into the signal region
341// calls SigmaLiMa to determine the significance of the gamma signal
342// in the range |alpha| < alphasig
343// calls FitGaussPoly to fit a (polynomial+Gauss) function in the
344// whole |alpha| region
345//
346//
347Bool_t MHFindSignificance::FindSigma(TH1 *fhist, Double_t alphamin,
348 Double_t alphamax, Int_t degree, Double_t alphasig,
349 Bool_t drawpoly, Bool_t fitgauss, Bool_t print)
350{
351 //*fLog << "MHFindSignificance::FindSigma;" << endl;
352
353 fHistOrig = fhist;
354
355 fHist = (TH1*)fHistOrig->Clone("|alpha| plot");
356 if ( !fHist )
357 {
358 *fLog << "MHFindSignificance::FindSigma; Clone of histogram could not be generated"
359 << endl;
360 return kFALSE;
361 }
362
363 fHist->Sumw2();
364 fHist->SetNameTitle("Alpha", "alpha plot");
365 fHist->SetXTitle("|alpha| [\\circ]");
366 fHist->SetYTitle("Counts");
367
368 fAlphamin = alphamin;
369 fAlphamax = alphamax;
370 fAlphammm = (alphamin+alphamax)/2.0;
371 fDegree = degree;
372 fAlphasig = alphasig;
373
374 fDraw = drawpoly;
375 fFitGauss = fitgauss;
376
377
378 //--------------------------------------------
379 // fit a polynomial in the background region
380
381 //*fLog << "MHFindSignificance::FindSigma; calling FitPolynomial()" << endl;
382 if ( !FitPolynomial() )
383 {
384 *fLog << "MHFindSignificance::FindSigma; FitPolynomial failed"
385 << endl;
386 return kFALSE;
387 }
388
389
390 //--------------------------------------------
391 // calculate the number of excess events in the signal region
392
393 //*fLog << "MHFindSignificance::FindSigma; calling DetExcess()" << endl;
394 if ( !DetExcess() )
395 {
396 *fLog << "MHFindSignificance::FindSigma; DetExcess failed"
397 << endl;
398 return kFALSE;
399 }
400
401
402 //--------------------------------------------
403 // calculate the significance of the excess
404
405 //*fLog << "MHFindSignificance::FindSigma; calling SigmaLiMa()" << endl;
406 Double_t siglima = 0.0;
407 if ( !SigmaLiMa(fNon, fNoff, fGamma, &siglima) )
408 {
409 *fLog << "MHFindSignificance::FindSigma; SigmaLiMa failed"
410 << endl;
411 return kFALSE;
412 }
413 fSigLiMa = siglima;
414
415 //--------------------------------------------
416 // calculate the error of the number of excess events
417
418 fdNex = fNex / fSigLiMa;
419
420
421 //--------------------------------------------
422
423 //*fLog << "MHFindSignificance::FindSigma; calling PrintPoly()" << endl;
424 if (print)
425 PrintPoly();
426
427
428 //--------------------------------------------
429 // fit a (polynomial + Gauss) function
430
431 if (fFitGauss)
432 {
433 //--------------------------------------------------
434 // delete objects from this fit
435 // in order to have independent starting conditions for the next fit
436
437 delete gMinuit;
438 gMinuit = NULL;
439 //--------------------------------------------------
440
441 //*fLog << "MHFindSignificance::FindSigma; calling FitGaussPoly()"
442 // << endl;
443 if ( !FitGaussPoly() )
444 {
445 *fLog << "MHFindSignificance::FindSigma; FitGaussPoly failed"
446 << endl;
447 return kFALSE;
448 }
449
450 if (print)
451 {
452 //*fLog << "MHFindSignificance::FindSigma; calling PrintPolyGauss()"
453 // << endl;
454 PrintPolyGauss();
455 }
456 }
457
458 //--------------------------------------------------
459 // draw the histogram if requested
460
461 if (fDraw)
462 {
463 //*fLog << "MHFindSignificance::FindSigma; calling DrawFit()" << endl;
464 if ( !DrawFit() )
465 {
466 *fLog << "MHFindSignificance::FindSigma; DrawFit failed"
467 << endl;
468 return kFALSE;
469 }
470 }
471
472
473 //--------------------------------------------------
474 // delete objects from this fit
475 // in order to have independent starting conditions for the next fit
476
477 delete gMinuit;
478 gMinuit = NULL;
479 //--------------------------------------------------
480
481 return kTRUE;
482}
483
484// --------------------------------------------------------------------------
485//
486// SigmaVsAlpha (like FindSigma. However, alphasig is scanned and
487// the significance is plotted versus alphasig)
488//
489// calls FitPolynomial to fit the background in the background region
490//
491// scan alphasig; for a given alphasig :
492// calls DetExcess to determine the number of excess events
493// calls SigmaLiMa to determine the significance of the gamma signal
494// in the range fAlphalow < |alpha| < alphasig
495//
496
497Bool_t MHFindSignificance::SigmaVsAlpha(TH1 *fhist, Double_t alphamin,
498 Double_t alphamax, Int_t degree, Bool_t print)
499{
500 //*fLog << "MHFindSignificance::SigmaVsAlpha;" << endl;
501
502 fHistOrig = fhist;
503
504 fHist = (TH1*)fHistOrig->Clone("|alpha| plot");
505 fHist->Sumw2();
506 fHist->SetNameTitle("alpha", "alpha plot");
507 fHist->SetXTitle("|alpha| [\\circ]");
508 fHist->SetYTitle("Counts");
509
510
511 fAlphamin = alphamin;
512 fAlphamax = alphamax;
513 fAlphammm = (alphamin+alphamax)/2.0;
514 fDegree = degree;
515
516
517 //--------------------------------------------
518 // fit a polynomial in the background region
519
520 //*fLog << "MHFindSignificance::SigmaVsAlpha calling FitPolynomial()"
521 // << endl;
522 if ( !FitPolynomial() )
523 {
524 *fLog << "MHFindSignificance::SigmaVsAlpha; FitPolynomial() failed"
525 << endl;
526 return kFALSE;
527 }
528
529
530 //--------------------------------------------
531 // loop over different signal regions
532
533 Int_t nsteps = 15;
534
535 fSigVsAlpha = new TH1D("SigVsAlpha","Sigma vs Alpha", nsteps, 0.0, alphamin);
536 fSigVsAlpha->SetXTitle("upper edge of signal region in |alpha| [\\circ]");
537 fSigVsAlpha->SetYTitle("Significance of gamma signal");
538
539 for (Int_t i=1; i<=nsteps; i++)
540 {
541 fAlphasig = fSigVsAlpha->GetBinCenter(i);
542
543 if ( !DetExcess() )
544 {
545 *fLog << "MHFindSignificance::SigmaVsAlpha; DetExcess() failed" << endl;
546 continue;
547 }
548
549 Double_t siglima = 0.0;
550 if ( !SigmaLiMa(fNon, fNoff, fGamma, &siglima) )
551 {
552 *fLog << "MHFindSignificance::SigmaVsAlpha; SigmaLiMa() failed" << endl;
553 continue;
554 }
555
556 fdNex = fNex / siglima;
557 fSigVsAlpha->SetBinContent(i, siglima);
558
559 if (print)
560 PrintPoly();
561 }
562
563 //--------------------------------------------
564 // plot significance versus alphasig
565
566 TCanvas *ccc = new TCanvas("SigVsAlpha", "Sigma vs Alpha", 600, 600);
567
568 gROOT->SetSelectedPad(NULL);
569 gStyle->SetPadLeftMargin(0.05);
570
571 ccc->cd();
572 fSigVsAlpha->DrawCopy();
573
574 ccc->Modified();
575 ccc->Update();
576
577 return kTRUE;
578}
579
580// --------------------------------------------------------------------------
581//
582// FitPolynomial
583//
584// - create a clone 'fHist' of the |alpha| distribution 'fHistOrig'
585// - fit a polynomial of degree 'fDegree' to the alpha distribution
586// 'fHist' in the region alphamin < |alpha| < alphamax
587//
588// in pathological cases the histogram is rebinned before fitting
589// (this is done only if fRebin is kTRUE)
590//
591// if the highest coefficient of the polynomial is compatible with zero
592// the fit is repeated with a polynomial of lower degree
593// (this is done only if fReduceDegree is kTRUE)
594//
595//
596Bool_t MHFindSignificance::FitPolynomial()
597{
598 //--------------------------------------------------
599 // check the histogram :
600 // - calculate initial values of the parameters
601 // - check for bins with zero entries
602 // - set minimum errors
603 // - save the original errors
604 // - set errors huge outside the fit range
605 // (in 'fcnpoly' points with huge errors will be ignored)
606
607
608 Double_t dummy = 1.e20;
609
610 Double_t mean;
611 Double_t rms;
612 Double_t nclose;
613 Double_t nfar;
614 Double_t a2init = 0.0;
615 TArrayD saveError;
616
617 Int_t nbins;
618 Int_t nrebin = 1;
619
620 //---------------- start while loop for rebinning -----------------
621 while(1)
622 {
623
624 fNzero = 0;
625 fMbins = 0;
626 fMlow = 0;
627 fNbgtot = 0.0;
628
629 fAlphami = 10000.0;
630 fAlphamm = 10000.0;
631 fAlphama = -10000.0;
632
633 mean = 0.0;
634 rms = 0.0;
635 nclose = 0.0;
636 nfar = 0.0;
637
638 nbins = fHist->GetNbinsX();
639 saveError.Set(nbins);
640
641 for (Int_t i=1; i<=nbins; i++)
642 {
643 saveError[i-1] = fHist->GetBinError(i);
644
645 // bin should be completely contained in the fit range
646 // (fAlphamin, fAlphamax)
647 Double_t xlo = fHist->GetBinLowEdge(i);
648 Double_t xup = fHist->GetBinLowEdge(i+1);
649
650 if ( xlo >= fAlphamin-fEps && xlo <= fAlphamax+fEps &&
651 xup >= fAlphamin-fEps && xup <= fAlphamax+fEps )
652 {
653 fMbins++;
654
655 if ( xlo < fAlphami )
656 fAlphami = xlo;
657
658 if ( xup > fAlphama )
659 fAlphama = xup;
660
661 Double_t content = fHist->GetBinContent(i);
662 fNbgtot += content;
663
664 mean += content;
665 rms += content*content;
666
667 // count events in low-alpha and high-alpha region
668 if ( xlo >= fAlphammm-fEps && xup >= fAlphammm-fEps)
669 {
670 nfar += content;
671 if ( xlo < fAlphamm )
672 fAlphamm = xlo;
673 if ( xup < fAlphamm )
674 fAlphamm = xup;
675 }
676 else
677 {
678 nclose += content;
679 if ( xlo > fAlphamm )
680 fAlphamm = xlo;
681 if ( xup > fAlphamm )
682 fAlphamm = xup;
683 }
684
685 // count bins with zero entry
686 if (content <= 0.0)
687 fNzero++;
688
689 // set minimum error
690 if (content < 9.0)
691 {
692 fMlow += 1;
693 fHist->SetBinError(i, 3.0);
694 }
695
696 //*fLog << "Take : i, content, error = " << i << ", "
697 // << fHist->GetBinContent(i) << ", "
698 // << fHist->GetBinError(i) << endl;
699
700 continue;
701 }
702 // bin is not completely contained in the fit range : set error huge
703
704 fHist->SetBinError(i, dummy);
705
706 //*fLog << "Omit : i, content, error = " << i << ", "
707 // << fHist->GetBinContent(i) << ", " << fHist->GetBinError(i)
708 // << endl;
709
710 }
711
712 // mean of entries/bin in the fit range
713 if (fMbins > 0)
714 {
715 mean /= ((Double_t) fMbins);
716 rms /= ((Double_t) fMbins);
717 }
718
719 rms = sqrt( rms - mean*mean );
720
721 // if there are no events in the background region
722 // there is no reason for rebinning
723 // and this is the condition for assuming a constant background (= 0)
724 if (mean <= 0.0)
725 break;
726
727 Double_t helpmi = fAlphami*fAlphami*fAlphami;
728 Double_t helpmm = fAlphamm*fAlphamm*fAlphamm;
729 Double_t helpma = fAlphama*fAlphama*fAlphama;
730 Double_t help = (helpma-helpmm) * (fAlphamm-fAlphami)
731 - (helpmm-helpmi) * (fAlphama-fAlphamm);
732 if (help != 0.0)
733 a2init = ( (fAlphamm-fAlphami)*nfar - (fAlphama-fAlphamm)*nclose )
734 * 1.5 * fHist->GetBinWidth(1) / help;
735 else
736 a2init = 0.0;
737
738
739 //--------------------------------------------
740 // rebin the histogram
741 // - if a bin has no entries
742 // - or if there are too many bins with too few entries
743 // - or if the new bin width would exceed half the size of the
744 // signal region
745
746 if ( !fRebin ||
747 ( fNzero <= 0 && (Double_t)fMlow<0.05*(Double_t)fMbins ) ||
748 (Double_t)(nrebin+1)/(Double_t)nrebin * fHist->GetBinWidth(1)
749 > fAlphasig/2.0 )
750 {
751 //*fLog << "before break" << endl;
752 break;
753 }
754
755 nrebin += 1;
756 //if (fHist)
757 delete fHist;
758 fHist = NULL;
759
760 *fLog << "MHFindSignificance::FitPolynomial; rebin the |alpha| plot, grouping "
761 << nrebin << " bins together" << endl;
762
763 // TH1::Rebin doesn't work properly
764 //fHist = fHistOrig->Rebin(nrebin, "Rebinned");
765 // use private routine RebinHistogram()
766 fHist = new TH1F;
767 fHist->Sumw2();
768 fHist->SetNameTitle("Rebinned", "Rebinned alpha plot");
769
770
771 // do rebinning such that x0 remains a lower bin edge
772 Double_t x0 = 0.0;
773 if ( !RebinHistogram(x0, nrebin) )
774 {
775 *fLog << "MHFindSignificance::FitPolynomial; RebinHistgram() failed"
776 << endl;
777 return kFALSE;
778 }
779
780 fHist->SetXTitle("|alpha| [\\circ]");
781 fHist->SetYTitle("Counts");
782
783 }
784 //---------------- end of while loop for rebinning -----------------
785
786
787 // if there are still too many bins with too few entries don't fit
788 // and assume a constant background
789
790 fConstantBackg = kFALSE;
791 if ( fNzero > 0 || (Double_t)fMlow>0.05*(Double_t)fMbins )
792 {
793 *fLog << "MHFindSignificance::FitPolynomial; polynomial fit not possible, fNzero, fMlow, fMbins = "
794 << fNzero << ", " << fMlow << ", " << fMbins << endl;
795 *fLog << " assume a constant background" << endl;
796
797 fConstantBackg = kTRUE;
798 fDegree = 0;
799
800 TString funcname = "Poly";
801 Double_t xmin = 0.0;
802 Double_t xmax = 90.0;
803
804 TString formula = "[0]";
805
806 fPoly = new TF1(funcname, formula, xmin, xmax);
807 TList *funclist = fHist->GetListOfFunctions();
808 funclist->Add(fPoly);
809
810 //--------------------
811 Int_t nparfree = 1;
812 fChisq = 0.0;
813 fNdf = fMbins - nparfree;
814 fProb = 0.0;
815 fIstat = 0;
816
817 fValues.Set(1);
818 fErrors.Set(1);
819
820 Double_t val, err;
821 val = mean;
822 err = sqrt( mean / (Double_t)fMbins );
823
824 fPoly->SetParameter(0, val);
825 fPoly->SetParError (0, err);
826
827 fValues[0] = val;
828 fErrors[0] = err;
829
830 fEma[0][0] = err*err;
831 fCorr[0][0] = 1.0;
832 //--------------------
833
834 //--------------------------------------------------
835 // reset the errors of the points in the histogram
836 for (Int_t i=1; i<=nbins; i++)
837 {
838 fHist->SetBinError(i, saveError[i-1]);
839 }
840
841
842 return kTRUE;
843 }
844
845
846 //=========== start loop for reducing the degree ==================
847 // of the polynomial
848 while (1)
849 {
850 //--------------------------------------------------
851 // prepare fit of a polynomial : (a0 + a1*x + a2*x**2 + a3*x**3 + ...)
852
853 TString funcname = "Poly";
854 Double_t xmin = 0.0;
855 Double_t xmax = 90.0;
856
857 TString formula = "[0]";
858 TString bra1 = "+[";
859 TString bra2 = "]";
860 TString xpower = "*x";
861 TString newpower = "*x";
862 for (Int_t i=1; i<=fDegree; i++)
863 {
864 formula += bra1;
865 formula += i;
866 formula += bra2;
867 formula += xpower;
868
869 xpower += newpower;
870 }
871
872 //*fLog << "FitPolynomial : formula = " << formula << endl;
873
874 fPoly = new TF1(funcname, formula, xmin, xmax);
875 TList *funclist = fHist->GetListOfFunctions();
876 funclist->Add(fPoly);
877
878 //------------------------
879 // attention : the dimensions must agree with those in CallMinuit()
880 const UInt_t npar = fDegree+1;
881
882 TString parname[npar];
883 TArrayD vinit(npar);
884 TArrayD step(npar);
885 TArrayD limlo(npar);
886 TArrayD limup(npar);
887 TArrayI fix(npar);
888
889 vinit[0] = mean;
890 vinit[2] = a2init;
891
892 for (UInt_t j=0; j<npar; j++)
893 {
894 parname[j] = "p";
895 parname[j] += j+1;
896
897 step[j] = vinit[j] != 0.0 ? TMath::Abs(vinit[j]) / 10.0 : 0.000001;
898 }
899
900 // limit the first coefficient of the polynomial to positive values
901 // because the background must not be negative
902 limup[0] = fHist->GetEntries();
903
904 // use the subsequernt loop if you want to apply the
905 // constraint : uneven derivatives (at alpha=0) = zero
906 for (UInt_t j=1; j<npar; j+=2)
907 {
908 vinit[j] = 0;
909 step[j] = 0;
910 fix[j] = 1;
911 }
912
913 *fLog << "FitPolynomial : before CallMinuit()" << endl;
914
915 MMinuitInterface inter;
916 const Bool_t rc = inter.CallMinuit(fcnpoly, parname, vinit, step,
917 limlo, limup, fix, fHist, "Migrad",
918 kFALSE);
919
920 *fLog << "FitPolynomial : after CallMinuit()" << endl;
921
922 if (rc != 0)
923 {
924 // *fLog << "MHFindSignificance::FitPolynomial; polynomial fit failed"
925 // << endl;
926 // return kFALSE;
927 }
928
929
930 //-------------------
931 // get status of minimization
932 Double_t fmin = 0;
933 Double_t fedm = 0;
934 Double_t errdef = 0;
935 Int_t npari = 0;
936 Int_t nparx = 0;
937
938 if (gMinuit)
939 gMinuit->mnstat(fmin, fedm, errdef, npari, nparx, fIstat);
940
941 *fLog << "MHFindSignificance::FitPolynomial; fmin, fedm, errdef, npari, nparx, fIstat = "
942 << fmin << ", " << fedm << ", " << errdef << ", " << npari
943 << ", " << nparx << ", " << fIstat << endl;
944
945
946 //-------------------
947 // store the results
948
949 Int_t nparfree = gMinuit!=NULL ? gMinuit->GetNumFreePars() : 0;
950 fChisq = fmin;
951 fNdf = fMbins - nparfree;
952 fProb = TMath::Prob(fChisq, fNdf);
953
954
955 // get fitted parameter values and errors
956 fValues.Set(npar);
957 fErrors.Set(npar);
958
959 for (Int_t j=0; j<=fDegree; j++)
960 {
961 Double_t val, err;
962 if (gMinuit)
963 gMinuit->GetParameter(j, val, err);
964
965 fPoly->SetParameter(j, val);
966 fPoly->SetParError(j, err);
967
968 fValues[j] = val;
969 fErrors[j] = err;
970 }
971
972
973 //--------------------------------------------------
974 // if the highest coefficient (j0) of the polynomial
975 // is consistent with zero reduce the degree of the polynomial
976
977 Int_t j0 = 0;
978 for (Int_t j=fDegree; j>1; j--)
979 {
980 // ignore fixed parameters
981 if (fErrors[j] == 0)
982 continue;
983
984 // this is the highest coefficient
985 j0 = j;
986 break;
987 }
988
989 if (!fReduceDegree || j0==0 || TMath::Abs(fValues[j0]) > fErrors[j0])
990 break;
991
992 // reduce the degree of the polynomial
993 *fLog << "MHFindSignificance::FitPolynomial; reduce the degree of the polynomial from "
994 << fDegree << " to " << (j0-2) << endl;
995 fDegree = j0 - 2;
996
997 funclist->Remove(fPoly);
998 //if (fPoly)
999 delete fPoly;
1000 fPoly = NULL;
1001
1002 // delete the Minuit object in order to have independent starting
1003 // conditions for the next minimization
1004 //if (gMinuit)
1005 delete gMinuit;
1006 gMinuit = NULL;
1007 }
1008 //=========== end of loop for reducing the degree ==================
1009 // of the polynomial
1010
1011
1012 //--------------------------------------------------
1013 // get the error matrix of the fitted parameters
1014
1015
1016 if (fIstat >= 1)
1017 {
1018 // error matrix was calculated
1019 if (gMinuit)
1020 gMinuit->mnemat(&fEmat[0][0], fNdim);
1021
1022 // copy covariance matrix into a matrix which includes also the fixed
1023 // parameters
1024 TString name;
1025 Double_t bnd1, bnd2, val, err;
1026 Int_t jvarbl;
1027 Int_t kvarbl;
1028 for (Int_t j=0; j<=fDegree; j++)
1029 {
1030 if (gMinuit)
1031 gMinuit->mnpout(j, name, val, err, bnd1, bnd2, jvarbl);
1032
1033 for (Int_t k=0; k<=fDegree; k++)
1034 {
1035 if (gMinuit)
1036 gMinuit->mnpout(k, name, val, err, bnd1, bnd2, kvarbl);
1037
1038 fEma[j][k] = jvarbl==0 || kvarbl==0 ? 0 : fEmat[jvarbl-1][kvarbl-1];
1039 }
1040 }
1041 }
1042 else
1043 {
1044 // error matrix was not calculated, construct it
1045 *fLog << "MHFindSignificance::FitPolynomial; error matrix not defined"
1046 << endl;
1047 for (Int_t j=0; j<=fDegree; j++)
1048 {
1049 for (Int_t k=0; k<=fDegree; k++)
1050 fEma[j][k] = 0;
1051
1052 fEma[j][j] = fErrors[j]*fErrors[j];
1053 }
1054 }
1055
1056
1057 //--------------------------------------------------
1058 // calculate correlation matrix
1059 for (Int_t j=0; j<=fDegree; j++)
1060 for (Int_t k=0; k<=fDegree; k++)
1061 {
1062 const Double_t sq = fEma[j][j]*fEma[k][k];
1063 fCorr[j][k] = sq==0 ? 0 : fEma[j][k] / TMath::Sqrt(fEma[j][j]*fEma[k][k]);
1064 }
1065
1066
1067 //--------------------------------------------------
1068 // reset the errors of the points in the histogram
1069 for (Int_t i=1; i<=nbins; i++)
1070 fHist->SetBinError(i, saveError[i-1]);
1071
1072
1073 return kTRUE;
1074}
1075
1076// --------------------------------------------------------------------------
1077//
1078// ReBinHistogram
1079//
1080// rebin the histogram 'fHistOrig' by grouping 'nrebin' bins together
1081// put the result into the histogram 'fHist'
1082// the rebinning is made such that 'x0' remains a lower bound of a bin
1083//
1084
1085Bool_t MHFindSignificance::RebinHistogram(Double_t x0, Int_t nrebin)
1086{
1087 //-----------------------------------------
1088 // search bin i0 which has x0 as lower edge
1089
1090 Int_t i0 = -1;
1091 Int_t nbold = fHistOrig->GetNbinsX();
1092 for (Int_t i=1; i<=nbold; i++)
1093 {
1094 if (TMath::Abs(fHistOrig->GetBinLowEdge(i) - x0) < 1.e-4 )
1095 {
1096 i0 = i;
1097 break;
1098 }
1099 }
1100
1101 if (i0 == -1)
1102 {
1103 i0 = 1;
1104 *fLog << "MHFindsignificance::Rebin; no bin found with " << x0
1105 << " as lower edge, start rebinning with bin 1" << endl;
1106 }
1107
1108 Int_t istart = i0 - nrebin * ( (i0-1)/nrebin );
1109
1110 //-----------------------------------------
1111 // get new bin edges
1112
1113 const Int_t nbnew = (nbold-istart+1) / nrebin;
1114 const Double_t xmin = fHistOrig->GetBinLowEdge(istart);
1115 const Double_t xmax = xmin + (Double_t)nbnew * nrebin * fHistOrig->GetBinWidth(1);
1116 fHist->SetBins(nbnew, xmin, xmax);
1117
1118 *fLog << "MHFindSignificance::ReBin; x0, i0, nbold, nbnew, xmin, xmax = "
1119 << x0 << ", " << i0 << ", " << nbold << ", " << nbnew << ", "
1120 << xmin << ", " << xmax << endl;
1121
1122 //-----------------------------------------
1123 // get new bin entries
1124
1125 for (Int_t i=1; i<=nbnew; i++)
1126 {
1127 Int_t j = nrebin*(i-1) + istart;
1128
1129 Double_t content = 0;
1130 Double_t error2 = 0;
1131 for (Int_t k=0; k<nrebin; k++)
1132 {
1133 content += fHistOrig->GetBinContent(j+k);
1134 error2 += fHistOrig->GetBinError(j+k) * fHistOrig->GetBinError(j+k);
1135 }
1136 fHist->SetBinContent(i, content);
1137 fHist->SetBinError (i, sqrt(error2));
1138 }
1139 fHist->SetEntries( fHistOrig->GetEntries() );
1140
1141 return kTRUE;
1142}
1143
1144// --------------------------------------------------------------------------
1145//
1146// FitGaussPoly
1147//
1148// fits a (Gauss + polynomial function) to the alpha distribution 'fhist'
1149//
1150//
1151Bool_t MHFindSignificance::FitGaussPoly()
1152{
1153 *fLog << "Entry FitGaussPoly" << endl;
1154
1155 //--------------------------------------------------
1156 // check the histogram :
1157 // - calculate initial values of the parameters
1158 // - check for bins with zero entries
1159 // - set minimum errors
1160 // - save the original errors
1161 // - set errors huge outside the fit range
1162 // (in 'fcnpoly' points with huge errors will be ignored)
1163
1164
1165 Double_t dummy = 1.e20;
1166
1167 fGNzero = 0;
1168 fGMbins = 0;
1169
1170 //------------------------------------------
1171 // if a constant background has been assumed (due to low statistics)
1172 // fit only in the signal region
1173 if ( !fConstantBackg )
1174 {
1175 fAlphalow = 0.0;
1176 fAlphahig = fAlphamax;
1177 }
1178 else
1179 {
1180 fAlphalow = 0.0;
1181 fAlphahig = 2.0*fAlphasig>25.0 ? 25.0 : 2.0*fAlphasig;
1182 }
1183 //------------------------------------------
1184
1185
1186 fAlphalo = 10000.0;
1187 fAlphahi = -10000.0;
1188
1189
1190 Int_t nbins = fHist->GetNbinsX();
1191 TArrayD saveError(nbins);
1192
1193 for (Int_t i=1; i<=nbins; i++)
1194 {
1195 saveError[i-1] = fHist->GetBinError(i);
1196
1197 // bin should be completely contained in the fit range
1198 // (fAlphalow, fAlphahig)
1199 Double_t xlo = fHist->GetBinLowEdge(i);
1200 Double_t xup = fHist->GetBinLowEdge(i+1);
1201
1202 if ( xlo >= fAlphalow-fEps && xlo <= fAlphahig+fEps &&
1203 xup >= fAlphalow-fEps && xup <= fAlphahig+fEps )
1204 {
1205 fGMbins++;
1206
1207 if ( xlo < fAlphalo )
1208 fAlphalo = xlo;
1209
1210 if ( xup > fAlphahi )
1211 fAlphahi = xup;
1212
1213 Double_t content = fHist->GetBinContent(i);
1214
1215
1216 // count bins with zero entry
1217 if (content <= 0.0)
1218 fGNzero++;
1219
1220 // set minimum error
1221 if (content < 9.0)
1222 fHist->SetBinError(i, 3.0);
1223
1224 //*fLog << "Take : i, content, error = " << i << ", "
1225 // << fHist->GetBinContent(i) << ", "
1226 // << fHist->GetBinError(i) << endl;
1227
1228 continue;
1229 }
1230 // bin is not completely contained in the fit range : set error huge
1231
1232 fHist->SetBinError(i, dummy);
1233
1234 //*fLog << "Omit : i, content, error = " << i << ", "
1235 // << fHist->GetBinContent(i) << ", " << fHist->GetBinError(i)
1236 // << endl;
1237
1238 }
1239
1240
1241 // if a bin has no entries don't fit
1242 if (fGNzero > 0)
1243 {
1244 *fLog << "MHFindSignificance::FitGaussPoly; out of " << fGMbins
1245 << " bins there are " << fGNzero
1246 << " bins with zero entry" << endl;
1247
1248 fGPoly = NULL;
1249 return kFALSE;
1250 }
1251
1252
1253 //--------------------------------------------------
1254 // prepare fit of a (polynomial+Gauss) :
1255 // (a0 + a1*x + a2*x**2 + a3*x**3 + ...) + A*exp( -0.5*((x-x0)/sigma)**2 )
1256
1257 TString funcname = "PolyGauss";
1258 Double_t xmin = 0.0;
1259 Double_t xmax = 90.0;
1260
1261 TString xpower = "*x";
1262 TString newpower = "*x";
1263
1264 TString formulaBackg = "[0]";
1265 for (Int_t i=1; i<=fDegree; i++)
1266 formulaBackg += Form("+[%d]*x^%d", i, i);
1267
1268 const TString formulaGauss =
1269 Form("[%d]/[%d]*exp(-0.5*((x-[%d])/[%d])^2)",
1270 fDegree+1, fDegree+3, fDegree+2, fDegree+3);
1271
1272 TString formula = formulaBackg;
1273 formula += "+";
1274 formula += formulaGauss;
1275
1276 *fLog << "FitGaussPoly : formulaBackg = " << formulaBackg << endl;
1277 *fLog << "FitGaussPoly : formulaGauss = " << formulaGauss << endl;
1278 *fLog << "FitGaussPoly : formula = " << formula << endl;
1279
1280 fGPoly = new TF1(funcname, formula, xmin, xmax);
1281 TList *funclist = fHist->GetListOfFunctions();
1282 funclist->Add(fGPoly);
1283
1284 fGBackg = new TF1("Backg", formulaBackg, xmin, xmax);
1285 funclist->Add(fGBackg);
1286
1287 //------------------------
1288 // attention : the dimensions must agree with those in CallMinuit()
1289 Int_t npar = fDegree+1 + 3;
1290
1291 TString parname[npar];
1292 TArrayD vinit(npar);
1293 TArrayD step(npar);
1294 TArrayD limlo(npar);
1295 TArrayD limup(npar);
1296 TArrayI fix(npar);
1297
1298
1299 // take as initial values for the polynomial
1300 // the result from the polynomial fit
1301 for (Int_t j=0; j<=fDegree; j++)
1302 vinit[j] = fPoly->GetParameter(j);
1303
1304 Double_t sigma = 8;
1305 vinit[fDegree+1] = 2.0 * fNex * fHist->GetBinWidth(1) / TMath::Sqrt(TMath::Pi()*2);
1306 vinit[fDegree+2] = 0;
1307 vinit[fDegree+3] = sigma;
1308
1309 *fLog << "FitGaussPoly : starting value for Gauss-amplitude = "
1310 << vinit[fDegree+1] << endl;
1311
1312 for (Int_t j=0; j<npar; j++)
1313 {
1314 parname[j] = "p";
1315 parname[j] += j+1;
1316
1317 step[j] = vinit[j]!=0 ? TMath::Abs(vinit[j]) / 10.0 : 0.000001;
1318 }
1319
1320 // limit the first coefficient of the polynomial to positive values
1321 // because the background must not be negative
1322 limup[0] = fHist->GetEntries()*10;
1323
1324 // limit the sigma of the Gauss function
1325 limup[fDegree+3] = 20;
1326
1327
1328 // use the subsequernt loop if you want to apply the
1329 // constraint : uneven derivatives (at alpha=0) = zero
1330 for (Int_t j=1; j<=fDegree; j+=2)
1331 {
1332 vinit[j] = 0;
1333 step[j] = 0;
1334 fix[j] = 1;
1335 }
1336
1337 // fix position of Gauss function
1338 vinit[fDegree+2] = 0;
1339 step[fDegree+2] = 0;
1340 fix[fDegree+2] = 1;
1341
1342 // if a constant background has been assumed (due to low statistics)
1343 // fix the background
1344 if (fConstantBackg)
1345 {
1346 step[0] = 0;
1347 fix[0] = 1;
1348 }
1349
1350 MMinuitInterface inter;
1351 const Bool_t rc = inter.CallMinuit(fcnpolygauss, parname, vinit, step,
1352 limlo, limup, fix, fHist, "Migrad",
1353 kFALSE);
1354
1355 if (rc != 0)
1356 {
1357 // *fLog << "MHFindSignificance::FitGaussPoly; (polynomial+Gauss) fit failed"
1358 // << endl;
1359 // return kFALSE;
1360 }
1361
1362
1363 //-------------------
1364 // get status of the minimization
1365 Double_t fmin;
1366 Double_t fedm;
1367 Double_t errdef;
1368 Int_t npari;
1369 Int_t nparx;
1370
1371 if (gMinuit)
1372 gMinuit->mnstat(fmin, fedm, errdef, npari, nparx, fGIstat);
1373
1374 *fLog << "MHFindSignificance::FitGaussPoly; fmin, fedm, errdef, npari, nparx, fGIstat = "
1375 << fmin << ", " << fedm << ", " << errdef << ", " << npari
1376 << ", " << nparx << ", " << fGIstat << endl;
1377
1378
1379 //-------------------
1380 // store the results
1381
1382 Int_t nparfree = gMinuit!=NULL ? gMinuit->GetNumFreePars() : 0;
1383 fGChisq = fmin;
1384 fGNdf = fGMbins - nparfree;
1385 fGProb = TMath::Prob(fGChisq, fGNdf);
1386
1387
1388 // get fitted parameter values and errors
1389 fGValues.Set(npar);
1390 fGErrors.Set(npar);
1391
1392 for (Int_t j=0; j<npar; j++)
1393 {
1394 Double_t val, err;
1395 if (gMinuit)
1396 gMinuit->GetParameter(j, val, err);
1397
1398 fGPoly->SetParameter(j, val);
1399 fGPoly->SetParError(j, err);
1400
1401 fGValues[j] = val;
1402 fGErrors[j] = err;
1403
1404 if (j <=fDegree)
1405 {
1406 fGBackg->SetParameter(j, val);
1407 fGBackg->SetParError(j, err);
1408 }
1409 }
1410
1411 fSigmaGauss = fGValues[fDegree+3];
1412 fdSigmaGauss = fGErrors[fDegree+3];
1413 // fitted total number of excess events
1414 fNexGauss = fGValues[fDegree+1] * TMath::Sqrt(TMath::Pi()*2) /
1415 (fHist->GetBinWidth(1)*2 );
1416 fdNexGauss = fNexGauss * fGErrors[fDegree+1]/fGValues[fDegree+1];
1417
1418 //--------------------------------------------------
1419 // get the error matrix of the fitted parameters
1420
1421
1422 if (fGIstat >= 1)
1423 {
1424 // error matrix was calculated
1425 if (gMinuit)
1426 gMinuit->mnemat(&fGEmat[0][0], fGNdim);
1427
1428 // copy covariance matrix into a matrix which includes also the fixed
1429 // parameters
1430 TString name;
1431 Double_t bnd1, bnd2, val, err;
1432 Int_t jvarbl;
1433 Int_t kvarbl;
1434 for (Int_t j=0; j<npar; j++)
1435 {
1436 if (gMinuit)
1437 gMinuit->mnpout(j, name, val, err, bnd1, bnd2, jvarbl);
1438
1439 for (Int_t k=0; k<npar; k++)
1440 {
1441 if (gMinuit)
1442 gMinuit->mnpout(k, name, val, err, bnd1, bnd2, kvarbl);
1443
1444 fGEma[j][k] = jvarbl==0 || kvarbl==0 ? 0 : fGEmat[jvarbl-1][kvarbl-1];
1445 }
1446 }
1447 }
1448 else
1449 {
1450 // error matrix was not calculated, construct it
1451 *fLog << "MHFindSignificance::FitPolynomial; error matrix not defined"
1452 << endl;
1453 for (Int_t j=0; j<npar; j++)
1454 {
1455 for (Int_t k=0; k<npar; k++)
1456 fGEma[j][k] = 0;
1457
1458 fGEma[j][j] = fGErrors[j]*fGErrors[j];
1459 }
1460 }
1461
1462
1463 //--------------------------------------------------
1464 // calculate correlation matrix
1465 for (Int_t j=0; j<npar; j++)
1466 {
1467 for (Int_t k=0; k<npar; k++)
1468 {
1469 const Double_t sq = fGEma[j][j]*fGEma[k][k];
1470 fGCorr[j][k] = sq==0 ? 0 : fGEma[j][k] / sqrt( fGEma[j][j]*fGEma[k][k] );
1471 }
1472 }
1473
1474
1475 //--------------------------------------------------
1476 // reset the errors of the points in the histogram
1477 for (Int_t i=1; i<=nbins; i++)
1478 fHist->SetBinError(i, saveError[i-1]);
1479
1480 return kTRUE;
1481
1482}
1483
1484// --------------------------------------------------------------------------
1485//
1486// DetExcess
1487//
1488// using the result of the polynomial fit (fValues), DetExcess determines
1489//
1490// - the total number of events in the signal region (fNon)
1491// - the number of backgound events in the signal region (fNbg)
1492// - the number of excess events (fNex)
1493// - the effective number of background events (fNoff), and fGamma :
1494// fNbg = fGamma * fNoff; fdNbg = fGamma * sqrt(fNoff);
1495//
1496// It assumed that the polynomial is defined as
1497// a0 + a1*x + a2*x**2 + a3*x**3 + ..
1498//
1499// and that the alpha distribution has the range 0 < alpha < 90 degrees
1500//
1501
1502Bool_t MHFindSignificance::DetExcess()
1503{
1504 //*fLog << "MHFindSignificance::DetExcess;" << endl;
1505
1506 //--------------------------------------------
1507 // calculate the total number of events (fNon) in the signal region
1508
1509 fNon = 0.0;
1510 fdNon = 0.0;
1511
1512 Double_t alphaup = -1000.0;
1513 Double_t binwidth = fHist->GetBinWidth(1);
1514
1515 Int_t nbins = fHist->GetNbinsX();
1516 for (Int_t i=1; i<=nbins; i++)
1517 {
1518 Double_t xlo = fHist->GetBinLowEdge(i);
1519 Double_t xup = fHist->GetBinLowEdge(i+1);
1520
1521 // bin must be completely contained in the signal region
1522 if ( xlo <= (fAlphasig+fEps) && xup <= (fAlphasig+fEps) )
1523 {
1524 Double_t width = fabs(xup-xlo);
1525 if (fabs(width-binwidth) > fEps)
1526 {
1527 *fLog << "MHFindSignificance::DetExcess; alpha plot has variable binning, which is not allowed"
1528 << endl;
1529 return kFALSE;
1530 }
1531
1532 if (xup > alphaup)
1533 alphaup = xup;
1534
1535 fNon += fHist->GetBinContent(i);
1536 fdNon += fHist->GetBinError(i) * fHist->GetBinError(i);
1537 }
1538 }
1539 fdNon = sqrt(fdNon);
1540
1541 // the actual signal range is :
1542 if (alphaup == -1000.0)
1543 return kFALSE;
1544
1545 fAlphasi = alphaup;
1546
1547 //*fLog << "fAlphasi, fNon, fdNon, binwidth, fDegree = " << fAlphasi << ", "
1548 // << fNon << ", " << fdNon << ", " << binwidth << ", "
1549 // << fDegree << endl;
1550
1551 //--------------------------------------------
1552 // calculate the number of background events (fNbg) in the signal region
1553 // and its error (fdNbg)
1554
1555 Double_t fac = 1.0/binwidth;
1556
1557 fNbg = 0.0;
1558 Double_t altothejplus1 = fAlphasi;
1559 for (Int_t j=0; j<=fDegree; j++)
1560 {
1561 fNbg += fValues[j] * altothejplus1 / ((Double_t)(j+1));
1562 altothejplus1 *= fAlphasi;
1563 }
1564 fNbg *= fac;
1565
1566 // derivative of Nbg
1567 Double_t facj;
1568 Double_t fack;
1569
1570 Double_t sum = 0.0;
1571 altothejplus1 = fAlphasi;
1572 for (Int_t j=0; j<=fDegree; j++)
1573 {
1574 facj = altothejplus1 / ((Double_t)(j+1));
1575
1576 Double_t altothekplus1 = fAlphasi;
1577 for (Int_t k=0; k<=fDegree; k++)
1578 {
1579 fack = altothekplus1 / ((Double_t)(k+1));
1580
1581 sum += facj * fack * fEma[j][k];
1582 altothekplus1 *= fAlphasi;
1583 }
1584 altothejplus1 *= fAlphasi;
1585 }
1586 sum *= fac*fac;
1587
1588 if (sum < 0.0)
1589 {
1590 *fLog << "MHFindsignificance::DetExcess; error squared is negative"
1591 << endl;
1592 return kFALSE;
1593 }
1594
1595 fdNbg = sqrt(sum);
1596
1597
1598 //--------------------------------------------
1599 // AS A CHECK :
1600 // calculate the number of background events (fNbgtotFitted) in the
1601 // background region, and its error (fdNbgtotFitted)
1602 // expect fdnbg to be approximately equal to sqrt(fNbgtotFitted)
1603
1604 Double_t fNmi = 0.0;
1605 altothejplus1 = fAlphami;
1606 for (Int_t j=0; j<=fDegree; j++)
1607 {
1608 fNmi += fValues[j] * altothejplus1 / ((Double_t)(j+1));
1609 altothejplus1 *= fAlphami;
1610 }
1611 fNmi *= fac;
1612
1613 Double_t fNma = 0.0;
1614 altothejplus1 = fAlphama;
1615 for (Int_t j=0; j<=fDegree; j++)
1616 {
1617 fNma += fValues[j] * altothejplus1 / ((Double_t)(j+1));
1618 altothejplus1 *= fAlphama;
1619 }
1620 fNma *= fac;
1621
1622 fNbgtotFitted = fNma - fNmi;
1623
1624 //----------------------
1625
1626 sum = 0.0;
1627 Double_t altothejma = fAlphama;
1628 Double_t altothejmi = fAlphami;
1629 for (Int_t j=0; j<=fDegree; j++)
1630 {
1631 facj = (altothejma-altothejmi) / ((Double_t)(j+1));
1632
1633 Double_t altothekma = fAlphama;
1634 Double_t altothekmi = fAlphami;
1635 for (Int_t k=0; k<=fDegree; k++)
1636 {
1637 fack = (altothekma-altothekmi) / ((Double_t)(k+1));
1638
1639 sum += facj * fack * fEma[j][k];
1640 altothekma *= fAlphama;
1641 altothekmi *= fAlphami;
1642 }
1643 altothejma *= fAlphama;
1644 altothejmi *= fAlphami;
1645 }
1646 sum *= fac*fac;
1647
1648 fdNbgtotFitted = sqrt(sum);
1649 if ( fabs(fdNbgtotFitted - sqrt(fNbgtotFitted)) > 0.2 * sqrt(fNbgtotFitted) )
1650 {
1651 *fLog << "MHFindSignificance::DetExcess; error of calculated number of background events (in the background region) does not agree with the expectation :"
1652 << endl;
1653 *fLog << " fNbgtotFitted, fdNbgtotFitted = "
1654 << fNbgtotFitted << ", " << fdNbgtotFitted
1655 << ", expected : " << sqrt(fNbgtotFitted) << endl;
1656 }
1657
1658
1659 //--------------------------------------------
1660 // calculate the number of excess events in the signal region
1661
1662 fNex = fNon - fNbg;
1663
1664 //--------------------------------------------
1665 // calculate the effective number of background events (fNoff) , and fGamma :
1666 // fNbg = fGamma * fNoff; dfNbg = fGamma * sqrt(fNoff);
1667
1668 if (fNbg < 0.0)
1669 {
1670 *fLog << "MHFindSignificamce::DetExcess; number of background events is negative, fNbg, fdNbg = "
1671 << fNbg << ", " << fdNbg << endl;
1672
1673 fGamma = 1.0;
1674 fNoff = 0.0;
1675 return kFALSE;
1676 }
1677
1678 if (fNbg > 0.0)
1679 {
1680 fGamma = fdNbg*fdNbg / fNbg;
1681 fNoff = fNbg*fNbg / (fdNbg*fdNbg);
1682 }
1683 else
1684 {
1685 fGamma = 1.0;
1686 fNoff = 0.0;
1687 }
1688
1689 //*fLog << "Exit DetExcess()" << endl;
1690
1691 return kTRUE;
1692}
1693
1694// --------------------------------------------------------------------------
1695//
1696// SigmaLiMa
1697//
1698// calculates the significance according to Li & Ma
1699// ApJ 272 (1983) 317
1700//
1701Bool_t MHFindSignificance::SigmaLiMa(Double_t non, Double_t noff,
1702 Double_t gamma, Double_t *siglima)
1703{
1704 if (gamma <= 0.0 || non <= 0.0 || noff <= 0.0)
1705 {
1706 *siglima = 0.0;
1707 return kFALSE;
1708 }
1709
1710 Double_t help1 = non * log( (1.0+gamma)*non / (gamma*(non+noff)) );
1711 Double_t help2 = noff * log( (1.0+gamma)*noff / ( non+noff ) );
1712 *siglima = sqrt( 2.0 * (help1+help2) );
1713
1714 Double_t nex = non - gamma*noff;
1715 if (nex < 0.0)
1716 *siglima = - *siglima;
1717
1718 //*fLog << "MHFindSignificance::SigmaLiMa; non, noff, gamma, *siglima = "
1719 // << non << ", " << noff << ", " << gamma << ", " << *siglima << endl;
1720
1721 return kTRUE;
1722}
1723
1724// --------------------------------------------------------------------------
1725//
1726Bool_t MHFindSignificance::DrawFit(const Option_t *opt)
1727{
1728 if (fHist == NULL)
1729 *fLog << "MHFindSignificance::DrawFit; fHist = NULL" << endl;
1730
1731
1732 //TCanvas *fCanvas = new TCanvas("Alpha", "Alpha plot", 600, 600);
1733 fCanvas = new TCanvas(fHist->GetName(), "Alpha plot", 600, 600);
1734
1735 //gStyle->SetOptFit(1011);
1736
1737 gROOT->SetSelectedPad(NULL);
1738 gStyle->SetPadLeftMargin(0.1);
1739
1740 fCanvas->cd();
1741
1742
1743 if (fHist)
1744 {
1745 fHist->DrawCopy();
1746 }
1747
1748 TF1 *fpoly = fHist->GetFunction("Poly");
1749 if (fpoly == NULL)
1750 *fLog << "MHFindSignificance::DrawFit; fpoly = NULL" << endl;
1751
1752 if (fpoly)
1753 {
1754 // 2, 1 is red and solid
1755 fpoly->SetLineColor(2);
1756 fpoly->SetLineStyle(1);
1757 fpoly->SetLineWidth(2);
1758 fpoly->DrawCopy("same");
1759 }
1760
1761 if (fFitGauss)
1762 {
1763 TF1 *fpolygauss = fHist->GetFunction("PolyGauss");
1764 if (fpolygauss == NULL)
1765 *fLog << "MHFindSignificance::DrawFit; fpolygauss = NULL" << endl;
1766
1767 if (fpolygauss)
1768 {
1769 // 4, 1 is blue and solid
1770 fpolygauss->SetLineColor(4);
1771 fpolygauss->SetLineStyle(1);
1772 fpolygauss->SetLineWidth(4);
1773 fpolygauss->DrawCopy("same");
1774 }
1775
1776 TF1 *fbackg = fHist->GetFunction("Backg");
1777 if (fbackg == NULL)
1778 *fLog << "MHFindSignificance::DrawFit; fbackg = NULL" << endl;
1779
1780 if (fbackg)
1781 {
1782 // 6, 4 is pink and dotted
1783 fbackg->SetLineColor(4);
1784 fbackg->SetLineStyle(4);
1785 fbackg->SetLineWidth(4);
1786 fbackg->DrawCopy("same");
1787 }
1788 }
1789
1790
1791 //-------------------------------
1792 // print results onto the figure
1793 TPaveText *pt = new TPaveText(0.30, 0.35, 0.70, 0.90, "NDC");
1794 char tx[100];
1795
1796 sprintf(tx, "Results of polynomial fit (order %2d) :", fDegree);
1797 TText *t1 = pt->AddText(tx);
1798 t1->SetTextSize(0.03);
1799 t1->SetTextColor(2);
1800
1801 sprintf(tx, " (%6.2f< |alpha| <%6.2f [\\circ])", fAlphami, fAlphama);
1802 pt->AddText(tx);
1803
1804 sprintf(tx, " chi2 = %8.2f, Ndof = %4d, Prob = %6.2f",
1805 fChisq, fNdf, fProb);
1806 pt->AddText(tx);
1807
1808 sprintf(tx, " Nbgtot(fit) = %8.1f #pm %8.1f",
1809 fNbgtotFitted, fdNbgtotFitted);
1810 pt->AddText(tx);
1811
1812 sprintf(tx, " Nbgtot(meas) = %8.1f", fNbgtot);
1813 pt->AddText(tx);
1814
1815
1816 //sprintf(tx, " ");
1817 //pt->AddText(tx);
1818
1819 //--------------
1820 sprintf(tx, "Results for |alpha|< %6.2f [\\circ] :", fAlphasi);
1821 TText *t6 = pt->AddText(tx);
1822 t6->SetTextSize(0.03);
1823 t6->SetTextColor(8);
1824
1825 sprintf(tx, " Non = %8.1f #pm %8.1f", fNon, fdNon);
1826 pt->AddText(tx);
1827
1828 sprintf(tx, " Nex = %8.1f #pm %8.1f", fNex, fdNex);
1829 pt->AddText(tx);
1830
1831 sprintf(tx, " Nbg = %8.1f #pm %8.1f, gamma = %6.1f",
1832 fNbg, fdNbg, fGamma);
1833 pt->AddText(tx);
1834
1835 Double_t ratio = fNbg>0.0 ? fNex/fNbg : 0.0;
1836 sprintf(tx, " Significance = %6.2f, Nex/Nbg = %6.2f",
1837 fSigLiMa, ratio);
1838 pt->AddText(tx);
1839
1840 //sprintf(tx, " ");
1841 //pt->AddText(tx);
1842
1843 //--------------
1844 if (fFitGauss)
1845 {
1846 sprintf(tx, "Results of (polynomial+Gauss) fit :");
1847 TText *t7 = pt->AddText(tx);
1848 t7->SetTextSize(0.03);
1849 t7->SetTextColor(4);
1850
1851 sprintf(tx, " chi2 = %8.2f, Ndof = %4d, Prob = %6.2f",
1852 fGChisq, fGNdf, fGProb);
1853 pt->AddText(tx);
1854
1855 sprintf(tx, " Sigma of Gauss = %8.1f #pm %8.1f [\\circ]",
1856 fSigmaGauss, fdSigmaGauss);
1857 pt->AddText(tx);
1858
1859 sprintf(tx, " total no.of excess events = %8.1f #pm %8.1f",
1860 fNexGauss, fdNexGauss);
1861 pt->AddText(tx);
1862 }
1863 //--------------
1864
1865 pt->SetFillStyle(0);
1866 pt->SetBorderSize(0);
1867 pt->SetTextAlign(12);
1868
1869 pt->Draw();
1870
1871 fCanvas->Modified();
1872 fCanvas->Update();
1873
1874 return kTRUE;
1875}
1876
1877
1878
1879// --------------------------------------------------------------------------
1880//
1881// Print the results of the polynomial fit to the alpha distribution
1882//
1883//
1884void MHFindSignificance::PrintPoly(Option_t *o)
1885{
1886 *fLog << "---------------------------" << endl;
1887 *fLog << "MHFindSignificance::PrintPoly :" << endl;
1888
1889 *fLog << "fAlphami, fAlphama, fDegree, fAlphasi = "
1890 << fAlphami << ", " << fAlphama << ", " << fDegree << ", "
1891 << fAlphasi << endl;
1892
1893 *fLog << "fMbins, fNzero, fIstat = " << fMbins << ", "
1894 << fNzero << ", " << fIstat << endl;
1895
1896 *fLog << "fChisq, fNdf, fProb = " << fChisq << ", "
1897 << fNdf << ", " << fProb << endl;
1898
1899 *fLog << "fNon, fNbg, fdNbg, fNbgtot, fNbgtotFitted, fdNbgtotFitted = "
1900 << fNon << ", " << fNbg << ", " << fdNbg << ", " << fNbgtot
1901 << ", " << fNbgtotFitted << ", " << fdNbgtotFitted << endl;
1902
1903 Double_t sigtoback = fNbg>0.0 ? fNex/fNbg : 0.0;
1904 *fLog << "fNex, fdNex, fGamma, fNoff, fSigLiMa, sigtoback = "
1905 << fNex << ", " << fdNex << ", " << fGamma << ", " << fNoff
1906 << ", " << fSigLiMa << ", " << sigtoback << endl;
1907
1908 //------------------------------------
1909 // get errors
1910
1911 Double_t eplus;
1912 Double_t eminus;
1913 Double_t eparab;
1914 Double_t gcc;
1915 Double_t errdiag;
1916
1917
1918 if ( !fConstantBackg )
1919 {
1920 *fLog << "parameter value error eplus eminus eparab errdiag gcc"
1921 << endl;
1922
1923 for (Int_t j=0; j<=fDegree; j++)
1924 {
1925 if (gMinuit)
1926 gMinuit->mnerrs(j, eplus, eminus, eparab, gcc);
1927 errdiag = sqrt(fEma[j][j]);
1928 *fLog << j << " " << fValues[j] << " " << fErrors[j] << " "
1929 << eplus << " " << eminus << " " << eparab << " "
1930 << errdiag << " " << gcc << endl;
1931 }
1932 }
1933 else
1934 {
1935 *fLog << "parameter value error errdiag "
1936 << endl;
1937
1938 for (Int_t j=0; j<=fDegree; j++)
1939 {
1940 errdiag = sqrt(fEma[j][j]);
1941 *fLog << j << " " << fValues[j] << " " << fErrors[j] << " "
1942 << errdiag << endl;
1943 }
1944 }
1945
1946 //----------------------------------------
1947 *fLog << "Covariance matrix :" << endl;
1948 for (Int_t j=0; j<=fDegree; j++)
1949 {
1950 *fLog << "j = " << j << " : ";
1951 for (Int_t k=0; k<=fDegree; k++)
1952 {
1953 *fLog << fEma[j][k] << " ";
1954 }
1955 *fLog << endl;
1956 }
1957
1958 *fLog << "Correlation matrix :" << endl;
1959 for (Int_t j=0; j<=fDegree; j++)
1960 {
1961 *fLog << "j = " << j << " : ";
1962 for (Int_t k=0; k<=fDegree; k++)
1963 {
1964 *fLog << fCorr[j][k] << " ";
1965 }
1966 *fLog << endl;
1967 }
1968
1969 *fLog << "---------------------------" << endl;
1970}
1971
1972// --------------------------------------------------------------------------
1973//
1974// Print the results of the (polynomial+Gauss) fit to the alpha distribution
1975//
1976//
1977void MHFindSignificance::PrintPolyGauss(Option_t *o)
1978{
1979 *fLog << "---------------------------" << endl;
1980 *fLog << "MHFindSignificance::PrintPolyGauss :" << endl;
1981
1982 *fLog << "fAlphalo, fAlphahi = "
1983 << fAlphalo << ", " << fAlphahi << endl;
1984
1985 *fLog << "fGMbins, fGNzero, fGIstat = " << fGMbins << ", "
1986 << fGNzero << ", " << fGIstat << endl;
1987
1988 *fLog << "fGChisq, fGNdf, fGProb = " << fGChisq << ", "
1989 << fGNdf << ", " << fGProb << endl;
1990
1991
1992 //------------------------------------
1993 // get errors
1994
1995 Double_t eplus;
1996 Double_t eminus;
1997 Double_t eparab;
1998 Double_t gcc;
1999 Double_t errdiag;
2000
2001 *fLog << "parameter value error eplus eminus eparab errdiag gcc"
2002 << endl;
2003 for (Int_t j=0; j<=(fDegree+3); j++)
2004 {
2005 if (gMinuit)
2006 gMinuit->mnerrs(j, eplus, eminus, eparab, gcc);
2007 errdiag = sqrt(fGEma[j][j]);
2008 *fLog << j << " " << fGValues[j] << " " << fGErrors[j] << " "
2009 << eplus << " " << eminus << " " << eparab << " "
2010 << errdiag << " " << gcc << endl;
2011 }
2012
2013
2014 *fLog << "Covariance matrix :" << endl;
2015 for (Int_t j=0; j<=(fDegree+3); j++)
2016 {
2017 *fLog << "j = " << j << " : ";
2018 for (Int_t k=0; k<=(fDegree+3); k++)
2019 {
2020 *fLog << fGEma[j][k] << " ";
2021 }
2022 *fLog << endl;
2023 }
2024
2025 *fLog << "Correlation matrix :" << endl;
2026 for (Int_t j=0; j<=(fDegree+3); j++)
2027 {
2028 *fLog << "j = " << j << " : ";
2029 for (Int_t k=0; k<=(fDegree+3); k++)
2030 {
2031 *fLog << fGCorr[j][k] << " ";
2032 }
2033 *fLog << endl;
2034 }
2035
2036 *fLog << "---------------------------" << endl;
2037}
2038
2039//============================================================================
2040
Note: See TracBrowser for help on using the repository browser.