source: trunk/MagicSoft/Mars/mhflux/MAlphaFitter.cc@ 7110

Last change on this file since 7110 was 7093, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 19.3 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): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MAlphaFitter
28//
29// Create a single Alpha-Plot. The alpha-plot is fitted online. You can
30// check the result when it is filles in the MStatusDisplay
31// For more information see MHFalseSource::FitSignificance
32//
33// For convinience (fit) the output significance is stored in a
34// container in the parlisrt
35//
36// PRELIMINARY!
37//
38//////////////////////////////////////////////////////////////////////////////
39#include "MAlphaFitter.h"
40
41#include <TF1.h>
42#include <TH1.h>
43#include <TH3.h>
44
45#include <TRandom.h>
46
47#include <TLatex.h>
48
49#include "MMath.h"
50
51#include "MLogManip.h"
52
53ClassImp(MAlphaFitter);
54
55using namespace std;
56
57void MAlphaFitter::Clear(Option_t *o)
58{
59 fSignificance=0;
60 fEventsExcess=0;
61 fEventsSignal=0;
62 fEventsBackground=0;
63
64 fChiSqSignal=0;
65 fChiSqBg=0;
66 fIntegralMax=0;
67 fScaleFactor=1;
68
69 fCoefficients.Reset();
70}
71
72// --------------------------------------------------------------------------
73//
74// This is a preliminary implementation of a alpha-fit procedure for
75// all possible source positions. It will be moved into its own
76// more powerfull class soon.
77//
78// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
79// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
80// or
81// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
82//
83// Parameter [1] is fixed to 0 while the alpha peak should be
84// symmetric around alpha=0.
85//
86// Parameter [4] is fixed to 0 because the first derivative at
87// alpha=0 should be 0, too.
88//
89// In a first step the background is fitted between bgmin and bgmax,
90// while the parameters [0]=0 and [2]=1 are fixed.
91//
92// In a second step the signal region (alpha<sigmax) is fittet using
93// the whole function with parameters [1], [3], [4] and [5] fixed.
94//
95// The number of excess and background events are calculated as
96// s = int(hist, 0, 1.25*sigint)
97// b = int(pol2(3), 0, 1.25*sigint)
98//
99// The Significance is calculated using the Significance() member
100// function.
101//
102Bool_t MAlphaFitter::Fit(TH1D &h, Bool_t paint)
103{
104 Clear();
105 if (h.GetEntries()==0)
106 return kFALSE;
107
108 Double_t sigmax=fSigMax;
109 Double_t bgmin =fBgMin;
110 Double_t bgmax =fBgMax;
111
112 //*fLog << inf << "Fit: " << sigmax << " " << fSigInt << " " << bgmin << " " << bgmax << endl;
113
114 //TF1 fFunc("", Form("gaus(0) + pol%d(3)", fPolynomOrder), 0, 90);
115
116 //fFunc->SetParameters(fCoefficients.GetArray());
117
118 fFunc->FixParameter(1, 0);
119 fFunc->FixParameter(4, 0);
120 fFunc->SetParLimits(2, 0, 90);
121 fFunc->SetParLimits(3, -1, 1);
122
123 const Double_t alpha0 = h.GetBinContent(1);
124 const Double_t alphaw = h.GetXaxis()->GetBinWidth(1);
125
126 // Check for the regios which is not filled...
127 if (alpha0==0)
128 return kFALSE; //*fLog << warn << "Histogram empty." << endl;
129
130 // First fit a polynom in the off region
131 fFunc->FixParameter(0, 0);
132 fFunc->FixParameter(2, 1);
133 fFunc->ReleaseParameter(3);
134
135 for (int i=5; i<fFunc->GetNpar(); i++)
136 fFunc->ReleaseParameter(i);
137
138 // options : N do not store the function, do not draw
139 // I use integral of function in bin rather than value at bin center
140 // R use the range specified in the function range
141 // Q quiet mode
142 // E Perform better Errors estimation using Minos technique
143 h.Fit(fFunc, "NQI", "", bgmin, bgmax);
144
145 fChiSqBg = fFunc->GetChisquare()/fFunc->GetNDF();
146
147 // ------------------------------------
148 if (paint)
149 {
150 fFunc->SetRange(0, 90);
151 fFunc->SetLineColor(kRed);
152 fFunc->SetLineWidth(2);
153 fFunc->Paint("same");
154 }
155 // ------------------------------------
156
157 fFunc->ReleaseParameter(0);
158 //func.ReleaseParameter(1);
159 fFunc->ReleaseParameter(2);
160 fFunc->FixParameter(3, fFunc->GetParameter(3));
161 for (int i=5; i<fFunc->GetNpar(); i++)
162 fFunc->FixParameter(i, fFunc->GetParameter(i));
163
164 // Do not allow signals smaller than the background
165 const Double_t A = alpha0-fFunc->GetParameter(3);
166 const Double_t dA = TMath::Abs(A);
167 fFunc->SetParLimits(0, -dA*4, dA*4);
168 fFunc->SetParLimits(2, 0, 90);
169
170 // Now fit a gaus in the on region on top of the polynom
171 fFunc->SetParameter(0, A);
172 fFunc->SetParameter(2, sigmax*0.75);
173
174 // options : N do not store the function, do not draw
175 // I use integral of function in bin rather than value at bin center
176 // R use the range specified in the function range
177 // Q quiet mode
178 // E Perform better Errors estimation using Minos technique
179 h.Fit(fFunc, "NQI", "", 0, sigmax);
180
181 fChiSqSignal = fFunc->GetChisquare()/fFunc->GetNDF();
182 fCoefficients.Set(fFunc->GetNpar(), fFunc->GetParameters());
183
184 //const Bool_t ok = NDF>0 && chi2<2.5*NDF;
185
186 // ------------------------------------
187 if (paint)
188 {
189 fFunc->SetLineColor(kGreen);
190 fFunc->SetLineWidth(2);
191 fFunc->Paint("same");
192 }
193 // ------------------------------------
194
195 //const Double_t s = fFunc->Integral(0, fSigInt)/alphaw;
196 fFunc->SetParameter(0, 0);
197 fFunc->SetParameter(2, 1);
198 //const Double_t b = fFunc->Integral(0, fSigInt)/alphaw;
199 //fSignificance = MMath::SignificanceLiMaSigned(s, b);
200
201 const Int_t bin = h.GetXaxis()->FindFixBin(fSigInt*0.999);
202
203 fIntegralMax = h.GetBinLowEdge(bin+1);
204 fEventsBackground = fFunc->Integral(0, fIntegralMax)/alphaw;
205 fEventsSignal = h.Integral(1, bin);
206 fEventsExcess = fEventsSignal-fEventsBackground;
207 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground);
208
209 if (fEventsExcess<0)
210 fEventsExcess=0;
211
212 return kTRUE;
213}
214
215Bool_t MAlphaFitter::Fit(const TH1D &hon, const TH1D &hof, Double_t alpha, Bool_t paint)
216{
217 TH1D h(hon);
218 h.Add(&hof, -1); // substracts also number of entries!
219 h.SetEntries(hon.GetEntries());
220
221 MAlphaFitter fit(*this);
222 fit.SetPolynomOrder(1);
223
224 if (alpha<=0 || !fit.Fit(h, paint))
225 return kFALSE;
226
227 fChiSqSignal = fit.GetChiSqSignal();
228 fChiSqBg = fit.GetChiSqBg();
229 fCoefficients = fit.GetCoefficients();
230
231 const Int_t bin = hon.GetXaxis()->FindFixBin(fSigInt*0.999);
232
233 fIntegralMax = hon.GetBinLowEdge(bin+1);
234 fEventsBackground = hof.Integral(1, bin);
235 fEventsSignal = hon.Integral(1, bin);
236 fEventsExcess = fEventsSignal-fEventsBackground;
237 fScaleFactor = alpha;
238 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground/alpha, alpha);
239
240 if (fEventsExcess<0)
241 fEventsExcess=0;
242/*
243 TF1 func("", "gaus(0)+pol0(3)", 0., 90.);
244
245 const Double_t A = fEventsSignal/bin;
246 const Double_t dA = TMath::Abs(A);
247 func.SetParLimits(0, -dA*4, dA*4);
248 func.SetParLimits(2, 0, 90);
249 func.SetParLimits(3, -dA, dA);
250
251 func.SetParameter(0, A);
252 func.FixParameter(1, 0);
253 func.SetParameter(2, fSigMax*0.75);
254 func.SetParameter(3, 0);
255
256 // options : N do not store the function, do not draw
257 // I use integral of function in bin rather than value at bin center
258 // R use the range specified in the function range
259 // Q quiet mode
260 // E Perform better Errors estimation using Minos technique
261 TH1D h(hon);
262 h.Add(&hof, -1);
263 h.Fit(&func, "NQI", "", 0, 90);
264
265 fChiSqSignal = func.GetChisquare()/func.GetNDF();
266
267 const Int_t bin1 = h.GetXaxis()->FindFixBin(func.GetParameter(2)*2);
268
269 fChiSqBg = 0;
270 for (int i=bin1; i<=h.GetNbinsX(); i++)
271 {
272 const Float_t val = h.GetBinContent(i);
273 fChiSqBg = val*val;
274 }
275 if (fChiSqBg>0)
276 fChiSqBg /= h.GetNbinsX()+1-bin1;
277
278 fCoefficients.Set(func.GetNpar(), func.GetParameters());
279
280 // ------------------------------------
281 if (paint)
282 {
283 func.SetLineColor(kBlue);
284 func.SetLineWidth(2);
285 func.Paint("same");
286 }
287 // ------------------------------------
288 */
289 return kTRUE;
290}
291
292void MAlphaFitter::PaintResult(Float_t x, Float_t y, Float_t size) const
293{
294 const Double_t w = GetGausSigma();
295 const Double_t m = fIntegralMax;
296
297 const Int_t l1 = w<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(w));
298 const Int_t l2 = m<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(m));
299 const TString fmt = Form("\\sigma_{Li/Ma}=%%.1f \\omega=%%.%df\\circ E=%%d B=%%d (x<%%.%df) (\\chi_{b}^{2}/ndf=%%.1f \\chi_{s}^{2}/ndf=%%.1f c_{0}=%%.1f)",
300 l1<1?1:l1+1, l2<1?1:l2+1);
301
302 TLatex text(x, y, Form(fmt.Data(), fSignificance, w, (int)fEventsExcess,
303 (int)fEventsBackground, m, fChiSqBg, fChiSqSignal,
304 fCoefficients[3]));
305
306 text.SetBit(TLatex::kTextNDC);
307 text.SetTextSize(size);
308 text.Paint();
309}
310
311void MAlphaFitter::Copy(TObject &o) const
312{
313 MAlphaFitter &f = static_cast<MAlphaFitter&>(o);
314
315 // Setup
316 f.fSigInt = fSigInt;
317 f.fSigMax = fSigMax;
318 f.fBgMin = fBgMin;
319 f.fBgMax = fBgMax;
320 f.fScaleMin = fScaleMin;
321 f.fScaleMax = fScaleMax;
322 f.fPolynomOrder = fPolynomOrder;
323 f.fScaleMode = fScaleMode;
324 f.fScaleUser = fScaleUser;
325 f.fStrategy = fStrategy;
326 f.fCoefficients.Set(fCoefficients.GetSize());
327 f.fCoefficients.Reset();
328
329 // Result
330 f.fSignificance = fSignificance;
331 f.fEventsExcess = fEventsExcess;
332 f.fEventsSignal = fEventsSignal;
333 f.fEventsBackground = fEventsBackground;
334 f.fChiSqSignal = fChiSqSignal;
335 f.fChiSqBg = fChiSqBg;
336 f.fIntegralMax = fIntegralMax;
337 f.fScaleFactor = fScaleFactor;
338
339 // Function
340 TF1 *fcn = f.fFunc;
341 f.fFunc = new TF1(*fFunc);
342 gROOT->GetListOfFunctions()->Remove(f.fFunc);
343 f.fFunc->SetName("Dummy");
344 delete fcn;
345}
346
347void MAlphaFitter::Print(Option_t *o) const
348{
349 *fLog << GetDescriptor() << ": Fitting..." << endl;
350 *fLog << " ...background from " << fBgMin << " to " << fBgMax << endl;
351 *fLog << " ...signal to " << fSigMax << " (integrate into bin at " << fSigInt << ")" << endl;
352 *fLog << " ...polynom order " << fPolynomOrder << endl;
353 *fLog << " ...scale mode: ";
354 switch (fScaleMode)
355 {
356 case kNone: *fLog << "none."; break;
357 case kEntries: *fLog << "entries."; break;
358 case kIntegral: *fLog << "integral."; break;
359 case kOffRegion: *fLog << "off region (integral between " << fScaleMin << " and " << fScaleMax << ")"; break;
360 case kBackground: *fLog << "background (integral between " << fBgMin << " and " << fBgMax << ")"; break;
361 case kLeastSquare: *fLog << "least square (N/A)"; break;
362 case kUserScale: *fLog << "user def (" << fScaleUser << ")"; break;
363 }
364 *fLog << endl;
365
366 if (TString(o).Contains("result"))
367 {
368 *fLog << "Result:" << endl;
369 *fLog << " - Significance (Li/Ma) " << fSignificance << endl;
370 *fLog << " - Excess Events " << fEventsExcess << endl;
371 *fLog << " - Signal Events " << fEventsSignal << endl;
372 *fLog << " - Background Events " << fEventsBackground << endl;
373 *fLog << " - Chi^2/ndf (Signal) " << fChiSqSignal << endl;
374 *fLog << " - Chi^2/ndf (Background) " << fChiSqBg << endl;
375 *fLog << " - Signal integrated up to " << fIntegralMax << "°" << endl;
376 *fLog << " - Scale Factor (Off) " << fScaleFactor << endl;
377 }
378}
379
380Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, UInt_t bin, Bool_t paint)
381{
382 const TString name(Form("TempAlphaEnergy%06d", gRandom->Integer(1000000)));
383 TH1D *h = hon.ProjectionZ(name, -1, 9999, bin, bin, "E");
384 h->SetDirectory(0);
385
386 const Bool_t rc = Fit(*h, paint);
387 delete h;
388 return rc;
389}
390
391Bool_t MAlphaFitter::FitTheta(const TH3D &hon, UInt_t bin, Bool_t paint)
392{
393 const TString name(Form("TempAlphaTheta%06d", gRandom->Integer(1000000)));
394 TH1D *h = hon.ProjectionZ(name, bin, bin, -1, 9999, "E");
395 h->SetDirectory(0);
396
397 const Bool_t rc = Fit(*h, paint);
398 delete h;
399 return rc;
400}
401
402Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, Bool_t paint)
403{
404 const TString name(Form("TempAlpha%06d", gRandom->Integer(1000000)));
405 TH1D *h = hon.ProjectionZ(name, -1, 9999, -1, 9999, "E");
406 h->SetDirectory(0);
407
408 const Bool_t rc = Fit(*h, paint);
409 delete h;
410 return rc;
411}
412
413Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
414{
415 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
416 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
417
418 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, bin, bin, "E");
419 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, bin, bin, "E");
420 h1->SetDirectory(0);
421 h0->SetDirectory(0);
422
423 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
424
425 delete h0;
426 delete h1;
427
428 return rc;
429}
430
431Bool_t MAlphaFitter::FitTheta(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
432{
433 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
434 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
435
436 TH1D *h1 = hon.ProjectionZ(name1, bin, bin, -1, 9999, "E");
437 TH1D *h0 = hof.ProjectionZ(name0, bin, bin, -1, 9999, "E");
438 h1->SetDirectory(0);
439 h0->SetDirectory(0);
440
441 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
442
443 delete h0;
444 delete h1;
445
446 return rc;
447}
448
449Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, const TH3D &hof, Bool_t paint)
450{
451 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
452 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
453
454 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, -1, 9999, "E");
455 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, -1, 9999, "E");
456 h1->SetDirectory(0);
457 h0->SetDirectory(0);
458
459 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
460
461 delete h0;
462 delete h1;
463
464 return rc;
465}
466
467Double_t MAlphaFitter::Scale(TH1D &of, const TH1D &on) const
468{
469 Float_t scaleon = 1;
470 Float_t scaleof = 1;
471 switch (fScaleMode)
472 {
473 case kNone:
474 return 1;
475
476 case kEntries:
477 scaleon = on.GetEntries();
478 scaleof = of.GetEntries();
479 break;
480
481 case kIntegral:
482 scaleon = on.Integral();
483 scaleof = of.Integral();
484 break;
485
486 case kOffRegion:
487 {
488 const Int_t min = on.GetXaxis()->FindFixBin(fScaleMin);
489 const Int_t max = on.GetXaxis()->FindFixBin(fScaleMax);
490 scaleon = on.Integral(min, max);
491 scaleof = of.Integral(min, max);
492 }
493 break;
494
495 case kBackground:
496 {
497 const Int_t min = on.GetXaxis()->FindFixBin(fBgMin);
498 const Int_t max = on.GetXaxis()->FindFixBin(fBgMax);
499 scaleon = on.Integral(min, max);
500 scaleof = of.Integral(min, max);
501 }
502 break;
503
504 case kUserScale:
505 scaleon = fScaleUser;
506 break;
507
508 // This is just to make some compiler happy
509 default:
510 return 1;
511 }
512
513 if (scaleof!=0)
514 {
515 of.Scale(scaleon/scaleof);
516 return scaleon/scaleof;
517 }
518 else
519 {
520 of.Reset();
521 return 0;
522 }
523}
524
525Double_t MAlphaFitter::GetMinimizationValue() const
526{
527 switch (fStrategy)
528 {
529 case kSignificance:
530 return -GetSignificance();
531 case kSignificanceChi2:
532 return -GetSignificance()/GetChiSqSignal();
533 case kSignificanceLogExcess:
534 if (GetEventsExcess()<1)
535 return 0;
536 return -GetSignificance()*TMath::Log10(GetEventsExcess());
537 case kSignificanceExcess:
538 return -GetSignificance()*GetEventsExcess();
539 case kExcess:
540 return -GetEventsExcess();
541 }
542 return 0;
543}
544
545Int_t MAlphaFitter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
546{
547 Bool_t rc = kFALSE;
548
549 //void SetScaleUser(Float_t scale) { fScaleUser = scale; fScaleMode=kUserScale; }
550 //void SetScaleMode(ScaleMode_t mode) { fScaleMode = mode; }
551
552 if (IsEnvDefined(env, prefix, "SignalIntegralMax", print))
553 {
554 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalIntegralMax", fSigInt));
555 rc = kTRUE;
556 }
557 if (IsEnvDefined(env, prefix, "SignalFitMax", print))
558 {
559 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalFitMax", fSigMax));
560 rc = kTRUE;
561 }
562 if (IsEnvDefined(env, prefix, "BackgroundFitMax", print))
563 {
564 SetBackgroundFitMax(GetEnvValue(env, prefix, "BackgroundFitMax", fBgMax));
565 rc = kTRUE;
566 }
567 if (IsEnvDefined(env, prefix, "BackgroundFitMin", print))
568 {
569 SetBackgroundFitMin(GetEnvValue(env, prefix, "BackgroundFitMin", fBgMin));
570 rc = kTRUE;
571 }
572 if (IsEnvDefined(env, prefix, "ScaleMin", print))
573 {
574 SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
575 rc = kTRUE;
576 }
577 if (IsEnvDefined(env, prefix, "ScaleMax", print))
578 {
579 SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
580 rc = kTRUE;
581 }
582 if (IsEnvDefined(env, prefix, "PolynomOrder", print))
583 {
584 SetPolynomOrder(GetEnvValue(env, prefix, "PolynomOrder", fPolynomOrder));
585 rc = kTRUE;
586 }
587
588 if (IsEnvDefined(env, prefix, "MinimizationStrategy", print))
589 {
590 TString txt = GetEnvValue(env, prefix, "MinimizationStrategy", "");
591 txt = txt.Strip(TString::kBoth);
592 txt.ToLower();
593 if (txt==(TString)"significance")
594 fStrategy = kSignificance;
595 if (txt==(TString)"significancechi2")
596 fStrategy = kSignificanceChi2;
597 if (txt==(TString)"significanceexcess")
598 fStrategy = kSignificanceExcess;
599 if (txt==(TString)"excess")
600 fStrategy = kExcess;
601 rc = kTRUE;
602 }
603
604 if (IsEnvDefined(env, prefix, "ScaleMode", print))
605 {
606 TString txt = GetEnvValue(env, prefix, "ScaleMode", "");
607 txt = txt.Strip(TString::kBoth);
608 txt.ToLower();
609 if (txt==(TString)"none")
610 fScaleMode = kNone;
611 if (txt==(TString)"entries")
612 fScaleMode = kEntries;
613 if (txt==(TString)"integral")
614 fScaleMode = kIntegral;
615 if (txt==(TString)"offregion")
616 fScaleMode = kOffRegion;
617 if (txt==(TString)"background")
618 fScaleMode = kBackground;
619 if (txt==(TString)"leastsquare")
620 fScaleMode = kLeastSquare;
621 if (txt==(TString)"userscale")
622 fScaleMode = kUserScale;
623 rc = kTRUE;
624 }
625 if (IsEnvDefined(env, prefix, "Scale", print))
626 {
627 fScaleUser = GetEnvValue(env, prefix, "Scale", fScaleUser);
628 rc = kTRUE;
629 }
630
631 return rc;
632}
Note: See TracBrowser for help on using the repository browser.