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

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