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

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