source: tags/Mars-V0.9.1/mhflux/MAlphaFitter.cc

Last change on this file was 6926, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 18.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): 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(0, 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(0, bin);
235 fEventsSignal = hon.Integral(0, 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 TLatex text(x, y, Form("\\sigma_{Li/Ma}=%.1f \\omega=%.1f\\circ E=%d (\\alpha<%.1f\\circ) (\\chi_{b}^{2}/ndf=%.1f \\chi_{s}^{2}/ndf=%.1f c_{0}=%.1f)",
295 fSignificance, GetGausSigma(),
296 (int)fEventsExcess, fIntegralMax,
297 fChiSqBg, fChiSqSignal, fCoefficients[3]));
298
299 text.SetBit(TLatex::kTextNDC);
300 text.SetTextSize(size);
301 text.Paint();
302}
303
304void MAlphaFitter::Copy(TObject &o) const
305{
306 MAlphaFitter &f = static_cast<MAlphaFitter&>(o);
307
308 f.fSigInt = fSigInt;
309 f.fSigMax = fSigMax;
310 f.fBgMin = fBgMin;
311 f.fBgMax = fBgMax;
312 f.fScaleMin = fScaleMin;
313 f.fScaleMax = fScaleMax;
314 f.fPolynomOrder = fPolynomOrder;
315 f.fScaleMode = fScaleMode;
316 f.fScaleUser = fScaleUser;
317 f.fStrategy = fStrategy;
318 f.fCoefficients.Set(fCoefficients.GetSize());
319 f.fCoefficients.Reset();
320
321 TF1 *fcn = f.fFunc;
322 f.fFunc = new TF1(*fFunc);
323 gROOT->GetListOfFunctions()->Remove(f.fFunc);
324 f.fFunc->SetName("Dummy");
325 delete fcn;
326}
327
328void MAlphaFitter::Print(Option_t *o) const
329{
330 *fLog << inf << GetDescriptor() << ": Fitting..." << endl;
331 *fLog << " ...background from " << fBgMin << " to " << fBgMax << endl;
332 *fLog << " ...signal to " << fSigMax << " (integrate into bin at " << fSigInt << ")" << endl;
333 *fLog << " ...polynom order " << fPolynomOrder << endl;
334 *fLog << " ...scale mode: ";
335 switch (fScaleMode)
336 {
337 case kNone: *fLog << "none."; break;
338 case kEntries: *fLog << "entries."; break;
339 case kIntegral: *fLog << "integral."; break;
340 case kOffRegion: *fLog << "off region."; break;
341 case kLeastSquare: *fLog << "least square."; break;
342 case kUserScale: *fLog << "user def (" << fScaleUser << ")"; break;
343 }
344 *fLog << endl;
345
346 if (TString(o).Contains("result"))
347 {
348 *fLog << "Result:" << endl;
349 *fLog << " - Significance " << fSignificance << endl;
350 *fLog << " - Excess Events " << fEventsExcess << endl;
351 *fLog << " - Signal Events " << fEventsSignal << endl;
352 *fLog << " - Background Events " << fEventsBackground << endl;
353 *fLog << " - Chi^2/ndf (Signal) " << fChiSqSignal << endl;
354 *fLog << " - Chi^2/ndf (Background) " << fChiSqBg << endl;
355 *fLog << " - Signal integrated " << fIntegralMax << "°" << endl;
356 *fLog << " - Scale Factor (Off) " << fScaleFactor << endl;
357 }
358}
359
360Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, UInt_t bin, Bool_t paint)
361{
362 const TString name(Form("TempAlphaEnergy%06d", gRandom->Integer(1000000)));
363 TH1D *h = hon.ProjectionZ(name, -1, 9999, bin, bin, "E");
364 h->SetDirectory(0);
365
366 const Bool_t rc = Fit(*h, paint);
367 delete h;
368 return rc;
369}
370
371Bool_t MAlphaFitter::FitTheta(const TH3D &hon, UInt_t bin, Bool_t paint)
372{
373 const TString name(Form("TempAlphaTheta%06d", gRandom->Integer(1000000)));
374 TH1D *h = hon.ProjectionZ(name, bin, bin, -1, 9999, "E");
375 h->SetDirectory(0);
376
377 const Bool_t rc = Fit(*h, paint);
378 delete h;
379 return rc;
380}
381
382Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, Bool_t paint)
383{
384 const TString name(Form("TempAlpha%06d", gRandom->Integer(1000000)));
385 TH1D *h = hon.ProjectionZ(name, -1, 9999, -1, 9999, "E");
386 h->SetDirectory(0);
387
388 const Bool_t rc = Fit(*h, paint);
389 delete h;
390 return rc;
391}
392
393Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
394{
395 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
396 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
397
398 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, bin, bin, "E");
399 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, bin, bin, "E");
400 h1->SetDirectory(0);
401 h0->SetDirectory(0);
402
403 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
404
405 delete h0;
406 delete h1;
407
408 return rc;
409}
410
411Bool_t MAlphaFitter::FitTheta(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
412{
413 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
414 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
415
416 TH1D *h1 = hon.ProjectionZ(name1, bin, bin, -1, 9999, "E");
417 TH1D *h0 = hof.ProjectionZ(name0, bin, bin, -1, 9999, "E");
418 h1->SetDirectory(0);
419 h0->SetDirectory(0);
420
421 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
422
423 delete h0;
424 delete h1;
425
426 return rc;
427}
428
429Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, const TH3D &hof, Bool_t paint)
430{
431 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
432 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
433
434 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, -1, 9999, "E");
435 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, -1, 9999, "E");
436 h1->SetDirectory(0);
437 h0->SetDirectory(0);
438
439 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
440
441 delete h0;
442 delete h1;
443
444 return rc;
445}
446
447Double_t MAlphaFitter::Scale(TH1D &of, const TH1D &on) const
448{
449 Float_t scaleon = 1;
450 Float_t scaleof = 1;
451 switch (fScaleMode)
452 {
453 case kNone:
454 return 1;
455
456 case kEntries:
457 scaleon = on.GetEntries();
458 scaleof = of.GetEntries();
459 break;
460
461 case kIntegral:
462 scaleon = on.Integral();
463 scaleof = of.Integral();
464 break;
465
466 case kOffRegion:
467 {
468 const Int_t min = on.GetXaxis()->FindFixBin(fScaleMin);
469 const Int_t max = on.GetXaxis()->FindFixBin(fScaleMax);
470 scaleon = on.Integral(min, max);
471 scaleof = of.Integral(min, max);
472 }
473 break;
474
475 case kUserScale:
476 scaleon = fScaleUser;
477 break;
478
479 // This is just to make some compiler happy
480 default:
481 return 1;
482 }
483
484 if (scaleof!=0)
485 {
486 of.Scale(scaleon/scaleof);
487 return scaleon/scaleof;
488 }
489 else
490 {
491 of.Reset();
492 return 0;
493 }
494}
495
496Double_t MAlphaFitter::GetMinimizationValue() const
497{
498 switch (fStrategy)
499 {
500 case kSignificance:
501 return -GetSignificance();
502 case kSignificanceChi2:
503 return -GetSignificance()/GetChiSqSignal();
504 case kSignificanceLogExcess:
505 if (GetEventsExcess()<1)
506 return 0;
507 return -GetSignificance()*TMath::Log10(GetEventsExcess());
508 case kSignificanceExcess:
509 return -GetSignificance()*GetEventsExcess();
510 case kExcess:
511 return -GetEventsExcess();
512 }
513 return 0;
514}
515
516Int_t MAlphaFitter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
517{
518 Bool_t rc = kFALSE;
519
520 //void SetScaleUser(Float_t scale) { fScaleUser = scale; fScaleMode=kUserScale; }
521 //void SetScaleMode(ScaleMode_t mode) { fScaleMode = mode; }
522
523 if (IsEnvDefined(env, prefix, "SignalIntegralMax", print))
524 {
525 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalIntegralMax", fSigInt));
526 rc = kTRUE;
527 }
528 if (IsEnvDefined(env, prefix, "SignalFitMax", print))
529 {
530 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalFitMax", fSigMax));
531 rc = kTRUE;
532 }
533 if (IsEnvDefined(env, prefix, "BackgroundFitMax", print))
534 {
535 SetBackgroundFitMax(GetEnvValue(env, prefix, "BackgroundFitMax", fBgMax));
536 rc = kTRUE;
537 }
538 if (IsEnvDefined(env, prefix, "BackgroundFitMin", print))
539 {
540 SetBackgroundFitMin(GetEnvValue(env, prefix, "BackgroundFitMin", fBgMin));
541 rc = kTRUE;
542 }
543 if (IsEnvDefined(env, prefix, "ScaleMin", print))
544 {
545 SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
546 rc = kTRUE;
547 }
548 if (IsEnvDefined(env, prefix, "ScaleMax", print))
549 {
550 SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
551 rc = kTRUE;
552 }
553 if (IsEnvDefined(env, prefix, "PolynomOrder", print))
554 {
555 SetPolynomOrder(GetEnvValue(env, prefix, "PolynomOrder", fPolynomOrder));
556 rc = kTRUE;
557 }
558
559 if (IsEnvDefined(env, prefix, "MinimizationStrategy", print))
560 {
561 TString txt = GetEnvValue(env, prefix, "MinimizationStrategy", "");
562 txt = txt.Strip(TString::kBoth);
563 txt.ToLower();
564 if (txt==(TString)"significance")
565 fStrategy = kSignificance;
566 if (txt==(TString)"significancechi2")
567 fStrategy = kSignificanceChi2;
568 if (txt==(TString)"significanceexcess")
569 fStrategy = kSignificanceExcess;
570 if (txt==(TString)"excess")
571 fStrategy = kExcess;
572 rc = kTRUE;
573 }
574
575 if (IsEnvDefined(env, prefix, "ScaleMode", print))
576 {
577 TString txt = GetEnvValue(env, prefix, "ScaleMode", "");
578 txt = txt.Strip(TString::kBoth);
579 txt.ToLower();
580 if (txt==(TString)"none")
581 fScaleMode = kNone;
582 if (txt==(TString)"entries")
583 fScaleMode = kEntries;
584 if (txt==(TString)"integral")
585 fScaleMode = kIntegral;
586 if (txt==(TString)"offregion")
587 fScaleMode = kOffRegion;
588 if (txt==(TString)"leastsquare")
589 fScaleMode = kLeastSquare;
590 if (txt==(TString)"userscale")
591 fScaleMode = kUserScale;
592 rc = kTRUE;
593 }
594 if (IsEnvDefined(env, prefix, "Scale", print))
595 {
596 fScaleUser = GetEnvValue(env, prefix, "Scale", fScaleUser);
597 rc = kTRUE;
598 }
599
600 return rc;
601}
Note: See TracBrowser for help on using the repository browser.