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

Last change on this file since 8785 was 8780, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 30.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-2005
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 "MSrcPosCam.h"
57#include "MHillas.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);
74
75using namespace std;
76
77// --------------------------------------------------------------------------
78//
79// Default Constructor
80//
81MHAlpha::MHAlpha(const char *name, const char *title)
82 : fNameParameter("MHillasSrc"), fParameter(0),
83 fOffData(0), fResult(0), fEnergy(0), fBin(0),
84 fPointPos(0), fTimeEffOn(0), fTime(0), fNumTimeBins(10),
85 fHillas(0), fMatrix(0), fSkipHistTime(kFALSE), fSkipHistTheta(kFALSE),
86 fSkipHistEnergy(kFALSE), fForceUsingSize(kFALSE)
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 fHist.SetName("Alpha");
95 fHist.SetTitle("Alpha");
96 fHist.SetXTitle("\\Theta [deg]");
97 //fHist.SetYTitle("E_{est} [GeV]");
98 fHist.SetZTitle("|\\alpha| [\\circ]");
99 fHist.SetDirectory(NULL);
100 fHist.UseCurrentStyle();
101
102 // Main histogram
103 fHistTime.SetName("Alpha");
104 fHistTime.SetXTitle("|\\alpha| [\\circ]");
105 fHistTime.SetYTitle("Counts");
106 fHistTime.UseCurrentStyle();
107 fHistTime.SetDirectory(NULL);
108
109
110 fHEnergy.SetName("Excess");
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.SetEdgesASin(67, -0.005, 0.665);
141 binse.Apply(fHEnergy);
142 binst.Apply(fHTheta);
143 binsa.Apply(fHistTime);
144
145 MH::SetBinning(&fHist, &binst, &binse, &binsa);
146}
147
148Float_t MHAlpha::FitEnergyBins(Bool_t paint)
149{
150 // Do not store the 'final' result in fFit
151 MAlphaFitter fit(fFit);
152
153 const Int_t n = fHist.GetNbinsY();
154
155 fHEnergy.SetEntries(0);
156
157 Float_t mean = 0;
158 Int_t num = 0;
159 for (int i=1; i<=n; i++)
160 {
161 if (!fit.FitEnergy(fHist, fOffData, i))
162 continue;
163
164 if (fit.GetSignificanceExc()<=0)
165 continue;
166
167 fHEnergy.SetBinContent(i, fit.GetEventsExcess());
168 fHEnergy.SetBinError(i, fit.GetEventsExcess()/fit.GetSignificanceExc());
169
170 mean += fit.GetSignificanceExc()*fit.GetSignificanceExc();
171 num++;
172 }
173 return TMath::Sqrt(mean)/num;
174}
175
176void MHAlpha::FitThetaBins(Bool_t paint)
177{
178 // Do not store the 'final' result in fFit
179 MAlphaFitter fit(fFit);
180
181 const Int_t n = fHist.GetNbinsX();
182
183 fHTheta.SetEntries(0);
184
185 for (int i=1; i<=n; i++)
186 {
187 if (!fit.FitTheta(fHist, fOffData, i))
188 continue;
189
190 if (fit.GetSignificanceExc()<=0)
191 continue;
192
193 fHTheta.SetBinContent(i, fit.GetEventsExcess());
194 fHTheta.SetBinError(i, fit.GetEventsExcess()/fit.GetSignificanceExc());
195 }
196}
197
198// --------------------------------------------------------------------------
199//
200// Return the value from fParemeter which should be filled into the plots
201//
202Double_t MHAlpha::GetVal() const
203{
204 return static_cast<const MHillasSrc*>(fParameter)->GetAlpha();
205}
206
207
208// --------------------------------------------------------------------------
209//
210// Store the pointer to the parameter container storing the plotted value
211// (here MHillasSrc) in fParameter.
212//
213// return whether it was found or not.
214//
215Bool_t MHAlpha::GetParameter(const MParList &pl)
216{
217 fParameter = (MParContainer*)pl.FindObject(fNameParameter, "MHillasSrc");
218 if (fParameter)
219 return kTRUE;
220
221 *fLog << err << fNameParameter << " [MHillasSrc] not found... abort." << endl;
222 return kFALSE;
223}
224
225Bool_t MHAlpha::SetupFill(const MParList *pl)
226{
227 fHist.Reset();
228 fHEnergy.Reset();
229 fHTheta.Reset();
230 fHTime.Reset();
231
232 const TString off(Form("%sOff", fName.Data()));
233 if (fName!=off && fOffData==NULL)
234 {
235 const TString desc(Form("%s [%s] found... using ", off.Data(), ClassName()));
236 MHAlpha *hoff = (MHAlpha*)pl->FindObject(off, ClassName());
237 if (!hoff)
238 *fLog << inf << "No " << desc << "current data only!" << endl;
239 else
240 {
241 *fLog << inf << desc << "on-off mode!" << endl;
242 SetOffData(*hoff);
243 }
244 }
245
246 if (!GetParameter(*pl))
247 return kFALSE;
248
249 fHillas = 0;
250 fEnergy = fForceUsingSize ? 0 : (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
251 if (!fEnergy)
252 {
253 *fLog << warn << "MEnergyEst [MParameterD] not found... " << flush;
254
255 if (!fHillas)
256 fHillas = (MHillas*)pl->FindObject("MHillas");
257 if (fHillas)
258 *fLog << "using SIZE instead!" << endl;
259 else
260 *fLog << "ignored." << endl;
261
262 fHEnergy.SetTitle(" N_{exc} vs. Size ");
263 fHEnergy.SetXTitle("Size [phe]");
264 fHist.SetYTitle("Size [phe]");
265 }
266 else
267 {
268 fHEnergy.SetTitle(" N_{exc} vs. E_{est} ");
269 fHEnergy.SetXTitle("E_{est} [GeV]");
270 fHist.SetYTitle("E_{est} [GeV]");
271 }
272
273 fPointPos = (MPointingPos*)pl->FindObject("MPointingPos");
274 if (!fPointPos)
275 *fLog << warn << "MPointingPos not found... ignored." << endl;
276
277 fTimeEffOn = (MTime*)pl->FindObject("MTimeEffectiveOnTime", "MTime");
278 if (!fTimeEffOn)
279 *fLog << warn << "MTimeEffectiveOnTime [MTime] not found... ignored." << endl;
280 else
281 *fTimeEffOn = MTime(); // FIXME: How to do it different?
282
283 fTime = (MTime*)pl->FindObject("MTime");
284 if (!fTime)
285 *fLog << warn << "MTime not found... ignored." << endl;
286
287 fResult = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MinimizationValue");
288 if (!fResult)
289 return kFALSE;
290 fSigma = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "GaussSigma");
291 if (!fSigma)
292 return kFALSE;
293 fBin = (MParameterI*)const_cast<MParList*>(pl)->FindCreateObj("MParameterI", "Bin");
294 if (!fBin)
295 return kFALSE;
296
297 //fExcess = (MParameterD*)const_cast<MParList*>(pl)->FindCreateObj("MParameterD", "MExcess");
298 //if (!fExcess)
299 // return kFALSE;
300
301 fLastTime = MTime();
302 fNumRebin = fNumTimeBins;
303
304 MBinning binst, binse, binsa;
305 binst.SetEdges(fOffData ? *fOffData : fHist, 'x');
306 binse.SetEdges(fOffData ? *fOffData : fHist, 'y');
307 binsa.SetEdges(fOffData ? *fOffData : fHist, 'z');
308 if (!fOffData)
309 {
310 if (!fPointPos)
311 binst.SetEdges(1, 0, 60);
312 else
313 binst.SetEdges(*pl, "BinningTheta");
314
315 if (!fEnergy && !fHillas)
316 binse.SetEdges(1, 10, 100000);
317 else
318 if (fEnergy)
319 binse.SetEdges(*pl, "BinningEnergyEst");
320 else
321 binse.SetEdges(*pl, "BinningSize");
322
323 binsa.SetEdges(*pl, Form("Binning%s", ClassName()+2));
324 }
325 else
326 {
327 fHEnergy.SetTitle(fOffData->GetTitle());
328 fHEnergy.SetXTitle(fOffData->GetYaxis()->GetTitle());
329 fHist.SetYTitle(fOffData->GetYaxis()->GetTitle());
330 }
331
332 binse.Apply(fHEnergy);
333 binst.Apply(fHTheta);
334 binsa.Apply(fHistTime);
335 MH::SetBinning(&fHist, &binst, &binse, &binsa);
336
337 MAlphaFitter *fit = (MAlphaFitter*)pl->FindObject("MAlphaFitter");
338 if (!fit)
339 *fLog << warn << "MAlphaFitter not found... ignored." << endl;
340 else
341 fit->Copy(fFit);
342
343 *fLog << inf;
344 fFit.Print();
345
346 return kTRUE;
347}
348
349void MHAlpha::InitAlphaTime(const MTime &t)
350{
351 //
352 // If this is the first call we have to initialize the time-histogram
353 //
354 MBinning bins;
355 bins.SetEdges(1, t.GetAxisTime()-60, t.GetAxisTime());
356 bins.Apply(fHTime);
357
358 fLastTime=t;
359}
360
361void MHAlpha::UpdateAlphaTime(Bool_t final)
362{
363 if (!fTimeEffOn)
364 return;
365
366 if (!final)
367 {
368 if (fLastTime==MTime() && *fTimeEffOn!=MTime())
369 {
370 InitAlphaTime(*fTimeEffOn);
371 return;
372 }
373
374 if (fLastTime!=*fTimeEffOn)
375 {
376 fLastTime=*fTimeEffOn;
377 fNumRebin--;
378 }
379
380 if (fNumRebin>0)
381 return;
382 }
383
384 // Check new 'last time'
385 MTime *time = final ? fTime : fTimeEffOn;
386
387 if (time->GetAxisTime()<=fHTime.GetXaxis()->GetXmax())
388 {
389 *fLog << warn << "WARNING - New time-stamp " << *time << " lower" << endl;
390 *fLog << "than upper edge of histogram... skipped." << endl;
391 *fLog << "This should not happen. Maybe you started you eventloop with" << endl;
392 *fLog << "an already initialized time stamp MTimeEffectiveOnTime?" << endl;
393 fNumRebin++;
394 return;
395 }
396
397 // Fit time histogram
398 MAlphaFitter fit(fFit);
399
400 TH1D *h = fOffData ? fOffData->ProjectionZ("ProjTimeTemp", -1, -1, -1, -1, "E") : 0;
401 const Bool_t rc = fit.ScaleAndFit(fHistTime, h);
402 if (h)
403 delete h;
404
405 if (!rc)
406 return;
407
408 // Reset Histogram
409 fHistTime.Reset();
410 fHistTime.SetEntries(0);
411
412 //
413 // Prepare Histogram
414 //
415 if (final)
416 time->Plus1ns();
417
418 // Get number of bins
419 const Int_t n = fHTime.GetNbinsX();
420
421 // Enhance binning
422 MBinning bins;
423 bins.SetEdges(fHTime, 'x');
424 bins.AddEdge(time->GetAxisTime());
425 bins.Apply(fHTime);
426
427 //
428 // Fill histogram
429 //
430 if (fit.GetSignificanceExc()>0)
431 {
432 fHTime.SetBinContent(n+1, fit.GetEventsExcess());
433 fHTime.SetBinError(n+1, fit.GetEventsExcess()/fit.GetSignificanceExc());
434 }
435
436 *fLog << all << *fTimeEffOn << ": " << fit.GetEventsExcess() << endl;
437
438 fNumRebin = fNumTimeBins;
439}
440
441void MHAlpha::SetBin(Int_t ibin)
442{
443 // Is this necessary?
444 // Could be speed up up searching for it only once.
445 const Float_t max = fFit.GetSignalIntegralMax();
446 const Int_t bin0 = fHist.GetZaxis()->FindFixBin(max);
447
448 const Int_t nbinsx = fHist.GetNbinsX();
449 const Int_t nbinsy = fHist.GetNbinsY();
450 const Int_t nxy = (nbinsx+2)*(nbinsy+2);
451
452 const Int_t binz = ibin/nxy;
453
454 const Bool_t issignal = binz>0 && binz<bin0;
455
456 fBin->SetVal(issignal ? binz : -binz);
457}
458
459// --------------------------------------------------------------------------
460//
461// Fill the histogram. For details see the code or the class description
462//
463Bool_t MHAlpha::Fill(const MParContainer *par, const Stat_t w)
464{
465 Double_t alpha, energy, theta;
466 Double_t size=-1;
467
468 if (fMatrix)
469 {
470 alpha = fMap[0]<0 ? GetVal() : (*fMatrix)[fMap[0]];
471 energy = fMap[1]<0 ? -1 : (*fMatrix)[fMap[1]];
472 size = fMap[2]<0 ? -1 : (*fMatrix)[fMap[2]];
473 //<0 ? 1000 : (*fMatrix)[fMap[1]];
474 theta = 0;
475
476 if (energy<0)
477 energy=size;
478 if (size<0)
479 size=energy;
480
481 if (energy<0 && size<0)
482 energy = size = 1000;
483 }
484 else
485 {
486 alpha = GetVal();
487
488 if (fHillas)
489 size = fHillas->GetSize();
490 energy = fEnergy ? fEnergy->GetVal() : (fHillas?fHillas->GetSize():1000);
491 theta = fPointPos ? fPointPos->GetZd() : 0;
492 }
493
494 // enhance histogram if necessary
495 UpdateAlphaTime();
496
497 // Fill histograms
498 const Int_t ibin = fHist.Fill(theta, energy, TMath::Abs(alpha), w);
499 SetBin(ibin);
500
501 if (!fSkipHistTime)
502 fHistTime.Fill(TMath::Abs(alpha), w);
503
504 return kTRUE;
505}
506
507// --------------------------------------------------------------------------
508//
509// Paint the integral and the error on top of the histogram
510//
511void MHAlpha::PaintText(Double_t val, Double_t error) const
512{
513 TLatex text(0.45, 0.94, Form("N_{exc} = %.1f \\pm %.1f", val, error));
514 text.SetBit(TLatex::kTextNDC);
515 text.SetTextSize(0.04);
516 text.Paint();
517}
518
519// --------------------------------------------------------------------------
520//
521// Paint the integral and the error on top of the histogram
522//
523void MHAlpha::PaintText(const TH1D &h) const
524{
525 Double_t sumv = 0;
526 Double_t sume = 0;
527
528 for (int i=0; i<h.GetNbinsX(); i++)
529 {
530 sumv += h.GetBinContent(i+1);
531 sume += h.GetBinError(i+1);
532 }
533 PaintText(sumv, sume);
534}
535// --------------------------------------------------------------------------
536//
537// Update the projections
538//
539void MHAlpha::Paint(Option_t *opt)
540{
541 // Note: Any object cannot be added to a pad twice!
542 // The object is searched by gROOT->FindObject only in
543 // gPad itself!
544 TVirtualPad *padsave = gPad;
545
546 TH1D *h0=0;
547
548 TString o(opt);
549 if (o==(TString)"proj")
550 {
551 TPaveStats *st=0;
552 for (int x=0; x<4; x++)
553 {
554 TVirtualPad *p=padsave->GetPad(x+1);
555 if (!p || !(st = (TPaveStats*)p->GetPrimitive("stats")))
556 continue;
557
558 if (st->GetOptStat()==11)
559 continue;
560
561 const Double_t y1 = st->GetY1NDC();
562 const Double_t y2 = st->GetY2NDC();
563 const Double_t x1 = st->GetX1NDC();
564 const Double_t x2 = st->GetX2NDC();
565
566 st->SetY1NDC((y2-y1)/3+y1);
567 st->SetX1NDC((x2-x1)/3+x1);
568 st->SetOptStat(11);
569 }
570
571 padsave->cd(1);
572 TH1D *hon = fHist.ProjectionZ("Proj");
573 if (fOffData)
574 {
575 TH1D *hoff = fOffData->ProjectionZ("ProjOff");
576 const Double_t alpha = fFit.Scale(*hoff, *hon);
577
578 hon->SetMaximum();
579 hon->SetMaximum(TMath::Max(hon->GetMaximum(), hoff->GetMaximum())*1.05);
580
581 // BE CARFEULL: This is a really weird workaround!
582 hoff->SetMaximum(alpha);
583
584 // For some reason the line-color is resetted
585 hoff->SetLineColor(kRed);
586
587 if ((h0=(TH1D*)gPad->FindObject("ProjOnOff")))
588 {
589 h0->Reset();
590 h0->Add(hoff, hon, -1);
591 const Float_t min = h0->GetMinimum()*1.05;
592 hon->SetMinimum(min<0 ? min : 0);
593 }
594 }
595 else
596 hon->SetMinimum(0);
597 FitEnergyBins();
598 FitThetaBins();
599 }
600
601 if (o==(TString)"variable")
602 if ((h0 = (TH1D*)gPad->FindObject("Proj")))
603 {
604 // Check whether we are dealing with on-off analysis
605 TH1D *hoff = (TH1D*)gPad->FindObject("ProjOff");
606
607 // BE CARFEULL: This is a really weird workaround!
608 const Double_t scale = !hoff || hoff->GetMaximum()<0 ? 1 : hoff->GetMaximum();
609
610 // Do not store the 'final' result in fFit
611 MAlphaFitter fit(fFit);
612 fit.Fit(*h0, hoff, scale, kTRUE);
613 fit.PaintResult();
614 }
615
616 if (o==(TString)"time")
617 PaintText(fHTime);//.Integral(), 0);
618
619 if (o==(TString)"theta")
620 {
621 TH1 *h = (TH1*)gPad->FindObject(Form("%s_x", fHist.GetName()));
622 if (h)
623 {
624 TH1D *h2 = (TH1D*)fHist.Project3D("dum_x");
625 h2->SetDirectory(0);
626 h2->Scale(fHTheta.Integral()/h2->Integral());
627 h->Reset();
628 h->Add(h2);
629 delete h2;
630 }
631 PaintText(fHTheta);//.Integral(), 0);
632 }
633
634 if (o==(TString)"energy")
635 {
636 TH1 *h = (TH1*)gPad->FindObject(Form("%s_y", fHist.GetName()));
637 if (h)
638 {
639 TH1D *h2 = (TH1D*)fHist.Project3D("dum_y");
640 h2->SetDirectory(0);
641 h2->Scale(fHEnergy.Integral()/h2->Integral());
642 h->Reset();
643 h->Add(h2);
644 delete h2;
645 }
646 PaintText(fHEnergy);//.Integral(), 0);
647
648 if (fHEnergy.GetMaximum()>1)
649 {
650 gPad->SetLogx();
651 gPad->SetLogy();
652 }
653 }
654
655 gPad = padsave;
656}
657
658// --------------------------------------------------------------------------
659//
660// Draw the histogram
661//
662void MHAlpha::Draw(Option_t *opt)
663{
664 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
665
666 /*
667 if (TString(opt).Contains("sizebins", TString::kIgnoreCase))
668 {
669 AppendPad("sizebins");
670 return;
671 }
672 */
673
674 // Do the projection before painting the histograms into
675 // the individual pads
676 AppendPad("proj");
677
678 pad->SetBorderMode(0);
679 pad->Divide(2,2);
680
681 TH1D *h=0;
682
683 pad->cd(1);
684 gPad->SetBorderMode(0);
685
686 h = fHist.ProjectionZ("Proj", -1, 9999, -1, 9999, "E");
687 h->SetBit(TH1::kNoTitle);
688 h->SetStats(kTRUE);
689 h->SetXTitle(fHist.GetZaxis()->GetTitle());
690 h->SetYTitle("Counts");
691 h->SetDirectory(NULL);
692 h->SetMarkerStyle(kFullDotMedium);
693 h->SetBit(kCanDelete);
694 h->Draw("");
695
696 if (fOffData)
697 {
698 // To get green on-data
699 //h->SetMarkerColor(kGreen);
700 //h->SetLineColor(kGreen);
701
702 h = fOffData->ProjectionZ("ProjOff", -1, 9999, -1, 9999, "E");
703 h->SetBit(TH1::kNoTitle);
704 h->SetXTitle(fHist.GetZaxis()->GetTitle());
705 h->SetYTitle("Counts");
706 h->SetDirectory(NULL);
707 h->SetMarkerStyle(kFullDotMedium);
708 h->SetBit(kCanDelete);
709 h->SetMarkerColor(kRed);
710 h->SetLineColor(kRed);
711 //h->SetFillColor(18);
712 h->Draw("same"/*"bar same"*/);
713
714 h = (TH1D*)h->Clone("ProjOnOff");
715 h->SetBit(TH1::kNoTitle);
716 h->SetXTitle(fHist.GetZaxis()->GetTitle());
717 h->SetYTitle("Counts");
718 h->SetDirectory(NULL);
719 h->SetMarkerStyle(kFullDotMedium);
720 h->SetBit(kCanDelete);
721 h->SetMarkerColor(kBlue);
722 h->SetLineColor(kBlue);
723 h->Draw("same");
724
725 TLine lin;
726 lin.SetLineStyle(kDashed);
727 lin.DrawLine(h->GetXaxis()->GetXmin(), 0, h->GetXaxis()->GetXmax(), 0);
728 }
729
730 // After the Alpha-projection has been drawn. Fit the histogram
731 // and paint the result into this pad
732 AppendPad("variable");
733
734 if (fHEnergy.GetNbinsX()>1 || fHEnergy.GetBinContent(1)>0)
735 {
736 pad->cd(2);
737 gPad->SetBorderMode(0);
738 gPad->SetGridx();
739 gPad->SetGridy();
740 fHEnergy.Draw();
741
742 AppendPad("energy");
743
744 h = (TH1D*)fHist.Project3D("y");
745 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
746 h->SetXTitle("E [GeV]");
747 h->SetYTitle("Counts");
748 h->SetDirectory(NULL);
749 h->SetMarkerStyle(kFullDotMedium);
750 h->SetBit(kCanDelete);
751 h->SetMarkerColor(kCyan);
752 h->SetLineColor(kCyan);
753 h->Draw("Psame");
754 }
755 else
756 delete pad->GetPad(2);
757
758 if (fTimeEffOn && fTime || fHTime.GetNbinsX()>1 || fHTime.GetBinContent(1)>0)
759 {
760 pad->cd(3);
761 gPad->SetBorderMode(0);
762 gPad->SetGridx();
763 gPad->SetGridy();
764 fHTime.Draw();
765 AppendPad("time");
766 }
767 else
768 delete pad->GetPad(3);
769
770 if (fHTheta.GetNbinsX()>1 || fHTheta.GetBinContent(1)>0)
771 {
772 pad->cd(4);
773 gPad->SetGridx();
774 gPad->SetGridy();
775 gPad->SetBorderMode(0);
776 fHTheta.Draw();
777
778 AppendPad("theta");
779
780 h = (TH1D*)fHist.Project3D("x");
781 h->SetBit(TH1::kNoTitle|TH1::kNoStats);
782 h->SetXTitle("\\theta [\\circ]");
783 h->SetYTitle("Counts");
784 h->SetDirectory(NULL);
785 h->SetMarkerStyle(kFullDotMedium);
786 h->SetBit(kCanDelete);
787 h->SetMarkerColor(kCyan);
788 h->SetLineColor(kCyan);
789 h->Draw("Psame");
790 }
791 else
792 delete pad->GetPad(4);
793}
794
795void MHAlpha::DrawAll(Bool_t newc)
796{
797 if (!newc && !fDisplay)
798 return;
799
800 // FIXME: Do in Paint if special option given!
801 TCanvas &c = newc ? *new TCanvas : fDisplay->AddTab("SizeBins");
802 Int_t n = fHist.GetNbinsY();
803 Int_t nc = (Int_t)(TMath::Sqrt((Float_t)n-1)+1);
804 c.Divide(nc, nc, 1e-10, 1e-10);
805 gPad->SetBorderMode(0);
806 gPad->SetFrameBorderMode(0);
807
808 // Do not store the 'final' result in fFit
809 MAlphaFitter fit(fFit);
810
811 for (int i=1; i<=fHist.GetNbinsY(); i++)
812 {
813 c.cd(i);
814 gPad->SetBorderMode(0);
815 gPad->SetFrameBorderMode(0);
816
817 TH1D *hon = fHist.ProjectionZ("Proj", -1, 9999, i, i, "E");
818 hon->SetBit(TH1::kNoTitle);
819 hon->SetStats(kTRUE);
820 hon->SetXTitle(fHist.GetZaxis()->GetTitle());
821 hon->SetYTitle("Counts");
822 hon->SetDirectory(NULL);
823 hon->SetMarkerStyle(kFullDotMedium);
824 hon->SetBit(kCanDelete);
825 hon->Draw("");
826
827 TH1D *hof = 0;
828 Double_t alpha = 1;
829
830 if (fOffData)
831 {
832 hon->SetMarkerColor(kGreen);
833
834 hof = fOffData->ProjectionZ("ProjOff", -1, 9999, i, i, "E");
835 hof->SetBit(TH1::kNoTitle|TH1::kNoStats);
836 hof->SetXTitle(fHist.GetZaxis()->GetTitle());
837 hof->SetYTitle("Counts");
838 hof->SetDirectory(NULL);
839 hof->SetMarkerStyle(kFullDotMedium);
840 hof->SetBit(kCanDelete);
841 hof->SetMarkerColor(kRed);
842 hof->SetLineColor(kRed);
843 hof->Draw("same");
844
845 alpha = fit.Scale(*hof, *hon);
846
847 hon->SetMaximum();
848 hon->SetMaximum(TMath::Max(hon->GetMaximum(), hof->GetMaximum())*1.05);
849
850 TH1D *diff = new TH1D(*hon);
851 diff->Add(hof, -1);
852 diff->SetBit(TH1::kNoTitle);
853 diff->SetXTitle(fHist.GetZaxis()->GetTitle());
854 diff->SetYTitle("Counts");
855 diff->SetDirectory(NULL);
856 diff->SetMarkerStyle(kFullDotMedium);
857 diff->SetBit(kCanDelete);
858 diff->SetMarkerColor(kBlue);
859 diff->Draw("same");
860
861 TLine lin;
862 lin.SetLineStyle(kDashed);
863 lin.DrawLine(diff->GetXaxis()->GetXmin(), 0, diff->GetXaxis()->GetXmax(), 0);
864
865 const Float_t min = diff->GetMinimum()*1.05;
866 hon->SetMinimum(min<0 ? min : 0);
867 }
868
869 if (hof ? fit.Fit(*hon, *hof, alpha) : fit.Fit(*hon))
870 {
871 *fLog << dbg << "Bin " << i << ": sigmaexc=" << fit.GetSignificanceExc() << " omega=" << fit.GetGausSigma() << " events=" << fit.GetEventsExcess() << " scale=" << fit.GetScaleFactor() << endl;
872 fit.DrawResult();
873 }
874 /*
875 if (fit.FitEnergy(fHist, fOffData, i, kTRUE))
876 {
877 fHEnergy.SetBinContent(i, fit.GetEventsExcess());
878 fHEnergy.SetBinError(i, fit.GetEventsExcess()*0.2);
879 }*/
880 }
881
882}
883
884void MHAlpha::DrawNicePlot(Bool_t newc, const char *title, const char *watermark)
885{
886 if (!newc && !fDisplay)
887 return;
888
889 if (!fOffData)
890 return;
891
892 // Open and setup canvas/pad
893 TCanvas &c = newc ? *new TCanvas : fDisplay->AddTab("ThetsSq");
894
895 c.SetBorderMode(0);
896 c.SetFrameBorderMode(0);
897 c.SetFillColor(kWhite);
898
899 c.SetLeftMargin(0.12);
900 c.SetRightMargin(0.01);
901 c.SetBottomMargin(0.16);
902 c.SetTopMargin(0.18);
903
904 c.SetGridy();
905
906 gStyle->SetOptStat(0);
907
908 // Get on-data
909 TH1D *hon = (TH1D*)fHist.ProjectionZ("Proj", -1, 9999, -1, 9999, "E");
910 hon->SetDirectory(NULL);
911 hon->SetBit(kCanDelete);
912 hon->SetMarkerSize(0);
913 hon->SetLineWidth(2);
914 hon->SetLineColor(kBlack);
915 hon->SetMarkerColor(kBlack);
916
917 // Get off-data
918 TH1D *hoff = 0;
919 if (fOffData)
920 {
921 hoff = (TH1D*)fOffData->ProjectionZ("ProjOff", -1, 9999, -1, 9999, "E");
922 hoff->SetDirectory(NULL);
923 hoff->SetBit(kCanDelete);
924 hoff->SetFillColor(17);
925 hoff->SetMarkerSize(0);
926 hoff->SetLineColor(kBlack);
927 hoff->SetMarkerColor(kBlack);
928 }
929
930 // Setup plot which is drawn first
931 TH1D *h = hoff ? hoff : hon;
932 h->GetXaxis()->SetLabelSize(0.06);
933 h->GetXaxis()->SetTitleSize(0.06);
934 h->GetXaxis()->SetTitleOffset(0.95);
935 h->GetYaxis()->SetLabelSize(0.06);
936 h->GetYaxis()->SetTitleSize(0.06);
937 h->GetYaxis()->CenterTitle();
938 h->SetYTitle("Counts");
939 h->SetTitleSize(0.07);
940 h->SetTitle("");
941
942 const Double_t imax = fFit.GetSignalIntegralMax();
943 if (imax<1)
944 h->GetXaxis()->SetRangeUser(0, 0.6*0.6);
945
946 // scale off-data
947 fFit.Scale(*hoff, *hon);
948
949 hon->SetMinimum(0);
950 hoff->SetMinimum(0);
951
952 // draw data
953 if (hoff)
954 {
955 hoff->SetMaximum(TMath::Max(hon->GetMaximum(),hoff->GetMaximum())*1.1);
956 hoff->Draw("bar");
957 hon->Draw("same");
958 }
959 else
960 {
961 hon->SetMaximum();
962 hon->Draw();
963 }
964
965 // draw a preliminary tag
966 TLatex text;
967 text.SetTextColor(kWhite);
968 text.SetBit(TLatex::kTextNDC);
969 text.SetTextSize(0.07);
970 text.SetTextAngle(2.5);
971 if (watermark)
972 text.DrawLatex(0.45, 0.2, watermark);
973 //enum { kTextNDC = BIT(14) };
974
975 // draw line showing cut
976 TLine line;
977 line.SetLineColor(14);
978 line.SetLineStyle(7);
979 line.DrawLine(imax, 0, imax, h->GetMaximum());
980
981 // Add a title above the plot
982 TPaveText *pave=new TPaveText(0.12, 0.83, 0.99, 0.975, "blNDC");
983 pave->SetBorderSize(1);
984 pave->SetLabel(title);
985
986 char txt[1000];
987 TText *ptxt;
988 sprintf(txt, " ");
989 ptxt = pave->AddText(txt);
990 ptxt->SetTextAlign(23);
991
992 //sprintf(txt, "Significance %.1f\\sigma, off-scale %.2f (\\omega=%.2f\\circ)",
993 // fFit.GetSignificance(), fFit.GetScaleFactor(), fFit.GetGausSigma());
994 sprintf(txt, "Significance %.1f\\sigma, off-scale %.2f",
995 fFit.GetSignificance(), fFit.GetScaleFactor());
996 ptxt = pave->AddText(txt);
997 ptxt->SetTextAlign(23);
998
999 sprintf(txt, "%.1f excess events, %.1f background events",
1000 fFit.GetEventsExcess(), fFit.GetEventsBackground());
1001 ptxt = pave->AddText(txt);
1002 ptxt->SetTextAlign(23);
1003 pave->SetBit(kCanDelete);
1004 pave->Draw();
1005}
1006
1007Bool_t MHAlpha::Finalize()
1008{
1009 if (!FitAlpha())
1010 {
1011 *fLog << warn << "MAlphaFitter - Fit failed..." << endl;
1012 return kTRUE;
1013 }
1014
1015 // Store the final result in fFit
1016 *fLog << all;
1017 fFit.Print("result");
1018
1019 fResult->SetVal(fFit.GetMinimizationValue());
1020 fSigma->SetVal(fFit.GetGausSigma());
1021
1022 if (!fSkipHistEnergy)
1023 {
1024 *fLog << inf << "Processing energy bins..." << endl;
1025 FitEnergyBins();
1026 }
1027 if (!fSkipHistTheta)
1028 {
1029 *fLog << inf << "Processing theta bins..." << endl;
1030 FitThetaBins();
1031 }
1032 if (!fSkipHistTime)
1033 {
1034 *fLog << inf << "Processing time bins..." << endl;
1035 UpdateAlphaTime(kTRUE);
1036 MH::RemoveFirstBin(fHTime);
1037 }
1038
1039 if (fOffData)
1040 DrawAll(kFALSE);
1041
1042 return kTRUE;
1043}
1044
1045// --------------------------------------------------------------------------
1046//
1047// You can use this function if you want to use a MHMatrix instead of
1048// MMcEvt. This function adds all necessary columns to the
1049// given matrix. Afterward you should fill the matrix with the corresponding
1050// data (eg from a file by using MHMatrix::Fill). If you now loop
1051// through the matrix (eg using MMatrixLoop) MHHadronness::Fill
1052// will take the values from the matrix instead of the containers.
1053//
1054// It takes fSkipHist* into account!
1055//
1056void MHAlpha::InitMapping(MHMatrix *mat, Int_t type)
1057{
1058 if (fMatrix)
1059 return;
1060
1061 fMatrix = mat;
1062
1063 fMap[0] = fMatrix->AddColumn(GetParameterRule());
1064 fMap[1] = -1;
1065 fMap[2] = -1;
1066 fMap[3] = -1;
1067 fMap[4] = -1;
1068
1069 if (!fSkipHistEnergy)
1070 if (type==0)
1071 {
1072 fMap[1] = fMatrix->AddColumn("MEnergyEst.fVal");
1073 fMap[2] = -1;
1074 }
1075 else
1076 {
1077 fMap[1] = -1;
1078 fMap[2] = fMatrix->AddColumn("MHillas.fSize");
1079 }
1080
1081 if (!fSkipHistTheta)
1082 fMap[3] = fMatrix->AddColumn("MPointingPos.fZd");
1083
1084 // if (!fSkipHistTime)
1085 // fMap[4] = fMatrix->AddColumn("MTime.GetAxisTime");
1086}
1087
1088void MHAlpha::StopMapping()
1089{
1090 fMatrix = NULL;
1091}
1092
1093void MHAlpha::ApplyScaling()
1094{
1095 if (!fOffData)
1096 return;
1097
1098 fFit.ApplyScaling(fHist, *const_cast<TH3D*>(fOffData));
1099}
1100
1101Int_t MHAlpha::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1102{
1103 Bool_t rc = kFALSE;
1104 if (IsEnvDefined(env, prefix, "NumTimeBins", print))
1105 {
1106 SetNumTimeBins(GetEnvValue(env, prefix, "NumTimeBins", fNumTimeBins));
1107 rc = kTRUE;
1108 }
1109 if (IsEnvDefined(env, prefix, "ForceUsingSize", print))
1110 {
1111 fForceUsingSize = GetEnvValue(env, prefix, "ForceUsingSize", fForceUsingSize);
1112 rc = kTRUE;
1113 }
1114 return rc;
1115}
1116
1117Int_t MHAlpha::DistancetoPrimitive(Int_t px, Int_t py)
1118{
1119 // If pad has no subpad return (we are in one of the subpads)
1120 if (gPad->GetPad(1)==NULL)
1121 return 9999;
1122
1123 // If pad has a subpad we are in the main pad. Check its value.
1124 return gPad->GetPad(1)->DistancetoPrimitive(px,py)==0 ? 0 : 9999;
1125}
Note: See TracBrowser for help on using the repository browser.