source: tags/Mars-V0.8.4/mhist/MHFalseSource.cc

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