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

Last change on this file since 4693 was 3779, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 16.1 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MHAlpha
28//
29// Create a single Alpha-Plot. The alpha-plot is fitted online. You can
30// check the result when it is filles in the MStatusDisplay
31// For more information see MHFalseSource::FitSignificance
32//
33// For convinience (fit) the output significance is stored in a
34// container in the parlisrt
35//
36// PRELIMINARY!
37//
38//////////////////////////////////////////////////////////////////////////////
39#include "MHAlpha.h"
40
41#include <TF1.h>
42#include <TH2.h>
43#include <TGraph.h>
44#include <TStyle.h>
45#include <TCanvas.h>
46#include <TPaveText.h>
47#include <TStopwatch.h>
48
49#include "MGeomCam.h"
50#include "MSrcPosCam.h"
51#include "MHillasSrc.h"
52#include "MTime.h"
53#include "MObservatory.h"
54#include "MPointingPos.h"
55#include "MAstroSky2Local.h"
56#include "MStatusDisplay.h"
57#include "MParameters.h"
58#include "MHMatrix.h"
59
60#include "MBinning.h"
61#include "MParList.h"
62
63#include "MLog.h"
64#include "MLogManip.h"
65
66ClassImp(MHAlpha);
67
68using namespace std;
69
70// --------------------------------------------------------------------------
71//
72// Default Constructor
73//
74MHAlpha::MHAlpha(const char *name, const char *title)
75 : fAlphaCut(12.5), fBgMean(55), fResult(0), fMatrix(0)
76{
77 //
78 // set the name and title of this object
79 //
80 fName = name ? name : "MHAlpha";
81 fTitle = title ? title : "Alpha plot";
82
83 fHist.SetDirectory(NULL);
84
85 fHist.SetName("Alpha");
86 fHist.SetTitle("Alpha");
87 fHist.SetXTitle("|\\alpha| [\\circ]");
88 fHist.SetYTitle("Counts");
89 fHist.UseCurrentStyle();
90
91 MBinning binsa;
92 binsa.SetEdges(18, 0, 90);
93 binsa.Apply(fHist);
94
95 fSigInt=15;
96 fSigMax=75;
97 fBgMin=45;
98 fBgMax=90;
99 fPolynom=1;
100}
101
102// --------------------------------------------------------------------------
103//
104// Calculate Significance as
105// significance = (s-b)/sqrt(s+k*k*b) mit k=s/b
106//
107Double_t MHAlpha::Significance(Double_t s, Double_t b)
108{
109 const Double_t k = b==0 ? 0 : s/b;
110 const Double_t f = s+k*k*b;
111
112 return f==0 ? 0 : (s-b)/TMath::Sqrt(f);
113}
114
115Bool_t MHAlpha::SetupFill(const MParList *pl)
116{
117 fHist.Reset();
118
119 fResult = (MParameterD*)pl->FindObject("Significance", "MParameterD");
120 if (!fResult)
121 *fLog << warn << "Significance [MParameterD] not found... ignored." << endl;
122
123 return kTRUE;
124}
125
126// --------------------------------------------------------------------------
127//
128// Fill the histogram. For details see the code or the class description
129//
130Bool_t MHAlpha::Fill(const MParContainer *par, const Stat_t w)
131{
132 Double_t alpha;
133
134 if (fMatrix)
135 alpha = (*fMatrix)[fMap];
136 else
137 {
138 const MHillasSrc *hil = dynamic_cast<const MHillasSrc*>(par);
139 if (!par)
140 return kFALSE;
141
142 alpha = hil->GetAlpha();
143 }
144
145 fHist.Fill(TMath::Abs(alpha), w);
146
147 return kTRUE;
148}
149
150// --------------------------------------------------------------------------
151//
152// Update the projections
153//
154void MHAlpha::Paint(Option_t *opt)
155{
156 Float_t sigint=fSigInt;
157 Float_t sigmax=fSigMax;
158 Float_t bgmin =fBgMin;
159 Float_t bgmax =fBgMax;
160 Byte_t polynom=fPolynom;
161
162 // Implementing the function yourself is only about 5% faster
163 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
164 //TF1 func("", Form("[0]*(TMath::Gaus(x, [1], [2])+TMath::Gaus(x, -[1], [2]))+pol%d(3)", polynom), 0, 90);
165 TArrayD maxpar(func.GetNpar());
166
167 func.FixParameter(1, 0);
168 func.FixParameter(4, 0);
169 func.SetParLimits(2, 0, 90);
170 func.SetParLimits(3, -1, 1);
171
172 TH1 *h=&fHist;
173 //const Double_t alpha0 = h->GetBinContent(1);
174 //const Double_t alphaw = h->GetXaxis()->GetBinWidth(1);
175
176 // Check for the regios which is not filled...
177 const Double_t alpha0 = h->GetBinContent(1);
178 if (alpha0==0)
179 return;
180
181 // First fit a polynom in the off region
182 func.FixParameter(0, 0);
183 func.FixParameter(2, 1);
184 func.ReleaseParameter(3);
185
186 for (int i=5; i<func.GetNpar(); i++)
187 func.ReleaseParameter(i);
188
189 //h->Fit(&func, "N0Q", "", bgmin, bgmax);
190 h->Fit(&func, "N0Q", "", bgmin, bgmax);
191 func.SetRange(0, 90);
192 func.SetLineColor(kRed);
193 func.Paint("same");
194
195 //h4a.Fill(func.GetChisquare());
196 //h5a.Fill(func.GetProb());
197
198 func.ReleaseParameter(0);
199 //func.ReleaseParameter(1);
200 func.ReleaseParameter(2);
201 func.FixParameter(3, func.GetParameter(3));
202 for (int i=5; i<func.GetNpar(); i++)
203 func.FixParameter(i, func.GetParameter(i));
204
205 // Do not allow signals smaller than the background
206 const Double_t A = alpha0-func.GetParameter(3);
207 const Double_t dA = TMath::Abs(A);
208 func.SetParLimits(0, -dA*4, dA*4);
209 func.SetParLimits(2, 0, 90);
210
211 // Now fit a gaus in the on region on top of the polynom
212 func.SetParameter(0, A);
213 func.SetParameter(2, sigmax*0.75);
214
215 //h->Fit(&func, "N0Q", "", 0, sigmax);
216 h->Fit(&func, "N0Q", "", 0, sigmax);
217 func.SetLineColor(kGreen);
218 func.Paint("same");
219}
220
221// --------------------------------------------------------------------------
222//
223// Draw the histogram
224//
225void MHAlpha::Draw(Option_t *opt)
226{
227 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
228 pad->SetBorderMode(0);
229
230 fHist.Draw();
231
232 AppendPad("");
233}
234
235// --------------------------------------------------------------------------
236//
237// This is a preliminary implementation of a alpha-fit procedure for
238// all possible source positions. It will be moved into its own
239// more powerfull class soon.
240//
241// The fit function is "gaus(0)+pol2(3)" which is equivalent to:
242// [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2
243// or
244// A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2
245//
246// Parameter [1] is fixed to 0 while the alpha peak should be
247// symmetric around alpha=0.
248//
249// Parameter [4] is fixed to 0 because the first derivative at
250// alpha=0 should be 0, too.
251//
252// In a first step the background is fitted between bgmin and bgmax,
253// while the parameters [0]=0 and [2]=1 are fixed.
254//
255// In a second step the signal region (alpha<sigmax) is fittet using
256// the whole function with parameters [1], [3], [4] and [5] fixed.
257//
258// The number of excess and background events are calculated as
259// s = int(0, sigint, gaus(0)+pol2(3))
260// b = int(0, sigint, pol2(3))
261//
262// The Significance is calculated using the Significance() member
263// function.
264//
265Bool_t MHAlpha::Finalize()
266{
267 //*fLog << dbg << "Fitting..." << endl;
268
269 Float_t sigint=fSigInt;
270 Float_t sigmax=fSigMax;
271 Float_t bgmin=fBgMin;
272 Float_t bgmax=fBgMax;
273 Byte_t polynom=fPolynom;
274
275 // xmin, xmax, npar
276 //TF1 func("MyFunc", fcn, 0, 90, 6);
277 // Implementing the function yourself is only about 5% faster
278 TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90);
279 TArrayD maxpar(func.GetNpar());
280
281 //*fLog << "Setup pars..." << endl;
282
283 func.FixParameter(1, 0);
284 func.FixParameter(4, 0);
285 func.SetParLimits(2, 0, 90);
286 func.SetParLimits(3, -1, 1);
287
288 // TStopwatch clk;
289 // clk.Start();
290 /*
291 *fLog << inf;
292 *fLog << "Signal fit: alpha < " << sigmax << endl;
293 *fLog << "Integration: alpha < " << sigint << endl;
294 *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl;
295 *fLog << "Polynom order: " << (int)polynom << endl;
296 *fLog << "Fitting False Source Plot..." << flush;
297 */
298 TH1 *h=&fHist;
299 const Double_t alpha0 = h->GetBinContent(1);
300 const Double_t alphaw = h->GetXaxis()->GetBinWidth(1);
301
302 // Check for the regios which is not filled...
303 if (alpha0==0)
304 {
305 *fLog << warn << "Histogram empty." << endl;
306 return kTRUE;
307 }
308
309 // First fit a polynom in the off region
310 func.FixParameter(0, 0);
311 func.FixParameter(2, 1);
312 func.ReleaseParameter(3);
313
314 //*fLog << "Release pars..." << endl;
315
316 for (int i=5; i<func.GetNpar(); i++)
317 func.ReleaseParameter(i);
318
319 //*fLog << dbg << "Fit 1" << endl;
320 //h->Fit(&func, "N0Q", "", bgmin, bgmax);
321 h->Fit(&func, "N0Q", "", bgmin, bgmax);
322 //*fLog << dbg << ix << "/" << iy << ": " << func.GetParameter(3) << " " << func.GetParameter(4) << endl;
323
324 //h4a.Fill(func.GetChisquare());
325 //h5a.Fill(func.GetProb());
326
327 //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1));
328 //func.SetParLimits(2, 5, 90);
329
330 //*fLog << dbg << "Change params..." << endl;
331
332 func.ReleaseParameter(0);
333 //func.ReleaseParameter(1);
334 func.ReleaseParameter(2);
335 func.FixParameter(3, func.GetParameter(3));
336 for (int i=5; i<func.GetNpar(); i++)
337 func.FixParameter(i, func.GetParameter(i));
338
339 // Do not allow signals smaller than the background
340 const Double_t A = alpha0-func.GetParameter(3);
341 const Double_t dA = TMath::Abs(A);
342 func.SetParLimits(0, -dA*4, dA*4);
343 func.SetParLimits(2, 0, 90);
344
345 // Now fit a gaus in the on region on top of the polynom
346 func.SetParameter(0, A);
347 func.SetParameter(2, sigmax*0.75);
348
349 //*fLog << dbg << "Fit 2..." << endl;
350 //h->Fit(&func, "N0Q", "", 0, sigmax);
351 h->Fit(&func, "N0Q", "", 0, sigmax);
352 //*fLog << dbg << " " << func.GetParameter(0) << " " << func.GetParameter(1) << " " << func.GetParameter(2) << endl;
353
354 //*fLog << dbg << "Copy par..." << endl;
355 TArrayD p(func.GetNpar(), func.GetParameters());
356
357 // Fill results into some histograms
358 /*
359 h0a.Fill(p[0]);
360 h0b.Fill(p[3]);
361 h1.Fill(p[1]);
362 h2.Fill(p[2]);
363 if (polynom>1)
364 h3.Fill(p[5]);
365 h4b.Fill(func.GetChisquare());
366 //h5b.Fill(func.GetProb());
367 */
368 // Implementing the integral as analytical function
369 // gives the same result in the order of 10e-5
370 // and it is not faster at all...
371
372 //*fLog << dbg << "Int1 par..." << endl;
373 const Double_t s = func.Integral(0, sigint)/alphaw;
374
375 func.SetParameter(0, 0);
376 //func.SetParameter(2, 1);
377
378 //*fLog << dbg << "Int2 par..." << endl;
379 const Double_t b = func.Integral(0, sigint)/alphaw;
380 const Double_t sig = Significance(s, b);
381
382 //*fLog << dbg << "Set Val "<<sig<<" ..." << fResult << endl;
383 if (fResult)
384 fResult->SetVal(sig);
385
386 *fLog << inf << "Found Significance: " << sig << " (width=";
387 *fLog << func.GetParameter(2) << "°, excess=" << s-b << ")" << endl;
388
389 //*fLog << dbg << "Done." << endl;
390 return kTRUE;
391 /*
392 if (sig!=0)
393 h6.Fill(sig);
394
395 if (sig<maxs)
396 {
397 maxs = sig;
398 maxx = ix;
399 maxy = iy;
400 maxpar = p;
401 }
402 */
403/*
404
405 *fLog << inf << "done." << endl;
406
407 h0a.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
408 h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5);
409
410 //hists->SetMinimum(GetMinimumGT(*hists));
411 histb->SetMinimum(GetMinimumGT(*histb));
412
413// clk.Stop();
414// clk.Print();
415 TCanvas *c=new TCanvas;
416
417 c->Divide(3,2, 0, 0);
418 c->cd(1);
419 gPad->SetBorderMode(0);
420 hists->Draw("colz");
421 hists->SetBit(kCanDelete);
422 c->cd(2);
423 gPad->SetBorderMode(0);
424 hist->Draw("colz");
425 hist->SetBit(kCanDelete);
426 c->cd(3);
427 gPad->SetBorderMode(0);
428 histb->Draw("colz");
429 histb->SetBit(kCanDelete);
430 c->cd(4);
431 gPad->Divide(1,3, 0, 0);
432 TVirtualPad *p=gPad;
433 p->SetBorderMode(0);
434 p->cd(1);
435 gPad->SetBorderMode(0);
436 h0b.DrawCopy();
437 h0a.DrawCopy("same");
438 p->cd(2);
439 gPad->SetBorderMode(0);
440 h3.DrawCopy();
441 p->cd(3);
442 gPad->SetBorderMode(0);
443 h2.DrawCopy();
444 c->cd(6);
445 gPad->Divide(1,2, 0, 0);
446 TVirtualPad *q=gPad;
447 q->SetBorderMode(0);
448 q->cd(1);
449 gPad->SetBorderMode(0);
450 gPad->SetBorderMode(0);
451 h4b.DrawCopy();
452 h4a.DrawCopy("same");
453 h6.DrawCopy("same");
454 q->cd(2);
455 gPad->SetBorderMode(0);
456 //h5b.DrawCopy();
457 //h5a.DrawCopy("same");
458 c->cd(5);
459 gPad->SetBorderMode(0);
460 if (maxx>0 && maxy>0)
461 {
462 const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ",
463 fHist.GetXaxis()->GetBinCenter(maxx),
464 fHist.GetYaxis()->GetBinCenter(maxy), maxs);
465
466 TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy);
467 result->SetDirectory(NULL);
468 result->SetNameTitle("Result \\alpha", title);
469 result->SetBit(kCanDelete);
470 result->SetXTitle("\\alpha [\\circ]");
471 result->SetYTitle("Counts");
472 result->Draw();
473
474 TF1 f1("", func.GetExpFormula(), 0, 90);
475 TF1 f2("", func.GetExpFormula(), 0, 90);
476 f1.SetParameters(maxpar.GetArray());
477 f2.SetParameters(maxpar.GetArray());
478 f2.FixParameter(0, 0);
479 f2.FixParameter(1, 0);
480 f2.FixParameter(2, 1);
481 f1.SetLineColor(kGreen);
482 f2.SetLineColor(kRed);
483
484 f2.DrawCopy("same");
485 f1.DrawCopy("same");
486
487 TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC");
488 leg->SetBorderSize(1);
489 leg->SetTextSize(0.04);
490 leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22);
491 //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22);
492 leg->AddLine(0, 0.65, 0, 0.65);
493 leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11);
494 leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11);
495 leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11);
496 leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11);
497 leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11);
498 if (polynom>1)
499 leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11);
500 leg->SetBit(kCanDelete);
501 leg->Draw();
502
503 q->cd(2);
504
505 TGraph *g = new TGraph;
506 g->SetPoint(0, 0, 0);
507
508 Int_t max=0;
509 Float_t maxsig=0;
510 for (int i=1; i<89; i++)
511 {
512 const Double_t s = f1.Integral(0, (float)i);
513 const Double_t b = f2.Integral(0, (float)i);
514
515 const Double_t sig = Significance(s, b);
516
517 g->SetPoint(g->GetN(), i, sig);
518
519 if (sig>maxsig)
520 {
521 max = i;
522 maxsig = sig;
523 }
524 }
525
526 g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha");
527 g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]");
528 g->GetHistogram()->SetYTitle("Significance");
529 g->SetBit(kCanDelete);
530 g->Draw("AP");
531
532 leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC");
533 leg->SetBorderSize(1);
534 leg->SetTextSize(0.1);
535 leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max));
536 leg->SetBit(kCanDelete);
537 leg->Draw();
538 }*/
539}
540
541// --------------------------------------------------------------------------
542//
543// You can use this function if you want to use a MHMatrix instead of
544// MMcEvt. This function adds all necessary columns to the
545// given matrix. Afterward you should fill the matrix with the corresponding
546// data (eg from a file by using MHMatrix::Fill). If you now loop
547// through the matrix (eg using MMatrixLoop) MHHadronness::Fill
548// will take the values from the matrix instead of the containers.
549//
550void MHAlpha::InitMapping(MHMatrix *mat)
551{
552 if (fMatrix)
553 return;
554
555 fMatrix = mat;
556
557 fMap = fMatrix->AddColumn("MHillasSrc.fAlpha");
558}
559
560void MHAlpha::StopMapping()
561{
562 fMatrix = NULL;
563}
Note: See TracBrowser for help on using the repository browser.