source: trunk/MagicSoft/Mars/mhist/MHAlpha.cc@ 3778

Last change on this file since 3778 was 3778, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 18.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// MHAlpha
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// MHAlpha 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("MHAlpha", "MHillas");
77// or
78// MHAlpha 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 "MHAlpha.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 "MParameters.h"
122#include "MHMatrix.h"
123
124#include "MBinning.h"
125#include "MParList.h"
126
127#include "MLog.h"
128#include "MLogManip.h"
129
130ClassImp(MHAlpha);
131
132using namespace std;
133
134// --------------------------------------------------------------------------
135//
136// Default Constructor
137//
138MHAlpha::MHAlpha(const char *name, const char *title)
139 : fAlphaCut(12.5), fBgMean(55), fResult(0), fMatrix(0)
140{
141 //
142 // set the name and title of this object
143 //
144 fName = name ? name : "MHAlpha";
145 fTitle = title ? title : "Alpha plot";
146
147 fHist.SetDirectory(NULL);
148
149 fHist.SetName("Alpha");
150 fHist.SetTitle("Alpha");
151 fHist.SetXTitle("|\\alpha| [\\circ]");
152 fHist.SetYTitle("Counts");
153 fHist.UseCurrentStyle();
154
155 MBinning binsa;
156 binsa.SetEdges(18, 0, 90);
157 binsa.Apply(fHist);
158
159 fSigInt=15;
160 fSigMax=75;
161 fBgMin=45;
162 fBgMax=90;
163 fPolynom=1;
164}
165
166// --------------------------------------------------------------------------
167//
168// Calculate Significance as
169// significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
170//
171Double_t MHAlpha::Significance(Double_t s, Double_t b)
172{
173 const Double_t k = b==0 ? 0 : s/b;
174 const Double_t f = s+k*k*b;
175
176 return f==0 ? 0 : (s-b)/TMath::Sqrt(f);
177}
178
179Bool_t MHAlpha::SetupFill(const MParList *pl)
180{
181 fHist.Reset();
182
183 fResult = (MParameterD*)pl->FindObject("Significance", "MParameterD");
184 if (!fResult)
185 *fLog << warn << "Significance [MParameterD] not found... ignored." << endl;
186
187 return kTRUE;
188}
189
190// --------------------------------------------------------------------------
191//
192// Fill the histogram. For details see the code or the class description
193//
194Bool_t MHAlpha::Fill(const MParContainer *par, const Stat_t w)
195{
196 Double_t alpha;
197
198 if (fMatrix)
199 alpha = (*fMatrix)[fMap];
200 else
201 {
202 const MHillasSrc *hil = dynamic_cast<const MHillasSrc*>(par);
203 if (!par)
204 return kFALSE;
205
206 alpha = hil->GetAlpha();
207 }
208
209 fHist.Fill(TMath::Abs(alpha), w);
210
211 return kTRUE;
212}
213
214// --------------------------------------------------------------------------
215//
216// Update the projections
217//
218void MHAlpha::Paint(Option_t *opt)
219{
220 Float_t sigint=fSigInt;
221 Float_t sigmax=fSigMax;
222 Float_t bgmin =fBgMin;
223 Float_t bgmax =fBgMax;
224 Byte_t polynom=fPolynom;
225
226 // Implementing the function yourself is only about 5% faster
227 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
228 //TF1 func("", Form("[0]*(TMath::Gaus(x, [1], [2])+TMath::Gaus(x, -[1], [2]))+pol%d(3)", polynom), 0, 90);
229 TArrayD maxpar(func.GetNpar());
230
231 func.FixParameter(1, 0);
232 func.FixParameter(4, 0);
233 func.SetParLimits(2, 0, 90);
234 func.SetParLimits(3, -1, 1);
235
236 TH1 *h=&fHist;
237 //const Double_t alpha0 = h->GetBinContent(1);
238 //const Double_t alphaw = h->GetXaxis()->GetBinWidth(1);
239
240 // Check for the regios which is not filled...
241 const Double_t alpha0 = h->GetBinContent(1);
242 if (alpha0==0)
243 return;
244
245 // First fit a polynom in the off region
246 func.FixParameter(0, 0);
247 func.FixParameter(2, 1);
248 func.ReleaseParameter(3);
249
250 for (int i=5; i<func.GetNpar(); i++)
251 func.ReleaseParameter(i);
252
253 //h->Fit(&func, "N0Q", "", bgmin, bgmax);
254 h->Fit(&func, "N0Q", "", bgmin, bgmax);
255 func.SetRange(0, 90);
256 func.SetLineColor(kRed);
257 func.Paint("same");
258
259 //h4a.Fill(func.GetChisquare());
260 //h5a.Fill(func.GetProb());
261
262 func.ReleaseParameter(0);
263 //func.ReleaseParameter(1);
264 func.ReleaseParameter(2);
265 func.FixParameter(3, func.GetParameter(3));
266 for (int i=5; i<func.GetNpar(); i++)
267 func.FixParameter(i, func.GetParameter(i));
268
269 // Do not allow signals smaller than the background
270 const Double_t A = alpha0-func.GetParameter(3);
271 const Double_t dA = TMath::Abs(A);
272 func.SetParLimits(0, -dA*4, dA*4);
273 func.SetParLimits(2, 0, 90);
274
275 // Now fit a gaus in the on region on top of the polynom
276 func.SetParameter(0, A);
277 func.SetParameter(2, sigmax*0.75);
278
279 //h->Fit(&func, "N0Q", "", 0, sigmax);
280 h->Fit(&func, "N0Q", "", 0, sigmax);
281 func.SetLineColor(kGreen);
282 func.Paint("same");
283}
284
285// --------------------------------------------------------------------------
286//
287// Draw the histogram
288//
289void MHAlpha::Draw(Option_t *opt)
290{
291 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
292 pad->SetBorderMode(0);
293
294 fHist.Draw();
295
296 AppendPad("");
297}
298
299// --------------------------------------------------------------------------
300//
301// This is a preliminary implementation of a alpha-fit procedure for
302// all possible source positions. It will be moved into its own
303// more powerfull class soon.
304//
305// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
306// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
307// or
308// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
309//
310// Parameter [1] is fixed to 0 while the alpha peak should be
311// symmetric around alpha=0.
312//
313// Parameter [4] is fixed to 0 because the first derivative at
314// alpha=0 should be 0, too.
315//
316// In a first step the background is fitted between bgmin and bgmax,
317// while the parameters [0]=0 and [2]=1 are fixed.
318//
319// In a second step the signal region (alpha<sigmax) is fittet using
320// the whole function with parameters [1], [3], [4] and [5] fixed.
321//
322// The number of excess and background events are calculated as
323// s = int(0, sigint, gaus(0)+pol2(3))
324// b = int(0, sigint, pol2(3))
325//
326// The Significance is calculated using the Significance() member
327// function.
328//
329Bool_t MHAlpha::Finalize()
330{
331 //*fLog << dbg << "Fitting..." << endl;
332
333 Float_t sigint=fSigInt;
334 Float_t sigmax=fSigMax;
335 Float_t bgmin=fBgMin;
336 Float_t bgmax=fBgMax;
337 Byte_t polynom=fPolynom;
338
339 // xmin, xmax, npar
340 //TF1 func("MyFunc", fcn, 0, 90, 6);
341 // Implementing the function yourself is only about 5% faster
342 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
343 TArrayD maxpar(func.GetNpar());
344
345 //*fLog << "Setup pars..." << endl;
346
347 func.FixParameter(1, 0);
348 func.FixParameter(4, 0);
349 func.SetParLimits(2, 0, 90);
350 func.SetParLimits(3, -1, 1);
351
352 // TStopwatch clk;
353 // clk.Start();
354 /*
355 *fLog << inf;
356 *fLog << "Signal fit: alpha < " << sigmax << endl;
357 *fLog << "Integration: alpha < " << sigint << endl;
358 *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
359 *fLog << "Polynom order: " << (int)polynom << endl;
360 *fLog << "Fitting False Source Plot..." << flush;
361 */
362 TH1 *h=&fHist;
363 const Double_t alpha0 = h->GetBinContent(1);
364 const Double_t alphaw = h->GetXaxis()->GetBinWidth(1);
365
366 // Check for the regios which is not filled...
367 if (alpha0==0)
368 {
369 *fLog << warn << "Histogram empty." << endl;
370 return kTRUE;
371 }
372
373 // First fit a polynom in the off region
374 func.FixParameter(0, 0);
375 func.FixParameter(2, 1);
376 func.ReleaseParameter(3);
377
378 //*fLog << "Release pars..." << endl;
379
380 for (int i=5; i<func.GetNpar(); i++)
381 func.ReleaseParameter(i);
382
383 //*fLog << dbg << "Fit 1" << endl;
384 //h->Fit(&func, "N0Q", "", bgmin, bgmax);
385 h->Fit(&func, "N0Q", "", bgmin, bgmax);
386 //*fLog << dbg << ix << "/" << iy << ": " << func.GetParameter(3) << " " << func.GetParameter(4) << endl;
387
388 //h4a.Fill(func.GetChisquare());
389 //h5a.Fill(func.GetProb());
390
391 //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1));
392 //func.SetParLimits(2, 5, 90);
393
394 //*fLog << dbg << "Change params..." << endl;
395
396 func.ReleaseParameter(0);
397 //func.ReleaseParameter(1);
398 func.ReleaseParameter(2);
399 func.FixParameter(3, func.GetParameter(3));
400 for (int i=5; i<func.GetNpar(); i++)
401 func.FixParameter(i, func.GetParameter(i));
402
403 // Do not allow signals smaller than the background
404 const Double_t A = alpha0-func.GetParameter(3);
405 const Double_t dA = TMath::Abs(A);
406 func.SetParLimits(0, -dA*4, dA*4);
407 func.SetParLimits(2, 0, 90);
408
409 // Now fit a gaus in the on region on top of the polynom
410 func.SetParameter(0, A);
411 func.SetParameter(2, sigmax*0.75);
412
413 //*fLog << dbg << "Fit 2..." << endl;
414 //h->Fit(&func, "N0Q", "", 0, sigmax);
415 h->Fit(&func, "N0Q", "", 0, sigmax);
416 //*fLog << dbg << " " << func.GetParameter(0) << " " << func.GetParameter(1) << " " << func.GetParameter(2) << endl;
417
418 //*fLog << dbg << "Copy par..." << endl;
419 TArrayD p(func.GetNpar(), func.GetParameters());
420
421 // Fill results into some histograms
422 /*
423 h0a.Fill(p[0]);
424 h0b.Fill(p[3]);
425 h1.Fill(p[1]);
426 h2.Fill(p[2]);
427 if (polynom>1)
428 h3.Fill(p[5]);
429 h4b.Fill(func.GetChisquare());
430 //h5b.Fill(func.GetProb());
431 */
432 // Implementing the integral as analytical function
433 // gives the same result in the order of 10e-5
434 // and it is not faster at all...
435
436 //*fLog << dbg << "Int1 par..." << endl;
437 const Double_t s = func.Integral(0, sigint)/alphaw;
438
439 func.SetParameter(0, 0);
440 //func.SetParameter(2, 1);
441
442 //*fLog << dbg << "Int2 par..." << endl;
443 const Double_t b = func.Integral(0, sigint)/alphaw;
444 const Double_t sig = Significance(s, b);
445
446 //*fLog << dbg << "Set Val "<<sig<<" ..." << fResult << endl;
447 if (fResult)
448 fResult->SetVal(sig);
449
450 *fLog << inf << "Found Significance: " << sig << " (width=";
451 *fLog << func.GetParameter(2) << "°, excess=" << s-b << ")" << endl;
452
453 //*fLog << dbg << "Done." << endl;
454 return kTRUE;
455 /*
456 if (sig!=0)
457 h6.Fill(sig);
458
459 if (sig<maxs)
460 {
461 maxs = sig;
462 maxx = ix;
463 maxy = iy;
464 maxpar = p;
465 }
466 */
467/*
468
469 *fLog << inf << "done." << endl;
470
471 h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
472 h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
473
474 //hists->SetMinimum(GetMinimumGT(*hists));
475 histb->SetMinimum(GetMinimumGT(*histb));
476
477// clk.Stop();
478// clk.Print();
479 TCanvas *c=new TCanvas;
480
481 c->Divide(3,2, 0, 0);
482 c->cd(1);
483 gPad->SetBorderMode(0);
484 hists->Draw("colz");
485 hists->SetBit(kCanDelete);
486 c->cd(2);
487 gPad->SetBorderMode(0);
488 hist->Draw("colz");
489 hist->SetBit(kCanDelete);
490 c->cd(3);
491 gPad->SetBorderMode(0);
492 histb->Draw("colz");
493 histb->SetBit(kCanDelete);
494 c->cd(4);
495 gPad->Divide(1,3, 0, 0);
496 TVirtualPad *p=gPad;
497 p->SetBorderMode(0);
498 p->cd(1);
499 gPad->SetBorderMode(0);
500 h0b.DrawCopy();
501 h0a.DrawCopy("same");
502 p->cd(2);
503 gPad->SetBorderMode(0);
504 h3.DrawCopy();
505 p->cd(3);
506 gPad->SetBorderMode(0);
507 h2.DrawCopy();
508 c->cd(6);
509 gPad->Divide(1,2, 0, 0);
510 TVirtualPad *q=gPad;
511 q->SetBorderMode(0);
512 q->cd(1);
513 gPad->SetBorderMode(0);
514 gPad->SetBorderMode(0);
515 h4b.DrawCopy();
516 h4a.DrawCopy("same");
517 h6.DrawCopy("same");
518 q->cd(2);
519 gPad->SetBorderMode(0);
520 //h5b.DrawCopy();
521 //h5a.DrawCopy("same");
522 c->cd(5);
523 gPad->SetBorderMode(0);
524 if (maxx>0 && maxy>0)
525 {
526 const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ",
527 fHist.GetXaxis()->GetBinCenter(maxx),
528 fHist.GetYaxis()->GetBinCenter(maxy), maxs);
529
530 TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
531 result->SetDirectory(NULL);
532 result->SetNameTitle("Result \\alpha", title);
533 result->SetBit(kCanDelete);
534 result->SetXTitle("\\alpha [\\circ]");
535 result->SetYTitle("Counts");
536 result->Draw();
537
538 TF1 f1("", func.GetExpFormula(), 0, 90);
539 TF1 f2("", func.GetExpFormula(), 0, 90);
540 f1.SetParameters(maxpar.GetArray());
541 f2.SetParameters(maxpar.GetArray());
542 f2.FixParameter(0, 0);
543 f2.FixParameter(1, 0);
544 f2.FixParameter(2, 1);
545 f1.SetLineColor(kGreen);
546 f2.SetLineColor(kRed);
547
548 f2.DrawCopy("same");
549 f1.DrawCopy("same");
550
551 TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
552 leg->SetBorderSize(1);
553 leg->SetTextSize(0.04);
554 leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
555 //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
556 leg->AddLine(0, 0.65, 0, 0.65);
557 leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
558 leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
559 leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
560 leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
561 leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
562 if (polynom>1)
563 leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11);
564 leg->SetBit(kCanDelete);
565 leg->Draw();
566
567 q->cd(2);
568
569 TGraph *g = new TGraph;
570 g->SetPoint(0, 0, 0);
571
572 Int_t max=0;
573 Float_t maxsig=0;
574 for (int i=1; i<89; i++)
575 {
576 const Double_t s = f1.Integral(0, (float)i);
577 const Double_t b = f2.Integral(0, (float)i);
578
579 const Double_t sig = Significance(s, b);
580
581 g->SetPoint(g->GetN(), i, sig);
582
583 if (sig>maxsig)
584 {
585 max = i;
586 maxsig = sig;
587 }
588 }
589
590 g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
591 g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
592 g->GetHistogram()->SetYTitle("Significance");
593 g->SetBit(kCanDelete);
594 g->Draw("AP");
595
596 leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
597 leg->SetBorderSize(1);
598 leg->SetTextSize(0.1);
599 leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max));
600 leg->SetBit(kCanDelete);
601 leg->Draw();
602 }*/
603}
604
605// --------------------------------------------------------------------------
606//
607// You can use this function if you want to use a MHMatrix instead of
608// MMcEvt. This function adds all necessary columns to the
609// given matrix. Afterward you should fill the matrix with the corresponding
610// data (eg from a file by using MHMatrix::Fill). If you now loop
611// through the matrix (eg using MMatrixLoop) MHHadronness::Fill
612// will take the values from the matrix instead of the containers.
613//
614void MHAlpha::InitMapping(MHMatrix *mat)
615{
616 if (fMatrix)
617 return;
618
619 fMatrix = mat;
620
621 fMap = fMatrix->AddColumn("MHillasSrc.fAlpha");
622}
623
624void MHAlpha::StopMapping()
625{
626 fMatrix = NULL;
627}
Note: See TracBrowser for help on using the repository browser.