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

Last change on this file since 9344 was 9302, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 29.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-2008
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// Version 2:
37// ----------
38// + Double_t fSignificanceExc; // significance of a known excess
39//
40// Version 3:
41// ----------
42// + TArrayD fErrors; // errors of coefficients
43//
44// Version 4:
45// ----------
46// + Double_t fErrorExcess;
47// - Double_t fSignificanceExc;
48//
49//
50//////////////////////////////////////////////////////////////////////////////
51#include "MAlphaFitter.h"
52
53#include <TF1.h>
54#include <TH1.h>
55#include <TH3.h>
56
57#include <TRandom.h>
58#include <TFeldmanCousins.h>
59
60#include <TLine.h>
61#include <TLatex.h>
62#include <TVirtualPad.h>
63
64#include "MMath.h"
65#include "MString.h"
66
67#include "MLogManip.h"
68
69ClassImp(MAlphaFitter);
70
71using namespace std;
72
73// --------------------------------------------------------------------------
74//
75// Default constructor.
76//
77MAlphaFitter::MAlphaFitter(const char *name, const char *title) : fSigInt(15),
78 fSigMax(75), fBgMin(45), fBgMax(85), fScaleMin(40), fScaleMax(80),
79 fPolynomOrder(2), fFitBackground(kTRUE), fFunc(0),
80 fScaleMode(kOffRegion), fScaleUser(1), fStrategy(kSignificance)
81{
82 fName = name ? name : "MAlphaFitter";
83 fTitle = title ? title : "Fit alpha";
84
85 SetSignalFunction(kGauss);
86
87 Clear();
88}
89
90// --------------------------------------------------------------------------
91//
92// Destructor
93//
94MAlphaFitter::~MAlphaFitter()
95{
96 delete fFunc;
97}
98
99// --------------------------------------------------------------------------
100//
101// Re-initializes fFunc either according to SignalFunc_t
102//
103void MAlphaFitter::SetSignalFunction(SignalFunc_t func)
104{
105 if (gROOT->GetListOfFunctions()->FindObject(""))
106 {
107 gLog << err << "MAlphaFitter::SetSignalFunction -- '' found!" << endl;
108 return;
109 }
110
111 delete fFunc;
112 fFunc = 0;
113
114 switch (func)
115 {
116 case kGauss:
117 fFunc=new TF1("", MString::Format("gaus(0) + pol%d(3)", fPolynomOrder));
118 break;
119 case kThetaSq:
120 if (fPolynomOrder>0)
121 fPolynomOrder = 1;
122 fFunc=new TF1("", "[0]*exp(-0.5*((sqrt(x)-[1])/[2])^2) + expo(3)");
123 break;
124 }
125
126 fSignalFunc=func;
127
128 fFunc->SetName("Dummy");
129 gROOT->GetListOfFunctions()->Remove(fFunc);
130
131 fCoefficients.Set(3+fPolynomOrder+1);
132 fCoefficients.Reset();
133
134 fErrors.Set(3+fPolynomOrder+1);
135 fErrors.Reset();
136}
137
138// --------------------------------------------------------------------------
139//
140// Reset variables which belong to results. Reset the arrays.
141//
142void MAlphaFitter::Clear(Option_t *o)
143{
144 fSignificance=0;
145 fErrorExcess=0;
146 fEventsExcess=0;
147 fEventsSignal=0;
148 fEventsBackground=0;
149
150 fChiSqSignal=0;
151 fChiSqBg=0;
152 fIntegralMax=0;
153 fScaleFactor=1;
154
155 fCoefficients.Reset();
156 fErrors.Reset();
157}
158
159// --------------------------------------------------------------------------
160//
161// Returns fFunc->Eval(d) or 0 if fFunc==NULL
162//
163Double_t MAlphaFitter::Eval(Double_t d) const
164{
165 return fFunc ? fFunc->Eval(d) : 0;
166}
167
168// --------------------------------------------------------------------------
169//
170// This function implementes the fit to the off-data as used in Fit()
171//
172Bool_t MAlphaFitter::FitOff(TH1D &h, Int_t paint)
173{
174 if (h.GetEntries()==0)
175 return kFALSE;
176
177 // First fit a polynom in the off region
178 fFunc->FixParameter(0, 0);
179 fFunc->FixParameter(1, 0);
180 fFunc->FixParameter(2, 1);
181 fFunc->ReleaseParameter(3);
182 if (fPolynomOrder!=1)
183 fFunc->FixParameter(4, 0);
184
185 for (int i=5; i<fFunc->GetNpar(); i++)
186 if (fFitBackground)
187 fFunc->ReleaseParameter(i);
188 else
189 fFunc->SetParameter(i, 0);
190
191 if (!fFitBackground)
192 return kTRUE;
193
194 if (fSignalFunc==kThetaSq)
195 {
196 const Double_t sum = h.Integral(1, 3)/3;
197 const Double_t a = sum<=1 ? 0 : TMath::Log(sum);
198 const Double_t b = -1.7;
199
200 // Do a best-guess
201 fFunc->SetParameter(3, a);
202 fFunc->SetParameter(4, b);
203 }
204
205 // options : N do not store the function, do not draw
206 // I use integral of function in bin rather than value at bin center
207 // R use the range specified in the function range
208 // Q quiet mode
209 // E Perform better Errors estimation using Minos technique
210 h.Fit(fFunc, "NQI", "", fBgMin, fBgMax);
211 fChiSqBg = fFunc->GetChisquare()/fFunc->GetNDF();
212
213 fCoefficients.Set(fFunc->GetNpar(), fFunc->GetParameters());
214 fErrors.Set(fFunc->GetNpar());
215 for (int i=3; i<fFunc->GetNpar(); i++)
216 fErrors[i] = fFunc->GetParError(i);
217
218 // ------------------------------------
219
220 if (paint)
221 {
222 if (paint==2)
223 {
224 fFunc->SetLineColor(kBlack);
225 fFunc->SetLineWidth(1);
226 }
227 else
228 {
229 fFunc->SetRange(0, 90);
230 fFunc->SetLineColor(kRed);
231 fFunc->SetLineWidth(2);
232 }
233 fFunc->Paint("same");
234 }
235
236 return kTRUE;
237}
238
239// --------------------------------------------------------------------------
240//
241// Calculate the result of the fit and set the corresponding data members
242//
243void MAlphaFitter::FitResult(const TH1D &h)
244{
245 const Double_t alphaw = h.GetXaxis()->GetBinWidth(1);
246
247 const Int_t bin = h.GetXaxis()->FindFixBin(fSigInt*0.999);
248
249 fIntegralMax = h.GetBinLowEdge(bin+1);
250 fEventsBackground = fFunc->Integral(0, fIntegralMax)/alphaw;
251 fEventsSignal = h.Integral(1, bin);
252 fEventsExcess = fEventsSignal-fEventsBackground;
253 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground);
254 fErrorExcess = MMath::ErrorExc(fEventsSignal, fEventsBackground);
255
256 // !Finitite includes IsNaN
257 if (!TMath::Finite(fSignificance))
258 fSignificance=0;
259}
260
261// --------------------------------------------------------------------------
262//
263// This is a preliminary implementation of a alpha-fit procedure for
264// all possible source positions. It will be moved into its own
265// more powerfull class soon.
266//
267// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
268// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
269// or
270// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
271//
272// Parameter [1] is fixed to 0 while the alpha peak should be
273// symmetric around alpha=0.
274//
275// Parameter [4] is fixed to 0 because the first derivative at
276// alpha=0 should be 0, too.
277//
278// In a first step the background is fitted between bgmin and bgmax,
279// while the parameters [0]=0 and [2]=1 are fixed.
280//
281// In a second step the signal region (alpha<sigmax) is fittet using
282// the whole function with parameters [1], [3], [4] and [5] fixed.
283//
284// The number of excess and background events are calculated as
285// s = int(hist, 0, 1.25*sigint)
286// b = int(pol2(3), 0, 1.25*sigint)
287//
288// The Significance is calculated using the Significance() member
289// function.
290//
291Bool_t MAlphaFitter::Fit(TH1D &h, Bool_t paint)
292{
293 Clear();
294
295 // Check for the region which is not filled...
296 // if (alpha0==0)
297 // return kFALSE;
298
299 // Perform fit to the off-data
300 if (!FitOff(h, paint))
301 return kFALSE;
302
303 fFunc->ReleaseParameter(0); // It is also released by SetParLimits later on
304 //func.ReleaseParameter(1); // It is also released by SetParLimits later on
305 fFunc->ReleaseParameter(2);
306 for (int i=3; i<fFunc->GetNpar(); i++)
307 fFunc->FixParameter(i, fFunc->GetParameter(i));
308
309
310 // Do not allow signals smaller than the background
311 const Double_t alpha0 = h.GetBinContent(1);
312 const Double_t s = fSignalFunc==kGauss ? fFunc->GetParameter(3) : TMath::Exp(fFunc->GetParameter(3));
313 const Double_t A = alpha0-s;
314 //const Double_t dA = TMath::Abs(A);
315 //fFunc->SetParLimits(0, -dA*4, dA*4); // SetParLimits also releases the parameter
316 fFunc->SetParLimits(2, 0, 90); // SetParLimits also releases the parameter
317
318 // Now fit a gaus in the on region on top of the polynom
319 fFunc->SetParameter(0, A);
320 fFunc->SetParameter(2, fSigMax*0.75);
321
322 // options : N do not store the function, do not draw
323 // I use integral of function in bin rather than value at bin center
324 // R use the range specified in the function range
325 // Q quiet mode
326 // E Perform better Errors estimation using Minos technique
327 h.Fit(fFunc, "NQI", "", 0, fSigMax);
328
329 fChiSqSignal = fFunc->GetChisquare()/fFunc->GetNDF();
330 fCoefficients.Set(fFunc->GetNpar(), fFunc->GetParameters());
331 for (int i=0; i<3; i++)
332 fErrors[i] = fFunc->GetParError(i);
333 //const Bool_t ok = NDF>0 && chi2<2.5*NDF;
334
335 // ------------------------------------
336 if (paint)
337 {
338 fFunc->SetLineColor(kGreen);
339 fFunc->SetLineWidth(2);
340 fFunc->Paint("same");
341 }
342 // ------------------------------------
343
344 //const Double_t s = fFunc->Integral(0, fSigInt)/alphaw;
345 fFunc->SetParameter(0, 0);
346 fFunc->SetParameter(2, 1);
347 //const Double_t b = fFunc->Integral(0, fSigInt)/alphaw;
348 //fSignificance = MMath::SignificanceLiMaSigned(s, b);
349
350 // Calculate the fit result and set the corresponding data members
351 FitResult(h);
352
353 return kTRUE;
354}
355
356Double_t MAlphaFitter::DoOffFit(const TH1D &hon, const TH1D &hof, Bool_t paint)
357{
358 if (fSignalFunc!=kThetaSq)
359 return 0;
360
361 // ----------------------------------------------------------------------------
362
363 const Int_t bin = hon.GetXaxis()->FindFixBin(fSigInt*0.999);
364
365 MAlphaFitter fit(*this);
366 fit.EnableBackgroundFit();
367 fit.SetBackgroundFitMin(0);
368
369 // produce a histogram containing the off-samples from on-source and
370 // off-source in the off-source region and the on-data in the source-region
371 TH1D h(hof);
372 h.SetDirectory(0);
373 h.Add(&hon);
374
375 h.Scale(0.5);
376 for (int i=1; i<=bin+3; i++)
377 {
378 h.SetBinContent(i, hof.GetBinContent(i));
379 h.SetBinError( i, hof.GetBinError(i));
380 }
381
382 // Now fit the off-data
383 if (!fit.FitOff(h, paint?2:0)) // FIXME: Show fit!
384 return -1;
385
386 // Calculate fit-result
387 fit.FitResult(h);
388
389 // Do a gaussian error propagation to calculated the error of
390 // the background estimated from the fit
391 const Double_t ea = fit.fErrors[3];
392 const Double_t eb = fit.fErrors[4];
393 const Double_t a = fit.fCoefficients[3];
394 const Double_t b = fit.fCoefficients[4];
395
396 const Double_t t = fIntegralMax;
397
398 const Double_t ex = TMath::Exp(t*b);
399 const Double_t eab = TMath::Exp(a)/b;
400
401 const Double_t eA = ex-1;
402 const Double_t eB = t*ex - eA/b;
403
404 const Double_t w = h.GetXaxis()->GetBinWidth(1);
405
406 // Error of estimated background
407 const Double_t er = TMath::Abs(eab)*TMath::Hypot(eA*ea, eB*eb)/w;
408
409 // Calculate arbitrary scale factor from propagated error from the
410 // condition: sqrt(alpha*background) = est.background/est.error
411 // const Double_t bg = hof.Integral(1, bin);
412 // const Double_t sc = bg * er*er / (fit2.GetEventsBackground()*fit2.GetEventsBackground());
413 // Assuming that bg and fit2.GetEventsBackground() are rather identical:
414 const Double_t sc = er*er / fit.fEventsBackground;
415
416
417 /*
418 cout << MMath::SignificanceLiMaSigned(hon.Integral(1, bin), fit.GetEventsBackground()/sc, sc) << " ";
419 cout << sc << " ";
420 cout << fit.fChiSqBg << endl;
421 */
422 return sc;
423}
424
425Bool_t MAlphaFitter::Fit(const TH1D &hon, const TH1D &hof, Double_t alpha, Bool_t paint)
426{
427 TH1D h(hon);
428 h.SetDirectory(0);
429 h.Add(&hof, -1); // substracts also number of entries!
430 h.SetEntries(hon.GetEntries());
431
432 MAlphaFitter fit(*this);
433 fit.SetPolynomOrder(0);
434 if (alpha<=0 || !fit.Fit(h, paint))
435 return kFALSE;
436
437 fChiSqSignal = fit.fChiSqSignal;
438 fChiSqBg = fit.fChiSqBg;
439 fCoefficients = fit.fCoefficients;
440 fErrors = fit.fErrors;
441
442 // ----------------------------------------------------------------------------
443
444 const Double_t scale = DoOffFit(hon, hof, paint);
445 if (scale<0)
446 return kFALSE;
447
448 // ----------------------------------------------------------------------------
449
450 const Int_t bin = hon.GetXaxis()->FindFixBin(fSigInt*0.999);
451
452 fIntegralMax = hon.GetBinLowEdge(bin+1);
453 fEventsBackground = hof.Integral(1, bin);
454 fEventsSignal = hon.Integral(1, bin);
455 fEventsExcess = fEventsSignal-fEventsBackground;
456 fScaleFactor = alpha;
457 fSignificance = MMath::SignificanceLiMaSigned(fEventsSignal, fEventsBackground/alpha, alpha);
458 fErrorExcess = MMath::ErrorExc(fEventsSignal, fEventsBackground/alpha, alpha);
459
460 // !Finitite includes IsNaN
461 if (!TMath::Finite(fSignificance))
462 fSignificance=0;
463
464 return kTRUE;
465}
466
467// --------------------------------------------------------------------------
468//
469// Calculate the upper limit for fEventsSignal number of observed events
470// and fEventsBackground number of background events.
471//
472// Therefor TFeldmanCousin is used.
473//
474// The Feldman-Cousins method as described in PRD V57 #7, p3873-3889
475//
476Double_t MAlphaFitter::CalcUpperLimit() const
477{
478 // get a FeldmanCousins calculation object with the default limits
479 // of calculating a 90% CL with the minimum signal value scanned
480 // = 0.0 and the maximum signal value scanned of 50.0
481 TFeldmanCousins f;
482 f.SetMuStep(0.05);
483 f.SetMuMax(100);
484 f.SetMuMin(0);
485 f.SetCL(90);
486
487 return f.CalculateUpperLimit(fEventsSignal, fEventsBackground);
488}
489
490void MAlphaFitter::PaintResult(Float_t x, Float_t y, Float_t size, Bool_t draw) const
491{
492 const Double_t w = GetGausSigma();
493 const Double_t m = fIntegralMax;
494
495 const Int_t l1 = w<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(w));
496 const Int_t l2 = m<=0 ? 0 : (Int_t)TMath::Ceil(-TMath::Log10(m));
497 const TString fmt = MString::Format("\\sigma_{L/M}=%%.1f \\omega=%%.%df\\circ E=%%d B=%%d x<%%.%df \\tilde\\chi_{b}=%%.1f \\tilde\\chi_{s}=%%.1f c=%%.1f f=%%.2f",
498 l1<1?1:l1+1, l2<1?1:l2+1);
499 const TString txt = MString::Format(fmt.Data(), fSignificance, w, (int)fEventsExcess,
500 (int)fEventsBackground, m, fChiSqBg, fChiSqSignal,
501 fCoefficients[3], fScaleFactor);
502
503 // This is totaly weired but the only way to get both options
504 // working with this nonsense implementation of TLatex
505 TLatex text(x, y, txt);
506 text.SetBit(TLatex::kTextNDC);
507 text.SetTextSize(size);
508 if (draw)
509 text.DrawLatex(x, y, txt);
510 else
511 text.Paint();
512
513 TLine line;
514 line.SetLineColor(14);
515 if (draw)
516 line.DrawLine(m, gPad->GetUymin(), m, gPad->GetUymax());
517 else
518 line.PaintLine(m, gPad->GetUymin(), m, gPad->GetUymax());
519}
520
521void MAlphaFitter::Copy(TObject &o) const
522{
523 MAlphaFitter &f = static_cast<MAlphaFitter&>(o);
524
525 // Setup
526 f.fSigInt = fSigInt;
527 f.fSigMax = fSigMax;
528 f.fBgMin = fBgMin;
529 f.fBgMax = fBgMax;
530 f.fScaleMin = fScaleMin;
531 f.fScaleMax = fScaleMax;
532 f.fPolynomOrder = fPolynomOrder;
533 f.fFitBackground= fFitBackground;
534 f.fSignalFunc = fSignalFunc;
535 f.fScaleMode = fScaleMode;
536 f.fScaleUser = fScaleUser;
537 f.fStrategy = fStrategy;
538 f.fCoefficients.Set(fCoefficients.GetSize());
539 f.fCoefficients.Reset();
540 f.fErrors.Set(fCoefficients.GetSize());
541 f.fErrors.Reset();
542
543 // Result
544 f.fSignificance = fSignificance;
545 f.fErrorExcess = fErrorExcess;
546 f.fEventsExcess = fEventsExcess;
547 f.fEventsSignal = fEventsSignal;
548 f.fEventsBackground = fEventsBackground;
549 f.fChiSqSignal = fChiSqSignal;
550 f.fChiSqBg = fChiSqBg;
551 f.fIntegralMax = fIntegralMax;
552 f.fScaleFactor = fScaleFactor;
553
554 // Function
555 delete f.fFunc;
556
557 f.fFunc = new TF1(*fFunc);
558 f.fFunc->SetName("Dummy");
559 gROOT->GetListOfFunctions()->Remove(f.fFunc);
560}
561
562void MAlphaFitter::Print(Option_t *o) const
563{
564 *fLog << GetDescriptor() << ": Fitting..." << endl;
565 *fLog << " ...signal to " << fSigMax << " (integrate into bin at " << fSigInt << ")" << endl;
566 *fLog << " ...signal function: ";
567 switch (fSignalFunc)
568 {
569 case kGauss: *fLog << "gauss(x)/pol" << fPolynomOrder; break;
570 case kThetaSq: *fLog << "gauss(sqrt(x))/expo"; break;
571 }
572 *fLog << endl;
573 if (!fFitBackground)
574 *fLog << " ...no background." << endl;
575 else
576 {
577 *fLog << " ...background from " << fBgMin << " to " << fBgMax << endl;
578 *fLog << " ...polynom order " << fPolynomOrder << endl;
579 *fLog << " ...scale mode: ";
580 switch (fScaleMode)
581 {
582 case kNone: *fLog << "none."; break;
583 case kEntries: *fLog << "entries."; break;
584 case kIntegral: *fLog << "integral."; break;
585 case kOffRegion: *fLog << "off region (integral between " << fScaleMin << " and " << fScaleMax << ")"; break;
586 case kBackground: *fLog << "background (integral between " << fBgMin << " and " << fBgMax << ")"; break;
587 case kLeastSquare: *fLog << "least square (N/A)"; break;
588 case kUserScale: *fLog << "user def (" << fScaleUser << ")"; break;
589 }
590 *fLog << endl;
591 }
592
593 if (TString(o).Contains("result"))
594 {
595 *fLog << "Result:" << endl;
596 *fLog << " - Significance (Li/Ma) " << fSignificance << endl;
597 *fLog << " - Excess Events " << fEventsExcess << endl;
598 *fLog << " - Signal Events " << fEventsSignal << endl;
599 *fLog << " - Background Events " << fEventsBackground << endl;
600 *fLog << " - Chi^2/ndf (Signal) " << fChiSqSignal << endl;
601 *fLog << " - Chi^2/ndf (Background) " << fChiSqBg << endl;
602 *fLog << " - Signal integrated up to " << fIntegralMax << "°" << endl;
603 *fLog << " - Scale Factor (Off) " << fScaleFactor << endl;
604 }
605}
606
607Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, UInt_t bin, Bool_t paint)
608{
609 const TString name(MString::Format("TempAlphaEnergy%06d", gRandom->Integer(1000000)));
610
611 TH1D *h = hon.ProjectionZ(name, 0, hon.GetNbinsX()+1, bin, bin, "E");
612 h->SetDirectory(0);
613
614 const Bool_t rc = Fit(*h, paint);
615
616 delete h;
617
618 return rc;
619}
620
621Bool_t MAlphaFitter::FitTheta(const TH3D &hon, UInt_t bin, Bool_t paint)
622{
623 const TString name(MString::Format("TempAlphaTheta%06d", gRandom->Integer(1000000)));
624
625 TH1D *h = hon.ProjectionZ(name, bin, bin, 0, hon.GetNbinsY()+1, "E");
626 h->SetDirectory(0);
627
628 const Bool_t rc = Fit(*h, paint);
629
630 delete h;
631
632 return rc;
633}
634/*
635Bool_t MAlphaFitter::FitTime(const TH3D &hon, UInt_t bin, Bool_t paint)
636{
637 const TString name(Form("TempAlphaTime%06d", gRandom->Integer(1000000)));
638
639 hon.GetZaxis()->SetRange(bin,bin);
640 TH1D *h = (TH1D*)hon.Project3D("ye");
641 hon.GetZaxis()->SetRange(-1,-1);
642
643 h->SetDirectory(0);
644
645 const Bool_t rc = Fit(*h, paint);
646 delete h;
647 return rc;
648}
649*/
650Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, Bool_t paint)
651{
652 const TString name(MString::Format("TempAlpha%06d", gRandom->Integer(1000000)));
653
654 TH1D *h = hon.ProjectionZ(name, 0, hon.GetNbinsX()+1, 0, hon.GetNbinsY()+1, "E");
655 h->SetDirectory(0);
656
657 const Bool_t rc = Fit(*h, paint);
658
659 delete h;
660
661 return rc;
662}
663
664Bool_t MAlphaFitter::FitEnergy(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
665{
666 const TString name1(MString::Format("TempAlpha%06d_on", gRandom->Integer(1000000)));
667 const TString name0(MString::Format("TempAlpha%06d_off", gRandom->Integer(1000000)));
668
669 TH1D *h1 = hon.ProjectionZ(name1, 0, hon.GetNbinsX()+1, bin, bin, "E");
670 h1->SetDirectory(0);
671
672 TH1D *h0 = hof.ProjectionZ(name0, 0, hof.GetNbinsX()+1, bin, bin, "E");
673 h0->SetDirectory(0);
674
675 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
676
677 delete h0;
678 delete h1;
679
680 return rc;
681}
682
683Bool_t MAlphaFitter::FitTheta(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
684{
685 const TString name1(MString::Format("TempAlpha%06d_on", gRandom->Integer(1000000)));
686 const TString name0(MString::Format("TempAlpha%06d_off", gRandom->Integer(1000000)));
687
688 TH1D *h1 = hon.ProjectionZ(name1, bin, bin, 0, hon.GetNbinsY()+1, "E");
689 h1->SetDirectory(0);
690
691 TH1D *h0 = hof.ProjectionZ(name0, bin, bin, 0, hof.GetNbinsY()+1, "E");
692 h0->SetDirectory(0);
693
694 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
695
696 delete h0;
697 delete h1;
698
699 return rc;
700}
701/*
702Bool_t MAlphaFitter::FitTime(const TH3D &hon, const TH3D &hof, UInt_t bin, Bool_t paint)
703{
704 const TString name1(Form("TempAlphaTime%06d_on", gRandom->Integer(1000000)));
705 const TString name0(Form("TempAlphaTime%06d_off", gRandom->Integer(1000000)));
706
707 hon.GetZaxis()->SetRange(bin,bin);
708 TH1D *h1 = (TH1D*)hon.Project3D("ye");
709 hon.GetZaxis()->SetRange(-1,-1);
710 h1->SetDirectory(0);
711
712 hof.GetZaxis()->SetRange(bin,bin);
713 TH1D *h0 = (TH1D*)hof.Project3D("ye");
714 hof.GetZaxis()->SetRange(-1,-1);
715 h0->SetDirectory(0);
716
717 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
718
719 delete h0;
720 delete h1;
721
722 return rc;
723}
724*/
725
726Bool_t MAlphaFitter::FitAlpha(const TH3D &hon, const TH3D &hof, Bool_t paint)
727{
728 const TString name1(MString::Format("TempAlpha%06d_on", gRandom->Integer(1000000)));
729 const TString name0(MString::Format("TempAlpha%06d_off", gRandom->Integer(1000000)));
730
731 TH1D *h1 = hon.ProjectionZ(name1, 0, hon.GetNbinsX()+1, 0, hon.GetNbinsY()+1, "E");
732 h1->SetDirectory(0);
733
734 TH1D *h0 = hof.ProjectionZ(name0, 0, hof.GetNbinsX()+1, 0, hof.GetNbinsY()+1, "E");
735 h0->SetDirectory(0);
736
737 const Bool_t rc = ScaleAndFit(*h1, h0, paint);
738
739 delete h0;
740 delete h1;
741
742 return rc;
743}
744
745Bool_t MAlphaFitter::ApplyScaling(const TH3D &hon, TH3D &hof, UInt_t bin) const
746{
747 const TString name1(MString::Format("TempAlpha%06d_on", gRandom->Integer(1000000)));
748 const TString name0(MString::Format("TempAlpha%06d_off", gRandom->Integer(1000000)));
749
750 TH1D *h1 = hon.ProjectionZ(name1, 0, hon.GetNbinsX()+1, bin, bin, "E");
751 h1->SetDirectory(0);
752
753 TH1D *h0 = hof.ProjectionZ(name0, 0, hof.GetNbinsX()+1, bin, bin, "E");
754 h0->SetDirectory(0);
755
756 const Double_t scale = Scale(*h0, *h1);
757
758 delete h0;
759 delete h1;
760
761 for (int x=0; x<=hof.GetNbinsX()+1; x++)
762 for (int z=0; z<=hof.GetNbinsZ()+1; z++)
763 {
764 hof.SetBinContent(x, bin, z, hof.GetBinContent(x, bin, z)*scale);
765 hof.SetBinError( x, bin, z, hof.GetBinError( x, bin, z)*scale);
766 }
767
768 return scale>0;
769}
770
771Bool_t MAlphaFitter::ApplyScaling(const TH3D &hon, TH3D &hof) const
772{
773 for (int y=0; y<=hof.GetNbinsY()+1; y++)
774 ApplyScaling(hon, hof, y);
775
776 return kTRUE;
777}
778
779Double_t MAlphaFitter::Scale(TH1D &of, const TH1D &on) const
780{
781 Float_t scaleon = 1;
782 Float_t scaleof = 1;
783 switch (fScaleMode)
784 {
785 case kNone:
786 return 1;
787
788 case kEntries:
789 scaleon = on.GetEntries();
790 scaleof = of.GetEntries();
791 break;
792
793 case kIntegral:
794 scaleon = on.Integral();
795 scaleof = of.Integral();
796 break;
797
798 case kOffRegion:
799 {
800 const Int_t min = on.GetXaxis()->FindFixBin(fScaleMin);
801 const Int_t max = on.GetXaxis()->FindFixBin(fScaleMax);
802 scaleon = on.Integral(min, max);
803 scaleof = of.Integral(min, max);
804 }
805 break;
806
807 case kBackground:
808 {
809 const Int_t min = on.GetXaxis()->FindFixBin(fBgMin);
810 const Int_t max = on.GetXaxis()->FindFixBin(fBgMax);
811 scaleon = on.Integral(min, max);
812 scaleof = of.Integral(min, max);
813 }
814 break;
815
816 case kUserScale:
817 scaleon = fScaleUser;
818 break;
819
820 // This is just to make some compiler happy
821 default:
822 return 1;
823 }
824
825 if (scaleof!=0)
826 {
827 of.Scale(scaleon/scaleof);
828 return scaleon/scaleof;
829 }
830 else
831 {
832 of.Reset();
833 return 0;
834 }
835}
836
837Double_t MAlphaFitter::GetMinimizationValue() const
838{
839 switch (fStrategy)
840 {
841 case kSignificance:
842 return -GetSignificance();
843 case kSignificanceChi2:
844 return -GetSignificance()/GetChiSqSignal();
845 case kSignificanceLogExcess:
846 if (GetEventsExcess()<1)
847 return 0;
848 return -GetSignificance()*TMath::Log10(GetEventsExcess());
849 case kSignificanceExcess:
850 return -GetSignificance()*GetEventsExcess();
851 case kExcess:
852 return -GetEventsExcess();
853 case kGaussSigma:
854 return GetGausSigma();
855 case kWeakSource:
856 return GetEventsBackground()<1 ? -GetEventsExcess() : -GetEventsExcess()/TMath::Sqrt(GetEventsBackground());
857 case kWeakSourceLogExcess:
858 return GetEventsBackground()<1 ? -GetEventsExcess() : -GetEventsExcess()/TMath::Sqrt(GetEventsBackground())*TMath::Log10(GetEventsExcess());
859 }
860 return 0;
861}
862
863Int_t MAlphaFitter::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
864{
865 Bool_t rc = kFALSE;
866
867 //void SetScaleUser(Float_t scale) { fScaleUser = scale; fScaleMode=kUserScale; }
868 //void SetScaleMode(ScaleMode_t mode) { fScaleMode = mode; }
869
870 if (IsEnvDefined(env, prefix, "SignalIntegralMax", print))
871 {
872 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalIntegralMax", fSigInt));
873 rc = kTRUE;
874 }
875 if (IsEnvDefined(env, prefix, "SignalFitMax", print))
876 {
877 SetSignalIntegralMax(GetEnvValue(env, prefix, "SignalFitMax", fSigMax));
878 rc = kTRUE;
879 }
880 if (IsEnvDefined(env, prefix, "BackgroundFitMax", print))
881 {
882 SetBackgroundFitMax(GetEnvValue(env, prefix, "BackgroundFitMax", fBgMax));
883 rc = kTRUE;
884 }
885 if (IsEnvDefined(env, prefix, "BackgroundFitMin", print))
886 {
887 SetBackgroundFitMin(GetEnvValue(env, prefix, "BackgroundFitMin", fBgMin));
888 rc = kTRUE;
889 }
890 if (IsEnvDefined(env, prefix, "ScaleMin", print))
891 {
892 SetScaleMin(GetEnvValue(env, prefix, "ScaleMin", fScaleMin));
893 rc = kTRUE;
894 }
895 if (IsEnvDefined(env, prefix, "ScaleMax", print))
896 {
897 SetScaleMax(GetEnvValue(env, prefix, "ScaleMax", fScaleMax));
898 rc = kTRUE;
899 }
900 if (IsEnvDefined(env, prefix, "PolynomOrder", print))
901 {
902 SetPolynomOrder(GetEnvValue(env, prefix, "PolynomOrder", fPolynomOrder));
903 rc = kTRUE;
904 }
905
906 if (IsEnvDefined(env, prefix, "MinimizationStrategy", print))
907 {
908 TString txt = GetEnvValue(env, prefix, "MinimizationStrategy", "");
909 txt = txt.Strip(TString::kBoth);
910 txt.ToLower();
911 if (txt==(TString)"significance")
912 fStrategy = kSignificance;
913 if (txt==(TString)"significancechi2")
914 fStrategy = kSignificanceChi2;
915 if (txt==(TString)"significanceexcess")
916 fStrategy = kSignificanceExcess;
917 if (txt==(TString)"excess")
918 fStrategy = kExcess;
919 if (txt==(TString)"gausssigma" || txt==(TString)"gaussigma")
920 fStrategy = kGaussSigma;
921 if (txt==(TString)"weaksource")
922 fStrategy = kWeakSource;
923 rc = kTRUE;
924 }
925 if (IsEnvDefined(env, prefix, "Scale", print))
926 {
927 fScaleUser = GetEnvValue(env, prefix, "Scale", fScaleUser);
928 rc = kTRUE;
929 }
930 if (IsEnvDefined(env, prefix, "ScaleMode", print))
931 {
932 TString txt = GetEnvValue(env, prefix, "ScaleMode", "");
933 txt = txt.Strip(TString::kBoth);
934 txt.ToLower();
935 if (txt==(TString)"none")
936 fScaleMode = kNone;
937 if (txt==(TString)"entries")
938 fScaleMode = kEntries;
939 if (txt==(TString)"integral")
940 fScaleMode = kIntegral;
941 if (txt==(TString)"offregion")
942 fScaleMode = kOffRegion;
943 if (txt==(TString)"background")
944 fScaleMode = kBackground;
945 if (txt==(TString)"leastsquare")
946 fScaleMode = kLeastSquare;
947 if (txt==(TString)"userscale")
948 fScaleMode = kUserScale;
949 if (txt==(TString)"fixed")
950 FixScale();
951 rc = kTRUE;
952 }
953 if (IsEnvDefined(env, prefix, "SignalFunction", print))
954 {
955 TString txt = GetEnvValue(env, prefix, "SignalFunction", "");
956 txt = txt.Strip(TString::kBoth);
957 txt.ToLower();
958 if (txt==(TString)"gauss" || txt==(TString)"gaus")
959 SetSignalFunction(kGauss);
960 if (txt==(TString)"thetasq")
961 SetSignalFunction(kThetaSq);
962 rc = kTRUE;
963 }
964
965 return rc;
966}
Note: See TracBrowser for help on using the repository browser.