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

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