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

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