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

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