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

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