source: trunk/MagicSoft/Mars/mhflux/MHFalseSource.cc@ 5002

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