source: trunk/MagicSoft/Mars/mhflux/MHAlpha.cc@ 6951

Last change on this file since 6951 was 6932, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 25.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// 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// ToDo:
40// =====
41//
42// - Replace time-binning by histogram (which is enhanced bin-wise)
43//
44//
45//////////////////////////////////////////////////////////////////////////////
46#include "MHAlpha.h"
47
48#include <TF1.h>
49#include <TGraph.h>
50#include <TStyle.h>
51#include <TLatex.h>
52#include <TCanvas.h>
53#include <TPaveStats.h>
54#include <TStopwatch.h>
55
56#include "MGeomCam.h"
57#include "MSrcPosCam.h"
58#include "MHillasSrc.h"
59#include "MTime.h"
60#include "MObservatory.h"
61#include "MPointingPos.h"
62#include "MAstroSky2Local.h"
63#include "MStatusDisplay.h"
64#include "MParameters.h"
65#include "MHMatrix.h"
66
67#include "MBinning.h"
68#include "MParList.h"
69
70#include "MLog.h"
71#include "MLogManip.h"
72
73ClassImp(MHAlpha);
74ClassImp(MAlphaFitter);
75
76using namespace std;
77
78// --------------------------------------------------------------------------
79//
80// Default Constructor
81//
82MHAlpha::MHAlpha(const char *name, const char *title)
83 : fOffData(0), fResult(0), /*fExcess(0),*/ fEnergy(0), fHillas(0),
84 fPointPos(0), fTimeEffOn(0), fTime(0),
85 fSkipHistTime(kFALSE), fSkipHistTheta(kFALSE), fSkipHistEnergy(kFALSE),
86 fNumTimeBins(10), fMatrix(0)
87{
88 //
89 // set the name and title of this object
90 //
91 fName = name ? name : "MHAlpha";
92 fTitle = title ? title : "Alpha plot";
93
94 fHAlpha.SetName("Alpha");
95 fHAlpha.SetTitle("Alpha");
96 fHAlpha.SetXTitle("\\Theta [deg]");
97 fHAlpha.SetYTitle("E_{est} [GeV]");
98 fHAlpha.SetZTitle("|\\alpha| [\\circ]");
99 fHAlpha.SetDirectory(NULL);
100 fHAlpha.UseCurrentStyle();
101
102 // Main histogram
103 fHAlphaTime.SetName("Alpha");
104 fHAlphaTime.SetXTitle("|\\alpha| [\\circ]");
105 fHAlphaTime.SetYTitle("Counts");
106 fHAlphaTime.UseCurrentStyle();
107 fHAlphaTime.SetDirectory(NULL);
108
109
110 //fHEnergy.SetName("Energy");
111 //fHEnergy.SetTitle(" N_{exc} vs. E_{est} ");
112 //fHEnergy.SetXTitle("E_{est} [GeV]");
113 fHEnergy.SetYTitle("N_{exc}");
114 fHEnergy.SetDirectory(NULL);
115 fHEnergy.UseCurrentStyle();
116
117 fHTheta.SetName("ExcessTheta");
118 fHTheta.SetTitle(" N_{exc} vs. \\Theta ");
119 fHTheta.SetXTitle("\\Theta [\\circ]");
120 fHTheta.SetYTitle("N_{exc}");
121 fHTheta.SetDirectory(NULL);
122 fHTheta.UseCurrentStyle();
123
124 // effective on time versus time
125 fHTime.SetName("ExcessTime");
126 fHTime.SetTitle(" N_{exc} vs. Time ");
127 fHTime.SetXTitle("Time");
128 fHTime.SetYTitle("N_{exc} [s]");
129 fHTime.UseCurrentStyle();
130 fHTime.SetDirectory(NULL);
131 fHTime.GetYaxis()->SetTitleOffset(1.2);
132 fHTime.GetXaxis()->SetLabelSize(0.033);
133 fHTime.GetXaxis()->SetTimeFormat("%H:%M:%S %F1995-01-01 00:00:00 GMT");
134 fHTime.GetXaxis()->SetTimeDisplay(1);
135 fHTime.Sumw2();
136
137 MBinning binsa, binse, binst;
138 binsa.SetEdges(18, 0, 90);
139 binse.SetEdgesLog(15, 10, 100000);
140 binst.SetEdgesCos(50, 0, 60);
141 binse.Apply(fHEnergy);
142 binst.Apply(fHTheta);
143 binsa.Apply(fHAlphaTime);
144
145 MH::SetBinning(&fHAlpha, &binst, &binse, &binsa);
146}
147
148void MHAlpha::FitEnergySpec(Bool_t paint)
149{
150 /*
151 TF1 f("Spectrum", "pow(x, [0])*exp(-x/[1])*[2]*x");//
152 f.SetParameter(0, -2.1);
153 f.SetParameter(1, 1400); // 50
154 f.SetParameter(2, fHEnergy.Integral()); // 50
155 f.SetLineWidth(2);
156
157 fHEnergy.Fit(&f, "0QI", "", 100, 2400); // Integral? 90, 2800
158
159 const Double_t chi2 = f.GetChisquare();
160 const Int_t NDF = f.GetNDF();
161
162 // was fit successful ?
163 const Bool_t ok = NDF>0 && chi2<2.5*NDF;
164 f.SetLineColor(ok ? kGreen : kRed);
165
166 if (paint)
167 {
168 f.Paint("same");
169
170 if (ok)
171 {
172 TString s;
173 s += Form("\\alpha=%.1f ", f.GetParameter(0));
174 s += Form("E_{c}=%.1f TeV ", f.GetParameter(1)/1000);
175 s += Form("N=%.1e", f.GetParameter(2));
176
177 TLatex text(0.25, 0.94, s.Data());
178 text.SetBit(TLatex::kTextNDC);
179 text.SetTextSize(0.04);
180 text.Paint();
181 }
182 }*/
183}
184
185Float_t MHAlpha::FitEnergyBins(Bool_t paint)
186{
187 // Do not store the 'final' result in fFit
188 MAlphaFitter fit(fFit);
189
190 const Int_t n = fHAlpha.GetNbinsY();
191
192 fHEnergy.SetEntries(0);
193
194 Float_t mean = 0;
195 Int_t num = 0;
196 for (int i=1; i<=n; i++)
197 {
198 if (fit.FitEnergy(fHAlpha, fOffData, i))
199 {
200 fHEnergy.SetBinContent(i, fit.GetEventsExcess());
201 if (fit.GetSignificance()>1)
202 fHEnergy.SetBinError(i, fit.GetEventsExcess()/fit.GetSignificance());
203 else
204 fHEnergy.SetBinError(i, fit.GetEventsExcess());
205
206 if (fit.GetSignificance()>1)
207 {
208 mean += fit.GetSignificance();
209 num++;
210 }
211 }
212 }
213 return mean/num;
214}
215
216void MHAlpha::FitThetaBins(Bool_t paint)
217{
218 // Do not store the 'final' result in fFit
219 MAlphaFitter fit(fFit);
220
221 const Int_t n = fHAlpha.GetNbinsX();
222
223 fHTheta.SetEntries(0);
224
225 for (int i=1; i<=n; i++)
226 {
227 if (fit.FitTheta(fHAlpha, fOffData, i))
228 {
229 fHTheta.SetBinContent(i, fit.GetEventsExcess());
230 if (fit.GetSignificance()>1)
231 fHTheta.SetBinError(i, fit.GetEventsExcess()/fit.GetSignificance());
232 else
233 fHTheta.SetBinError(i, fit.GetEventsExcess());
234 }
235 }
236}
237
238Bool_t MHAlpha::SetupFill(const MParList *pl)
239{
240 fHAlpha.Reset();
241 fHEnergy.Reset();
242 fHTheta.Reset();
243 fHTime.Reset();
244
245 if (fName!=(TString)"MHAlphaOff" && fOffData==NULL)
246 {
247 MHAlpha *hoff = (MHAlpha*)pl->FindObject("MHAlphaOff", "MHAlpha");
248 if (!hoff)
249 *fLog << inf << "No MHAlphaOff [MHAlpha] found... using current data only!" << endl;
250 else
251 {
252 *fLog << inf << "MHAlphaOff [MHAlpha] found... using on-off mode!" << endl;
253 SetOffData(*hoff);
254 }
255 }
256
257 fHillas = 0;
258 /*
259 if (fSizeMin>=0 || fSizeMax>=0)
260 {
261 fHillas = (MHillas*)pl->FindObject("MHillas");
262 if (!fHillas)
263 {
264 *fLog << warn << "Size cut set, but MHillas not found... abort." << endl;
265 return kFALSE;
266 }
267 }
268 */
269 fEnergy = (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
270 if (!fEnergy)
271 { /*
272 if (fEnergyMin>=0 || fEnergyMax>=0)
273 {
274 *fLog << warn << "Energy cut set, but MEnergyEst not found... abort." << endl;
275 return kFALSE;
276 } */
277
278 *fLog << warn << "MEnergyEst [MParameterD] not found... " << flush;
279
280 if (!fHillas)
281 fHillas = (MHillas*)pl->FindObject("MHillas");
282 if (fHillas)
283 *fLog << "using SIZE instead!" << endl;
284 else
285 *fLog << "ignored." << endl;
286
287 fHEnergy.SetName("ExcessSize");
288 fHEnergy.SetTitle(" N_{exc} vs. Size ");
289 fHEnergy.SetXTitle("Size [\\gamma]");
290 }
291 else
292 {
293 fHEnergy.SetName("ExcessEnergy");
294 fHEnergy.SetTitle(" N_{exc} vs. E_{est} ");
295 fHEnergy.SetXTitle("E_{est} [GeV]");
296 }
297
298 fPointPos = (MPointingPos*)pl->FindObject("MPointingPos");
299 if (!fPointPos)
300 *fLog << warn << "MPointingPos not found... ignored." << endl;
301
302 fTimeEffOn = (MTime*)pl->FindObject("MTimeEffectiveOnTime", "MTime");
303 if (!fTimeEffOn)
304 *fLog << warn << "MTimeEffectiveOnTime [MTime] not found... ignored." << endl;
305 else
306 *fTimeEffOn = MTime(); // FIXME: How to do it different?
307
308 fTime = (MTime*)pl->FindObject("MTime");
309 if (!fTime)
310 *fLog << warn << "MTime not found... ignored." << endl;
311
312 fResult = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MinimizationValue");
313 if (!fResult)
314 return kFALSE;
315
316 //fExcess = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MExcess");
317 //if (!fExcess)
318 // return kFALSE;
319
320 fLastTime = MTime();
321 fNumRebin = fNumTimeBins;
322
323 MBinning binst, binse, binsa;
324 binst.SetEdges(fOffData ? *fOffData : fHAlpha, 'x');
325 binse.SetEdges(fOffData ? *fOffData : fHAlpha, 'y');
326 binsa.SetEdges(fOffData ? *fOffData : fHAlpha, 'z');
327
328 if (!fOffData)
329 {
330 if (!fPointPos)
331 binst.SetEdges(1, 0, 60);
332 else
333 binst.SetEdges(*pl, "BinningTheta");
334
335 if (!fEnergy && !fHillas)
336 binse.SetEdges(1, 10, 100000);
337 else
338 binse.SetEdges(*pl, "BinningEnergyEst");
339
340 binsa.SetEdges(*pl, "BinningAlpha");
341 }
342
343 binse.Apply(fHEnergy);
344 binst.Apply(fHTheta);
345 binsa.Apply(fHAlphaTime);
346 MH::SetBinning(&fHAlpha, &binst, &binse, &binsa);
347
348 MAlphaFitter *fit = (MAlphaFitter*)pl->FindObject("MAlphaFitter");
349 if (!fit)
350 *fLog << warn << "MAlphaFitter not found... ignored." << endl;
351 else
352 fit->Copy(fFit);
353
354 fFit.Print();
355
356 return kTRUE;
357}
358
359void MHAlpha::InitAlphaTime(const MTime &t)
360{
361 //
362 // If this is the first call we have to initialize the time-histogram
363 //
364 MBinning bins;
365 bins.SetEdges(1, t.GetAxisTime()-60, t.GetAxisTime());
366 bins.Apply(fHTime);
367
368 fLastTime=t;
369}
370
371void MHAlpha::UpdateAlphaTime(Bool_t final)
372{
373 if (!fTimeEffOn)
374 return;
375
376 if (!final)
377 {
378 if (fLastTime==MTime() && *fTimeEffOn!=MTime())
379 {
380 InitAlphaTime(*fTimeEffOn);
381 return;
382 }
383
384 if (fLastTime!=*fTimeEffOn)
385 {
386 fLastTime=*fTimeEffOn;
387 fNumRebin--;
388 }
389
390 if (fNumRebin>0)
391 return;
392 }
393
394 // Check new 'last time'
395 MTime *time = final ? fTime : fTimeEffOn;
396
397 if (time->GetAxisTime()<=fHTime.GetXaxis()->GetXmax())
398 {
399 *fLog << warn << "WARNING - New time-stamp " << *time << " lower" << endl;
400 *fLog << "than upper edge of histogram... skipped." << endl;
401 *fLog << "This should not happen. Maybe you started you eventloop with" << endl;
402 *fLog << "an already initialized time stamp MTimeEffectiveOnTime?" << endl;
403 fNumRebin++;
404 return;
405 }
406
407 // Fit time histogram
408 MAlphaFitter fit(fFit);
409
410 TH1D *h = fOffData ? fOffData->ProjectionZ("ProjTimeTemp", -1, 9999, -1, 9999, "E") : 0;
411 const Bool_t rc = fit.ScaleAndFit(fHAlphaTime, h);
412 if (h)
413 delete h;
414
415 if (!rc)
416 return;
417
418 // Reset Histogram
419 fHAlphaTime.Reset();
420 fHAlphaTime.SetEntries(0);
421
422 //
423 // Prepare Histogram
424 //
425 if (final)
426 time->Plus1ns();
427
428 // Get number of bins
429 const Int_t n = fHTime.GetNbinsX();
430
431 // Enhance binning
432 MBinning bins;
433 bins.SetEdges(fHTime, 'x');
434 bins.AddEdge(time->GetAxisTime());
435 bins.Apply(fHTime);
436
437 //
438 // Fill histogram
439 //
440 fHTime.SetBinContent(n+1, fit.GetEventsExcess());
441 if (fit.GetSignificance()>1)
442 fHTime.SetBinError(n+1, fit.GetEventsExcess()/fit.GetSignificance());
443 else
444 fHTime.SetBinError(n+1, fit.GetEventsExcess());
445
446 *fLog << all << *fTimeEffOn << ": " << fit.GetEventsExcess() << endl;
447
448 fNumRebin = fNumTimeBins;
449}
450
451// --------------------------------------------------------------------------
452//
453// Fill the histogram. For details see the code or the class description
454//
455Bool_t MHAlpha::Fill(const MParContainer *par, const Stat_t w)
456{
457 Double_t alpha, energy, theta;
458 Double_t size=-1;
459
460 if (fMatrix)
461 {
462 alpha = (*fMatrix)[fMap[0]];
463 energy = fMap[1]<0 ? -1 : (*fMatrix)[fMap[1]];
464 size = fMap[2]<0 ? -1 : (*fMatrix)[fMap[2]];
465 //<0 ? 1000 : (*fMatrix)[fMap[1]];
466 theta = 0;
467
468 if (energy<0)
469 energy=size;
470 if (size<0)
471 size=energy;
472
473 if (energy<0 && size<0)
474 energy = size = 1000;
475 }
476 else
477 {
478 const MHillasSrc *hil = dynamic_cast<const MHillasSrc*>(par);
479 if (!par)
480 {
481 *fLog << err << dbginf << "MHillasSrc not found... abort." << endl;
482 return kFALSE;
483 }
484
485 alpha = hil->GetAlpha();
486
487 if (fHillas)
488 size = fHillas->GetSize();
489 energy = fEnergy ? fEnergy->GetVal() : (fHillas?fHillas->GetSize():1000);
490 theta = fPointPos ? fPointPos->GetZd() : 0;
491 }
492
493 //if (size>0)
494 // alpha /= (2.4 + 1.13*(log10((energy-14)/0.37)-5)*(log10((energy-14)/0.37)-5))/15;
495
496 // enhance histogram if necessary
497 UpdateAlphaTime();
498
499 // Fill histograms
500 fHAlpha.Fill(theta, energy, TMath::Abs(alpha), w);
501
502 // Check cuts
503 /*
504 if ( (fEnergyMin>=0 && energy<fEnergyMin) ||
505 (fEnergyMax>=0 && energy>fEnergyMax) ||
506 (fSizeMin>=0 && size <fSizeMin) ||
507 (fSizeMax>=0 && size >fSizeMin) )
508 return kTRUE;
509 */
510
511 if (!fSkipHistTime)
512 fHAlphaTime.Fill(TMath::Abs(alpha), w);
513
514 return kTRUE;
515}
516
517// --------------------------------------------------------------------------
518//
519// Paint the integral and the error on top of the histogram
520//
521void MHAlpha::PaintText(Double_t val, Double_t error) const
522{
523 TLatex text(0.45, 0.94, Form("N_{exc} = %.1f \\pm %.1f", val, error));
524 text.SetBit(TLatex::kTextNDC);
525 text.SetTextSize(0.04);
526 text.Paint();
527}
528
529// --------------------------------------------------------------------------
530//
531// Update the projections
532//
533void MHAlpha::Paint(Option_t *opt)
534{
535 // Note: Any object cannot be added to a pad twice!
536 // The object is searched by gROOT->FindObject only in
537 // gPad itself!
538 TVirtualPad *padsave = gPad;
539
540 TH1D *h0=0;
541
542 TString o(opt);
543 if (o==(TString)"proj")
544 {
545 TPaveStats *st=0;
546 for (int x=0; x<4; x++)
547 {
548 TVirtualPad *p=padsave->GetPad(x+1);
549 if (!p || !(st = (TPaveStats*)p->GetPrimitive("stats")))
550 continue;
551
552 if (st->GetOptStat()==11)
553 continue;
554
555 const Double_t y1 = st->GetY1NDC();
556 const Double_t y2 = st->GetY2NDC();
557 const Double_t x1 = st->GetX1NDC();
558 const Double_t x2 = st->GetX2NDC();
559
560 st->SetY1NDC((y2-y1)/3+y1);
561 st->SetX1NDC((x2-x1)/3+x1);
562 st->SetOptStat(11);
563 }
564
565 padsave->cd(1);
566 TH1D *hon = fHAlpha.ProjectionZ("ProjAlpha");
567 if (fOffData)
568 {
569 TH1D *hoff = fOffData->ProjectionZ("ProjAlphaOff");
570 const Double_t alpha = fFit.Scale(*hoff, *hon);
571
572 hon->SetMaximum();
573 hon->SetMaximum(TMath::Max(hon->GetMaximum(), hoff->GetMaximum())*1.05);
574
575 // BE CARFEULL: This is a really weird workaround!
576 hoff->SetMaximum(alpha);
577
578 if ((h0=(TH1D*)gPad->FindObject("ProjAlphaOnOff")))
579 {
580 h0->Reset();
581 h0->Add(hoff, hon, -1);
582 const Float_t min = h0->GetMinimum()*1.05;
583 hon->SetMinimum(min<0 ? min : 0);
584 }
585 }
586 else
587 hon->SetMinimum(0);
588 FitEnergyBins();
589 FitThetaBins();
590 }
591
592 if (o==(TString)"alpha")
593 if ((h0 = (TH1D*)gPad->FindObject("ProjAlpha")))
594 {
595 // Check whether we are dealing with on-off analysis
596 TH1D *hoff = (TH1D*)gPad->FindObject("ProjAlphaOff");
597
598 // BE CARFEULL: This is a really weird workaround!
599 const Double_t scale = !hoff || hoff->GetMaximum()<0 ? 1 : hoff->GetMaximum();
600
601 // Do not store the 'final' result in fFit
602 MAlphaFitter fit(fFit);
603 fit.Fit(*h0, hoff, scale, kTRUE);
604 fit.PaintResult();
605 }
606
607 if (o==(TString)"time")
608 PaintText(fHTime.Integral(), 0);
609
610 if (o==(TString)"theta")
611 {
612 TH1 *h = (TH1*)gPad->FindObject("Alpha_x");
613 if (h)
614 {
615 TH1D *h2 = (TH1D*)fHAlpha.Project3D("dum_x");
616 h2->SetDirectory(0);
617 h2->Scale(fHTheta.Integral()/h2->Integral());
618 h->Reset();
619 h->Add(h2);
620 delete h2;
621 }
622 PaintText(fHTheta.Integral(), 0);
623 }
624
625 if (o==(TString)"energy")
626 {
627 TH1 *h = (TH1*)gPad->FindObject("Alpha_y");
628 if (h)
629 {
630 TH1D *h2 = (TH1D*)fHAlpha.Project3D("dum_y");
631 h2->SetDirectory(0);
632 h2->Scale(fHEnergy.Integral()/h2->Integral());
633 h->Reset();
634 h->Add(h2);
635 delete h2;
636 }
637
638 if (fHEnergy.GetMaximum()>1)
639 {
640 gPad->SetLogx();
641 gPad->SetLogy();
642 }
643 FitEnergySpec(kTRUE);
644 }
645
646 gPad = padsave;
647}
648
649// --------------------------------------------------------------------------
650//
651// Draw the histogram
652//
653void MHAlpha::Draw(Option_t *opt)
654{
655 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
656
657 // Do the projection before painting the histograms into
658 // the individual pads
659 AppendPad("proj");
660
661 pad->SetBorderMode(0);
662 pad->Divide(2,2);
663
664 TH1D *h=0;
665
666 pad->cd(1);
667 gPad->SetBorderMode(0);
668
669 h = fHAlpha.ProjectionZ("ProjAlpha", -1, 9999, -1, 9999, "E");
670 h->SetBit(TH1::kNoTitle);
671 h->SetStats(kTRUE);
672 h->SetXTitle("\\alpha [\\circ]");
673 h->SetYTitle("Counts");
674 h->SetDirectory(NULL);
675 h->SetMarkerStyle(kFullDotMedium);
676 h->SetBit(kCanDelete);
677 h->Draw("");
678
679 if (fOffData)
680 {
681 h->SetMarkerColor(kGreen);
682
683 h = fOffData->ProjectionZ("ProjAlphaOff", -1, 9999, -1, 9999, "E");
684 h->SetBit(TH1::kNoTitle);
685 h->SetXTitle("\\alpha [\\circ]");
686 h->SetYTitle("Counts");
687 h->SetDirectory(NULL);
688 h->SetMarkerStyle(kFullDotMedium);
689 h->SetBit(kCanDelete);
690 h->SetMarkerColor(kRed);
691 h->Draw("same");
692
693 h = (TH1D*)h->Clone("ProjAlphaOnOff");
694 h->SetBit(TH1::kNoTitle);
695 h->SetXTitle("\\alpha [\\circ]");
696 h->SetYTitle("Counts");
697 h->SetDirectory(NULL);
698 h->SetMarkerStyle(kFullDotMedium);
699 h->SetBit(kCanDelete);
700 h->SetMarkerColor(kBlue);
701 h->Draw("same");
702
703 TLine lin;
704 lin.SetLineStyle(kDashed);
705 lin.DrawLine(0, 0, 90, 0);
706 }
707
708 // After the Alpha-projection has been drawn. Fit the histogram
709 // and paint the result into this pad
710 AppendPad("alpha");
711
712 if (fHEnergy.GetNbinsX()>1)
713 {
714 pad->cd(2);
715 gPad->SetBorderMode(0);
716 gPad->SetGridx();
717 gPad->SetGridy();
718 fHEnergy.Draw();
719
720 AppendPad("energy");
721
722 h = (TH1D*)fHAlpha.Project3D("y");
723 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
724 h->SetXTitle("E [GeV]");
725 h->SetYTitle("Counts");
726 h->SetDirectory(NULL);
727 h->SetMarkerStyle(kFullDotMedium);
728 h->SetBit(kCanDelete);
729 h->SetMarkerColor(kCyan);
730 h->SetLineColor(kCyan);
731 h->Draw("Psame");
732 }
733 else
734 delete pad->GetPad(2);
735
736 if (fTimeEffOn && fTime || fHTime.GetNbinsX()>1)
737 {
738 pad->cd(3);
739 gPad->SetBorderMode(0);
740 gPad->SetGridx();
741 gPad->SetGridy();
742 fHTime.Draw();
743 AppendPad("time");
744 }
745 else
746 delete pad->GetPad(3);
747
748 if (fHTheta.GetNbinsX()>1)
749 {
750 pad->cd(4);
751 gPad->SetGridx();
752 gPad->SetGridy();
753 gPad->SetBorderMode(0);
754 fHTheta.Draw();
755
756 AppendPad("theta");
757
758 h = (TH1D*)fHAlpha.Project3D("x");
759 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
760 h->SetXTitle("\\theta [\\circ]");
761 h->SetYTitle("Counts");
762 h->SetDirectory(NULL);
763 h->SetMarkerStyle(kFullDotMedium);
764 h->SetBit(kCanDelete);
765 h->SetMarkerColor(kCyan);
766 h->SetLineColor(kCyan);
767 h->Draw("Psame");
768 }
769 else
770 delete pad->GetPad(4);
771}
772
773void MHAlpha::DrawAll()
774{
775 // FIXME: Do in Paint if special option given!
776 TCanvas *c = new TCanvas;
777 Int_t n = fHAlpha.GetNbinsY();
778 Int_t nc = (Int_t)(TMath::Sqrt((Float_t)n-1)+1);
779 c->Divide(nc, nc, 1e-10, 1e-10);
780
781 // Do not store the 'final' result in fFit
782 MAlphaFitter fit(fFit);
783
784 for (int i=1; i<=fHAlpha.GetNbinsY(); i++)
785 {
786 c->cd(i);
787
788 TH1D *hon = fHAlpha.ProjectionZ("ProjAlpha", -1, 9999, i, i, "E");
789 hon->SetBit(TH1::kNoTitle);
790 hon->SetStats(kTRUE);
791 hon->SetXTitle("\\alpha [\\circ]");
792 hon->SetYTitle("Counts");
793 hon->SetDirectory(NULL);
794 hon->SetMarkerStyle(kFullDotMedium);
795 hon->SetBit(kCanDelete);
796 hon->Draw("");
797
798 TH1D *hof = 0;
799 Double_t alpha = 1;
800
801 if (fOffData)
802 {
803 hon->SetMarkerColor(kGreen);
804
805 hof = fOffData->ProjectionZ("ProjAlphaOff", -1, 9999, i, i, "E");
806 hof->SetBit(TH1::kNoTitle|TH1::kNoStats);
807 hof->SetXTitle("\\alpha [\\circ]");
808 hof->SetYTitle("Counts");
809 hof->SetDirectory(NULL);
810 hof->SetMarkerStyle(kFullDotMedium);
811 hof->SetBit(kCanDelete);
812 hof->SetMarkerColor(kRed);
813 hof->Draw("same");
814
815 alpha = fit.Scale(*hof, *hon);
816
817 hon->SetMaximum();
818 hon->SetMaximum(TMath::Max(hon->GetMaximum(), hof->GetMaximum())*1.05);
819
820 TH1D *diff = new TH1D(*hon);
821 diff->Add(hof, -1);
822 diff->SetBit(TH1::kNoTitle);
823 diff->SetXTitle("\\alpha [\\circ]");
824 diff->SetYTitle("Counts");
825 diff->SetDirectory(NULL);
826 diff->SetMarkerStyle(kFullDotMedium);
827 diff->SetBit(kCanDelete);
828 diff->SetMarkerColor(kBlue);
829 diff->Draw("same");
830
831 TLine lin;
832 lin.SetLineStyle(kDashed);
833 lin.DrawLine(0, 0, 90, 0);
834
835 const Float_t min = diff->GetMinimum()*1.05;
836 hon->SetMinimum(min<0 ? min : 0);
837 }
838
839 if (hof ? fit.Fit(*hon, *hof, alpha) : fit.Fit(*hon))
840 *fLog << dbg << "Bin " << i << ": sigma=" << fit.GetSignificance() << " omega=" << fit.GetGausSigma() << " events=" << fit.GetEventsExcess() << endl;
841 /*
842 if (fit.FitEnergy(fHAlpha, fOffData, i, kTRUE))
843 {
844 fHEnergy.SetBinContent(i, fit.GetEventsExcess());
845 fHEnergy.SetBinError(i, fit.GetEventsExcess()*0.2);
846 }*/
847 }
848
849}
850
851
852Bool_t MHAlpha::Finalize()
853{
854 //TH1D *h = fHAlpha.ProjectionZ("AlphaExc_px", -1, 9999, -1, 9999, "E");
855 //h->SetDirectory(0);
856 //Bool_t rc = fFit.Fit(*h);
857 //delete h;
858
859 if (!fFit.FitAlpha(fHAlpha, fOffData))
860 {
861 *fLog << warn << "MAlphaFitter - Fit failed..." << endl;
862 return kTRUE;
863 }
864
865 // Store the final result in fFit
866 fFit.Print("result");
867
868 fResult->SetVal(fFit.GetMinimizationValue());
869
870 if (!fSkipHistEnergy)
871 {
872 *fLog << inf << "Processing energy bins..." << endl;
873 FitEnergyBins();
874 }
875 if (!fSkipHistTheta)
876 {
877 *fLog << inf << "Processing theta bins..." << endl;
878 FitThetaBins();
879 }
880 if (!fSkipHistTime)
881 {
882 *fLog << inf << "Processing time bins..." << endl;
883 UpdateAlphaTime(kTRUE);
884 MH::RemoveFirstBin(fHTime);
885 }
886
887 return kTRUE;
888}
889
890// --------------------------------------------------------------------------
891//
892// You can use this function if you want to use a MHMatrix instead of
893// MMcEvt. This function adds all necessary columns to the
894// given matrix. Afterward you should fill the matrix with the corresponding
895// data (eg from a file by using MHMatrix::Fill). If you now loop
896// through the matrix (eg using MMatrixLoop) MHHadronness::Fill
897// will take the values from the matrix instead of the containers.
898//
899// It takes fSkipHist* into account!
900//
901void MHAlpha::InitMapping(MHMatrix *mat, Int_t type)
902{
903 if (fMatrix)
904 return;
905
906 fMatrix = mat;
907
908 fMap[0] = fMatrix->AddColumn("MHillasSrc.fAlpha");
909 fMap[1] = -1;
910 fMap[2] = -1;
911 fMap[3] = -1;
912 fMap[4] = -1;
913
914 if (!fSkipHistEnergy)
915 if (type==0)
916 {
917 fMap[1] = fMatrix->AddColumn("MEnergyEst.fVal");
918 fMap[2] = -1;
919 }
920 else
921 {
922 fMap[1] = -1;
923 fMap[2] = fMatrix->AddColumn("MHillas.fSize");
924 }
925
926 if (!fSkipHistTheta)
927 fMap[3] = fMatrix->AddColumn("MPointingPos.fZd");
928
929 // if (!fSkipHistTime)
930 // fMap[4] = fMatrix->AddColumn("MTime.GetAxisTime");
931}
932
933void MHAlpha::StopMapping()
934{
935 fMatrix = NULL;
936}
937
938Int_t MHAlpha::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
939{
940 Bool_t rc = kFALSE;
941 if (IsEnvDefined(env, prefix, "NumTimeBins", print))
942 {
943 SetNumTimeBins(GetEnvValue(env, prefix, "NumTimeBins", fNumTimeBins));
944 rc = kTRUE;
945 }
946 return rc;
947}
Note: See TracBrowser for help on using the repository browser.