source: trunk/MagicSoft/Mars/mhist/MHFalseSource.cc@ 3571

Last change on this file since 3571 was 3568, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 32.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// MHFalseSource
28//
29// Create a false source plot. For the Binning in x,y the object MBinning
30// "BinningFalseSource" is taken from the paremeter list.
31//
32// The binning in alpha is currently fixed as 18bins from 0 to 90deg.
33//
34// If MTime, MObservatory and MPointingPos is available in the paremeter
35// list each image is derotated.
36//
37// MHFalseSource fills a 3D histogram with alpha distribution for
38// each (derotated) x and y position.
39//
40// Only a radius of 1.2deg around the camera center is taken into account.
41//
42// The displayed histogram is the projection of alpha<fAlphaCut into
43// the x,y plain and the projection of alpha>90-fAlphaCut
44//
45// By using the context menu (find it in a small region between the single
46// pads) you can change the AlphaCut 'online'
47//
48// Here is a slightly simplified version of the algorithm:
49// ------------------------------------------------------
50// MHillas hil; // Taken as second argument in MFillH
51//
52// MSrcPosCam src;
53// MHillasSrc hsrc;
54// hsrc.SetSrcPos(&src);
55//
56// for (int ix=0; ix<nx; ix++)
57// for (int iy=0; iy<ny; iy++)
58// {
59// TVector2 v(cx[ix], cy[iy]); //cx center of x-bin ix
60// if (rho!=0) //cy center of y-bin iy
61// v.Rotate(-rho); //image rotation angle
62//
63// src.SetXY(v); //source position in the camera
64//
65// if (!hsrc.Calc(hil)) //calc source dependant hillas
66// return;
67//
68// //fill absolute alpha into histogram
69// const Double_t alpha = hsrc.GetAlpha();
70// fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
71// }
72// }
73//
74// Use MHFalse Source like this:
75// -----------------------------
76// MFillH fill("MHFalseSource", "MHillas");
77// or
78// MHFalseSource hist;
79// hist.SetAlphaCut(12.5); // The default value
80// hist.SetBgmean(55); // The default value
81// MFillH fill(&hist, "MHillas");
82//
83// To be implemented:
84// ------------------
85// - a switch to use alpha or |alpha|
86// - taking the binning for alpha from the parlist (binning is
87// currently fixed)
88// - a possibility to change the fit-function (eg using a different TF1)
89// - a possibility to change the fit algorithm (eg which paremeters
90// are fixed at which time)
91// - use a different varaible than alpha
92// - a possiblity to change the algorithm which is used to calculate
93// alpha (or another parameter) is missing (this could be done using
94// a tasklist somewhere...)
95// - a more clever (and faster) algorithm to fill the histogram, eg by
96// calculating alpha once and fill the two coils around the mean
97// - more drawing options...
98// - Move Significance() to a more 'general' place and implement
99// also different algorithms like (Li/Ma)
100// - implement fit for best alpha distribution -- online (Paint)
101//
102//////////////////////////////////////////////////////////////////////////////
103#include "MHFalseSource.h"
104
105#include <TF1.h>
106#include <TH2.h>
107#include <TGraph.h>
108#include <TStyle.h>
109#include <TCanvas.h>
110#include <TPaveText.h>
111#include <TStopwatch.h>
112
113#include "MGeomCam.h"
114#include "MSrcPosCam.h"
115#include "MHillasSrc.h"
116#include "MTime.h"
117#include "MObservatory.h"
118#include "MPointingPos.h"
119#include "MAstroSky2Local.h"
120#include "MStatusDisplay.h"
121
122#include "MBinning.h"
123#include "MParList.h"
124
125#include "MLog.h"
126#include "MLogManip.h"
127
128ClassImp(MHFalseSource);
129
130using namespace std;
131
132// --------------------------------------------------------------------------
133//
134// Default Constructor
135//
136MHFalseSource::MHFalseSource(const char *name, const char *title)
137 : fMm2Deg(-1), fUseMmScale(kTRUE), fAlphaCut(12.5), fBgMean(55)
138{
139 //
140 // set the name and title of this object
141 //
142 fName = name ? name : "MHFalseSource";
143 fTitle = title ? title : "3D-plot of Alpha vs x, y";
144
145 fHist.SetDirectory(NULL);
146
147 fHist.SetName("Alpha");
148 fHist.SetTitle("3D-plot of Alpha vs x, y");
149 fHist.SetXTitle("x [\\circ]");
150 fHist.SetYTitle("y [\\circ]");
151 fHist.SetZTitle("\\alpha [\\circ]");
152}
153
154// --------------------------------------------------------------------------
155//
156// Set the alpha cut (|alpha|<fAlphaCut) which is use to estimate the
157// number of excess events
158//
159void MHFalseSource::SetAlphaCut(Float_t alpha)
160{
161 if (alpha<0)
162 *fLog << warn << "Alpha<0... taking absolute value." << endl;
163
164 fAlphaCut = TMath::Abs(alpha);
165
166 Modified();
167}
168
169// --------------------------------------------------------------------------
170//
171// Set mean alpha around which the off sample is determined
172// (fBgMean-fAlphaCut/2<|fAlpha|<fBgMean+fAlphaCut/2) which is use
173// to estimate the number of off events
174//
175void MHFalseSource::SetBgMean(Float_t alpha)
176{
177 if (alpha<0)
178 *fLog << warn << "Alpha<0... taking absolute value." << endl;
179
180 fBgMean = TMath::Abs(alpha);
181
182 Modified();
183}
184
185// --------------------------------------------------------------------------
186//
187// Use this function to setup your own conversion factor between degrees
188// and millimeters. The conversion factor should be the one calculated in
189// MGeomCam. Use this function with Caution: You could create wrong values
190// by setting up your own scale factor.
191//
192void MHFalseSource::SetMm2Deg(Float_t mmdeg)
193{
194 if (mmdeg<0)
195 {
196 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
197 return;
198 }
199
200 if (fMm2Deg>=0)
201 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
202
203 fMm2Deg = mmdeg;
204}
205
206// --------------------------------------------------------------------------
207//
208// With this function you can convert the histogram ('on the fly') between
209// degrees and millimeters.
210//
211void MHFalseSource::SetMmScale(Bool_t mmscale)
212{
213 if (fUseMmScale == mmscale)
214 return;
215
216 if (fMm2Deg<0)
217 {
218 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
219 return;
220 }
221
222 if (fUseMmScale)
223 {
224 fHist.SetXTitle("x [mm]");
225 fHist.SetYTitle("y [mm]");
226
227 fHist.Scale(1./fMm2Deg);
228 }
229 else
230 {
231 fHist.SetXTitle("x [\\circ]");
232 fHist.SetYTitle("y [\\circ]");
233
234 fHist.Scale(1./fMm2Deg);
235 }
236
237 fUseMmScale = mmscale;
238}
239
240// --------------------------------------------------------------------------
241//
242// Calculate Significance as
243// significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
244//
245// s: total number of events in signal region
246// b: number of background events in signal region
247//
248Double_t MHFalseSource::Significance(Double_t s, Double_t b)
249{
250 const Double_t k = b==0 ? 0 : s/b;
251 const Double_t f = s+k*k*b;
252
253 return f==0 ? 0 : (s-b)/TMath::Sqrt(f);
254}
255
256// --------------------------------------------------------------------------
257//
258// calculates the significance according to Li & Ma
259// ApJ 272 (1983) 317
260//
261// s: total number of events in signal region
262//
263// Double_t fGamma; // Nbg = Non - gamma * Noff
264// - the effective number of background events (fNoff), and fGamma :
265//
266// fGamma = b / fNoff;
267// fGamma = fdNbg / sqrt(fNoff);
268// fGamma = fdNbg*fdNbg / fNbg;
269// fNoff = b*b / (fdNbg*fdNbg);
270// Double_t fNbg; // number of background events in signal region
271// b = fNOff *fGamma;
272/*
273Double_t MHFindSignificance::SignificanceLiMa(Double_t non, Double_t noff, Double_t gamma)
274{
275 if (gamma <= 0.0 || non <= 0.0 || noff <= 0.0)
276 {
277 *siglima = 0.0;
278 return kFALSE;
279 }
280
281 Double_t help1 = TMath::Log( (gamma+1)*s / (gamma*(s+noff)) );
282 Double_t help2 = TMath::Log( (gamma+1)*noff / ( s+noff ) );
283
284 Double_t siglima = TMath::Sqrt((s*help1+noff*help2)*2);
285
286 return non<gamma*noff ? -siglima : siglima;
287}
288*/
289// --------------------------------------------------------------------------
290//
291// Set binnings (takes BinningFalseSource) and prepare filling of the
292// histogram.
293//
294// Also search for MTime, MObservatory and MPointingPos
295//
296Bool_t MHFalseSource::SetupFill(const MParList *plist)
297{
298 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
299 if (geom)
300 {
301 fMm2Deg = geom->GetConvMm2Deg();
302 fUseMmScale = kFALSE;
303
304 fHist.SetXTitle("x [\\circ]");
305 fHist.SetYTitle("y [\\circ]");
306 }
307
308 MBinning binsa;
309 binsa.SetEdges(18, 0, 90);
310
311 const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource");
312 if (!bins)
313 {
314 Float_t r = geom ? geom->GetMaxRadius() : 600;
315 r /= 3;
316 if (!fUseMmScale)
317 r *= fMm2Deg;
318
319 MBinning b;
320 b.SetEdges(20, -r, r);
321 SetBinning(&fHist, &b, &b, &binsa);
322 }
323 else
324 SetBinning(&fHist, bins, bins, &binsa);
325
326 fPointPos = (MPointingPos*)plist->FindObject(AddSerialNumber("MPointingPos"));
327 if (!fPointPos)
328 *fLog << warn << "MPointingPos not found... no derotation." << endl;
329
330 fTime = (MTime*)plist->FindObject(AddSerialNumber("MTime"));
331 if (!fTime)
332 *fLog << warn << "MTime not found... no derotation." << endl;
333
334 fObservatory = (MObservatory*)plist->FindObject(AddSerialNumber("MObservatory"));
335 if (!fObservatory)
336 *fLog << err << "MObservatory not found... no derotation." << endl;
337
338
339 return kTRUE;
340}
341
342// --------------------------------------------------------------------------
343//
344// Fill the histogram. For details see the code or the class description
345//
346Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w)
347{
348 MHillas *hil = (MHillas*)par;
349
350 const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
351
352 // Get camera rotation angle
353 Double_t rho = 0;
354 if (fTime && fObservatory && fPointPos)
355 rho = fPointPos->RotationAngle(*fObservatory, *fTime);
356 //if (fPointPos)
357 // rho = fPointPos->RotationAngle(*fObservatory);
358
359 MSrcPosCam src;
360 MHillasSrc hsrc;
361 hsrc.SetSrcPos(&src);
362
363 const Int_t nx = fHist.GetNbinsX();
364 const Int_t ny = fHist.GetNbinsY();
365
366 Axis_t cx[nx];
367 Axis_t cy[ny];
368 fHist.GetXaxis()->GetCenter(cx);
369 fHist.GetYaxis()->GetCenter(cy);
370
371 for (int ix=0; ix<nx; ix++)
372 {
373 for (int iy=0; iy<ny; iy++)
374 {
375 if (TMath::Hypot(cx[ix], cy[iy])>maxr)
376 continue;
377
378 TVector2 v(cx[ix], cy[iy]);
379 if (rho!=0)
380 v.Rotate(-rho);
381
382 if (!fUseMmScale)
383 v *= 1./fMm2Deg;
384
385 src.SetXY(v);
386
387 if (!hsrc.Calc(hil))
388 {
389 *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl;
390 return kFALSE;
391 }
392
393 const Double_t alpha = hsrc.GetAlpha();
394 fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
395 }
396 }
397
398 return kTRUE;
399}
400
401// --------------------------------------------------------------------------
402//
403// Create projection for off data, taking sum of bin contents of
404// range (fBgMean-fAlphaCut/2, fBgMean+fAlphaCut) Making sure to take
405// the same number of bins than for on-data
406//
407void MHFalseSource::ProjectOff(TH2D *h2)
408{
409 TAxis &axe = *fHist.GetZaxis();
410
411 // Find range to cut (left edge and width)
412 const Int_t f = axe.FindBin(fBgMean-fAlphaCut/2);
413 const Int_t l = axe.FindBin(fAlphaCut)+f-1;
414
415 axe.SetRange(f, l);
416 const Float_t cut1 = axe.GetBinLowEdge(f);
417 const Float_t cut2 = axe.GetBinUpEdge(l);
418 h2->SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2));
419
420 // Get projection for range
421 TH2D *p = (TH2D*)fHist.Project3D("xy_off");
422
423 // Reset range
424 axe.SetRange(0,9999);
425
426 // Move contents from projection to h2
427 h2->Reset();
428 h2->Add(p);
429
430 // Delete p
431 delete p;
432
433 // Set Minimum as minimum value Greater Than 0
434 h2->SetMinimum(GetMinimumGT(*h2));
435}
436
437// --------------------------------------------------------------------------
438//
439// Create projection for on data, taking sum of bin contents of
440// range (0, fAlphaCut)
441//
442void MHFalseSource::ProjectOn(TH2D *h3)
443{
444 TAxis &axe = *fHist.GetZaxis();
445
446 // Find range to cut
447 axe.SetRangeUser(0, fAlphaCut);
448 const Float_t cut = axe.GetBinUpEdge(axe.GetLast());
449 h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut));
450
451 // Get projection for range
452 TH2D *p = (TH2D*)fHist.Project3D("xy_on");
453
454 // Reset range
455 axe.SetRange(0,9999);
456
457 // Move contents from projection to h3
458 h3->Reset();
459 h3->Add(p);
460 delete p;
461
462 // Set Minimum as minimum value Greater Than 0
463 h3->SetMinimum(GetMinimumGT(*h3));
464}
465
466// --------------------------------------------------------------------------
467//
468// Update the projections
469//
470void MHFalseSource::Paint(Option_t *opt)
471{
472 // sigma = (s-b)/sqrt(s+k*k*b) mit k=s/b
473
474 gStyle->SetPalette(1, 0);
475
476 TVirtualPad *padsave = gPad;
477
478 TH1D* h1;
479 TH2D* h2;
480 TH2D* h3;
481 TH2D* h4;
482
483 padsave->cd(3);
484 if ((h3 = (TH2D*)gPad->FindObject("Alpha_xy_on")))
485 ProjectOn(h3);
486
487 padsave->cd(4);
488 if ((h2 = (TH2D*)gPad->FindObject("Alpha_xy_off")))
489 ProjectOff(h2);
490
491 padsave->cd(2);
492 if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_xy_sig")))
493 {
494 const Int_t nx = h4->GetXaxis()->GetNbins();
495 const Int_t ny = h4->GetYaxis()->GetNbins();
496
497 Int_t maxx=nx/2;
498 Int_t maxy=ny/2;
499
500 Int_t max = h4->GetBin(maxx, maxy);
501
502 for (int ix=0; ix<nx; ix++)
503 for (int iy=0; iy<ny; iy++)
504 {
505 const Int_t n = h4->GetBin(ix+1, iy+1);
506
507 const Double_t s = h3->GetBinContent(n);
508 const Double_t b = h2->GetBinContent(n);
509
510 const Double_t sig = Significance(s, b);
511
512 h4->SetBinContent(n, sig);
513
514 if (sig>h4->GetBinContent(max) && sig!=0)
515 {
516 max = n;
517 maxx=ix;
518 maxy=iy;
519 }
520 }
521
522 padsave->cd(1);
523 if ((h1 = (TH1D*)gPad->FindObject("Alpha")))
524 {
525 //h1->Reset();
526
527 const Double_t x = fHist.GetXaxis()->GetBinCenter(maxx);
528 const Double_t y = fHist.GetYaxis()->GetBinCenter(maxy);
529 const Double_t s = h4->GetBinContent(max);
530
531 TH1 *h = fHist.ProjectionZ("Alpha", maxx, maxx, maxy, maxy);
532 h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f)", x, y, s));
533 }
534 }
535
536 gPad = padsave;
537}
538
539// --------------------------------------------------------------------------
540//
541// Draw the histogram
542//
543void MHFalseSource::Draw(Option_t *opt)
544{
545 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
546 pad->SetBorderMode(0);
547
548 AppendPad("");
549
550 pad->Divide(2, 2);
551
552 // draw the 2D histogram Sigmabar versus Theta
553 pad->cd(1);
554 gPad->SetBorderMode(0);
555 TH1 *h1 = fHist.ProjectionZ("Alpha");
556 h1->SetDirectory(NULL);
557 h1->SetTitle("Distribution of \\alpha");
558 h1->SetXTitle(fHist.GetZaxis()->GetTitle());
559 h1->SetYTitle("Counts");
560 h1->Draw(opt);
561 h1->SetBit(kCanDelete);
562
563 pad->cd(4);
564 gPad->SetBorderMode(0);
565 fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2);
566 TH1 *h2 = fHist.Project3D("xy_off");
567 h2->SetDirectory(NULL);
568 h2->SetXTitle(fHist.GetXaxis()->GetTitle());
569 h2->SetYTitle(fHist.GetYaxis()->GetTitle());
570 h2->Draw("colz");
571 h2->SetBit(kCanDelete);
572
573 pad->cd(3);
574 gPad->SetBorderMode(0);
575 fHist.GetZaxis()->SetRangeUser(0,fAlphaCut);
576 TH1 *h3 = fHist.Project3D("xy_on");
577 fHist.GetZaxis()->SetRange(0,9999);
578 h3->SetDirectory(NULL);
579 h3->SetXTitle(fHist.GetXaxis()->GetTitle());
580 h3->SetYTitle(fHist.GetYaxis()->GetTitle());
581 h3->Draw("colz");
582 h3->SetBit(kCanDelete);
583
584 pad->cd(2);
585 gPad->SetBorderMode(0);
586 fHist.GetZaxis()->SetRange(0,0);
587 TH1 *h4 = fHist.Project3D("xy_sig"); // Do this to get the correct binning....
588 fHist.GetZaxis()->SetRange(0,9999);
589 h4->SetTitle("Significance");
590 h4->SetDirectory(NULL);
591 h4->SetXTitle(fHist.GetXaxis()->GetTitle());
592 h4->SetYTitle(fHist.GetYaxis()->GetTitle());
593 h4->Draw("colz");
594 h4->SetBit(kCanDelete);
595}
596
597// --------------------------------------------------------------------------
598//
599// Everything which is in the main pad belongs to this class!
600//
601Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py)
602{
603 return 0;
604}
605
606// --------------------------------------------------------------------------
607//
608// Set all sub-pads to Modified()
609//
610void MHFalseSource::Modified()
611{
612 if (!gPad)
613 return;
614
615 TVirtualPad *padsave = gPad;
616 padsave->Modified();
617 padsave->cd(1);
618 gPad->Modified();
619 padsave->cd(2);
620 gPad->Modified();
621 padsave->cd(3);
622 gPad->Modified();
623 padsave->cd(4);
624 gPad->Modified();
625 gPad->cd();
626}
627/*
628Double_t fcn(Double_t *arg, Double_t *p)
629{
630 const Double_t x = arg[0];
631
632 const Double_t dx = (x-p[1])/p[2];
633
634 const Double_t f1 = p[0]*TMath::Exp(-0.5*dx*dx);
635 const Double_t f2 = p[3] + p[5]*x*x;
636
637 return f1 + f2;
638}
639
640Double_t FcnI1(Double_t x, Double_t *p)
641{
642 return (p[5]*x*x/3+p[3])*x;
643}
644Double_t FcnI2(Double_t x, Double_t *p)
645{
646 static const Double_t sqrt2 = TMath::Sqrt(2.);
647 static const Double_t sqrt2pi = TMath::Sqrt(TMath::TwoPi());
648
649 const Double_t dx = (x-p[1])/p[2];
650
651 const Double_t f2 = p[0]*p[2]*sqrt2pi*TMath::Erf(dx/sqrt2)/2;
652
653 return f2;
654}
655*/
656/*
657class MHSignificance : public MH
658{
659private:
660 TH1D fHist;
661
662 MParameterD *fParam;
663
664public:
665 MHSignificance() : fParam(0)
666 {
667 fHist.SetName("Alpha");
668 fHist.SetTitle("Distribution of \\alpha");
669 fHist.SetXTitle("\\alpha [\\circ]");
670 fHist.SetYTitle("Counts");
671 }
672 Int_t SetupFill(MParList *p)
673 {
674 fHist.Reset();
675
676 fParam = (MParameterD*)p->FindCreateObj("Significance", "MParameterD");
677 if (fParam)
678 return kFALSE;
679
680 return kTRUE;
681 }
682 Int_t Process(MParContainer *p, Double_t w=1)
683 {
684 MHillasSrc *hil = dynamic_cast<MHillasSrc*>(p);
685 if (!hil)
686 {
687 *fLog << err << dbginf << "Got no MHillasSrc as argument of Fill()..." << endl;
688 return kFALSE;
689 }
690
691 fHist->Fill(hil->GetAlpha(), w);
692
693 return kTRUE;
694 }
695 Int_t Finalize()
696 {
697 if (fHist.GetEntries()==0)
698 {
699 *fLog << err << "Histogram empty." << endl;
700 return kFALSE;
701 }
702
703 Float_t sigmax=15;
704 Float_t bgmin =45;
705 Float_t bgmax =80;
706
707 fHist.SetNameTitle("Significance",
708 Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ",
709 sigmax, bgmin, bgmax));
710
711 // Implementing the function yourself is not faster at all!
712 TF1 func("gaus(0) + pol2(3)", fcn, 0, 90, 6);
713 TArrayD maxpar(func.GetNpar());
714
715 func.FixParameter(1, 0);
716 func.FixParameter(4, 0);
717 func.SetParLimits(3, -1, 1);
718
719 const Double_t alpha0 = fHist.GetBinContent(1);
720
721 // First fit a polynom in the off region
722 func.FixParameter(0, 0);
723 func.FixParameter(2, 1);
724 func.ReleaseParameter(3);
725 func.ReleaseParameter(5);
726
727 h->Fit(&func, "N0Q", "", bgmin, bgmax);
728
729 // Now fit a gaus in the on region on top of the polynom
730 func.SetParameter(0, alpha0-func.GetParameter(4));
731 func.SetParameter(2, sigmax*0.75);
732
733 func.ReleaseParameter(0);
734 func.ReleaseParameter(2);
735 func.FixParameter(3, func.GetParameter(3));
736 func.FixParameter(5, func.GetParameter(5));
737
738 func.SetParLimits(2, 0, 80);
739 h->Fit(&func, "N0Q", "", 0, sigmax);
740
741 TArrayD p(func.GetNpar(), func.GetParameters());
742
743 Double_t sig=0;
744
745 const Int_t n = hist->GetBin(ix+1, iy+1);
746 if (!(func.GetParameter(0)>alpha0*2 ||
747 func.GetParameter(2)<2.5 ||
748 func.GetParameter(2)>70))
749 {
750 // Implementing the integral as analytical function
751 // gives the same result in the order of 10e-5
752 // and it is not faster at all...
753 const Double_t s = func.Integral(0, 15);
754
755 func.SetParameter(0, 0);
756 func.SetParameter(2, 1);
757
758 const Double_t b = func.Integral(0, 15);
759
760 sig = Significance(s, b);
761 }
762
763 fParam->SetValue(sig);
764 fParam->SetReadyToSave();
765 return kTRUE;
766 }
767 ClassDef(MHSignificance, 0)
768};
769*/
770
771
772// --------------------------------------------------------------------------
773//
774// This is a preliminary implementation of a alpha-fit procedure for
775// all possible source positions. It will be moved into its own
776// more powerfull class soon.
777//
778// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
779// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
780// or
781// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
782//
783// Parameter [1] is fixed to 0 while the alpha peak should be
784// symmetric around alpha=0.
785//
786// Parameter [4] is fixed to 0 because the first derivative at
787// alpha=0 should be 0, too.
788//
789// In a first step the background is fitted between bgmin and bgmax,
790// while the parameters [0]=0 and [2]=1 are fixed.
791//
792// In a second step the signal region (alpha<sigmax) is fittet using
793// the whole function with parameters [1], [3], [4] and [5] fixed.
794//
795// The number of excess and background events are calculated as
796// s = int(0, sigint, gaus(0)+pol2(3))
797// b = int(0, sigint, pol2(3))
798//
799// The Significance is calculated using the Significance() member
800// function.
801//
802void MHFalseSource::FitSignificance(Float_t sigint, Float_t sigmax, Float_t bgmin, Float_t bgmax, Byte_t polynom)
803{
804 TH1D h0a("A", "", 50, 0, 4000);
805 TH1D h4a("chisq1", "", 50, 0, 35);
806 //TH1D h5a("prob1", "", 50, 0, 1.1);
807 TH1D h6("signifcance", "", 50, -20, 20);
808
809 TH1D h1("mu", "Parameter \\mu", 50, -1, 1);
810 TH1D h2("sigma", "Parameter \\sigma", 50, 0, 90);
811 TH1D h3("b", "Parameter b", 50, -0.1, 0.1);
812
813 TH1D h0b("a", "Parameter a (red), A (blue)", 50, 0, 4000);
814 TH1D h4b("\\chi^{2}", "\\chi^{2} (red, green) / significance (black)", 50, 0, 35);
815 //TH1D h5b("prob", "Fit probability: Bg(red), F(blue)", 50, 0, 1.1);
816
817 h0a.SetLineColor(kBlue);
818 h4a.SetLineColor(kBlue);
819 //h5a.SetLineColor(kBlue);
820 h0b.SetLineColor(kRed);
821 h4b.SetLineColor(kRed);
822 //h5b.SetLineColor(kRed);
823
824 TH1 *hist = fHist.Project3D("xy_fit");
825 hist->SetDirectory(0);
826 TH1 *hists = fHist.Project3D("xy_fit");
827 hists->SetDirectory(0);
828 TH1 *histb = fHist.Project3D("xy_fit");
829 histb->SetDirectory(0);
830 hist->Reset();
831 hists->Reset();
832 histb->Reset();
833 hist->SetNameTitle("Significance",
834 Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ",
835 sigmax, bgmin, bgmax));
836 hists->SetNameTitle("Excess", Form("Number of excess events for \\alpha<%.0f\\circ", sigint));
837 histb->SetNameTitle("Background", Form("Number of background events for \\alpha<%.0f\\circ", sigint));
838 hist->SetXTitle(fHist.GetXaxis()->GetTitle());
839 hists->SetXTitle(fHist.GetXaxis()->GetTitle());
840 histb->SetXTitle(fHist.GetXaxis()->GetTitle());
841 hist->SetYTitle(fHist.GetXaxis()->GetTitle());
842 hists->SetYTitle(fHist.GetXaxis()->GetTitle());
843 histb->SetYTitle(fHist.GetXaxis()->GetTitle());
844
845 // xmin, xmax, npar
846 //TF1 func("MyFunc", fcn, 0, 90, 6);
847 // Implementing the function yourself is only about 5% faster
848 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
849 TArrayD maxpar(func.GetNpar());
850
851 /*
852 func.SetParName(0, "A");
853 func.SetParName(1, "mu");
854 func.SetParName(2, "sigma");
855 func.SetParName(3, "a");
856 func.SetParName(4, "b");
857 func.SetParName(5, "c");
858 */
859
860 func.FixParameter(1, 0);
861 func.FixParameter(4, 0);
862 func.SetParLimits(2, 0, 90);
863 func.SetParLimits(3, -1, 1);
864
865 const Int_t nx = fHist.GetXaxis()->GetNbins();
866 const Int_t ny = fHist.GetYaxis()->GetNbins();
867
868 Double_t maxalpha0=0;
869 Double_t maxs=3;
870
871 Int_t maxx=0;
872 Int_t maxy=0;
873
874 TStopwatch clk;
875 clk.Start();
876
877 *fLog << inf;
878 *fLog << "Signal fit: alpha < " << sigmax << endl;
879 *fLog << "Integration: alpha < " << sigint << endl;
880 *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
881 *fLog << "Polynom order: " << (int)polynom << endl;
882 *fLog << "Fitting False Source Plot..." << flush;
883
884 TH1 *h=0;
885 for (int ix=1; ix<=nx; ix++)
886 for (int iy=1; iy<=ny; iy++)
887 {
888 h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy);
889
890 const Double_t alpha0 = h->GetBinContent(1);
891 const Double_t alphaw = h->GetXaxis()->GetBinWidth(1);
892
893 // Check for the regios which is not filled...
894 if (alpha0==0)
895 continue;
896
897 if (alpha0>maxalpha0)
898 maxalpha0=alpha0;
899
900 // First fit a polynom in the off region
901 func.FixParameter(0, 0);
902 func.FixParameter(2, 1);
903 func.ReleaseParameter(3);
904
905 for (int i=5; i<func.GetNpar(); i++)
906 func.ReleaseParameter(i);
907
908 h->Fit(&func, "N0Q", "", bgmin, bgmax);
909 //*fLog << dbg << ix << "/" << iy << ": " << func.GetParameter(3) << " " << func.GetParameter(4) << endl;
910
911 h4a.Fill(func.GetChisquare());
912 //h5a.Fill(func.GetProb());
913
914 //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1));
915 //func.SetParLimits(2, 5, 90);
916
917 func.ReleaseParameter(0);
918 //func.ReleaseParameter(1);
919 func.ReleaseParameter(2);
920 func.FixParameter(3, func.GetParameter(3));
921 for (int i=5; i<func.GetNpar(); i++)
922 func.FixParameter(i, func.GetParameter(i));
923
924 // Do not allow signals smaller than the background
925 const Double_t A = alpha0-func.GetParameter(3);
926 const Double_t dA = TMath::Abs(A);
927 func.SetParLimits(0, -dA*4, dA*4);
928 func.SetParLimits(2, 0, 90);
929
930 // Now fit a gaus in the on region on top of the polynom
931 func.SetParameter(0, A);
932 func.SetParameter(2, sigmax*0.75);
933
934 h->Fit(&func, "N0Q", "", 0, sigmax);
935 //*fLog << dbg << " " << func.GetParameter(0) << " " << func.GetParameter(1) << " " << func.GetParameter(2) << endl;
936
937 TArrayD p(func.GetNpar(), func.GetParameters());
938
939 // Fill results into some histograms
940 h0a.Fill(p[0]);
941 h0b.Fill(p[3]);
942 h1.Fill(p[1]);
943 h2.Fill(p[2]);
944 if (polynom>1)
945 h3.Fill(p[5]);
946 h4b.Fill(func.GetChisquare());
947 //h5b.Fill(func.GetProb());
948
949 // Implementing the integral as analytical function
950 // gives the same result in the order of 10e-5
951 // and it is not faster at all...
952 //const Bool_t ok = /*p[0]>=0 && /*p[0]<alpha0*2 &&*/ p[2]>1.75;// && p[2]<88.5;
953 /*
954 if (p[0]<-fabs(alpha0-p[3])*7 && p[2]<alphaw/3)
955 {
956 func.SetParameter(0, alpha0-p[3]);
957 cout << p[2] << " " << p[0] << " " << alpha0-p[3] << endl;
958 }
959 */
960
961 const Double_t s = func.Integral(0, sigint);
962
963 func.SetParameter(0, 0);
964 func.SetParameter(2, 1);
965
966 const Double_t b = func.Integral(0, sigint);
967 const Double_t sig = Significance(s, b);
968
969 const Int_t n = hist->GetBin(ix, iy);
970 hists->SetBinContent(n, s-b);
971 histb->SetBinContent(n, b);
972
973 hist->SetBinContent(n, sig);
974 if (sig!=0)
975 h6.Fill(sig);
976
977 if (sig>maxs)
978 {
979 maxs = sig;
980 maxx = ix;
981 maxy = iy;
982 maxpar = p;
983 }
984 }
985
986 *fLog << inf << "done." << endl;
987
988 h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
989 h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
990
991 //hists->SetMinimum(GetMinimumGT(*hists));
992 histb->SetMinimum(GetMinimumGT(*histb));
993
994 clk.Stop();
995 clk.Print();
996
997 TCanvas *c=new TCanvas;
998
999 c->Divide(3,2, 0, 0);
1000 c->cd(1);
1001 gPad->SetBorderMode(0);
1002 hists->Draw("colz");
1003 hists->SetBit(kCanDelete);
1004 c->cd(2);
1005 gPad->SetBorderMode(0);
1006 hist->Draw("colz");
1007 hist->SetBit(kCanDelete);
1008 c->cd(3);
1009 gPad->SetBorderMode(0);
1010 histb->Draw("colz");
1011 histb->SetBit(kCanDelete);
1012 c->cd(4);
1013 gPad->Divide(1,3, 0, 0);
1014 TVirtualPad *p=gPad;
1015 p->SetBorderMode(0);
1016 p->cd(1);
1017 gPad->SetBorderMode(0);
1018 h0b.DrawCopy();
1019 h0a.DrawCopy("same");
1020 p->cd(2);
1021 gPad->SetBorderMode(0);
1022 h3.DrawCopy();
1023 p->cd(3);
1024 gPad->SetBorderMode(0);
1025 h2.DrawCopy();
1026 c->cd(6);
1027 gPad->Divide(1,2, 0, 0);
1028 TVirtualPad *q=gPad;
1029 q->SetBorderMode(0);
1030 q->cd(1);
1031 gPad->SetBorderMode(0);
1032 gPad->SetBorderMode(0);
1033 h4b.DrawCopy();
1034 h4a.DrawCopy("same");
1035 h6.DrawCopy("same");
1036 q->cd(2);
1037 gPad->SetBorderMode(0);
1038 //h5b.DrawCopy();
1039 //h5a.DrawCopy("same");
1040 c->cd(5);
1041 gPad->SetBorderMode(0);
1042 if (maxx>0 && maxy>0)
1043 {
1044 const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ",
1045 fHist.GetXaxis()->GetBinCenter(maxx),
1046 fHist.GetYaxis()->GetBinCenter(maxy), maxs);
1047
1048 TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
1049 result->SetDirectory(NULL);
1050 result->SetNameTitle("Result \\alpha", title);
1051 result->SetBit(kCanDelete);
1052 result->SetXTitle("\\alpha [\\circ]");
1053 result->SetYTitle("Counts");
1054 result->Draw();
1055
1056 TF1 f1("", func.GetExpFormula(), 0, 90);
1057 TF1 f2("", func.GetExpFormula(), 0, 90);
1058 f1.SetParameters(maxpar.GetArray());
1059 f2.SetParameters(maxpar.GetArray());
1060 f2.FixParameter(0, 0);
1061 f2.FixParameter(1, 0);
1062 f2.FixParameter(2, 1);
1063 f1.SetLineColor(kGreen);
1064 f2.SetLineColor(kRed);
1065
1066 f2.DrawCopy("same");
1067 f1.DrawCopy("same");
1068
1069 TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
1070 leg->SetBorderSize(1);
1071 leg->SetTextSize(0.04);
1072 leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
1073 //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
1074 leg->AddLine(0, 0.65, 0, 0.65);
1075 leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
1076 leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
1077 leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
1078 leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
1079 leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
1080 if (polynom>1)
1081 leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11);
1082 leg->SetBit(kCanDelete);
1083 leg->Draw();
1084
1085 q->cd(2);
1086
1087 TGraph *g = new TGraph;
1088 g->SetPoint(0, 0, 0);
1089
1090 Int_t max=0;
1091 Float_t maxsig=0;
1092 for (int i=1; i<89; i++)
1093 {
1094 const Double_t s = f1.Integral(0, (float)i);
1095 const Double_t b = f2.Integral(0, (float)i);
1096
1097 const Double_t sig = Significance(s, b);
1098
1099 g->SetPoint(g->GetN(), i, sig);
1100
1101 if (sig>maxsig)
1102 {
1103 max = i;
1104 maxsig = sig;
1105 }
1106 }
1107
1108 g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
1109 g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
1110 g->GetHistogram()->SetYTitle("Significance");
1111 g->SetBit(kCanDelete);
1112 g->Draw("AP");
1113
1114 leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
1115 leg->SetBorderSize(1);
1116 leg->SetTextSize(0.1);
1117 leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max));
1118 leg->SetBit(kCanDelete);
1119 leg->Draw();
1120 }
1121}
Note: See TracBrowser for help on using the repository browser.