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

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