source: tags/Mars-V0.8.4/mhist/MHCamera.cc

Last change on this file was 3667, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 41.8 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, 05/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Harald Kornmayer, 1/2001
20! Author(s): Markus Gaug, 03/2004 <mailto:markus@ifae.es>
21!
22! Copyright: MAGIC Software Development, 2000-2004
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28//
29// MHCamera
30//
31// Camera Display, based on a TH1D. Pleas be carefull using the
32// underlaying TH1D.
33//
34// To change the scale to a logarithmic scale SetLogy() of the Pad.
35//
36// Be carefull: Entries in this context means Entries/bin or Events
37//
38// FIXME? Maybe MHCamera can take the fLog object from MGeomCam?
39//
40////////////////////////////////////////////////////////////////////////////
41#include "MHCamera.h"
42
43#include <fstream>
44#include <iostream>
45
46#include <TBox.h>
47#include <TArrow.h>
48#include <TLatex.h>
49#include <TStyle.h>
50#include <TCanvas.h>
51#include <TArrayF.h>
52#include <TRandom.h>
53#include <TPaveText.h>
54#include <TPaveStats.h>
55#include <TClonesArray.h>
56#include <THistPainter.h>
57#include <THLimitsFinder.h>
58#include <TProfile.h>
59
60#include "MLog.h"
61#include "MLogManip.h"
62
63#include "MH.h"
64#include "MHexagon.h"
65
66#include "MGeomPix.h"
67#include "MGeomCam.h"
68
69#include "MCamEvent.h"
70
71
72#define kItemsLegend 48 // see SetPalette(1,0)
73
74ClassImp(MHCamera);
75
76using namespace std;
77
78void MHCamera::Init()
79{
80 UseCurrentStyle();
81
82 SetDirectory(NULL);
83
84 SetLineColor(kGreen);
85 SetMarkerStyle(kFullDotMedium);
86 SetXTitle("Pixel Index");
87
88 fNotify = new TList;
89 fNotify->SetBit(kMustCleanup);
90 gROOT->GetListOfCleanups()->Add(fNotify);
91
92#if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06)
93 SetPalette(1, 0);
94#else
95 SetInvDeepBlueSeaPalette();
96#endif
97}
98
99// ------------------------------------------------------------------------
100//
101// Default Constructor. To be used by the root system ONLY.
102//
103MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fColors(kItemsLegend)
104{
105 Init();
106}
107
108// ------------------------------------------------------------------------
109//
110// Constructor. Makes a clone of MGeomCam. Removed the TH1D from the
111// current directory. Calls Sumw2(). Set the histogram line color
112// (for error bars) to Green and the marker style to kFullDotMedium.
113//
114MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title)
115: fGeomCam(NULL), fColors(kItemsLegend)
116{
117 //fGeomCam = (MGeomCam*)geom.Clone();
118 SetGeometry(geom, name, title);
119 Init();
120
121 //
122 // root 3.02
123 // * base/inc/TObject.h:
124 // register BIT(8) as kNoContextMenu. If an object has this bit set it will
125 // not get an automatic context menu when clicked with the right mouse button.
126}
127
128void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title)
129{
130 SetNameTitle(name, title);
131
132 TAxis &x = *GetXaxis();
133
134 SetBins(geom.GetNumPixels(), 0, 1);
135 x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
136
137 //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
138 //Rebuild();
139
140 Sumw2(); // necessary?
141
142 if (fGeomCam)
143 delete fGeomCam;
144 fGeomCam = (MGeomCam*)geom.Clone();
145
146 fUsed.Set(geom.GetNumPixels());
147 for (Int_t i=0; i<fNcells-2; i++)
148 ResetUsed(i);
149}
150
151// ------------------------------------------------------------------------
152//
153// Destructor. Deletes the cloned fGeomCam and the notification list.
154//
155MHCamera::~MHCamera()
156{
157 if (fGeomCam)
158 delete fGeomCam;
159 if (fNotify)
160 delete fNotify;
161}
162
163// ------------------------------------------------------------------------
164//
165// Return kTRUE for sector<0. Otherwise return kTRUE only if the specified
166// sector idx matches the sector of the pixel with index idx.
167//
168Bool_t MHCamera::MatchSector(Int_t idx, const TArrayI &sector, const TArrayI &aidx) const
169{
170 const MGeomPix &pix = (*fGeomCam)[idx];
171 return FindVal(sector, pix.GetSector()) && FindVal(aidx, pix.GetAidx());
172}
173
174// ------------------------------------------------------------------------
175//
176// Taken from TH1D::Fill(). Uses the argument directly as bin index.
177// Doesn't increment the number of entries.
178//
179// -*-*-*-*-*-*-*-*Increment bin with abscissa X by 1*-*-*-*-*-*-*-*-*-*-*
180// ==================================
181//
182// if x is less than the low-edge of the first bin, the Underflow bin is incremented
183// if x is greater than the upper edge of last bin, the Overflow bin is incremented
184//
185// If the storage of the sum of squares of weights has been triggered,
186// via the function Sumw2, then the sum of the squares of weights is incremented
187// by 1 in the bin corresponding to x.
188//
189// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
190Int_t MHCamera::Fill(Axis_t x)
191{
192#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
193 if (fBuffer) return BufferFill(x,1);
194#endif
195 const Int_t bin = (Int_t)x+1;
196 AddBinContent(bin);
197 if (fSumw2.fN)
198 fSumw2.fArray[bin]++;
199
200 if (bin<=0 || bin>fNcells-2)
201 return -1;
202
203 fTsumw++;
204 fTsumw2++;
205 fTsumwx += x;
206 fTsumwx2 += x*x;
207 return bin;
208}
209
210// ------------------------------------------------------------------------
211//
212// Taken from TH1D::Fill(). Uses the argument directly as bin index.
213// Doesn't increment the number of entries.
214//
215// -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-*
216// =============================================
217//
218// if x is less than the low-edge of the first bin, the Underflow bin is incremented
219// if x is greater than the upper edge of last bin, the Overflow bin is incremented
220//
221// If the storage of the sum of squares of weights has been triggered,
222// via the function Sumw2, then the sum of the squares of weights is incremented
223// by w^2 in the bin corresponding to x.
224//
225// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
226Int_t MHCamera::Fill(Axis_t x, Stat_t w)
227{
228#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
229 if (fBuffer) return BufferFill(x,w);
230#endif
231 const Int_t bin = (Int_t)x+1;
232 AddBinContent(bin, w);
233 if (fSumw2.fN)
234 fSumw2.fArray[bin] += w*w;
235
236 if (bin<=0 || bin>fNcells-2)
237 return -1;
238
239 const Stat_t z = (w > 0 ? w : -w);
240 fTsumw += z;
241 fTsumw2 += z*z;
242 fTsumwx += z*x;
243 fTsumwx2 += z*x*x;
244 return bin;
245}
246
247// ------------------------------------------------------------------------
248//
249// Use x and y in millimeters
250//
251Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w)
252{
253 if (fNcells<=1 || IsFreezed())
254 return -1;
255
256 for (Int_t idx=0; idx<fNcells-2; idx++)
257 {
258 MHexagon hex((*fGeomCam)[idx]);
259 if (hex.DistanceToPrimitive(x, y)>0)
260 continue;
261
262 SetUsed(idx);
263 return Fill(idx, w);
264 }
265 return -1;
266}
267
268// ------------------------------------------------------------------------
269//
270// Return the mean value of all entries which are used if all=kFALSE and
271// of all entries if all=kTRUE if sector<0. If sector>=0 only
272// entries with match the given sector are taken into account.
273//
274Stat_t MHCamera::GetMeanSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
275{
276 if (fNcells<=1)
277 return 0;
278
279 Int_t n=0;
280
281 Stat_t mean = 0;
282 for (int i=0; i<fNcells-2; i++)
283 {
284 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
285 {
286 mean += fArray[i+1];
287 n++;
288 }
289 }
290
291 return n==0 ? 0 : Profile(mean/n);
292}
293
294// ------------------------------------------------------------------------
295//
296// Return the variance of all entries which are used if all=kFALSE and
297// of all entries if all=kTRUE if sector<0. If sector>=0 only
298// entries with match the given sector are taken into account.
299//
300Stat_t MHCamera::GetRmsSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
301{
302 if (fNcells<=1)
303 return -1;
304
305 Int_t n=0;
306
307 Stat_t sum = 0;
308 Stat_t sq = 0;
309 for (int i=0; i<fNcells-2; i++)
310 {
311 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
312 {
313 sum += fArray[i+1];
314 sq += fArray[i+1]*fArray[i+1];
315 n++;
316 }
317 }
318
319 if (n==0)
320 return 0;
321
322 sum /= n;
323 sq /= n;
324
325 return Profile(TMath::Sqrt(sq-sum*sum));
326}
327
328// ------------------------------------------------------------------------
329//
330// Return the minimum contents of all pixels (if all is set, otherwise
331// only of all 'used' pixels), fMinimum if fMinimum set. If sector>=0
332// only pixels with matching sector number are taken into account.
333//
334Double_t MHCamera::GetMinimumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
335{
336 if (fMinimum != -1111)
337 return fMinimum;
338
339 if (fNcells<=1)
340 return 0;
341
342 Double_t minimum=FLT_MAX;
343
344 for (Int_t idx=0; idx<fNcells-2; idx++)
345 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]<minimum)
346 minimum = fArray[idx+1];
347
348 return Profile(minimum);
349}
350
351// ------------------------------------------------------------------------
352//
353// Return the maximum contents of all pixels (if all is set, otherwise
354// only of all 'used' pixels), fMaximum if fMaximum set. If sector>=0
355// only pixels with matching sector number are taken into account.
356//
357Double_t MHCamera::GetMaximumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
358{
359 if (fMaximum!=-1111)
360 return fMaximum;
361
362 if (fNcells<=1)
363 return 1;
364
365 Double_t maximum=-FLT_MAX;
366 for (Int_t idx=0; idx<fNcells-2; idx++)
367 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]>maximum)
368 maximum = fArray[idx+1];
369
370 return Profile(maximum);
371}
372
373// ------------------------------------------------------------------------
374//
375// Call this function to draw the camera layout into your canvas.
376// Setup a drawing canvas. Add this object and all child objects
377// (hexagons, etc) to the current pad. If no pad exists a new one is
378// created. (To access the 'real' pad containing the camera you have
379// to do a cd(1) in the current layer.
380//
381// To draw a camera into its own pad do something like:
382//
383// MGeomCamMagic m;
384// MHCamera *d=new MHCamera(m);
385//
386// TCanvas *c = new TCanvas;
387// c->Divide(2,1);
388// c->cd(1);
389//
390// d->FillRandom();
391// d->Draw();
392// d->SetBit(kCanDelete);
393//
394// There are several drawing options:
395// 'hist' Draw as a standard TH1 histogram (value vs. pixel index)
396// 'box' Draw hexagons which size is in respect to its contents
397// 'nocol' Leave the 'boxed' hexagons empty
398// 'pixelindex' Display the pixel index in each pixel
399// 'sectorindex' Display the sector index in each pixel
400// 'content' Display the relative content aligned to GetMaximum() and
401// GeMinimum() ((val-min)/(max-min))
402// 'proj' Display the y-projection of the histogram
403// 'same' Draw trandparent pixels on top of an existing pad. This
404// makes it possible to draw the camera image on top of an
405// existing TH2, but also allows for distorted camera images
406//
407void MHCamera::Draw(Option_t *option)
408{
409 const Bool_t hassame = TString(option).Contains("same", TString::kIgnoreCase) && gPad;
410
411 // root 3.02:
412 // gPad->SetFixedAspectRatio()
413 const Color_t col = gPad ? gPad->GetFillColor() : 16;
414 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600);
415
416 if (!hassame)
417 {
418 pad->SetBorderMode(0);
419 pad->SetFillColor(col);
420
421 //
422 // Create an own pad for the MHCamera-Object which can be
423 // resized in paint to keep the correct aspect ratio
424 //
425 pad->Divide(1, 1, 0, 0, col);
426 pad->cd(1);
427 gPad->SetBorderMode(0);
428 }
429
430 AppendPad(option);
431 //fGeomCam->AppendPad();
432
433 //
434 // Do not change gPad. The user should not see, that Draw
435 // changes gPad...
436 //
437 if (!hassame)
438 pad->cd();
439}
440
441// ------------------------------------------------------------------------
442//
443// This is TObject::DrawClone but completely ignores
444// gROOT->GetSelectedPad(). tbretz had trouble with this in the past.
445// If this makes trouble please write a bug report.
446//
447TObject *MHCamera::DrawClone(Option_t *option) const
448{
449 // Draw a clone of this object in the current pad
450
451 //TVirtualPad *pad = gROOT->GetSelectedPad();
452 TVirtualPad *padsav = gPad;
453 //if (pad) pad->cd();
454
455 TObject *newobj = Clone();
456
457 if (!newobj)
458 return 0;
459
460 /*
461 if (pad) {
462 if (strlen(option)) pad->GetListOfPrimitives()->Add(newobj,option);
463 else pad->GetListOfPrimitives()->Add(newobj,GetDrawOption());
464 pad->Modified(kTRUE);
465 pad->Update();
466 if (padsav) padsav->cd();
467 return newobj;
468 }
469 */
470
471 TString opt(option);
472 opt.ToLower();
473
474 newobj->Draw(opt.IsNull() ? GetDrawOption() : option);
475
476 if (padsav)
477 padsav->cd();
478
479 return newobj;
480}
481
482// ------------------------------------------------------------------------
483//
484// Creates a TH1D which contains the projection of the contents of the
485// MHCamera onto the y-axis. The maximum and minimum are calculated
486// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
487// displayed using THLimitsFinder::OptimizeLimits.
488//
489// If no name is given the newly allocated histogram is removed from
490// the current directory calling SetDirectory(0) in any other case
491// the newly created histogram is removed from the current directory
492// and added to gROOT such the gROOT->FindObject can find the histogram.
493//
494// If the standard name "_py" is given "_py" is appended to the name
495// of the MHCamera and the corresponding histogram is searched using
496// gROOT->FindObject and updated with the present projection.
497//
498// It is the responsibility of the user to make sure, that the newly
499// created histogram is freed correctly.
500//
501// Currently the new histogram is restrictred to 50 bins.
502// Maybe a optimal number can be calulated from the number of
503// bins on the x-axis of the MHCamera?
504//
505// The code was taken mainly from TH2::ProjectX such the interface
506// is more or less the same than to TH2-projections.
507//
508// If sector>=0 only entries with matching sector index are taken
509// into account.
510//
511TH1D *MHCamera::ProjectionS(const TArrayI &sector, const TArrayI &aidx, const char *name) const
512{
513 Int_t nbins = 50;
514
515 // Create the projection histogram
516 TString pname(name);
517 if (name=="_py")
518 {
519 pname.Prepend(GetName());
520 if (sector.GetSize()>0)
521 {
522 pname += ";";
523 for (int i=0; i<sector.GetSize(); i++)
524 pname += sector[i];
525 }
526 if (aidx.GetSize()>0)
527 {
528 pname += ";";
529 for (int i=0; i<aidx.GetSize(); i++)
530 pname += aidx[i];
531 }
532 }
533
534 TH1D *h1=0;
535
536 //check if histogram with identical name exist
537 TObject *h1obj = gROOT->FindObject(pname);
538 if (h1obj && h1obj->InheritsFrom("TH1D")) {
539 h1 = (TH1D*)h1obj;
540 h1->Reset();
541 }
542
543 if (!h1)
544 {
545 Double_t min = GetMinimumSectors(sector, aidx);
546 Double_t max = GetMaximumSectors(sector, aidx);
547
548 Int_t newbins=0;
549
550 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
551
552 h1 = new TH1D(pname, GetTitle(), nbins, min, max);
553 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
554 h1->SetXTitle(GetYaxis()->GetTitle());
555 h1->SetYTitle("Counts");
556 //h1->Sumw2();
557 }
558
559 // Fill the projected histogram
560 for (Int_t idx=0; idx<fNcells-2; idx++)
561 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
562 h1->Fill(GetBinContent(idx+1));
563
564 return h1;
565}
566
567// ------------------------------------------------------------------------
568//
569// Creates a TH1D which contains the projection of the contents of the
570// MHCamera onto the radius from the camera center.
571// The maximum and minimum are calculated
572// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
573// displayed using THLimitsFinder::OptimizeLimits.
574//
575// If no name is given the newly allocated histogram is removed from
576// the current directory calling SetDirectory(0) in any other case
577// the newly created histogram is removed from the current directory
578// and added to gROOT such the gROOT->FindObject can find the histogram.
579//
580// If the standard name "_rad" is given "_rad" is appended to the name
581// of the MHCamera and the corresponding histogram is searched using
582// gROOT->FindObject and updated with the present projection.
583//
584// It is the responsibility of the user to make sure, that the newly
585// created histogram is freed correctly.
586//
587// Currently the new histogram is restrictred to 50 bins.
588// Maybe a optimal number can be calulated from the number of
589// bins on the x-axis of the MHCamera?
590//
591// The code was taken mainly from TH2::ProjectX such the interface
592// is more or less the same than to TH2-projections.
593//
594// If sector>=0 only entries with matching sector index are taken
595// into account.
596//
597TProfile *MHCamera::RadialProfileS(const TArrayI &sector, const TArrayI &aidx, const char *name, const Int_t nbins) const
598{
599
600 // Create the projection histogram
601 TString pname(name);
602 if (name=="_rad")
603 {
604 pname.Prepend(GetName());
605 if (sector.GetSize()>0)
606 {
607 pname += ";";
608 for (int i=0; i<sector.GetSize(); i++)
609 pname += sector[i];
610 }
611 if (aidx.GetSize()>0)
612 {
613 pname += ";";
614 for (int i=0; i<aidx.GetSize(); i++)
615 pname += aidx[i];
616 }
617 }
618
619 TProfile *h1=0;
620
621 //check if histogram with identical name exist
622 TObject *h1obj = gROOT->FindObject(pname);
623 if (h1obj && h1obj->InheritsFrom("TProfile")) {
624 h1 = (TProfile*)h1obj;
625 h1->Reset();
626 }
627
628 if (!h1)
629 {
630
631 Double_t min = 0.;
632 Double_t max = fGeomCam->GetMaxRadius();
633
634 Int_t newbins=0;
635
636 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
637
638 h1 = new TProfile(pname, GetTitle(), nbins, min, max);
639 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
640 h1->SetXTitle("Radius from camera center [mm]");
641 h1->SetYTitle(GetYaxis()->GetTitle());
642 }
643
644 // Fill the projected histogram
645 for (Int_t idx=0; idx<fNcells-2; idx++)
646 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
647 h1->Fill(TMath::Hypot((*fGeomCam)[idx].GetX(),(*fGeomCam)[idx].GetY()),
648 GetBinContent(idx+1));
649 return h1;
650}
651
652
653// ------------------------------------------------------------------------
654//
655// Resizes the current pad so that the camera is displayed in its
656// correct aspect ratio
657//
658void MHCamera::SetRange()
659{
660 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
661
662 //
663 // Maintain aspect ratio
664 //
665 const float ratio = TestBit(kNoLegend) ? 1 : 1.15;
666
667 //
668 // Calculate width and height of the current pad in pixels
669 //
670 Float_t w = gPad->GetWw();
671 Float_t h = gPad->GetWh()*ratio;
672
673 //
674 // This prevents the pad from resizing itself wrongly
675 //
676 if (gPad->GetMother() != gPad)
677 {
678 w *= gPad->GetMother()->GetAbsWNDC();
679 h *= gPad->GetMother()->GetAbsHNDC();
680 }
681
682 //
683 // Set Range (coordinate system) of pad
684 //
685 gPad->Range(-range, -range, (2*ratio-1)*range, range);
686
687 //
688 // Resize Pad to given ratio
689 //
690 if (h<w)
691 gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1);
692 else
693 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2);
694}
695
696// ------------------------------------------------------------------------
697//
698// Updates the pixel colors and paints the pixels
699//
700void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol, Bool_t issame)
701{
702 Double_t min = GetMinimum(kFALSE);
703 Double_t max = GetMaximum(kFALSE);
704 if (min==FLT_MAX)
705 {
706 min = 0;
707 max = 1;
708 }
709
710 if (min==max)
711 max += 1;
712
713 if (!issame)
714 UpdateLegend(min, max, islog);
715
716 // Try to estimate the units of the current display. This is only
717 // necessary for 'same' option and allows distorted images of the camera!
718 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
719 const Float_t conv = !issame ||
720 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
721 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
722
723 MHexagon hex;
724 for (Int_t i=0; i<fNcells-2; i++)
725 {
726 hex.SetFillStyle(issame ? 0 : 1001);
727
728 if (!issame)
729 {
730 if (IsUsed(i) && iscol)
731 {
732 if (TMath::IsNaN(fArray[i+1]))
733 gLog << warn << "MHCamera::Update: " << GetName() << " <" << GetTitle() << "> - Pixel Index #" << i << " contents is NaN (Not a Number)..." << endl;
734
735 hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog));
736 }
737 else
738 hex.SetFillColor(10);
739 }
740
741 const MGeomPix &pix = (*fGeomCam)[i];
742
743 const Float_t x = pix.GetX()*conv;
744 const Float_t y = pix.GetY()*conv;
745 const Float_t d = pix.GetD()*conv;
746
747 if (!isbox)
748 hex.PaintHexagon(x, y, d);
749 else
750 if (IsUsed(i) && !TMath::IsNaN(fArray[i+1]))
751 {
752 Float_t size = d*(GetBinContent(i+1)-min)/(max-min);
753 if (size>d)
754 size=d;
755 hex.PaintHexagon(x, y, size);
756 }
757 }
758}
759
760// ------------------------------------------------------------------------
761//
762// Print minimum and maximum
763//
764void MHCamera::Print(Option_t *) const
765{
766 gLog << all << "Minimum: " << GetMinimum();
767 if (fMinimum==-1111)
768 gLog << " <autoscaled>";
769 gLog << endl;
770 gLog << "Maximum: " << GetMaximum();
771 if (fMaximum==-1111)
772 gLog << " <autoscaled>";
773 gLog << endl;
774}
775
776// ------------------------------------------------------------------------
777//
778// Paint the y-axis title
779//
780void MHCamera::PaintAxisTitle()
781{
782 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
783 const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range;
784
785 TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle());
786
787 ptitle->SetTextSize(0.05);
788 ptitle->SetTextAlign(21);
789
790 // box with the histogram title
791 ptitle->SetTextColor(gStyle->GetTitleTextColor());
792#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01)
793 ptitle->SetTextFont(gStyle->GetTitleFont(""));
794#endif
795 ptitle->Paint();
796}
797
798// ------------------------------------------------------------------------
799//
800// Paints the camera.
801//
802void MHCamera::Paint(Option_t *o)
803{
804 if (fNcells<=1)
805 return;
806
807 TString opt(o);
808 opt.ToLower();
809
810 if (opt.Contains("hist"))
811 {
812 opt.ReplaceAll("hist", "");
813 TH1D::Paint(opt);
814 return;
815 }
816
817 if (opt.Contains("proj"))
818 {
819 opt.ReplaceAll("proj", "");
820 Projection(GetName())->Paint(opt);
821 return;
822 }
823
824 const Bool_t hassame = opt.Contains("same");
825 const Bool_t hasbox = opt.Contains("box");
826 const Bool_t hascol = hasbox ? !opt.Contains("nocol") : kTRUE;
827
828 if (!hassame)
829 {
830 gPad->Clear();
831
832 // Maintain aspect ratio
833 SetRange();
834
835 if (GetPainter())
836 {
837 // Paint statistics
838 if (!TestBit(TH1::kNoStats))
839 fPainter->PaintStat(gStyle->GetOptStat(), NULL);
840
841 // Paint primitives (pixels, color legend, photons, ...)
842 if (fPainter->InheritsFrom(THistPainter::Class()))
843 {
844 static_cast<THistPainter*>(fPainter)->MakeChopt("");
845 static_cast<THistPainter*>(fPainter)->PaintTitle();
846 }
847 }
848 }
849
850 // Update Contents of the pixels and paint legend
851 Update(gPad->GetLogy(), hasbox, hascol, hassame);
852
853 if (!hassame)
854 PaintAxisTitle();
855
856 if (opt.Contains("pixelindex"))
857 PaintIndices(0);
858 if (opt.Contains("sectorindex"))
859 PaintIndices(1);
860 if (opt.Contains("content"))
861 PaintIndices(2);
862}
863
864// ------------------------------------------------------------------------
865//
866// With this function you can change the color palette. For more
867// information see TStyle::SetPalette. Only palettes with 50 colors
868// are allowed.
869// In addition you can use SetPalette(52, 0) to create an inverse
870// deep blue sea palette
871//
872void MHCamera::SetPalette(Int_t ncolors, Int_t *colors)
873{
874 //
875 // If not enough colors are specified skip this.
876 //
877 if (ncolors>1 && ncolors<50)
878 {
879 gLog << err << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl;
880 return;
881 }
882
883 //
884 // If ncolors==52 create a reversed deep blue sea palette
885 //
886 if (ncolors==52)
887 {
888 gStyle->SetPalette(51, NULL);
889 TArrayI c(kItemsLegend);
890 for (int i=0; i<kItemsLegend; i++)
891 c[kItemsLegend-i-1] = gStyle->GetColorPalette(i);
892 gStyle->SetPalette(kItemsLegend, c.GetArray());
893 }
894 else
895 gStyle->SetPalette(ncolors, colors);
896
897 fColors.Set(kItemsLegend);
898 for (int i=0; i<kItemsLegend; i++)
899 fColors[i] = gStyle->GetColorPalette(i);
900}
901
902
903void MHCamera::SetPrettyPalette()
904{
905 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
906 SetPalette(1, 0);
907}
908
909void MHCamera::SetDeepBlueSeaPalette()
910{
911 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
912 SetPalette(51, 0);
913}
914
915void MHCamera::SetInvDeepBlueSeaPalette()
916{
917 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
918 SetPalette(52, 0);
919}
920
921void MHCamera::PaintIndices(Int_t type)
922{
923 if (fNcells<=1)
924 return;
925
926 const Double_t min = GetMinimum();
927 const Double_t max = GetMaximum();
928
929 if (type==2 && max==min)
930 return;
931
932 TText txt;
933 txt.SetTextFont(122);
934 txt.SetTextAlign(22); // centered/centered
935
936 for (Int_t i=0; i<fNcells-2; i++)
937 {
938 const MGeomPix &h = (*fGeomCam)[i];
939
940 TString num;
941 switch (type)
942 {
943 case 0: num += i; break;
944 case 1: num += h.GetSector(); break;
945 case 2: num += (Int_t)((fArray[i+1]-min)/(max-min)); break;
946 }
947
948 // FIXME: Should depend on the color of the pixel...
949 //(GetColor(GetBinContent(i+1), min, max, 0));
950 txt.SetTextColor(kRed);
951 txt.SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05);
952 txt.PaintText(h.GetX(), h.GetY(), num);
953 }
954}
955
956// ------------------------------------------------------------------------
957//
958// Call this function to add a MCamEvent on top of the present contents.
959//
960void MHCamera::AddCamContent(const MCamEvent &event, Int_t type)
961{
962 if (fNcells<=1 || IsFreezed())
963 return;
964
965 // FIXME: Security check missing!
966 for (Int_t idx=0; idx<fNcells-2; idx++)
967 {
968 Double_t val=0;
969 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
970 SetUsed(idx);
971
972 Fill(idx, val); // FIXME: Slow!
973 }
974 fEntries++;
975}
976
977// ------------------------------------------------------------------------
978//
979// Call this function to add a MCamEvent on top of the present contents.
980//
981void MHCamera::SetCamError(const MCamEvent &evt, Int_t type)
982{
983
984 if (fNcells<=1 || IsFreezed())
985 return;
986
987 // FIXME: Security check missing!
988 for (Int_t idx=0; idx<fNcells-2; idx++)
989 {
990 Double_t val=0;
991 if (evt.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
992 SetUsed(idx);
993
994 SetBinError(idx+1, val); // FIXME: Slow!
995 }
996}
997
998Stat_t MHCamera::GetBinError(Int_t bin) const
999{
1000 Double_t rc = 0;
1001 if (TestBit(kVariance) && GetEntries()>0) // error on the mean
1002 {
1003 const Double_t error = fSumw2.fArray[bin]/GetEntries();
1004 const Double_t val = fArray[bin]/GetEntries();
1005 rc = TMath::Sqrt(error - val*val);
1006 }
1007 else
1008 {
1009 rc = TH1D::GetBinError(bin);
1010 rc = Profile(rc);
1011 }
1012
1013 return rc;
1014}
1015
1016// ------------------------------------------------------------------------
1017//
1018// Call this function to add a MHCamera on top of the present contents.
1019// Type:
1020// 0) bin content
1021// 1) errors
1022// 2) rel. errors
1023//
1024void MHCamera::AddCamContent(const MHCamera &d, Int_t type)
1025{
1026 if (fNcells!=d.fNcells || IsFreezed())
1027 return;
1028
1029 // FIXME: Security check missing!
1030 for (Int_t idx=0; idx<fNcells-2; idx++)
1031 if (d.IsUsed(idx))
1032 SetUsed(idx);
1033
1034 switch (type)
1035 {
1036 case 1:
1037 for (Int_t idx=0; idx<fNcells-2; idx++)
1038 Fill(idx, d.GetBinError(idx+1));
1039 break;
1040 case 2:
1041 for (Int_t idx=0; idx<fNcells-2; idx++)
1042 if (d.GetBinContent(idx+1)!=0)
1043 Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1)));
1044 break;
1045 default:
1046 for (Int_t idx=0; idx<fNcells-2; idx++)
1047 Fill(idx, d.GetBinContent(idx+1));
1048 break;
1049 }
1050 fEntries++;
1051}
1052
1053// ------------------------------------------------------------------------
1054//
1055// Call this function to add a TArrayD on top of the present contents.
1056//
1057void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used)
1058{
1059 if (event.GetSize()!=fNcells-2 || IsFreezed())
1060 return;
1061
1062 if (used && used->GetSize()!=fNcells-2)
1063 return;
1064
1065 for (Int_t idx=0; idx<fNcells-2; idx++)
1066 {
1067 Fill(idx, const_cast<TArrayD&>(event)[idx]); // FIXME: Slow!
1068
1069 if (!used || (*const_cast<TArrayC*>(used))[idx])
1070 SetUsed(idx);
1071 }
1072 fEntries++;
1073}
1074
1075// ------------------------------------------------------------------------
1076//
1077// Call this function to add a MCamEvent on top of the present contents.
1078// 1 is added to each pixel if the contents of MCamEvent>threshold
1079//
1080void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type)
1081{
1082 if (fNcells<=1 || IsFreezed())
1083 return;
1084
1085 // FIXME: Security check missing!
1086 for (Int_t idx=0; idx<fNcells-2; idx++)
1087 {
1088 Double_t val=threshold;
1089 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1090 SetUsed(idx);
1091
1092 if (val>threshold)
1093 Fill(idx);
1094 }
1095 fEntries++;
1096}
1097
1098// ------------------------------------------------------------------------
1099//
1100// Call this function to add a TArrayD on top of the present contents.
1101// 1 is added to each pixel if the contents of MCamEvent>threshold
1102//
1103void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
1104{
1105 if (event.GetSize()!=fNcells-2 || IsFreezed())
1106 return;
1107
1108 for (Int_t idx=0; idx<fNcells-2; idx++)
1109 {
1110 if (const_cast<TArrayD&>(event)[idx]>threshold)
1111 Fill(idx);
1112
1113 if (!ispos || fArray[idx+1]>0)
1114 SetUsed(idx);
1115 }
1116 fEntries++;
1117}
1118
1119// ------------------------------------------------------------------------
1120//
1121// Fill the pixels with random contents.
1122//
1123void MHCamera::FillRandom()
1124{
1125 if (fNcells<=1 || IsFreezed())
1126 return;
1127
1128 Reset();
1129
1130 // FIXME: Security check missing!
1131 for (Int_t idx=0; idx<fNcells-2; idx++)
1132 {
1133 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
1134 SetUsed(idx);
1135 }
1136 fEntries=1;
1137}
1138
1139
1140// ------------------------------------------------------------------------
1141//
1142// The array must be in increasing order, eg: 2.5, 3.7, 4.9
1143// The values in each bin are replaced by the interval in which the value
1144// fits. In the example we have four intervals
1145// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
1146// accordingly.
1147//
1148void MHCamera::SetLevels(const TArrayF &arr)
1149{
1150 if (fNcells<=1)
1151 return;
1152
1153 for (Int_t i=0; i<fNcells-2; i++)
1154 {
1155 if (!IsUsed(i))
1156 continue;
1157
1158 Int_t j = arr.GetSize();
1159 while (j && fArray[i+1]<arr[j-1])
1160 j--;
1161
1162 fArray[i+1] = j;
1163 }
1164 SetMaximum(arr.GetSize());
1165 SetMinimum(0);
1166}
1167
1168// ------------------------------------------------------------------------
1169//
1170// Reset the all pixel colors to a default value
1171//
1172void MHCamera::Reset(Option_t *opt)
1173{
1174 if (fNcells<=1 || IsFreezed())
1175 return;
1176
1177 TH1::Reset(opt);
1178
1179 for (Int_t i=0; i<fNcells-2; i++)
1180 {
1181 fArray[i+1]=0;
1182 ResetUsed(i);
1183 }
1184 fArray[0] = 0;
1185 fArray[fNcells-1] = 0;
1186}
1187
1188// ------------------------------------------------------------------------
1189//
1190// Here we calculate the color index for the current value.
1191// The color index is defined with the class TStyle and the
1192// Color palette inside. We use the command gStyle->SetPalette(1,0)
1193// for the display. So we have to convert the value "wert" into
1194// a color index that fits the color palette.
1195// The range of the color palette is defined by the values fMinPhe
1196// and fMaxRange. Between this values we have 50 color index, starting
1197// with 0 up to 49.
1198//
1199Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
1200{
1201 if (TMath::IsNaN(val)) // FIXME: gLog!
1202 return 10;
1203
1204 //
1205 // first treat the over- and under-flows
1206 //
1207 const Int_t maxcolidx = kItemsLegend-1;
1208
1209 if (val >= max)
1210 return fColors[maxcolidx];
1211
1212 if (val <= min)
1213 return fColors[0];
1214
1215 //
1216 // calculate the color index
1217 //
1218 Float_t ratio;
1219 if (islog && min>0)
1220 ratio = log10(val/min) / log10(max/min);
1221 else
1222 ratio = (val-min) / (max-min);
1223
1224 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
1225 return fColors[colidx];
1226}
1227
1228TPaveStats *MHCamera::GetStatisticBox()
1229{
1230 TObject *obj = 0;
1231
1232 TIter Next(fFunctions);
1233 while ((obj = Next()))
1234 if (obj->InheritsFrom(TPaveStats::Class()))
1235 return static_cast<TPaveStats*>(obj);
1236
1237 return NULL;
1238}
1239
1240// ------------------------------------------------------------------------
1241//
1242// Change the text on the legend according to the range of the Display
1243//
1244void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
1245{
1246 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
1247
1248 TArrow arr;
1249 arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025);
1250 arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025);
1251
1252 TString text;
1253 text += (int)(range*.3);
1254 text += "mm";
1255
1256 TText newtxt2;
1257 newtxt2.SetTextSize(0.04);
1258 newtxt2.PaintText(-range*.85, -range*.85, text);
1259
1260 text = "";
1261 text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10;
1262 text += "\\circ";
1263 text = text.Strip(TString::kLeading);
1264
1265 TLatex latex;
1266 latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text);
1267
1268 if (TestBit(kNoLegend))
1269 return;
1270
1271 TPaveStats *stats = GetStatisticBox();
1272
1273 const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1);
1274 const Float_t H = (0.75-hndc)*range;
1275 const Float_t offset = hndc*range;
1276
1277 const Float_t h = 2./kItemsLegend;
1278 const Float_t w = range/sqrt((float)(fNcells-2));
1279
1280 TBox newbox;
1281 TText newtxt;
1282 newtxt.SetTextSize(0.03);
1283 newtxt.SetTextAlign(12);
1284#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
1285 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
1286 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
1287#endif
1288
1289 const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend;
1290 const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3));
1291 const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
1292
1293 for (Int_t i=0; i<kItemsLegend+1; i+=3)
1294 {
1295 Float_t val;
1296 if (islog && min>0)
1297 val = pow(10, step*i) * min;
1298 else
1299 val = min + step*i;
1300
1301 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
1302 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
1303 }
1304
1305 for (Int_t i=0; i<kItemsLegend; i++)
1306 {
1307 newbox.SetFillColor(fColors[i]);
1308 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
1309 }
1310}
1311
1312// ------------------------------------------------------------------------
1313//
1314// Save primitive as a C++ statement(s) on output stream out
1315//
1316void MHCamera::SavePrimitive(ofstream &out, Option_t *opt)
1317{
1318 gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
1319 /*
1320 if (!gROOT->ClassSaved(TCanvas::Class()))
1321 fDrawingPad->SavePrimitive(out, opt);
1322
1323 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
1324 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
1325 */
1326}
1327
1328// ------------------------------------------------------------------------
1329//
1330// compute the distance of a point (px,py) to the Camera
1331// this functions needed for graphical primitives, that
1332// means without this function you are not able to interact
1333// with the graphical primitive with the mouse!!!
1334//
1335// All calcutations are done in pixel coordinates
1336//
1337Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
1338{
1339 if (fNcells<=1)
1340 return 999999;
1341
1342 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
1343 if (box)
1344 {
1345 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
1346 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
1347 box->SetY1NDC(gStyle->GetStatY()-w);
1348 box->SetX2NDC(gStyle->GetStatX());
1349 box->SetY2NDC(gStyle->GetStatY());
1350 }
1351
1352 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1353 return TH1D::DistancetoPrimitive(px, py);
1354
1355 const Bool_t issame = TString(GetDrawOption()).Contains("same", TString::kIgnoreCase);
1356
1357 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
1358 const Float_t conv = !issame ||
1359 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
1360 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
1361
1362 if (GetPixelIndex(px, py, conv)>=0)
1363 return 0;
1364
1365 if (!box)
1366 return 999999;
1367
1368 const Int_t dist = box->DistancetoPrimitive(px, py);
1369 if (dist > TPad::GetMaxPickDistance())
1370 return 999999;
1371
1372 gPad->SetSelected(box);
1373 return dist;
1374}
1375
1376// ------------------------------------------------------------------------
1377//
1378//
1379Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py, Float_t conv) const
1380{
1381 if (fNcells<=1)
1382 return -1;
1383
1384 Int_t i;
1385 for (i=0; i<fNcells-2; i++)
1386 {
1387 MHexagon hex((*fGeomCam)[i]);
1388 if (hex.DistancetoPrimitive(px, py, conv)>0)
1389 continue;
1390
1391 return i;
1392 }
1393 return -1;
1394}
1395
1396// ------------------------------------------------------------------------
1397//
1398// Returns string containing info about the object at position (px,py).
1399// Returned string will be re-used (lock in MT environment).
1400//
1401char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
1402{
1403 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1404 return TH1D::GetObjectInfo(px, py);
1405
1406 static char info[128];
1407
1408 const Int_t idx=GetPixelIndex(px, py);
1409
1410 if (idx<0)
1411 return TObject::GetObjectInfo(px, py);
1412
1413 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
1414 return info;
1415}
1416
1417// ------------------------------------------------------------------------
1418//
1419// Add a MCamEvent which should be displayed when the user clicks on a
1420// pixel.
1421// Warning: The object MUST inherit from TObject AND MCamEvent
1422//
1423void MHCamera::AddNotify(TObject *obj)
1424{
1425 // Make sure, that the object derives from MCamEvent!
1426 MCamEvent *evt = dynamic_cast<MCamEvent*>(obj);
1427 if (!evt)
1428 {
1429 gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl;
1430 return;
1431 }
1432
1433 // Make sure, that it is deleted from the list too, if the obj is deleted
1434 obj->SetBit(kMustCleanup);
1435
1436 // Add object to list
1437 fNotify->Add(obj);
1438}
1439
1440// ------------------------------------------------------------------------
1441//
1442// Execute a mouse event on the camera
1443//
1444void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
1445{
1446 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1447 {
1448 TH1D::ExecuteEvent(event, px, py);
1449 return;
1450 }
1451 //if (event==kMouseMotion && fStatusBar)
1452 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
1453 if (event!=kButton1Down)
1454 return;
1455
1456 const Int_t idx = GetPixelIndex(px, py);
1457 if (idx<0)
1458 return;
1459
1460 gLog << all << GetTitle() << " <" << GetName() << ">" << endl;
1461 gLog << "Software Pixel Idx: " << idx << endl;
1462 gLog << "Hardware Pixel Id: " << idx+1 << endl;
1463 gLog << "Contents: " << GetBinContent(idx+1);
1464 if (GetBinError(idx+1)>0)
1465 gLog << " +/- " << GetBinError(idx+1);
1466 gLog << " <" << (IsUsed(idx)?"on":"off") << ">" << endl;
1467
1468 if (fNotify && fNotify->GetSize()>0)
1469 {
1470 // FIXME: Is there a simpler and more convinient way?
1471
1472 // The name which is created here depends on the instance of
1473 // MHCamera and on the pad on which it is drawn --> The name
1474 // is unique. For ExecuteEvent gPad is always correctly set.
1475 const TString name = Form("%p;%p;PixelContent", this, gPad);
1476
1477 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
1478 if (old)
1479 old->cd();
1480 else
1481 new TCanvas(name);
1482
1483 /*
1484 TIter Next(gPad->GetListOfPrimitives());
1485 TObject *o;
1486 while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl;
1487 */
1488
1489 // FIXME: Make sure, that the old histograms are really deleted.
1490 // Are they already deleted?
1491
1492 // The dynamic_cast is necessary here: We cannot use ForEach
1493 TIter Next(fNotify);
1494 MCamEvent *evt;
1495 while ((evt=dynamic_cast<MCamEvent*>(Next())))
1496 evt->DrawPixelContent(idx);
1497
1498 gPad->Modified();
1499 gPad->Update();
1500 }
1501}
1502
1503UInt_t MHCamera::GetNumPixels() const
1504{
1505 return fGeomCam->GetNumPixels();
1506}
1507
1508TH1 *MHCamera::DrawCopy() const
1509{
1510 gPad=NULL;
1511 return TH1D::DrawCopy(fName+";cpy");
1512}
Note: See TracBrowser for help on using the repository browser.