source: trunk/MagicSoft/Mars/mhist/MHCamera.cc@ 3528

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