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

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