source: trunk/MagicSoft/Mars/mgui/MCamDisplay.cc@ 2178

Last change on this file since 2178 was 2178, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 20.2 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 05/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Harald Kornmayer, 1/2001
20!
21! Copyright: MAGIC Software Development, 2000-2003
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MCamDisplay
29//
30// Camera Display. The Pixels are displayed in
31// contents/area [somthing/mm^2]
32//
33// To change the scale to a logarithmic scale SetLogz() of the Pad.
34//
35//
36////////////////////////////////////////////////////////////////////////////
37#include "MCamDisplay.h"
38
39#include <fstream>
40#include <iostream>
41
42#include <TBox.h>
43#include <TArrow.h>
44#include <TLatex.h>
45#include <TStyle.h>
46#include <TMarker.h>
47#include <TCanvas.h>
48#include <TArrayF.h>
49#include <TClonesArray.h>
50
51#include "MH.h"
52#include "MHexagon.h"
53
54#include "MGeomCam.h"
55
56#include "MRflEvtData.h"
57
58#include "MCerPhotPix.h"
59#include "MCerPhotEvt.h"
60
61#include "MPedestalPix.h"
62#include "MPedestalCam.h"
63
64#include "MCurrents.h"
65#include "MCamEvent.h"
66
67#include "MImgCleanStd.h"
68
69#define kItemsLegend 48 // see SetPalette(1,0)
70
71ClassImp(MCamDisplay);
72
73using namespace std;
74
75// ------------------------------------------------------------------------
76//
77// Default Constructor. To be used by the root system ONLY.
78//
79MCamDisplay::MCamDisplay() : fGeomCam(NULL), fAutoScale(kTRUE)
80{
81 fNumPixels = 0;
82 fRange = 0;
83
84 fPixels = NULL;
85 fLegend = NULL;
86 fLegText = NULL;
87 fArrowX = NULL;
88 fArrowY = NULL;
89 fLegRadius = NULL;
90 fLegDegree = NULL;
91
92 fMinimum = 0;
93 fMaximum = 1;
94
95 fNotify = new TList;
96}
97
98// ------------------------------------------------------------------------
99//
100// Constructor. Makes a clone of MGeomCam.
101//
102MCamDisplay::MCamDisplay(MGeomCam *geom)
103 : fGeomCam(NULL), fAutoScale(kTRUE), fColors(kItemsLegend), fData(geom->GetNumPixels()), fMinimum(0), fMaximum(1)
104{
105 fGeomCam = (MGeomCam*)geom->Clone();
106
107 fNotify = new TList;
108
109 //
110 // create the hexagons of the display
111 //
112 fNumPixels = fGeomCam->GetNumPixels();
113 fRange = fGeomCam->GetMaxRadius();
114
115 // root 3.02
116 // * base/inc/TObject.h:
117 // register BIT(8) as kNoContextMenu. If an object has this bit set it will
118 // not get an automatic context menu when clicked with the right mouse button.
119
120 fPhotons = new TClonesArray("TMarker", 0);
121
122 //
123 // Construct all hexagons. Use new-operator with placement
124 //
125 fPixels = new TClonesArray("MHexagon", fNumPixels);
126 for (UInt_t i=0; i<fNumPixels; i++)
127 {
128 MHexagon &pix = *new ((*fPixels)[i]) MHexagon((*fGeomCam)[i]);
129#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
130 pix.SetBit(/*kNoContextMenu|*/kCannotPick);
131#endif
132 pix.SetFillColor(16);
133 pix.ResetBit(kIsUsed);
134 }
135
136 //
137 // set up the Legend
138 //
139 const Float_t H = 0.9*fRange;
140 const Float_t h = 2./kItemsLegend;
141
142 const Float_t w = fRange/sqrt((float)fNumPixels);
143
144 fLegend = new TClonesArray("TBox", kItemsLegend);
145 fLegText = new TClonesArray("TText", kItemsLegend+1);
146
147 for (Int_t i=0; i<kItemsLegend; i++)
148 {
149 TBox &newbox = *new ((*fLegend)[i]) TBox;
150 newbox.SetX1(fRange);
151 newbox.SetX2(fRange+w);
152 newbox.SetY1(H*( i *h - 1.));
153 newbox.SetY2(H*((i+1)*h - 1.));
154 newbox.SetFillColor(16);
155#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
156 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
157#endif
158 }
159
160 for (Int_t i=0; i<kItemsLegend+1; i++)
161 {
162 TText &newtxt = *new ((*fLegText)[i]) TText;
163 newtxt.SetTextSize(0.025);
164 newtxt.SetTextAlign(12);
165 newtxt.SetX(fRange+1.5*w);
166 newtxt.SetY(H*((i+0.5)*h - 1.));
167#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
168 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
169#endif
170 }
171
172 fArrowX = new TArrow(-fRange*.9, -fRange*.9, -fRange*.6, -fRange*.9, 0.025);
173 fArrowY = new TArrow(-fRange*.9, -fRange*.9, -fRange*.9, -fRange*.6, 0.025);
174
175 TString text;
176 text += (int)(fRange*.3);
177 text += "mm";
178
179 fLegRadius = new TText(-fRange*.85, -fRange*.85, text);
180 text = "";
181 text += (float)((int)(fRange*.3*fGeomCam->GetConvMm2Deg()*10))/10;
182 text += "\\circ";
183 text = text.Strip(TString::kLeading);
184 fLegDegree = new TLatex(-fRange*.85, -fRange*.75, text);
185 fLegRadius->SetTextSize(0.04);
186 fLegDegree->SetTextSize(0.04);
187
188#if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06)
189 SetPalette(1, 0);
190#else
191 SetPalette(51, 0);
192#endif
193}
194
195// ------------------------------------------------------------------------
196//
197// Destructor. Deletes TClonesArrays for hexagons and legend elements.
198//
199MCamDisplay::~MCamDisplay()
200{
201 fPixels->Delete();
202 fLegend->Delete();
203 fLegText->Delete();
204 fPhotons->Delete();
205
206 delete fPixels;
207 delete fLegend;
208 delete fLegText;
209 delete fPhotons;
210
211 delete fArrowX;
212 delete fArrowY;
213
214 delete fLegRadius;
215 delete fLegDegree;
216
217 delete fGeomCam;
218
219 delete fNotify;
220}
221
222// ------------------------------------------------------------------------
223//
224// Call this function to draw the camera layout into your canvas.
225// Setup a drawing canvas. Add this object and all child objects
226// (hexagons, etc) to the current pad. If no pad exists a new one is
227// created.
228//
229void MCamDisplay::Draw(Option_t *option)
230{
231 // root 3.02:
232 // gPad->SetFixedAspectRatio()
233
234 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600);
235 pad->SetBorderMode(0);
236 pad->SetFillColor(16);
237 //pad->Modified();
238 /*
239 pad->Divide(1,1);
240 pad->cd(1);
241 gPad->SetBorderMode(0);
242 gPad->SetFillColor(16);
243*/
244 AppendPad("");
245}
246
247
248void MCamDisplay::SetRange()
249{
250 //
251 // Maintain aspect ratio
252 //
253 const float ratio = 1.15;
254
255 //
256 // Calculate width and height of the current pad in pixels
257 //
258 Float_t w = gPad->GetWw();
259 Float_t h = gPad->GetWh()*ratio;
260
261 //
262 // This prevents the pad from resizing itself wrongly
263 //
264 if (gPad->GetMother() != gPad)
265 {
266 w *= gPad->GetMother()->GetAbsWNDC();
267 h *= gPad->GetMother()->GetAbsHNDC();
268 }
269
270 //
271 // Set Range (coordinate system) of pad
272 //
273 gPad->Range(-fRange, -fRange, (2*ratio-1)*fRange, fRange);
274
275 //
276 // Resize Pad to given ratio
277 //
278 if (h<w)
279 gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1);
280 else
281 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2);
282}
283
284void MCamDisplay::Update(Bool_t islog)
285{
286 // FIXME: Don't do this if nothing changed!
287 if (fAutoScale)
288 {
289 fMinimum = FLT_MAX;
290 fMaximum = -FLT_MAX;
291
292 for (UInt_t i=0; i<fNumPixels; i++)
293 {
294 if (!(*this)[i].TestBit(kIsUsed))
295 continue;
296
297 if (fData[i]<fMinimum)
298 fMinimum = fData[i];
299 if (fData[i]>fMaximum)
300 fMaximum = fData[i];
301 }
302 }
303
304 if (fMinimum==FLT_MAX && fMaximum==-FLT_MAX)
305 fMinimum = fMaximum = 0;
306 if (fMinimum==fMaximum)
307 fMaximum = fMinimum + 1;
308
309 UpdateLegend(fMinimum, fMaximum, islog);
310
311 for (UInt_t i=0; i<fNumPixels; i++)
312 {
313 if ((*this)[i].TestBit(kIsUsed))
314 (*this)[i].SetFillColor(GetColor(fData[i], fMinimum, fMaximum, islog));
315 else
316 (*this)[i].SetFillColor(10);
317 }
318}
319
320// ------------------------------------------------------------------------
321//
322// This is called at any time the canvas should get repainted.
323// Here we maintain an aspect ratio of 1.15. This makes sure,
324// that the camera image doesn't get distorted by resizing the canvas.
325//
326void MCamDisplay::Paint(Option_t *opt)
327{
328 if (!fPixels)
329 return;
330
331 // Maintain aspect ratio
332 SetRange();
333
334 // Maintain colors
335 SetPalette();
336
337 // Update Contents of the pixels
338 Update(gPad->GetLogz());
339
340 // Paint Legend
341 fArrowX->Paint(">");
342 fArrowY->Paint(">");
343
344 fLegRadius->Paint();
345 fLegDegree->Paint();
346
347 // Paint primitives (pixels, color legend, photons, ...)
348 { fPixels->ForEach( TObject, Paint)(); }
349 { fLegend->ForEach( TObject, Paint)(); }
350 { fLegText->ForEach(TObject, Paint)(); }
351 { fPhotons->ForEach(TObject, Paint)(); }
352}
353
354// ------------------------------------------------------------------------
355//
356// With this function you can change the color palette. For more
357// information see TStyle::SetPalette. Only palettes with 50 colors
358// are allowed.
359// In addition you can use SetPalette(52, 0) to create an inverse
360// deep blue sea palette
361//
362void MCamDisplay::SetPalette(Int_t ncolors, Int_t *colors)
363{
364 //
365 // If not enough colors are specified skip this.
366 //
367 if (ncolors>1 && ncolors<50)
368 {
369 cout << "MCamDisplay::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl;
370 return;
371 }
372
373 //
374 // If ncolors==52 create a reversed deep blue sea palette
375 //
376 if (ncolors==52)
377 {
378 gStyle->SetPalette(51, NULL);
379 TArrayI c(kItemsLegend);
380 for (int i=0; i<kItemsLegend; i++)
381 c[kItemsLegend-i-1] = gStyle->GetColorPalette(i);
382 gStyle->SetPalette(kItemsLegend, c.GetArray());
383 }
384 else
385 gStyle->SetPalette(ncolors, colors);
386
387 //
388 // Change the colors of the pixels
389 //
390 for (unsigned int i=0; i<fNumPixels; i++)
391 {
392 //
393 // Get the old color index and check whether it is
394 // background or transparent
395 //
396 Int_t col = (*this)[i].GetFillColor();
397 if (col==10 || col==16)
398 continue;
399
400 //
401 // Search for the color index (level) in the color table
402 //
403 int idx;
404 for (idx=0; idx<kItemsLegend; idx++)
405 if (col==fColors[idx])
406 break;
407
408 //
409 // Should not happen
410 //
411 if (idx==kItemsLegend)
412 {
413 cout << "MCamDisplay::SetPalette: Strange... FIXME!" << endl;
414 continue;
415 }
416
417 //
418 // Convert the old color index (level) into the new one
419 //
420 (*this)[i].SetFillColor(gStyle->GetColorPalette(idx));
421 }
422
423 //
424 // Store the color palette used for a leter reverse lookup
425 //
426 for (int i=0; i<kItemsLegend; i++)
427 {
428 fColors[i] = gStyle->GetColorPalette(i);
429 GetBox(i)->SetFillColor(fColors[i]);
430 }
431}
432
433void MCamDisplay::SetPrettyPalette()
434{
435 SetPalette(1, 0);
436}
437
438void MCamDisplay::SetDeepBlueSeaPalette()
439{
440 SetPalette(51, 0);
441}
442
443void MCamDisplay::SetInvDeepBlueSeaPalette()
444{
445 SetPalette(52, 0);
446}
447
448void MCamDisplay::SetPalette()
449{
450 for (int i=0; i<kItemsLegend; i++)
451 GetBox(i)->SetFillColor(fColors[i]);
452}
453
454void MCamDisplay::DrawPixelNumbers()
455{
456 for (int i=0; i<kItemsLegend; i++)
457 fColors[i] = 16;
458
459 if (!gPad)
460 Draw();
461
462 TText txt;
463 txt.SetTextFont(122);
464 txt.SetTextAlign(22); // centered/centered
465
466 for (UInt_t i=0; i<fNumPixels; i++)
467 {
468 TString num;
469 num += i;
470
471 const MHexagon &h = (*this)[i];
472 TText *nt = txt.DrawText(h.GetX(), h.GetY(), num);
473 nt->SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius());
474 }
475}
476
477// ------------------------------------------------------------------------
478//
479// Call this function to fill the currents
480//
481void MCamDisplay::Fill(const MCamEvent &event, Int_t type)
482{
483 Reset();
484
485 // FIXME: Security check missing!
486 for (UInt_t idx=0; idx<fNumPixels; idx++)
487 {
488 fData[idx] = 0;
489 if (event.GetPixelContent(fData[idx], idx, fGeomCam->GetPixRatio(idx), type))
490 (*this)[idx].SetBit(kIsUsed);
491 }
492}
493
494// ------------------------------------------------------------------------
495//
496// Call this function to fill the currents
497//
498void MCamDisplay::Fill(const TArrayF &event, Bool_t ispos)
499{
500 Reset();
501
502 fData = event;
503
504 for (UInt_t idx=0; idx<fNumPixels; idx++)
505 if (!ispos || fData[idx]>0)
506 (*this)[idx].SetBit(kIsUsed);
507}
508
509// ------------------------------------------------------------------------
510//
511// Fill the colors in respect to the cleaning levels
512//
513void MCamDisplay::FillLevels(const MCerPhotEvt &event, Float_t lvl1, Float_t lvl2)
514{
515 Fill(event, 2);
516
517 for (UInt_t i=0; i<fNumPixels; i++)
518 {
519 if (!(*this)[i].TestBit(kIsUsed))
520 continue;
521
522 if (fData[i]>lvl1)
523 fData[i] = 0;
524 else
525 if (fData[i]>lvl2)
526 fData[i] = 1;
527 else
528 fData[i] = 2;
529 }
530}
531
532// ------------------------------------------------------------------------
533//
534// Fill the colors in respect to the cleaning levels
535//
536void MCamDisplay::FillLevels(const MCerPhotEvt &event, const MImgCleanStd &clean)
537{
538 FillLevels(event, clean.GetCleanLvl1(), clean.GetCleanLvl2());
539}
540
541// ------------------------------------------------------------------------
542//
543// Show a reflector event. EMarkerStyle is defined in root/include/Gtypes.h
544// To remove the photons from the display call FillRflEvent(NULL)
545//
546void MCamDisplay::ShowRflEvent(const MRflEvtData *event, EMarkerStyle ms)
547{
548 const Int_t num = event ? event->GetNumPhotons() : 0;
549
550 fPhotons->ExpandCreate(num);
551 if (num < 1)
552 return;
553
554 Int_t i=num-1;
555 do
556 {
557 const MRflSinglePhoton &ph = event->GetPhoton(i);
558 TMarker &m = (TMarker&)*fPhotons->UncheckedAt(i);
559 m.SetX(ph.GetX());
560 m.SetY(ph.GetY());
561 m.SetMarkerStyle(ms);
562 } while (i--);
563}
564
565// ------------------------------------------------------------------------
566//
567// Fill a reflector event. Sums all pixels in each pixel as the
568// pixel contents.
569//
570// WARNING: Due to the estimation in DistanceToPrimitive and the
571// calculation in pixels instead of x, y this is only a
572// rough estimate.
573//
574void MCamDisplay::FillRflEvent(const MRflEvtData &event)
575{
576 //
577 // reset pixel colors to background color
578 //
579 Reset();
580
581 //
582 // sum the photons content in each pixel
583 //
584 const Int_t entries = event.GetNumPhotons();
585
586 for (Int_t i=0; i<entries; i++)
587 {
588 const MRflSinglePhoton &ph = event.GetPhoton(i);
589
590 UInt_t id;
591 for (id=0; id<fNumPixels; id++)
592 {
593 if ((*this)[id].DistanceToPrimitive(ph.GetX(), ph.GetY())<0)
594 break;
595 }
596 if (id==fNumPixels)
597 continue;
598
599 fData[id] += fGeomCam->GetPixRatio(id);
600 }
601
602 //
603 // Set color of pixels
604 //
605 for (UInt_t id=0; id<fNumPixels; id++)
606 if (fData[id]>0)
607 (*this)[id].SetBit(kIsUsed);
608}
609
610// ------------------------------------------------------------------------
611//
612// Reset the all pixel colors to a default value
613//
614void MCamDisplay::Reset()
615{
616 for (UInt_t i=0; i<fNumPixels; i++)
617 (*this)[i].ResetBit(kIsUsed);
618
619 if (fAutoScale)
620 {
621 fMinimum = 0;
622 fMaximum = 0;
623 }
624}
625
626// ------------------------------------------------------------------------
627//
628// Here we calculate the color index for the current value.
629// The color index is defined with the class TStyle and the
630// Color palette inside. We use the command gStyle->SetPalette(1,0)
631// for the display. So we have to convert the value "wert" into
632// a color index that fits the color palette.
633// The range of the color palette is defined by the values fMinPhe
634// and fMaxRange. Between this values we have 50 color index, starting
635// with 0 up to 49.
636//
637Int_t MCamDisplay::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
638{
639 //
640 // first treat the over- and under-flows
641 //
642 const Int_t maxcolidx = kItemsLegend-1;
643
644 if (val >= max)
645 return fColors[maxcolidx];
646
647 if (val <= min)
648 return fColors[0];
649
650 //
651 // calculate the color index
652 //
653 Float_t ratio;
654 if (islog && min>0)
655 ratio = log10(val/min) / log10(max/min);
656 else
657 ratio = (val-min) / (max-min);
658 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
659 return fColors[colidx];
660}
661
662// ------------------------------------------------------------------------
663//
664// Change the text on the legend according to the range of the Display
665//
666void MCamDisplay::UpdateLegend(Float_t minphe, Float_t maxphe, Bool_t islog)
667{
668 for (Int_t i=0; i<kItemsLegend+1; i+=3)
669 {
670 const Float_t pos = (Float_t)i/kItemsLegend;
671
672 Float_t val;
673 if (islog && minphe>0)
674 val = pow(10, log10(maxphe/minphe)*pos) * minphe;
675 else
676 val = minphe + pos * (maxphe-minphe);
677
678 TText &txt = *(TText*)fLegText->At(i);
679 txt.SetText(txt.GetX(), txt.GetY(), Form(val<1e6?"%5.1f":"%5.1e", val));
680 }
681}
682
683// ------------------------------------------------------------------------
684//
685// Save primitive as a C++ statement(s) on output stream out
686//
687void MCamDisplay::SavePrimitive(ofstream &out, Option_t *opt)
688{
689 cout << "MCamDisplay::SavePrimitive: Must be rewritten!" << endl;
690 /*
691 if (!gROOT->ClassSaved(TCanvas::Class()))
692 fDrawingPad->SavePrimitive(out, opt);
693
694 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
695 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
696 */
697}
698
699// ------------------------------------------------------------------------
700//
701// compute the distance of a point (px,py) to the Camera
702// this functions needed for graphical primitives, that
703// means without this function you are not able to interact
704// with the graphical primitive with the mouse!!!
705//
706// All calcutations are running in pixel coordinates
707//
708Int_t MCamDisplay::DistancetoPrimitive(Int_t px, Int_t py)
709{
710 Int_t dist = 999999;
711
712 for (UInt_t i=0; i<fNumPixels; i++)
713 {
714 Int_t d = (*fPixels)[i]->DistancetoPrimitive(px, py);
715
716 if (d<dist)
717 dist=d;
718 }
719 return dist==0?0:999999;
720}
721
722// ------------------------------------------------------------------------
723//
724// Execute a mouse event on the camera
725//
726/*
727 void MCamDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
728 {
729 cout << "Execute Event Camera " << event << " @ " << px << " " << py << endl;
730 }
731 */
732
733
734// ------------------------------------------------------------------------
735//
736// Function introduced (31-01-03) WILL BE REMOVED IN THE FUTURE! DON'T
737// USE IT!
738//
739void MCamDisplay::SetPix(const Int_t pixnum, const Int_t color, Float_t min, Float_t max)
740{
741 fData[pixnum] = color;
742 (*this)[pixnum].SetBit(kIsUsed);
743 (*this)[pixnum].SetFillColor(GetColor(color, min, max, 0));
744}
745
746Int_t MCamDisplay::GetPixelIndex(Int_t px, Int_t py) const
747{
748 UInt_t i;
749 for (i=0; i<fNumPixels; i++)
750 {
751 if ((*fPixels)[i]->DistancetoPrimitive(px, py)>0)
752 continue;
753
754 return i;
755 }
756 return -1;
757}
758
759// ------------------------------------------------------------------------
760//
761// Returns string containing info about the object at position (px,py).
762// Returned string will be re-used (lock in MT environment).
763//
764char *MCamDisplay::GetObjectInfo(Int_t px, Int_t py) const
765{
766 static char info[128];
767
768 const Int_t idx=GetPixelIndex(px, py);
769
770 if (idx<0)
771 return TObject::GetObjectInfo(px, py);
772
773 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
774 return info;
775}
776
777// ------------------------------------------------------------------------
778//
779// Execute a mouse event on the camera
780//
781void MCamDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
782{
783 //if (event==kMouseMotion && fStatusBar)
784 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
785 if (event!=kButton1Down)
786 return;
787
788 const Int_t idx = GetPixelIndex(px, py);
789 if (idx<0)
790 return;
791
792 cout << "Software Pixel Index: " << idx << endl;
793 cout << "Hardware Pixel Id: " << idx+1 << endl;
794 cout << "Contents: " << fData[idx] << endl;
795
796 //fNotify->Print();
797 if (fNotify->GetSize()>0)
798 new TCanvas;
799 fNotify->ForEach(MCamEvent, DrawPixelContent)(idx);
800}
Note: See TracBrowser for help on using the repository browser.