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

Last change on this file since 8019 was 7287, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 38.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MHFalseSource
28//
29// Create a false source plot. For the Binning in x,y the object MBinning
30// "BinningFalseSource" is taken from the paremeter list.
31//
32// The binning in alpha is currently fixed as 18bins from 0 to 90deg.
33//
34// If MTime, MObservatory and MPointingPos is available in the paremeter
35// list each image is derotated.
36//
37// MHFalseSource fills a 3D histogram with alpha distribution for
38// each (derotated) x and y position.
39//
40// Only a radius of 1.2deg around the camera center is taken into account.
41//
42// The displayed histogram is the projection of alpha<fAlphaCut into
43// the x,y plain and the projection of alpha>90-fAlphaCut
44//
45// By using the context menu (find it in a small region between the single
46// pads) you can change the AlphaCut 'online'
47//
48// Each Z-Projection (Alpha-histogram) is scaled such, that the number
49// of entries fits the maximum number of entries in all Z-Projections.
50// This should correct for the different acceptance in the center
51// and at the vorder of the camera. This, however, produces more noise
52// at the border.
53//
54// Here is a slightly simplified version of the algorithm:
55// ------------------------------------------------------
56// MHillas hil; // Taken as second argument in MFillH
57//
58// MSrcPosCam src;
59// MHillasSrc hsrc;
60// hsrc.SetSrcPos(&src);
61//
62// for (int ix=0; ix<nx; ix++)
63// for (int iy=0; iy<ny; iy++)
64// {
65// TVector2 v(cx[ix], cy[iy]); //cx center of x-bin ix
66// if (rho!=0) //cy center of y-bin iy
67// v=v.Rotate(rho); //image rotation angle
68//
69// src.SetXY(v); //source position in the camera
70//
71// if (!hsrc.Calc(hil)) //calc source dependant hillas
72// return;
73//
74// //fill absolute alpha into histogram
75// const Double_t alpha = hsrc.GetAlpha();
76// fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
77// }
78// }
79//
80// Use MHFalse Source like this:
81// -----------------------------
82// MFillH fill("MHFalseSource", "MHillas");
83// or
84// MHFalseSource hist;
85// hist.SetAlphaCut(12.5); // The default value
86// hist.SetBgmean(55); // The default value
87// MFillH fill(&hist, "MHillas");
88//
89// To be implemented:
90// ------------------
91// - a switch to use alpha or |alpha|
92// - taking the binning for alpha from the parlist (binning is
93// currently fixed)
94// - a possibility to change the fit-function (eg using a different TF1)
95// - a possibility to change the fit algorithm (eg which paremeters
96// are fixed at which time)
97// - use a different varaible than alpha
98// - a possiblity to change the algorithm which is used to calculate
99// alpha (or another parameter) is missing (this could be done using
100// a tasklist somewhere...)
101// - a more clever (and faster) algorithm to fill the histogram, eg by
102// calculating alpha once and fill the two coils around the mean
103// - more drawing options...
104// - Move Significance() to a more 'general' place and implement
105// also different algorithms like (Li/Ma)
106// - implement fit for best alpha distribution -- online (Paint)
107// - currently a constant pointing position is assumed in Fill()
108// - the center of rotation need not to be the camera center
109// - 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// - use the g/h-separation filters from the tasklists ("Cut1") as filters
113// two
114//
115//////////////////////////////////////////////////////////////////////////////
116#include "MHFalseSource.h"
117
118#include <TF1.h>
119#include <TF2.h>
120#include <TH2.h>
121#include <TGraph.h>
122#include <TStyle.h>
123#include <TCanvas.h>
124#include <TRandom.h>
125#include <TPaveText.h>
126#include <TStopwatch.h>
127
128#include "MGeomCam.h"
129#include "MSrcPosCam.h"
130#include "MHillasSrc.h"
131#include "MTime.h"
132#include "MObservatory.h"
133#include "MPointingPos.h"
134#include "MAstroCatalog.h"
135#include "MAstroSky2Local.h"
136#include "MStatusDisplay.h"
137
138#include "MMath.h"
139#include "MAlphaFitter.h"
140
141#include "MBinning.h"
142#include "MParList.h"
143
144#include "MLog.h"
145#include "MLogManip.h"
146
147ClassImp(MHFalseSource);
148
149using namespace std;
150
151//class MHillasExt;
152
153// --------------------------------------------------------------------------
154//
155// Default Constructor
156//
157MHFalseSource::MHFalseSource(const char *name, const char *title)
158 : fTime(0), fPointPos(0), fObservatory(0), fMm2Deg(-1), fAlphaCut(12.5),
159 fBgMean(55), fMinDist(-1), fMaxDist(-1), fMinDW(-1), fMaxDW(-1),
160 fHistOff(0)
161{
162 //
163 // set the name and title of this object
164 //
165 fName = name ? name : "MHFalseSource";
166 fTitle = title ? title : "3D-plot of Alpha vs x, y";
167
168 fHist.SetDirectory(NULL);
169
170 fHist.SetName("Alpha");
171 fHist.SetTitle("3D-plot of Alpha vs x, y");
172 fHist.SetXTitle("x [\\circ]");
173 fHist.SetYTitle("y [\\circ]");
174 fHist.SetZTitle("\\alpha [\\circ]");
175}
176
177void MHFalseSource::MakeSymmetric(TH1 *h)
178{
179 h->SetMinimum();
180 h->SetMaximum();
181
182 const Float_t min = TMath::Abs(h->GetMinimum());
183 const Float_t max = TMath::Abs(h->GetMaximum());
184
185 const Float_t absmax = TMath::Max(min, max)*1.002;
186
187 h->SetMaximum( absmax);
188 h->SetMinimum(-absmax);
189}
190
191// --------------------------------------------------------------------------
192//
193// Set the alpha cut (|alpha|<fAlphaCut) which is use to estimate the
194// number of excess events
195//
196void MHFalseSource::SetAlphaCut(Float_t alpha)
197{
198 if (alpha<0)
199 *fLog << warn << "Alpha<0... taking absolute value." << endl;
200
201 fAlphaCut = TMath::Abs(alpha);
202
203 Modified();
204}
205
206// --------------------------------------------------------------------------
207//
208// Set mean alpha around which the off sample is determined
209// (fBgMean-fAlphaCut/2<|fAlpha|<fBgMean+fAlphaCut/2) which is use
210// to estimate the number of off events
211//
212void MHFalseSource::SetBgMean(Float_t alpha)
213{
214 if (alpha<0)
215 *fLog << warn << "Alpha<0... taking absolute value." << endl;
216
217 fBgMean = TMath::Abs(alpha);
218
219 Modified();
220}
221
222// --------------------------------------------------------------------------
223//
224// Set binnings (takes BinningFalseSource) and prepare filling of the
225// histogram.
226//
227// Also search for MTime, MObservatory and MPointingPos
228//
229Bool_t MHFalseSource::SetupFill(const MParList *plist)
230{
231 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
232 if (!geom)
233 {
234 *fLog << err << "MGeomCam not found... aborting." << endl;
235 return kFALSE;
236 }
237 fMm2Deg = geom->GetConvMm2Deg();
238
239 const TString off(Form("%sOff", fName.Data()));
240 if (fName!=off && fHistOff==NULL)
241 {
242 const TString desc(Form("%s [%s] found... using ", off.Data(), ClassName()));
243 MHFalseSource *hoff = (MHFalseSource*)plist->FindObject(off, ClassName());
244 if (!hoff)
245 *fLog << inf << "No " << desc << "current data only!" << endl;
246 else
247 {
248 *fLog << inf << desc << "on-off mode!" << endl;
249 SetOffData(*hoff);
250 }
251 }
252
253 if (fHistOff)
254 MH::SetBinning(&fHist, fHistOff);
255 else
256 {
257 MBinning binsa(18, 0, 90);
258 binsa.SetEdges(*plist, "BinningAlpha");
259
260 const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource");
261 if (!bins || bins->IsDefault())
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
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 MPointingPos *point = (MPointingPos*)plist->FindObject("MSourcePos", "MPointingPos");
290 if (!point)
291 point = fPointPos;
292
293 // FIXME: Because the pointing position could change we must check
294 // for the current pointing position and add a offset in the
295 // Fill function!
296 fRa = point ? point->GetRa() : 0;
297 fDec = point ? point->GetDec() : 90;
298
299 return kTRUE;
300}
301
302// --------------------------------------------------------------------------
303//
304// Fill the histogram. For details see the code or the class description
305//
306Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w)
307{
308 const MHillas *hil = dynamic_cast<const MHillas*>(par);
309 if (!hil)
310 {
311 *fLog << err << "MHFalseSource::Fill: No container specified!" << endl;
312 return kFALSE;
313 }
314
315 // Get max radius...
316 const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
317
318 // Get camera rotation angle
319 Double_t rho = 0;
320 if (fTime && fObservatory && fPointPos)
321 rho = fPointPos->RotationAngle(*fObservatory, *fTime);
322 //if (fPointPos)
323 // rho = fPointPos->RotationAngle(*fObservatory);
324
325 // Create necessary containers for calculation
326 MSrcPosCam src;
327 MHillasSrc hsrc;
328 hsrc.SetSrcPos(&src);
329
330 // Get number of bins and bin-centers
331 const Int_t nx = fHist.GetNbinsX();
332 const Int_t ny = fHist.GetNbinsY();
333 Axis_t cx[nx];
334 Axis_t cy[ny];
335 fHist.GetXaxis()->GetCenter(cx);
336 fHist.GetYaxis()->GetCenter(cy);
337
338 for (int ix=0; ix<nx; ix++)
339 {
340 for (int iy=0; iy<ny; iy++)
341 {
342 // check distance... to get a circle plot
343 if (TMath::Hypot(cx[ix], cy[iy])>maxr)
344 continue;
345
346 // rotate center of bin
347 // precalculation of sin/cos doesn't accelerate
348 TVector2 v(cx[ix], cy[iy]);
349 if (rho!=0)
350 v=v.Rotate(rho);
351
352 // convert degrees to millimeters
353 v *= 1./fMm2Deg;
354
355 if (fSrcPos)
356 v += fSrcPos->GetXY();
357
358 src.SetXY(v);
359
360 // Source dependant hillas parameters
361 if (hsrc.Calc(*hil/*, ext*/)>0)
362 {
363 *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl;
364 return kFALSE;
365 }
366
367 // FIXME: This should be replaced by an external MFilter
368 // and/or MTaskList
369 // Source dependant distance cut
370 if (fMinDist>0 && hsrc.GetDist()*fMm2Deg<fMinDist)
371 continue;
372 if (fMaxDist>0 && hsrc.GetDist()*fMm2Deg>fMaxDist)
373 continue;
374
375 if (fMaxDW>0 && hsrc.GetDist()>fMaxDW*hil->GetWidth())
376 continue;
377 if (fMinDW<0 && hsrc.GetDist()<fMinDW*hil->GetWidth())
378 continue;
379
380 // Fill histogram
381 const Double_t alpha = hsrc.GetAlpha();
382 fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
383 }
384 }
385
386 return kTRUE;
387}
388
389// --------------------------------------------------------------------------
390//
391// Create projection for off data, taking sum of bin contents of
392// range (fBgMean-fAlphaCut/2, fBgMean+fAlphaCut) Making sure to take
393// the same number of bins than for on-data
394//
395void MHFalseSource::ProjectOff(const TH3D &src, TH2D *h2, TH2D *all)
396{
397 TAxis &axe = *src.GetZaxis();
398
399 // Find range to cut (left edge and width)
400 const Int_t f = axe.FindBin(fBgMean-fAlphaCut/2);
401 const Int_t l = axe.FindBin(fAlphaCut)+f-1;
402
403 axe.SetRange(f, l);
404 const Float_t cut1 = axe.GetBinLowEdge(f);
405 const Float_t cut2 = axe.GetBinUpEdge(l);
406 h2->SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2));
407
408 // Get projection for range
409 TH2D *p = (TH2D*)src.Project3D("yx_off_NULL");
410
411 // Reset range
412 axe.SetRange(0,9999);
413
414//#if ROOT_VERSION_CODE < ROOT_VERSION(4,02,00)
415 // Move contents from projection to h2
416 h2->Reset();
417 h2->Add(p, all->GetMaximum());
418 h2->Divide(all);
419
420 // Delete p
421 delete p;
422/*#else
423 p->Scale(all->GetMaximum());
424 p->Divide(all);
425#endif*/
426
427 // Set Minimum as minimum value Greater Than 0
428 h2->SetMinimum(GetMinimumGT(*h2));
429}
430
431// --------------------------------------------------------------------------
432//
433// Create projection for on data, taking sum of bin contents of
434// range (0, fAlphaCut)
435//
436void MHFalseSource::ProjectOn(const TH3D &src, TH2D *h3, TH2D *all)
437{
438 TAxis &axe = *src.GetZaxis();
439
440 // Find range to cut
441 axe.SetRangeUser(0, fAlphaCut);
442 const Float_t cut = axe.GetBinUpEdge(axe.GetLast());
443 h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut));
444
445 // Get projection for range
446 TH2D *p = (TH2D*)src.Project3D("yx_on_NULL");
447
448 // Reset range
449 axe.SetRange(0,9999);
450
451//#if ROOT_VERSION_CODE < ROOT_VERSION(4,02,00)
452 // Move contents from projection to h3
453 h3->Reset();
454 h3->Add(p, all->GetMaximum());
455 h3->Divide(all);
456
457 // Delete p
458 delete p;
459/*#else
460 p->Scale(all->GetMaximum());
461 p->Divide(all);
462#endif*/
463
464 // Set Minimum as minimum value Greater Than 0
465 h3->SetMinimum(GetMinimumGT(*h3));
466}
467
468// --------------------------------------------------------------------------
469//
470// Create projection for all data, taking sum of bin contents of
471// range (0, 90) - corresponding to the number of entries in this slice.
472//
473void MHFalseSource::ProjectAll(TH2D *h3)
474{
475 h3->SetTitle("Number of entries");
476
477 // Get projection for range
478#if ROOT_VERSION_CODE < ROOT_VERSION(4,02,00)
479 TH2D *p = (TH2D*)fHist.Project3D("yx_all");
480
481 // Move contents from projection to h3
482 h3->Reset();
483 h3->Add(p);
484 delete p;
485#else
486 fHist.Project3D("yx_all");
487#endif
488
489 // Set Minimum as minimum value Greater Than 0
490 h3->SetMinimum(GetMinimumGT(*h3));
491}
492
493void MHFalseSource::ProjectOnOff(TH2D *h2, TH2D *h0)
494{
495 ProjectOn(*fHistOff, h2, h0);
496
497 TH2D h;
498 MH::SetBinning(&h, h2);
499
500 // Divide by number of entries in off region (of off-data)
501 ProjectOff(*fHistOff, &h, h0);
502 h2->Divide(&h);
503
504 // Multiply by number of entries in off region (of on-data)
505 ProjectOff(fHist, &h, h0);
506 h2->Multiply(&h);
507
508 // Recalculate Minimum
509 h2->SetMinimum(GetMinimumGT(*h2));
510}
511
512// --------------------------------------------------------------------------
513//
514// Update the projections and paint them
515//
516void MHFalseSource::Paint(Option_t *opt)
517{
518 // Set pretty color palette
519 gStyle->SetPalette(1, 0);
520
521 TVirtualPad *padsave = gPad;
522
523 TH1D* h1;
524 TH2D* h0;
525 TH2D* h2;
526 TH2D* h3;
527 TH2D* h4;
528 TH2D* h5;
529
530 // Update projection of all-events
531 padsave->GetPad(2)->cd(3);
532 if ((h0 = (TH2D*)gPad->FindObject("Alpha_yx_all")))
533 ProjectAll(h0);
534
535 // Update projection of on-events
536 padsave->GetPad(1)->cd(1);
537 if ((h3 = (TH2D*)gPad->FindObject("Alpha_yx_on")))
538 ProjectOn(fHist, h3, h0);
539
540 // Update projection of off-events
541 padsave->GetPad(1)->cd(3);
542 if ((h2 = (TH2D*)gPad->FindObject("Alpha_yx_off")))
543 {
544 if (!fHistOff)
545 ProjectOff(fHist, h2, h0);
546 else
547 ProjectOnOff(h2, h0);
548 }
549
550 padsave->GetPad(2)->cd(2);
551 if ((h5 = (TH2D*)gPad->FindObject("Alpha_yx_diff")))
552 {
553 h5->Add(h2, h3, -1);
554 MakeSymmetric(h5);
555 }
556
557 // Update projection of significance
558 padsave->GetPad(1)->cd(2);
559 if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_yx_sig")))
560 {
561 const Int_t nx = h4->GetXaxis()->GetNbins();
562 const Int_t ny = h4->GetYaxis()->GetNbins();
563
564 Int_t maxx=nx/2;
565 Int_t maxy=ny/2;
566
567 Int_t max = h4->GetBin(nx, ny);
568
569 h4->SetEntries(0);
570 for (int ix=1; ix<=nx; ix++)
571 for (int iy=1; iy<=ny; iy++)
572 {
573 const Int_t n = h4->GetBin(ix, iy);
574
575 const Double_t s = h3->GetBinContent(n);
576 const Double_t b = h2->GetBinContent(n);
577
578 const Double_t sig = MMath::SignificanceLiMaSigned(s, b);
579
580 h4->SetBinContent(n, sig);
581
582 if (sig>h4->GetBinContent(max) && sig>0/* && (ix-nx/2)*(ix-nx/2)+(iy-ny/2)*(iy-ny/2)<nr*nr/9*/)
583 {
584 max = n;
585 maxx=ix;
586 maxy=iy;
587 }
588 }
589
590 MakeSymmetric(h4);
591
592 // Update projection of 'the best alpha-plot'
593 padsave->GetPad(2)->cd(1);
594 if ((h1 = (TH1D*)gPad->FindObject("Alpha_z")) && max>0)
595 {
596 const Double_t x = h4->GetXaxis()->GetBinCenter(maxx);
597 const Double_t y = h4->GetYaxis()->GetBinCenter(maxy);
598 const Double_t s = h4->GetBinContent(max);
599
600 // This is because a 3D histogram x and y are vice versa
601 // Than for their projections
602 TH1 *h = fHist.ProjectionZ("Alpha_z", maxx, maxx, maxy, maxy);
603 h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (S_{max}=%.1f\\sigma)", x, y, s));
604
605 TH1D *h0=0;
606 if ((h0 = (TH1D*)gPad->FindObject("AlphaOff_z")))
607 {
608 fHistOff->ProjectionZ("AlphaOff_z", maxx, maxx, maxy, maxy);
609
610 /* ============= local scaling ================ */
611 const Int_t f = h0->GetXaxis()->FindFixBin(fBgMean-1.5*fAlphaCut);
612 const Int_t l = h0->GetXaxis()->FindFixBin(fAlphaCut*3)+f-1;
613 h0->Scale(h1->Integral(f, l)/h0->Integral(f, l));
614
615
616 //h0->Scale(h1->GetEntries()/h0->GetEntries());
617
618 /* ============= global scaling ================
619 const Int_t f = fHistOff->GetZaxis()->FindFixBin(fBgMean-1.5*fAlphaCut);
620 const Int_t l = fHistOff->GetZaxis()->FindFixBin(fAlphaCut*3)+f-1;
621
622 Double_t c0 = fHist.Integral(0, 9999, 0, 9999, f, l);
623 Double_t c1 = fHistOff->Integral(0, 9999, 0, 9999, f, l);
624
625 h0->Scale(c0/c1);
626 */
627 }
628 }
629 }
630
631 gPad = padsave;
632}
633
634// --------------------------------------------------------------------------
635//
636// Get the MAstroCatalog corresponding to fRa, fDec. The limiting magnitude
637// is set to 9, while the fov is equal to the current fov of the false
638// source plot.
639//
640TObject *MHFalseSource::GetCatalog()
641{
642 const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
643
644 // Create catalog...
645 MAstroCatalog *stars = new MAstroCatalog;
646 stars->SetMarkerColor(kWhite);
647 stars->SetLimMag(9);
648 stars->SetGuiActive(kFALSE);
649 stars->SetRadiusFOV(maxr);
650 stars->SetRaDec(fRa*TMath::DegToRad()*15, fDec*TMath::DegToRad());
651 stars->ReadBSC("bsc5.dat");
652
653 stars->SetBit(kCanDelete);
654 return stars;
655}
656
657// --------------------------------------------------------------------------
658//
659// Draw the histogram
660//
661void MHFalseSource::Draw(Option_t *opt)
662{
663 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
664 pad->SetBorderMode(0);
665
666 AppendPad("");
667
668 pad->Divide(1, 2, 1e-10, 0.03);
669
670// TObject *catalog = GetCatalog();
671
672 // Initialize upper part
673 pad->cd(1);
674 // Make sure that the catalog is deleted only once
675 // Normally this is not done by root, because it is not necessary...
676 // but in this case it is necessary, because the catalog is
677 // deleted by the first pad and the second one tries to do the same!
678// gROOT->GetListOfCleanups()->Add(gPad);
679 gPad->SetBorderMode(0);
680 gPad->Divide(3, 1);
681
682 // PAD #1
683 pad->GetPad(1)->cd(1);
684 gPad->SetBorderMode(0);
685 fHist.GetZaxis()->SetRangeUser(0,fAlphaCut);
686 TH1 *h3 = fHist.Project3D("yx_on");
687 fHist.GetZaxis()->SetRange(0,9999);
688 h3->SetDirectory(NULL);
689 h3->SetXTitle(fHist.GetXaxis()->GetTitle());
690 h3->SetYTitle(fHist.GetYaxis()->GetTitle());
691 h3->Draw("colz");
692 h3->SetBit(kCanDelete);
693// catalog->Draw("mirror same *");
694
695 // PAD #2
696 pad->GetPad(1)->cd(2);
697 gPad->SetBorderMode(0);
698 fHist.GetZaxis()->SetRange(0,0);
699 TH1 *h4 = fHist.Project3D("yx_sig"); // Do this to get the correct binning....
700 fHist.GetZaxis()->SetRange(0,9999);
701 h4->SetTitle("Significance");
702 h4->SetDirectory(NULL);
703 h4->SetXTitle(fHist.GetXaxis()->GetTitle());
704 h4->SetYTitle(fHist.GetYaxis()->GetTitle());
705 h4->Draw("colz");
706 h4->SetBit(kCanDelete);
707// catalog->Draw("mirror same *");
708
709 // PAD #3
710 pad->GetPad(1)->cd(3);
711 gPad->SetBorderMode(0);
712 TH1 *h2 = 0;
713 if (fHistOff)
714 {
715 fHistOff->GetZaxis()->SetRangeUser(0,fAlphaCut);
716 h2 = fHistOff->Project3D("yx_off");
717 fHistOff->GetZaxis()->SetRange(0,9999);
718 }
719 else
720 {
721 fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2);
722 h2 = fHist.Project3D("yx_off");
723 fHist.GetZaxis()->SetRange(0,9999);
724 }
725 h2->SetDirectory(NULL);
726 h2->SetXTitle(fHist.GetXaxis()->GetTitle());
727 h2->SetYTitle(fHist.GetYaxis()->GetTitle());
728 h2->Draw("colz");
729 h2->SetBit(kCanDelete);
730// catalog->Draw("mirror same *");
731
732 // Initialize lower part
733 pad->cd(2);
734 // Make sure that the catalog is deleted only once
735// gROOT->GetListOfCleanups()->Add(gPad);
736 gPad->SetBorderMode(0);
737 gPad->Divide(3, 1);
738
739 // PAD #4
740 pad->GetPad(2)->cd(1);
741 gPad->SetBorderMode(0);
742 TH1 *h1 = fHist.ProjectionZ("Alpha_z");
743 h1->SetDirectory(NULL);
744 h1->SetTitle("Distribution of \\alpha");
745 h1->SetXTitle(fHist.GetZaxis()->GetTitle());
746 h1->SetYTitle("Counts");
747 h1->Draw();
748 h1->SetBit(kCanDelete);
749 if (fHistOff)
750 {
751 h1->SetLineColor(kGreen);
752
753 h1 = fHistOff->ProjectionZ("AlphaOff_z");
754 h1->SetDirectory(NULL);
755 h1->SetTitle("Distribution of \\alpha");
756 h1->SetXTitle(fHistOff->GetZaxis()->GetTitle());
757 h1->SetYTitle("Counts");
758 h1->Draw("same");
759 h1->SetBit(kCanDelete);
760 h1->SetLineColor(kRed);
761 }
762
763 // PAD #5
764 pad->GetPad(2)->cd(2);
765 gPad->SetBorderMode(0);
766 TH1 *h5 = (TH1*)h3->Clone("Alpha_yx_diff");
767 h5->Add(h2, -1);
768 h5->SetTitle("Difference of on- and off-distribution");
769 h5->SetDirectory(NULL);
770 h5->Draw("colz");
771 h5->SetBit(kCanDelete);
772// catalog->Draw("mirror same *");
773
774 // PAD #6
775 pad->GetPad(2)->cd(3);
776 gPad->SetBorderMode(0);
777 TH1 *h0 = fHist.Project3D("yx_all");
778 h0->SetDirectory(NULL);
779 h0->SetXTitle(fHist.GetXaxis()->GetTitle());
780 h0->SetYTitle(fHist.GetYaxis()->GetTitle());
781 h0->Draw("colz");
782 h0->SetBit(kCanDelete);
783// catalog->Draw("mirror same *");
784}
785
786// --------------------------------------------------------------------------
787//
788// Everything which is in the main pad belongs to this class!
789//
790Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py)
791{
792 return 0;
793}
794
795// --------------------------------------------------------------------------
796//
797// Set all sub-pads to Modified()
798//
799void MHFalseSource::Modified()
800{
801 if (!gPad)
802 return;
803
804 TVirtualPad *padsave = gPad;
805 padsave->Modified();
806 padsave->GetPad(1)->cd(1);
807 gPad->Modified();
808 padsave->GetPad(1)->cd(3);
809 gPad->Modified();
810 padsave->GetPad(2)->cd(1);
811 gPad->Modified();
812 padsave->GetPad(2)->cd(2);
813 gPad->Modified();
814 padsave->GetPad(2)->cd(3);
815 gPad->Modified();
816 gPad->cd();
817}
818
819// --------------------------------------------------------------------------
820//
821// The following resources are available:
822//
823// MHFalseSource.DistMin: 0.5
824// MHFalseSource.DistMax: 1.4
825// MHFalseSource.DWMin: 0.1
826// MHFalseSource.DWMax: 0.3
827//
828Int_t MHFalseSource::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
829{
830 Bool_t rc = kFALSE;
831 if (IsEnvDefined(env, prefix, "DistMin", print))
832 {
833 rc = kTRUE;
834 SetMinDist(GetEnvValue(env, prefix, "DistMin", fMinDist));
835 }
836 if (IsEnvDefined(env, prefix, "DistMax", print))
837 {
838 rc = kTRUE;
839 SetMaxDist(GetEnvValue(env, prefix, "DistMax", fMaxDist));
840 }
841
842 if (IsEnvDefined(env, prefix, "DWMin", print))
843 {
844 rc = kTRUE;
845 SetMinDW(GetEnvValue(env, prefix, "DWMin", fMinDist));
846 }
847 if (IsEnvDefined(env, prefix, "DWMax", print))
848 {
849 rc = kTRUE;
850 SetMaxDW(GetEnvValue(env, prefix, "DWMax", fMaxDist));
851 }
852
853 return rc;
854}
855
856static Double_t FcnGauss2d(Double_t *x, Double_t *par)
857{
858 TVector2 v = TVector2(x[0], x[1]).Rotate(par[5]*TMath::DegToRad());
859
860 const Double_t g0 = TMath::Gaus(v.X(), par[1], par[2]);
861 const Double_t g1 = TMath::Gaus(v.Y(), par[3], par[4]);
862
863 //Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE);
864 return par[0]*g0*g1;
865}
866
867// --------------------------------------------------------------------------
868//
869// This is a preliminary implementation of a alpha-fit procedure for
870// all possible source positions. It will be moved into its own
871// more powerfull class soon.
872//
873// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
874// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
875// or
876// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
877//
878// Parameter [1] is fixed to 0 while the alpha peak should be
879// symmetric around alpha=0.
880//
881// Parameter [4] is fixed to 0 because the first derivative at
882// alpha=0 should be 0, too.
883//
884// In a first step the background is fitted between bgmin and bgmax,
885// while the parameters [0]=0 and [2]=1 are fixed.
886//
887// In a second step the signal region (alpha<sigmax) is fittet using
888// the whole function with parameters [1], [3], [4] and [5] fixed.
889//
890// The number of excess and background events are calculated as
891// s = int(0, sigint, gaus(0)+pol2(3))
892// b = int(0, sigint, pol2(3))
893//
894// The Significance is calculated using the Significance() member
895// function.
896//
897void MHFalseSource::FitSignificance(Float_t sigint, Float_t sigmax, Float_t bgmin, Float_t bgmax, Byte_t polynom)
898{
899// TObject *catalog = GetCatalog();
900
901 TH1D h0a("A", "", 50, 0, 4000);
902 TH1D h4a("chisq1", "", 50, 0, 35);
903 //TH1D h5a("prob1", "", 50, 0, 1.1);
904 TH1D h6("signifcance", "", 50, -20, 20);
905
906 TH1D h1("mu", "Parameter \\mu", 50, -1, 1);
907 TH1D h2("sigma", "Parameter \\sigma", 50, 0, 90);
908 TH1D h3("b", "Parameter b", 50, -0.1, 0.1);
909
910 TH1D h0b("a", "Parameter a (red), A (blue)", 50, 0, 4000);
911 TH1D h4b("\\chi^{2}", "\\chi^{2} (red, green) / significance (black)", 50, 0, 35);
912 //TH1D h5b("prob", "Fit probability: Bg(red), F(blue)", 50, 0, 1.1);
913
914 h0a.SetLineColor(kBlue);
915 h4a.SetLineColor(kBlue);
916 //h5a.SetLineColor(kBlue);
917 h0b.SetLineColor(kRed);
918 h4b.SetLineColor(kRed);
919 //h5b.SetLineColor(kRed);
920
921 TH1 *hist = fHist.Project3D("yx_fit");
922 hist->SetDirectory(0);
923 TH1 *hists = fHist.Project3D("yx_fit");
924 hists->SetDirectory(0);
925 TH1 *histb = fHist.Project3D("yx_fit");
926 histb->SetDirectory(0);
927 hist->Reset();
928 hists->Reset();
929 histb->Reset();
930 hist->SetNameTitle("Significance",
931 Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ",
932 sigmax, bgmin, bgmax));
933 hists->SetName("Excess");
934 histb->SetName("Background");
935 hist->SetXTitle(fHist.GetXaxis()->GetTitle());
936 hists->SetXTitle(fHist.GetXaxis()->GetTitle());
937 histb->SetXTitle(fHist.GetXaxis()->GetTitle());
938 hist->SetYTitle(fHist.GetYaxis()->GetTitle());
939 hists->SetYTitle(fHist.GetYaxis()->GetTitle());
940 histb->SetYTitle(fHist.GetYaxis()->GetTitle());
941
942 const Double_t w = fHist.GetZaxis()->GetBinWidth(1);
943
944 TArrayD maxpar;
945
946 /* func.SetParName(0, "A");
947 * func.SetParName(1, "mu");
948 * func.SetParName(2, "sigma");
949 */
950
951 const Int_t nx = hist->GetXaxis()->GetNbins();
952 const Int_t ny = hist->GetYaxis()->GetNbins();
953 //const Int_t nr = nx*nx+ny*ny;
954
955 Double_t maxalpha0=0;
956 Double_t maxs=3;
957
958 Int_t maxx=0;
959 Int_t maxy=0;
960
961 TStopwatch clk;
962 clk.Start();
963
964 *fLog << inf;
965 *fLog << "Signal fit: alpha < " << sigmax << endl;
966 *fLog << "Integration: alpha < " << sigint << endl;
967 *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
968 *fLog << "Polynom order: " << (int)polynom << endl;
969 *fLog << "Fitting False Source Plot..." << flush;
970
971 TH1 *h0 = fHist.Project3D("yx_entries");
972 Float_t entries = h0->GetMaximum();
973 delete h0;
974
975 MAlphaFitter fit;
976 fit.SetSignalIntegralMax(sigint);
977 fit.SetSignalFitMax(sigmax);
978 fit.SetBackgroundFitMin(bgmin);
979 fit.SetBackgroundFitMax(bgmax);
980 fit.SetPolynomOrder(polynom);
981
982 TH1D *h=0, *hoff=0;
983 Double_t scale = 1;
984 for (int ix=1; ix<=nx; ix++)
985 for (int iy=1; iy<=ny; iy++)
986 {
987 // This is because a 3D histogram x and y are vice versa
988 // Than for their projections
989 h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy);
990
991 if (h->GetEntries()==0)
992 continue;
993
994 // This is a quick hack to correct for the lower acceptance at
995 // the edge of the camera
996 h->Scale(entries/h->GetEntries());
997
998 if (fHistOff)
999 {
1000 hoff = fHistOff->ProjectionZ("AlphaFitOff", ix, ix, iy, iy);
1001 // This is a quick hack to correct for the lower acceptance at
1002 // the edge of the camera
1003 hoff->Scale(entries/h->GetEntries());
1004 scale = fit.Scale(*hoff, *h);
1005 }
1006
1007 if (!fit.Fit(*h, hoff, scale))
1008 continue;
1009
1010 const Double_t alpha0 = h->GetBinContent(1);
1011 if (alpha0>maxalpha0)
1012 maxalpha0=alpha0;
1013
1014 // Fill results into some histograms
1015 h0a.Fill(fit.GetGausA()); // gaus-A
1016 h0b.Fill(fit.GetCoefficient(3)); // 3
1017 h1.Fill(fit.GetGausMu()); // mu
1018 h2.Fill(fit.GetGausSigma()); // sigma-gaus
1019 if (polynom>1 && !fHistOff)
1020 h3.Fill(fit.GetCoefficient(5));
1021 h4b.Fill(fit.GetChiSqSignal());
1022
1023 const Double_t sig = fit.GetSignificance();
1024 const Double_t b = fit.GetEventsBackground();
1025 const Double_t s = fit.GetEventsSignal();
1026
1027 const Int_t n = hist->GetBin(ix, iy);
1028 hists->SetBinContent(n, s-b);
1029 histb->SetBinContent(n, b);
1030
1031 hist->SetBinContent(n, sig);
1032 if (sig!=0)
1033 h6.Fill(sig);
1034
1035 if (sig>maxs)
1036 {
1037 maxs = sig;
1038 maxx = ix;
1039 maxy = iy;
1040 maxpar = fit.GetCoefficients();
1041 }
1042 }
1043
1044 *fLog << "Done." << endl;
1045
1046 h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
1047 h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
1048
1049 hists->SetTitle(Form("Excess events for \\alpha<%.0f\\circ (N_{max}=%d)", sigint, (int)hists->GetMaximum()));
1050 histb->SetTitle(Form("Background events for \\alpha<%.0f\\circ", sigint));
1051
1052 //hists->SetMinimum(GetMinimumGT(*hists));
1053 histb->SetMinimum(GetMinimumGT(*histb));
1054
1055 MakeSymmetric(hists);
1056 MakeSymmetric(hist);
1057
1058 clk.Stop();
1059 clk.Print("m");
1060
1061 TCanvas *c=new TCanvas;
1062
1063 gStyle->SetPalette(1, 0);
1064
1065 c->Divide(3,2, 1e-10, 1e-10);
1066 c->cd(1);
1067 gPad->SetBorderMode(0);
1068 hists->Draw("colz");
1069 hists->SetBit(kCanDelete);
1070// catalog->Draw("mirror same *");
1071 c->cd(2);
1072 gPad->SetBorderMode(0);
1073 hist->Draw("colz");
1074 hist->SetBit(kCanDelete);
1075
1076
1077 const Double_t maxr = 0.9*TMath::Abs(fHist.GetBinCenter(1));
1078 TF2 f2d("Gaus-2D", FcnGauss2d, -maxr, maxr, -maxr, maxr, 6);
1079 f2d.SetLineWidth(1);
1080 f2d.SetParName(0, "Max sigma");
1081 f2d.SetParName(1, "Mean_1 deg");
1082 f2d.SetParName(2, "Sigma_1 deg");
1083 f2d.SetParName(3, "Mean_2 deg");
1084 f2d.SetParName(4, "Sigma_2 deg");
1085 f2d.SetParName(5, "Phi deg");
1086 f2d.SetParLimits(1, -maxr/2, maxr/2); // mu_1
1087 f2d.SetParLimits(3, -maxr/2, maxr/2); // mu_2
1088 f2d.SetParLimits(2, 0, maxr); // sigma_1
1089 f2d.SetParLimits(4, 0, maxr); // sigma_2
1090 f2d.SetParLimits(5, 0, 45); // phi
1091 f2d.SetParameter(0, maxs); // A
1092 f2d.SetParameter(1, hist->GetXaxis()->GetBinCenter(maxx)); // mu_1
1093 f2d.SetParameter(2, 0.1); // sigma_1
1094 f2d.SetParameter(3, hist->GetYaxis()->GetBinCenter(maxy)); // mu_2
1095 f2d.SetParameter(4, 0.1); // sigma_2
1096 f2d.FixParameter(5, 0); // phi
1097 hist->Fit(&f2d, "NI0R");
1098 f2d.DrawCopy("cont2same");
1099
1100
1101// catalog->Draw("mirror same *");
1102 c->cd(3);
1103 gPad->SetBorderMode(0);
1104 histb->Draw("colz");
1105 histb->SetBit(kCanDelete);
1106// catalog->Draw("mirror same *");
1107 c->cd(4);
1108 gPad->Divide(1,3, 0, 0);
1109 TVirtualPad *p=gPad;
1110 p->SetBorderMode(0);
1111 p->cd(1);
1112 gPad->SetBorderMode(0);
1113 h0b.DrawCopy();
1114 h0a.DrawCopy("same");
1115 p->cd(2);
1116 gPad->SetBorderMode(0);
1117 h3.DrawCopy();
1118 p->cd(3);
1119 gPad->SetBorderMode(0);
1120 h2.DrawCopy();
1121 c->cd(6);
1122 gPad->Divide(1,2, 0, 0);
1123 TVirtualPad *q=gPad;
1124 q->SetBorderMode(0);
1125 q->cd(1);
1126 gPad->SetBorderMode(0);
1127 gPad->SetBorderMode(0);
1128 h4b.DrawCopy();
1129 h4a.DrawCopy("same");
1130 h6.DrawCopy("same");
1131 q->cd(2);
1132 gPad->SetBorderMode(0);
1133 //h5b.DrawCopy();
1134 //h5a.DrawCopy("same");
1135 c->cd(5);
1136 gPad->SetBorderMode(0);
1137 if (maxx>0 && maxy>0)
1138 {
1139 const char *title = Form(" \\alpha for x=%.2f y=%.2f (S_{max}=%.1f\\sigma) ",
1140 hist->GetXaxis()->GetBinCenter(maxx),
1141 hist->GetYaxis()->GetBinCenter(maxy), maxs);
1142
1143 h = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
1144 h->Scale(entries/h->GetEntries());
1145
1146 h->SetDirectory(NULL);
1147 h->SetNameTitle("Result \\alpha", title);
1148 h->SetBit(kCanDelete);
1149 h->SetXTitle("\\alpha [\\circ]");
1150 h->SetYTitle("Counts");
1151 h->Draw();
1152
1153 if (fHistOff)
1154 {
1155 h->SetLineColor(kGreen);
1156
1157 TH1D *hof=fHistOff->ProjectionZ("AlphaFitOff", maxx, maxx, maxy, maxy);
1158 hof->Scale(entries/hof->GetEntries());
1159
1160 fit.Scale(*(TH1D*)hof, *(TH1D*)h);
1161
1162 hof->SetLineColor(kRed);
1163 hof->SetDirectory(NULL);
1164 hof->SetNameTitle("Result \\alpha", title);
1165 hof->SetBit(kCanDelete);
1166 hof->SetXTitle("\\alpha [\\circ]");
1167 hof->SetYTitle("Counts");
1168 hof->Draw("same");
1169
1170 TH1D *diff = (TH1D*)h->Clone("AlphaFitOnOff");
1171 diff->Add(hof, -1);
1172 diff->SetLineColor(kBlue);
1173 diff->SetNameTitle("Result On-Off \\alpha", title);
1174 diff->SetBit(kCanDelete);
1175 diff->SetXTitle("\\alpha [\\circ]");
1176 diff->SetYTitle("Counts");
1177 diff->Draw("same");
1178
1179 h->SetMinimum(diff->GetMinimum()<0 ? diff->GetMinimum()*1.2 : 0);
1180
1181 TLine l;
1182 l.DrawLine(0, 0, 90, 0);
1183 }
1184
1185 TF1 f1("f1", Form("gaus(0) + pol%d(3)", fHistOff ? 0 : polynom), 0, 90);
1186 TF1 f2("f2", Form("gaus(0) + pol%d(3)", fHistOff ? 0 : polynom), 0, 90);
1187 f1.SetParameters(maxpar.GetArray());
1188 f2.SetParameters(maxpar.GetArray());
1189 f2.FixParameter(0, 0);
1190 f2.FixParameter(1, 0);
1191 f2.FixParameter(2, 1);
1192 f1.SetLineColor(kGreen);
1193 f2.SetLineColor(kRed);
1194
1195 f2.DrawCopy("same");
1196 f1.DrawCopy("same");
1197
1198 TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
1199 leg->SetBorderSize(1);
1200 leg->SetTextSize(0.04);
1201 leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
1202 //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
1203 leg->AddLine(0, 0.65, 0, 0.65);
1204 leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
1205 leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
1206 leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
1207 leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
1208 leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
1209 if (polynom>1)
1210 leg->AddText(0.60, 0.14, Form("c=%.2f", !fHistOff?maxpar[5]:0))->SetTextAlign(11);
1211 leg->SetBit(kCanDelete);
1212 leg->Draw();
1213
1214 q->cd(2);
1215
1216 TGraph *g = new TGraph;
1217 g->SetPoint(0, 0, 0);
1218
1219 Int_t max=0;
1220 Float_t maxsig=0;
1221 for (int i=1; i<89; i++)
1222 {
1223 const Double_t s = f1.Integral(0, (float)i)/w;
1224 const Double_t b = f2.Integral(0, (float)i)/w;
1225
1226 const Double_t sig = MMath::SignificanceLiMaSigned(s, b);
1227
1228 g->SetPoint(g->GetN(), i, sig);
1229
1230 if (sig>maxsig)
1231 {
1232 max = i;
1233 maxsig = sig;
1234 }
1235 }
1236
1237 g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
1238 g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
1239 g->GetHistogram()->SetYTitle("Significance");
1240 g->SetBit(kCanDelete);
1241 g->Draw("AP");
1242
1243 leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
1244 leg->SetBorderSize(1);
1245 leg->SetTextSize(0.1);
1246 leg->AddText(Form("S_{max}=%.1f\\sigma at \\alpha_{max}=%d\\circ", maxsig, max));
1247 leg->SetBit(kCanDelete);
1248 leg->Draw();
1249 }
1250}
Note: See TracBrowser for help on using the repository browser.