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

Last change on this file since 6456 was 6282, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 16.8 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(TH1D &hon, TH1D &hof, Double_t alpha, Bool_t paint)
216{
217 /*
218 Clear();
219 if (hon.GetEntries()==0)
220 return kFALSE;
221 */
222
223 TH1D h(hon);
224 h.Add(&hof, -1);
225
226 MAlphaFitter fit(*this);
227 fit.SetPolynomOrder(1);
228
229 if (alpha<=0 || !fit.Fit(h, paint))
230 return kFALSE;
231
232 fChiSqSignal = fit.GetChiSqSignal();
233 fChiSqBg = fit.GetChiSqBg();
234 fCoefficients = fit.GetCoefficients();
235
236 const Int_t bin = hon.GetXaxis()->FindFixBin(fSigInt*0.999);
237
238 fIntegralMax = hon.GetBinLowEdge(bin+1);
239 fEventsBackground = hof.Integral(0, bin);
240 fEventsSignal = hon.Integral(0, bin);
241 fEventsExcess = fEventsSignal-fEventsBackground;
242 fScaleFactor = alpha;
243 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground/alpha, alpha);
244
245 if (fEventsExcess<0)
246 fEventsExcess=0;
247/*
248 TF1 func("", "gaus(0)+pol0(3)", 0., 90.);
249
250 const Double_t A = fEventsSignal/bin;
251 const Double_t dA = TMath::Abs(A);
252 func.SetParLimits(0, -dA*4, dA*4);
253 func.SetParLimits(2, 0, 90);
254 func.SetParLimits(3, -dA, dA);
255
256 func.SetParameter(0, A);
257 func.FixParameter(1, 0);
258 func.SetParameter(2, fSigMax*0.75);
259 func.SetParameter(3, 0);
260
261 // options : N do not store the function, do not draw
262 // I use integral of function in bin rather than value at bin center
263 // R use the range specified in the function range
264 // Q quiet mode
265 // E Perform better Errors estimation using Minos technique
266 TH1D h(hon);
267 h.Add(&hof, -1);
268 h.Fit(&func, "NQI", "", 0, 90);
269
270 fChiSqSignal = func.GetChisquare()/func.GetNDF();
271
272 const Int_t bin1 = h.GetXaxis()->FindFixBin(func.GetParameter(2)*2);
273
274 fChiSqBg = 0;
275 for (int i=bin1; i<=h.GetNbinsX(); i++)
276 {
277 const Float_t val = h.GetBinContent(i);
278 fChiSqBg = val*val;
279 }
280 if (fChiSqBg>0)
281 fChiSqBg /= h.GetNbinsX()+1-bin1;
282
283 fCoefficients.Set(func.GetNpar(), func.GetParameters());
284
285 // ------------------------------------
286 if (paint)
287 {
288 func.SetLineColor(kBlue);
289 func.SetLineWidth(2);
290 func.Paint("same");
291 }
292 // ------------------------------------
293 */
294 return kTRUE;
295}
296
297void MAlphaFitter::PaintResult(Float_t x, Float_t y, Float_t size) const
298{
299 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)",
300 fSignificance, GetGausSigma(),
301 (int)fEventsExcess, fIntegralMax,
302 fChiSqBg, fChiSqSignal, fCoefficients[3]));
303
304 text.SetBit(TLatex::kTextNDC);
305 text.SetTextSize(size);
306 text.Paint();
307}
308
309void MAlphaFitter::Copy(TObject &o) const
310{
311 MAlphaFitter &f = static_cast<MAlphaFitter&>(o);
312
313 f.fSigInt = fSigInt;
314 f.fSigMax = fSigMax;
315 f.fBgMin = fBgMin;
316 f.fBgMax = fBgMax;
317 f.fScaleMin = fScaleMin;
318 f.fScaleMax = fScaleMax;
319 f.fPolynomOrder = fPolynomOrder;
320 f.fScaleMode = fScaleMode;
321 f.fScaleUser = fScaleUser;
322 f.fStrategy = fStrategy;
323 f.fCoefficients.Set(fCoefficients.GetSize());
324 f.fCoefficients.Reset();
325
326 TF1 *fcn = f.fFunc;
327 f.fFunc = new TF1(*fFunc);
328 gROOT->GetListOfFunctions()->Remove(f.fFunc);
329 f.fFunc->SetName("Dummy");
330 delete fcn;
331}
332
333void MAlphaFitter::Print(Option_t *o) const
334{
335 *fLog << inf << GetDescriptor() << ": Fitting..." << endl;
336 *fLog << " ...background from " << fBgMin << " to " << fBgMax << endl;
337 *fLog << " ...signal to " << fSigMax << " (integrate into bin at " << fSigInt << ")" << endl;
338 *fLog << " ...polynom order " << fPolynomOrder << endl;
339 *fLog << " ...scale mode: ";
340 switch (fScaleMode)
341 {
342 case kNone: *fLog << "none."; break;
343 case kEntries: *fLog << "entries."; break;
344 case kIntegral: *fLog << "integral."; break;
345 case kOffRegion: *fLog << "off region."; break;
346 case kLeastSquare: *fLog << "least square."; break;
347 case kUserScale: *fLog << "user def (" << fScaleUser << ")"; break;
348 }
349 *fLog << endl;
350
351 if (TString(o).Contains("result"))
352 {
353 *fLog << "Result:" << endl;
354 *fLog << " - Significance " << fSignificance << endl;
355 *fLog << " - Excess Events " << fEventsExcess << endl;
356 *fLog << " - Signal Events " << fEventsSignal << endl;
357 *fLog << " - Background Events " << fEventsBackground << endl;
358 *fLog << " - Chi^2/ndf (Signal) " << fChiSqSignal << endl;
359 *fLog << " - Chi^2/ndf (Background) " << fChiSqBg << endl;
360 *fLog << " - Signal integrated " << fIntegralMax << "°" << endl;
361 }
362}
363
364Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, UInt_t bin, Bool_t paint)
365{
366 const TString name(Form("TempAlphaEnergy%06d", gRandom->Integer(1000000)));
367 TH1D *h = hon.ProjectionZ(name, -1, 9999, bin, bin, "E");
368 h->SetDirectory(0);
369
370 const Bool_t rc = Fit(*h, paint);
371 delete h;
372 return rc;
373}
374
375Bool_t MAlphaFitter::FitTheta(const TH3D &hon, UInt_t bin, Bool_t paint)
376{
377 const TString name(Form("TempAlphaTheta%06d", gRandom->Integer(1000000)));
378 TH1D *h = hon.ProjectionZ(name, bin, bin, -1, 9999, "E");
379 h->SetDirectory(0);
380
381 const Bool_t rc = Fit(*h, paint);
382 delete h;
383 return rc;
384}
385
386Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, Bool_t paint)
387{
388 const TString name(Form("TempAlpha%06d", gRandom->Integer(1000000)));
389 TH1D *h = hon.ProjectionZ(name, -1, 9999, -1, 9999, "E");
390 h->SetDirectory(0);
391
392 const Bool_t rc = Fit(*h, paint);
393 delete h;
394 return rc;
395}
396
397Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
398{
399 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
400 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
401
402 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, bin, bin, "E");
403 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, bin, bin, "E");
404 h1->SetDirectory(0);
405 h0->SetDirectory(0);
406
407 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
408
409 delete h0;
410 delete h1;
411
412 return rc;
413}
414
415Bool_t MAlphaFitter::FitTheta(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
416{
417 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
418 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
419
420 TH1D *h1 = hon.ProjectionZ(name1, bin, bin, -1, 9999, "E");
421 TH1D *h0 = hof.ProjectionZ(name0, bin, bin, -1, 9999, "E");
422 h1->SetDirectory(0);
423 h0->SetDirectory(0);
424
425 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
426
427 delete h0;
428 delete h1;
429
430 return rc;
431}
432
433Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, const TH3D &hof, Bool_t paint)
434{
435 const TString name1(Form("TempAlpha%06d_on", gRandom->Integer(1000000)));
436 const TString name0(Form("TempAlpha%06d_off", gRandom->Integer(1000000)));
437
438 TH1D *h1 = hon.ProjectionZ(name1, -1, 9999, -1, 9999, "E");
439 TH1D *h0 = hof.ProjectionZ(name0, -1, 9999, -1, 9999, "E");
440 h1->SetDirectory(0);
441 h0->SetDirectory(0);
442
443 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
444
445 delete h0;
446 delete h1;
447
448 return rc;
449}
450
451Double_t MAlphaFitter::Scale(TH1D &of, const TH1D &on) const
452{
453 Float_t scaleon = 1;
454 Float_t scaleof = 1;
455 switch (fScaleMode)
456 {
457 case kNone:
458 return 1;
459
460 case kEntries:
461 scaleon = on.GetEntries();
462 scaleof = of.GetEntries();
463 break;
464
465 case kIntegral:
466 scaleon = on.Integral();
467 scaleof = of.Integral();
468 break;
469
470 case kOffRegion:
471 {
472 const Int_t min = on.GetXaxis()->FindFixBin(fScaleMin);
473 const Int_t max = on.GetXaxis()->FindFixBin(fScaleMax);
474 scaleon = on.Integral(min, max);
475 scaleof = of.Integral(min, max);
476 }
477 break;
478
479 case kUserScale:
480 scaleon = fScaleUser;
481 break;
482
483 // This is just to make some compiler happy
484 default:
485 return 1;
486 }
487
488 if (scaleof!=0)
489 {
490 of.Scale(scaleon/scaleof);
491 return scaleon/scaleof;
492 }
493 else
494 {
495 of.Reset();
496 return 0;
497 }
498}
499
500Double_t MAlphaFitter::GetMinimizationValue() const
501{
502 switch (fStrategy)
503 {
504 case kSignificance:
505 return -GetSignificance();
506 case kSignificanceChi2:
507 return -GetSignificance()/GetChiSqSignal();
508 }
509 return 0;
510}
511
512Int_t MAlphaFitter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
513{
514 Bool_t rc = kFALSE;
515
516 //void SetScaleUser(Float_t scale) { fScaleUser = scale; fScaleMode=kUserScale; }
517 //void SetScaleMode(ScaleMode_t mode) { fScaleMode = mode; }
518
519 if (IsEnvDefined(env, prefix, "SignalIntegralMax", print))
520 {
521 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalIntegralMax", fSigInt));
522 rc = kTRUE;
523 }
524 if (IsEnvDefined(env, prefix, "SignalFitMax", print))
525 {
526 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalFitMax", fSigMax));
527 rc = kTRUE;
528 }
529 if (IsEnvDefined(env, prefix, "BackgroundFitMax", print))
530 {
531 SetBackgroundFitMax(GetEnvValue(env, prefix, "BackgroundFitMax", fBgMax));
532 rc = kTRUE;
533 }
534 if (IsEnvDefined(env, prefix, "BackgroundFitMin", print))
535 {
536 SetBackgroundFitMin(GetEnvValue(env, prefix, "BackgroundFitMin", fBgMin));
537 rc = kTRUE;
538 }
539 if (IsEnvDefined(env, prefix, "ScaleMin", print))
540 {
541 SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
542 rc = kTRUE;
543 }
544 if (IsEnvDefined(env, prefix, "ScaleMax", print))
545 {
546 SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
547 rc = kTRUE;
548 }
549 if (IsEnvDefined(env, prefix, "PolynomOrder", print))
550 {
551 SetPolynomOrder(GetEnvValue(env, prefix, "PolynomOrder", fPolynomOrder));
552 rc = kTRUE;
553 }
554
555 if (IsEnvDefined(env, prefix, "MinimizationStrategy", print))
556 {
557 TString txt = GetEnvValue(env, prefix, "MinimizationStrategy", "");
558 txt = txt.Strip(TString::kBoth);
559 txt.ToLower();
560 if (txt==(TString)"Significance")
561 fStrategy = kSignificance;
562 if (txt==(TString)"SignificanceChi2")
563 fStrategy = kSignificanceChi2;
564 rc = kTRUE;
565 }
566
567 return rc;
568}
Note: See TracBrowser for help on using the repository browser.