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

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