source: tags/Mars-V0.8.6/mhflux/MHFalseSource.cc

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