source: trunk/MagicSoft/Mars/mhist/MHFalseSource.cc@ 3580

Last change on this file since 3580 was 3578, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 29.0 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// Here is a slightly simplified version of the algorithm:
49// ------------------------------------------------------
50// MHillas hil; // Taken as second argument in MFillH
51//
52// MSrcPosCam src;
53// MHillasSrc hsrc;
54// hsrc.SetSrcPos(&src);
55//
56// for (int ix=0; ix<nx; ix++)
57// for (int iy=0; iy<ny; iy++)
58// {
59// TVector2 v(cx[ix], cy[iy]); //cx center of x-bin ix
60// if (rho!=0) //cy center of y-bin iy
61// v.Rotate(-rho); //image rotation angle
62//
63// src.SetXY(v); //source position in the camera
64//
65// if (!hsrc.Calc(hil)) //calc source dependant hillas
66// return;
67//
68// //fill absolute alpha into histogram
69// const Double_t alpha = hsrc.GetAlpha();
70// fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
71// }
72// }
73//
74// Use MHFalse Source like this:
75// -----------------------------
76// MFillH fill("MHFalseSource", "MHillas");
77// or
78// MHFalseSource hist;
79// hist.SetAlphaCut(12.5); // The default value
80// hist.SetBgmean(55); // The default value
81// MFillH fill(&hist, "MHillas");
82//
83// To be implemented:
84// ------------------
85// - a switch to use alpha or |alpha|
86// - taking the binning for alpha from the parlist (binning is
87// currently fixed)
88// - a possibility to change the fit-function (eg using a different TF1)
89// - a possibility to change the fit algorithm (eg which paremeters
90// are fixed at which time)
91// - use a different varaible than alpha
92// - a possiblity to change the algorithm which is used to calculate
93// alpha (or another parameter) is missing (this could be done using
94// a tasklist somewhere...)
95// - a more clever (and faster) algorithm to fill the histogram, eg by
96// calculating alpha once and fill the two coils around the mean
97// - more drawing options...
98// - Move Significance() to a more 'general' place and implement
99// also different algorithms like (Li/Ma)
100// - implement fit for best alpha distribution -- online (Paint)
101//
102//////////////////////////////////////////////////////////////////////////////
103#include "MHFalseSource.h"
104
105#include <TF1.h>
106#include <TH2.h>
107#include <TGraph.h>
108#include <TStyle.h>
109#include <TCanvas.h>
110#include <TPaveText.h>
111#include <TStopwatch.h>
112
113#include "MGeomCam.h"
114#include "MSrcPosCam.h"
115#include "MHillasSrc.h"
116#include "MTime.h"
117#include "MObservatory.h"
118#include "MPointingPos.h"
119#include "MAstroSky2Local.h"
120#include "MStatusDisplay.h"
121
122#include "MBinning.h"
123#include "MParList.h"
124
125#include "MLog.h"
126#include "MLogManip.h"
127
128ClassImp(MHFalseSource);
129
130using namespace std;
131
132// --------------------------------------------------------------------------
133//
134// Default Constructor
135//
136MHFalseSource::MHFalseSource(const char *name, const char *title)
137 : fMm2Deg(-1), fUseMmScale(kTRUE), fAlphaCut(12.5), fBgMean(55)
138{
139 //
140 // set the name and title of this object
141 //
142 fName = name ? name : "MHFalseSource";
143 fTitle = title ? title : "3D-plot of Alpha vs x, y";
144
145 fHist.SetDirectory(NULL);
146
147 fHist.SetName("Alpha");
148 fHist.SetTitle("3D-plot of Alpha vs x, y");
149 fHist.SetXTitle("x [\\circ]");
150 fHist.SetYTitle("y [\\circ]");
151 fHist.SetZTitle("\\alpha [\\circ]");
152}
153
154// --------------------------------------------------------------------------
155//
156// Set the alpha cut (|alpha|<fAlphaCut) which is use to estimate the
157// number of excess events
158//
159void MHFalseSource::SetAlphaCut(Float_t alpha)
160{
161 if (alpha<0)
162 *fLog << warn << "Alpha<0... taking absolute value." << endl;
163
164 fAlphaCut = TMath::Abs(alpha);
165
166 Modified();
167}
168
169// --------------------------------------------------------------------------
170//
171// Set mean alpha around which the off sample is determined
172// (fBgMean-fAlphaCut/2<|fAlpha|<fBgMean+fAlphaCut/2) which is use
173// to estimate the number of off events
174//
175void MHFalseSource::SetBgMean(Float_t alpha)
176{
177 if (alpha<0)
178 *fLog << warn << "Alpha<0... taking absolute value." << endl;
179
180 fBgMean = TMath::Abs(alpha);
181
182 Modified();
183}
184
185// --------------------------------------------------------------------------
186//
187// Use this function to setup your own conversion factor between degrees
188// and millimeters. The conversion factor should be the one calculated in
189// MGeomCam. Use this function with Caution: You could create wrong values
190// by setting up your own scale factor.
191//
192void MHFalseSource::SetMm2Deg(Float_t mmdeg)
193{
194 if (mmdeg<0)
195 {
196 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
197 return;
198 }
199
200 if (fMm2Deg>=0)
201 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
202
203 fMm2Deg = mmdeg;
204}
205
206// --------------------------------------------------------------------------
207//
208// With this function you can convert the histogram ('on the fly') between
209// degrees and millimeters.
210//
211void MHFalseSource::SetMmScale(Bool_t mmscale)
212{
213 if (fUseMmScale == mmscale)
214 return;
215
216 if (fMm2Deg<0)
217 {
218 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
219 return;
220 }
221
222 if (fUseMmScale)
223 {
224 fHist.SetXTitle("x [mm]");
225 fHist.SetYTitle("y [mm]");
226
227 fHist.Scale(1./fMm2Deg);
228 }
229 else
230 {
231 fHist.SetXTitle("x [\\circ]");
232 fHist.SetYTitle("y [\\circ]");
233
234 fHist.Scale(1./fMm2Deg);
235 }
236
237 fUseMmScale = mmscale;
238}
239
240// --------------------------------------------------------------------------
241//
242// Calculate Significance as
243// significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
244//
245// s: total number of events in signal region
246// b: number of background events in signal region
247//
248Double_t MHFalseSource::Significance(Double_t s, Double_t b)
249{
250 const Double_t k = b==0 ? 0 : s/b;
251 const Double_t f = s+k*k*b;
252
253 return f==0 ? 0 : (s-b)/TMath::Sqrt(f);
254}
255
256// --------------------------------------------------------------------------
257//
258// calculates the significance according to Li & Ma
259// ApJ 272 (1983) 317, Formula 17
260//
261// s // s: number of on events
262// b // b: number of off events
263// alpha = t_on/t_off; // t: observation time
264//
265Double_t MHFalseSource::SignificanceLiMa(Double_t s, Double_t b, Double_t alpha)
266{
267 const Double_t sum = s+b;
268
269 if (sum<=0 || alpha<=0)
270 return 0;
271
272 const Double_t l = s*TMath::Log(s/sum*(alpha+1)/alpha);
273 const Double_t m = b*TMath::Log(b/sum*(alpha+1) );
274
275 return TMath::Sqrt((l+m)*2);
276}
277
278// --------------------------------------------------------------------------
279//
280// Set binnings (takes BinningFalseSource) and prepare filling of the
281// histogram.
282//
283// Also search for MTime, MObservatory and MPointingPos
284//
285Bool_t MHFalseSource::SetupFill(const MParList *plist)
286{
287 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
288 if (geom)
289 {
290 fMm2Deg = geom->GetConvMm2Deg();
291 fUseMmScale = kFALSE;
292
293 fHist.SetXTitle("x [\\circ]");
294 fHist.SetYTitle("y [\\circ]");
295 }
296
297 MBinning binsa;
298 binsa.SetEdges(18, 0, 90);
299
300 const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource");
301 if (!bins)
302 {
303 Float_t r = geom ? geom->GetMaxRadius() : 600;
304 r /= 3;
305 if (!fUseMmScale)
306 r *= fMm2Deg;
307
308 MBinning b;
309 b.SetEdges(20, -r, r);
310 SetBinning(&fHist, &b, &b, &binsa);
311 }
312 else
313 SetBinning(&fHist, bins, bins, &binsa);
314
315 fPointPos = (MPointingPos*)plist->FindObject(AddSerialNumber("MPointingPos"));
316 if (!fPointPos)
317 *fLog << warn << "MPointingPos not found... no derotation." << endl;
318
319 fTime = (MTime*)plist->FindObject(AddSerialNumber("MTime"));
320 if (!fTime)
321 *fLog << warn << "MTime not found... no derotation." << endl;
322
323 fObservatory = (MObservatory*)plist->FindObject(AddSerialNumber("MObservatory"));
324 if (!fObservatory)
325 *fLog << err << "MObservatory not found... no derotation." << endl;
326
327
328 return kTRUE;
329}
330
331// --------------------------------------------------------------------------
332//
333// Fill the histogram. For details see the code or the class description
334//
335Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w)
336{
337 MHillas *hil = (MHillas*)par;
338
339 const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1));
340
341 // Get camera rotation angle
342 Double_t rho = 0;
343 if (fTime && fObservatory && fPointPos)
344 rho = fPointPos->RotationAngle(*fObservatory, *fTime);
345 //if (fPointPos)
346 // rho = fPointPos->RotationAngle(*fObservatory);
347
348 MSrcPosCam src;
349 MHillasSrc hsrc;
350 hsrc.SetSrcPos(&src);
351
352 const Int_t nx = fHist.GetNbinsX();
353 const Int_t ny = fHist.GetNbinsY();
354
355 Axis_t cx[nx];
356 Axis_t cy[ny];
357 fHist.GetXaxis()->GetCenter(cx);
358 fHist.GetYaxis()->GetCenter(cy);
359
360 for (int ix=0; ix<nx; ix++)
361 {
362 for (int iy=0; iy<ny; iy++)
363 {
364 if (TMath::Hypot(cx[ix], cy[iy])>maxr)
365 continue;
366
367 TVector2 v(cx[ix], cy[iy]);
368 if (rho!=0)
369 v.Rotate(-rho);
370
371 if (!fUseMmScale)
372 v *= 1./fMm2Deg;
373
374 src.SetXY(v);
375
376 if (!hsrc.Calc(hil))
377 {
378 *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl;
379 return kFALSE;
380 }
381
382 const Double_t alpha = hsrc.GetAlpha();
383 fHist.Fill(cx[ix], cy[iy], TMath::Abs(alpha), w);
384 }
385 }
386
387 return kTRUE;
388}
389
390// --------------------------------------------------------------------------
391//
392// Create projection for off data, taking sum of bin contents of
393// range (fBgMean-fAlphaCut/2, fBgMean+fAlphaCut) Making sure to take
394// the same number of bins than for on-data
395//
396void MHFalseSource::ProjectOff(TH2D *h2)
397{
398 TAxis &axe = *fHist.GetZaxis();
399
400 // Find range to cut (left edge and width)
401 const Int_t f = axe.FindBin(fBgMean-fAlphaCut/2);
402 const Int_t l = axe.FindBin(fAlphaCut)+f-1;
403
404 axe.SetRange(f, l);
405 const Float_t cut1 = axe.GetBinLowEdge(f);
406 const Float_t cut2 = axe.GetBinUpEdge(l);
407 h2->SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2));
408
409 // Get projection for range
410 TH2D *p = (TH2D*)fHist.Project3D("xy_off");
411
412 // Reset range
413 axe.SetRange(0,9999);
414
415 // Move contents from projection to h2
416 h2->Reset();
417 h2->Add(p);
418
419 // Delete p
420 delete p;
421
422 // Set Minimum as minimum value Greater Than 0
423 h2->SetMinimum(GetMinimumGT(*h2));
424}
425
426// --------------------------------------------------------------------------
427//
428// Create projection for on data, taking sum of bin contents of
429// range (0, fAlphaCut)
430//
431void MHFalseSource::ProjectOn(TH2D *h3)
432{
433 TAxis &axe = *fHist.GetZaxis();
434
435 // Find range to cut
436 axe.SetRangeUser(0, fAlphaCut);
437 const Float_t cut = axe.GetBinUpEdge(axe.GetLast());
438 h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut));
439
440 // Get projection for range
441 TH2D *p = (TH2D*)fHist.Project3D("xy_on");
442
443 // Reset range
444 axe.SetRange(0,9999);
445
446 // Move contents from projection to h3
447 h3->Reset();
448 h3->Add(p);
449 delete p;
450
451 // Set Minimum as minimum value Greater Than 0
452 h3->SetMinimum(GetMinimumGT(*h3));
453}
454
455// --------------------------------------------------------------------------
456//
457// Update the projections
458//
459void MHFalseSource::Paint(Option_t *opt)
460{
461 // sigma = (s-b)/sqrt(s+k*k*b) mit k=s/b
462
463 gStyle->SetPalette(1, 0);
464
465 TVirtualPad *padsave = gPad;
466
467 TH1D* h1;
468 TH2D* h2;
469 TH2D* h3;
470 TH2D* h4;
471
472 padsave->cd(3);
473 if ((h3 = (TH2D*)gPad->FindObject("Alpha_xy_on")))
474 ProjectOn(h3);
475
476 padsave->cd(4);
477 if ((h2 = (TH2D*)gPad->FindObject("Alpha_xy_off")))
478 ProjectOff(h2);
479
480 padsave->cd(2);
481 if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_xy_sig")))
482 {
483 const Int_t nx = h4->GetXaxis()->GetNbins();
484 const Int_t ny = h4->GetYaxis()->GetNbins();
485
486 Int_t maxx=nx/2;
487 Int_t maxy=ny/2;
488
489 Int_t max = h4->GetBin(maxx, maxy);
490
491 for (int ix=0; ix<nx; ix++)
492 for (int iy=0; iy<ny; iy++)
493 {
494 const Int_t n = h4->GetBin(ix+1, iy+1);
495
496 const Double_t s = h3->GetBinContent(n);
497 const Double_t b = h2->GetBinContent(n);
498
499 const Double_t sig = Significance(s, b);
500
501 h4->SetBinContent(n, sig);
502
503 if (sig>h4->GetBinContent(max) && sig!=0)
504 {
505 max = n;
506 maxx=ix;
507 maxy=iy;
508 }
509 }
510
511 padsave->cd(1);
512 if ((h1 = (TH1D*)gPad->FindObject("Alpha")))
513 {
514 //h1->Reset();
515
516 const Double_t x = fHist.GetXaxis()->GetBinCenter(maxx);
517 const Double_t y = fHist.GetYaxis()->GetBinCenter(maxy);
518 const Double_t s = h4->GetBinContent(max);
519
520 TH1 *h = fHist.ProjectionZ("Alpha", maxx, maxx, maxy, maxy);
521 h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f)", x, y, s));
522 }
523 }
524
525 gPad = padsave;
526}
527
528// --------------------------------------------------------------------------
529//
530// Draw the histogram
531//
532void MHFalseSource::Draw(Option_t *opt)
533{
534 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
535 pad->SetBorderMode(0);
536
537 AppendPad("");
538
539 pad->Divide(2, 2);
540
541 // draw the 2D histogram Sigmabar versus Theta
542 pad->cd(1);
543 gPad->SetBorderMode(0);
544 TH1 *h1 = fHist.ProjectionZ("Alpha");
545 h1->SetDirectory(NULL);
546 h1->SetTitle("Distribution of \\alpha");
547 h1->SetXTitle(fHist.GetZaxis()->GetTitle());
548 h1->SetYTitle("Counts");
549 h1->Draw(opt);
550 h1->SetBit(kCanDelete);
551
552 pad->cd(4);
553 gPad->SetBorderMode(0);
554 fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2);
555 TH1 *h2 = fHist.Project3D("xy_off");
556 h2->SetDirectory(NULL);
557 h2->SetXTitle(fHist.GetXaxis()->GetTitle());
558 h2->SetYTitle(fHist.GetYaxis()->GetTitle());
559 h2->Draw("colz");
560 h2->SetBit(kCanDelete);
561
562 pad->cd(3);
563 gPad->SetBorderMode(0);
564 fHist.GetZaxis()->SetRangeUser(0,fAlphaCut);
565 TH1 *h3 = fHist.Project3D("xy_on");
566 fHist.GetZaxis()->SetRange(0,9999);
567 h3->SetDirectory(NULL);
568 h3->SetXTitle(fHist.GetXaxis()->GetTitle());
569 h3->SetYTitle(fHist.GetYaxis()->GetTitle());
570 h3->Draw("colz");
571 h3->SetBit(kCanDelete);
572
573 pad->cd(2);
574 gPad->SetBorderMode(0);
575 fHist.GetZaxis()->SetRange(0,0);
576 TH1 *h4 = fHist.Project3D("xy_sig"); // Do this to get the correct binning....
577 fHist.GetZaxis()->SetRange(0,9999);
578 h4->SetTitle("Significance");
579 h4->SetDirectory(NULL);
580 h4->SetXTitle(fHist.GetXaxis()->GetTitle());
581 h4->SetYTitle(fHist.GetYaxis()->GetTitle());
582 h4->Draw("colz");
583 h4->SetBit(kCanDelete);
584}
585
586// --------------------------------------------------------------------------
587//
588// Everything which is in the main pad belongs to this class!
589//
590Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py)
591{
592 return 0;
593}
594
595// --------------------------------------------------------------------------
596//
597// Set all sub-pads to Modified()
598//
599void MHFalseSource::Modified()
600{
601 if (!gPad)
602 return;
603
604 TVirtualPad *padsave = gPad;
605 padsave->Modified();
606 padsave->cd(1);
607 gPad->Modified();
608 padsave->cd(2);
609 gPad->Modified();
610 padsave->cd(3);
611 gPad->Modified();
612 padsave->cd(4);
613 gPad->Modified();
614 gPad->cd();
615}
616
617// --------------------------------------------------------------------------
618//
619// This is a preliminary implementation of a alpha-fit procedure for
620// all possible source positions. It will be moved into its own
621// more powerfull class soon.
622//
623// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
624// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
625// or
626// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
627//
628// Parameter [1] is fixed to 0 while the alpha peak should be
629// symmetric around alpha=0.
630//
631// Parameter [4] is fixed to 0 because the first derivative at
632// alpha=0 should be 0, too.
633//
634// In a first step the background is fitted between bgmin and bgmax,
635// while the parameters [0]=0 and [2]=1 are fixed.
636//
637// In a second step the signal region (alpha<sigmax) is fittet using
638// the whole function with parameters [1], [3], [4] and [5] fixed.
639//
640// The number of excess and background events are calculated as
641// s = int(0, sigint, gaus(0)+pol2(3))
642// b = int(0, sigint, pol2(3))
643//
644// The Significance is calculated using the Significance() member
645// function.
646//
647void MHFalseSource::FitSignificance(Float_t sigint, Float_t sigmax, Float_t bgmin, Float_t bgmax, Byte_t polynom)
648{
649 TH1D h0a("A", "", 50, 0, 4000);
650 TH1D h4a("chisq1", "", 50, 0, 35);
651 //TH1D h5a("prob1", "", 50, 0, 1.1);
652 TH1D h6("signifcance", "", 50, -20, 20);
653
654 TH1D h1("mu", "Parameter \\mu", 50, -1, 1);
655 TH1D h2("sigma", "Parameter \\sigma", 50, 0, 90);
656 TH1D h3("b", "Parameter b", 50, -0.1, 0.1);
657
658 TH1D h0b("a", "Parameter a (red), A (blue)", 50, 0, 4000);
659 TH1D h4b("\\chi^{2}", "\\chi^{2} (red, green) / significance (black)", 50, 0, 35);
660 //TH1D h5b("prob", "Fit probability: Bg(red), F(blue)", 50, 0, 1.1);
661
662 h0a.SetLineColor(kBlue);
663 h4a.SetLineColor(kBlue);
664 //h5a.SetLineColor(kBlue);
665 h0b.SetLineColor(kRed);
666 h4b.SetLineColor(kRed);
667 //h5b.SetLineColor(kRed);
668
669 TH1 *hist = fHist.Project3D("xy_fit");
670 hist->SetDirectory(0);
671 TH1 *hists = fHist.Project3D("xy_fit");
672 hists->SetDirectory(0);
673 TH1 *histb = fHist.Project3D("xy_fit");
674 histb->SetDirectory(0);
675 hist->Reset();
676 hists->Reset();
677 histb->Reset();
678 hist->SetNameTitle("Significance",
679 Form("Fit Region: Signal<%.1f\\circ, %.1f\\circ<Bg<%.1f\\circ",
680 sigmax, bgmin, bgmax));
681 hists->SetNameTitle("Excess", Form("Number of excess events for \\alpha<%.0f\\circ", sigint));
682 histb->SetNameTitle("Background", Form("Number of background events for \\alpha<%.0f\\circ", sigint));
683 hist->SetXTitle(fHist.GetXaxis()->GetTitle());
684 hists->SetXTitle(fHist.GetXaxis()->GetTitle());
685 histb->SetXTitle(fHist.GetXaxis()->GetTitle());
686 hist->SetYTitle(fHist.GetXaxis()->GetTitle());
687 hists->SetYTitle(fHist.GetXaxis()->GetTitle());
688 histb->SetYTitle(fHist.GetXaxis()->GetTitle());
689
690 const Double_t w = fHist.GetZaxis()->GetBinWidth(1);
691
692 // xmin, xmax, npar
693 //TF1 func("MyFunc", fcn, 0, 90, 6);
694 // Implementing the function yourself is only about 5% faster
695 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
696 TArrayD maxpar(func.GetNpar());
697
698 /*
699 func.SetParName(0, "A");
700 func.SetParName(1, "mu");
701 func.SetParName(2, "sigma");
702 func.SetParName(3, "a");
703 func.SetParName(4, "b");
704 func.SetParName(5, "c");
705 */
706
707 func.FixParameter(1, 0);
708 func.FixParameter(4, 0);
709 func.SetParLimits(2, 0, 90);
710 func.SetParLimits(3, -1, 1);
711
712 const Int_t nx = fHist.GetXaxis()->GetNbins();
713 const Int_t ny = fHist.GetYaxis()->GetNbins();
714
715 Double_t maxalpha0=0;
716 Double_t maxs=3;
717
718 Int_t maxx=0;
719 Int_t maxy=0;
720
721 TStopwatch clk;
722 clk.Start();
723
724 *fLog << inf;
725 *fLog << "Signal fit: alpha < " << sigmax << endl;
726 *fLog << "Integration: alpha < " << sigint << endl;
727 *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
728 *fLog << "Polynom order: " << (int)polynom << endl;
729 *fLog << "Fitting False Source Plot..." << flush;
730
731 TH1 *h=0;
732 for (int ix=1; ix<=nx; ix++)
733 for (int iy=1; iy<=ny; iy++)
734 {
735 h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy);
736
737 const Double_t alpha0 = h->GetBinContent(1);
738
739 // Check for the regios which is not filled...
740 if (alpha0==0)
741 continue;
742
743 if (alpha0>maxalpha0)
744 maxalpha0=alpha0;
745
746 // First fit a polynom in the off region
747 func.FixParameter(0, 0);
748 func.FixParameter(2, 1);
749 func.ReleaseParameter(3);
750
751 for (int i=5; i<func.GetNpar(); i++)
752 func.ReleaseParameter(i);
753
754 h->Fit(&func, "N0Q", "", bgmin, bgmax);
755 //*fLog << dbg << ix << "/" << iy << ": " << func.GetParameter(3) << " " << func.GetParameter(4) << endl;
756
757 h4a.Fill(func.GetChisquare());
758 //h5a.Fill(func.GetProb());
759
760 //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1));
761 //func.SetParLimits(2, 5, 90);
762
763 func.ReleaseParameter(0);
764 //func.ReleaseParameter(1);
765 func.ReleaseParameter(2);
766 func.FixParameter(3, func.GetParameter(3));
767 for (int i=5; i<func.GetNpar(); i++)
768 func.FixParameter(i, func.GetParameter(i));
769
770 // Do not allow signals smaller than the background
771 const Double_t A = alpha0-func.GetParameter(3);
772 const Double_t dA = TMath::Abs(A);
773 func.SetParLimits(0, -dA*4, dA*4);
774 func.SetParLimits(2, 0, 90);
775
776 // Now fit a gaus in the on region on top of the polynom
777 func.SetParameter(0, A);
778 func.SetParameter(2, sigmax*0.75);
779
780 h->Fit(&func, "N0Q", "", 0, sigmax);
781 //*fLog << dbg << " " << func.GetParameter(0) << " " << func.GetParameter(1) << " " << func.GetParameter(2) << endl;
782
783 TArrayD p(func.GetNpar(), func.GetParameters());
784
785 // Fill results into some histograms
786 h0a.Fill(p[0]);
787 h0b.Fill(p[3]);
788 h1.Fill(p[1]);
789 h2.Fill(p[2]);
790 if (polynom>1)
791 h3.Fill(p[5]);
792 h4b.Fill(func.GetChisquare());
793 //h5b.Fill(func.GetProb());
794
795 // Implementing the integral as analytical function
796 // gives the same result in the order of 10e-5
797 // and it is not faster at all...
798 //const Bool_t ok = /*p[0]>=0 && /*p[0]<alpha0*2 &&*/ p[2]>1.75;// && p[2]<88.5;
799 /*
800 if (p[0]<-fabs(alpha0-p[3])*7 && p[2]<alphaw/3)
801 {
802 func.SetParameter(0, alpha0-p[3]);
803 cout << p[2] << " " << p[0] << " " << alpha0-p[3] << endl;
804 }
805 */
806
807 // The fitted function returned units of
808 // counts bin binwidth. To get the correct number
809 // of events we must adapt the functions by dividing
810 // the result of the integration by the bin-width
811 const Double_t s = func.Integral(0, sigint)/w;
812
813 func.SetParameter(0, 0);
814 func.SetParameter(2, 1);
815
816 const Double_t b = func.Integral(0, sigint)/w;
817 const Double_t sig = Significance(s, b);
818
819 const Int_t n = hist->GetBin(ix, iy);
820 hists->SetBinContent(n, s-b);
821 histb->SetBinContent(n, b);
822
823 hist->SetBinContent(n, sig);
824 if (sig!=0)
825 h6.Fill(sig);
826
827 if (sig>maxs)
828 {
829 maxs = sig;
830 maxx = ix;
831 maxy = iy;
832 maxpar = p;
833 }
834 }
835
836 *fLog << inf << "done." << endl;
837
838 h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
839 h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
840
841 //hists->SetMinimum(GetMinimumGT(*hists));
842 histb->SetMinimum(GetMinimumGT(*histb));
843
844 clk.Stop();
845 clk.Print();
846
847 TCanvas *c=new TCanvas;
848
849 c->Divide(3,2, 0, 0);
850 c->cd(1);
851 gPad->SetBorderMode(0);
852 hists->Draw("colz");
853 hists->SetBit(kCanDelete);
854 c->cd(2);
855 gPad->SetBorderMode(0);
856 hist->Draw("colz");
857 hist->SetBit(kCanDelete);
858 c->cd(3);
859 gPad->SetBorderMode(0);
860 histb->Draw("colz");
861 histb->SetBit(kCanDelete);
862 c->cd(4);
863 gPad->Divide(1,3, 0, 0);
864 TVirtualPad *p=gPad;
865 p->SetBorderMode(0);
866 p->cd(1);
867 gPad->SetBorderMode(0);
868 h0b.DrawCopy();
869 h0a.DrawCopy("same");
870 p->cd(2);
871 gPad->SetBorderMode(0);
872 h3.DrawCopy();
873 p->cd(3);
874 gPad->SetBorderMode(0);
875 h2.DrawCopy();
876 c->cd(6);
877 gPad->Divide(1,2, 0, 0);
878 TVirtualPad *q=gPad;
879 q->SetBorderMode(0);
880 q->cd(1);
881 gPad->SetBorderMode(0);
882 gPad->SetBorderMode(0);
883 h4b.DrawCopy();
884 h4a.DrawCopy("same");
885 h6.DrawCopy("same");
886 q->cd(2);
887 gPad->SetBorderMode(0);
888 //h5b.DrawCopy();
889 //h5a.DrawCopy("same");
890 c->cd(5);
891 gPad->SetBorderMode(0);
892 if (maxx>0 && maxy>0)
893 {
894 const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ",
895 fHist.GetXaxis()->GetBinCenter(maxx),
896 fHist.GetYaxis()->GetBinCenter(maxy), maxs);
897
898 TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
899 result->SetDirectory(NULL);
900 result->SetNameTitle("Result \\alpha", title);
901 result->SetBit(kCanDelete);
902 result->SetXTitle("\\alpha [\\circ]");
903 result->SetYTitle("Counts");
904 result->Draw();
905
906 TF1 f1("", func.GetExpFormula(), 0, 90);
907 TF1 f2("", func.GetExpFormula(), 0, 90);
908 f1.SetParameters(maxpar.GetArray());
909 f2.SetParameters(maxpar.GetArray());
910 f2.FixParameter(0, 0);
911 f2.FixParameter(1, 0);
912 f2.FixParameter(2, 1);
913 f1.SetLineColor(kGreen);
914 f2.SetLineColor(kRed);
915
916 f2.DrawCopy("same");
917 f1.DrawCopy("same");
918
919 TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
920 leg->SetBorderSize(1);
921 leg->SetTextSize(0.04);
922 leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
923 //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
924 leg->AddLine(0, 0.65, 0, 0.65);
925 leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
926 leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
927 leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
928 leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
929 leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
930 if (polynom>1)
931 leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11);
932 leg->SetBit(kCanDelete);
933 leg->Draw();
934
935 q->cd(2);
936
937 TGraph *g = new TGraph;
938 g->SetPoint(0, 0, 0);
939
940 Int_t max=0;
941 Float_t maxsig=0;
942 for (int i=1; i<89; i++)
943 {
944 const Double_t s = f1.Integral(0, (float)i);
945 const Double_t b = f2.Integral(0, (float)i);
946
947 const Double_t sig = Significance(s, b);
948
949 g->SetPoint(g->GetN(), i, sig);
950
951 if (sig>maxsig)
952 {
953 max = i;
954 maxsig = sig;
955 }
956 }
957
958 g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
959 g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
960 g->GetHistogram()->SetYTitle("Significance");
961 g->SetBit(kCanDelete);
962 g->Draw("AP");
963
964 leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
965 leg->SetBorderSize(1);
966 leg->SetTextSize(0.1);
967 leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max));
968 leg->SetBit(kCanDelete);
969 leg->Draw();
970 }
971}
Note: See TracBrowser for help on using the repository browser.