source: trunk/MagicSoft/Mars/mmain/MEventDisplay.cc@ 5198

Last change on this file since 5198 was 5145, checked in by tbretz, 20 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, 10/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24#include "MEventDisplay.h"
25
26//
27// C-lib
28//
29#include <stdlib.h> // atoi
30
31//
32// root
33//
34#include <TFile.h> // TFile
35#include <TTree.h> // TTree
36#include <TList.h> // TList::Add
37#include <TStyle.h> // gStyle->SetOptStat
38#include <TCanvas.h> // TCanvas::cd
39#include <TSystem.h> // TSystem::BaseName
40
41//
42// root GUI
43//
44#include <TGLabel.h> // TGLabel
45#include <TGButton.h> // TGPictureButton
46#include <TG3DLine.h> // TGHorizontal3DLine
47#include <TGTextEntry.h> // TGTextEntry
48#include <TGButtonGroup.h> // TGVButtonGroup
49#include <TRootEmbeddedCanvas.h> // TRootEmbeddedCanvas
50
51//
52// General
53//
54#include "MGList.h" // MGList
55
56#include "MLog.h"
57#include "MLogManip.h"
58
59#include "MParList.h" // MParList::AddToList
60#include "MEvtLoop.h" // MEvtLoop::GetParList
61#include "MTaskList.h" // MTaskList::AddToList
62
63//
64// Tasks
65//
66#include "MReadMarsFile.h" // MReadMarsFile
67#include "MGeomApply.h" // MGeomApply
68#include "MFDataMember.h" // MFDataMember
69#include "MMcPedestalCopy.h" // MMcPedestalCopy
70#include "MMcPedestalNSBAdd.h" // MMcPedestalNSBAdd
71
72#include "MCerPhotAnal2.h" // MCerPhotAnal2
73#include "MImgCleanStd.h" // MImgCleanStd
74#include "MHillasCalc.h" // MHillasCalc
75//#include "MHillasSrcCalc.h" // MHillasSrcCalc
76//#include "MBlindPixelCalc.h" // MBlindPixelCalc
77#include "MArrivalTimeCalc.h" // MArrivalTimeCalc
78#include "MFillH.h" // MFillH
79#include "MExtractSignal.h" // MExtractsignal
80#include "MMcCalibrationUpdate.h" // MMcCalibrationUpdate
81#include "MCalibrate.h" // MCalibrate
82#include "MMcTriggerLvl2Calc.h" // MMcTriggerLvl2Calc
83
84//
85// Container
86//
87#include "MHillas.h" // MHillas::Print(const MGeomCam&)
88#include "MHillasExt.h" // MHillasExt::Print(const MGeomCam&)
89#include "MHillasSrc.h" // MHillasSrc::Print(const MGeomCam&)
90#include "MImagePar.h" // MImagePar::Print(const MGeomCam&)
91#include "MNewImagePar.h" // MNewImagePar::Print(const MGeomCam&)
92#include "MHEvent.h" // MHEvent
93#include "MHCamera.h" // MHCamera
94#include "MRawEvtData.h" // MRawEvtData
95#include "MArrivalTimeCam.h" // MArrivalTimeCam
96#include "MBadPixelsCam.h" // MBadPixelsCam
97#include "MPedPhotCam.h" // MPedPhotCam
98#include "MCalibrationChargeCam.h" // MCalibrationChargeCam
99#include "MMcTriggerLvl2.h" // MMcTriggerLvl2
100
101using namespace std;
102
103ClassImp(MEventDisplay);
104
105// --------------------------------------------------------------------------
106//
107// Constructor.
108//
109MEventDisplay::MEventDisplay(const char *fname, const char *pname, const char *cname) : MStatusDisplay()
110{
111 //
112 // Setup Task list for hillas calculation
113 //
114 SetupTaskList("Events", fname, pname, cname);
115
116 //
117 // Add MEventDisplay GUI elements to the display
118 //
119 AddUserFrame(fname);
120
121 //
122 // Show new part of the window, resize to correct aspect ratio
123 //
124 // FIXME: This should be done by MStatusDisplay automatically
125 Resize(GetWidth(), GetHeight() + fUserFrame->GetDefaultHeight());
126 SetWindowName("Event Display");
127 MapSubwindows();
128
129 //
130 // Readin first event and display it
131 //
132 ReadFirstEvent();
133
134 SetNoContextMenu(kFALSE);
135}
136
137// --------------------------------------------------------------------------
138//
139// Destructor: PostProcess eventloop, delete eventloop with all containers
140//
141MEventDisplay::~MEventDisplay()
142{
143 fEvtLoop->PostProcess();
144 delete fEvtLoop;
145}
146
147Int_t MEventDisplay::GetFileType(const char *tree, const char *fname) const
148{
149 TFile f(fname, "READ");
150 TTree *t = (TTree*)f.Get(tree);
151 if (!t)
152 return -1;
153
154 if (t->FindBranch("MCerPhotEvt."))
155 return 1;
156
157 if (t->FindBranch("MRawEvtData."))
158 return 2;
159
160 return -2;
161}
162
163// --------------------------------------------------------------------------
164//
165// Setup Task and parameter list for hillas calculation,
166// preprocess tasks and read in first event (process)
167//
168void MEventDisplay::SetupTaskList(const char *tname, const char *fname,
169 const char *pname, const char *cname)
170{
171 MCalibrationChargeCam *ccam=0;
172 MPedPhotCam *pcam=0;
173
174 MBadPixelsCam *badpix = new MBadPixelsCam;
175
176 const Int_t type = GetFileType(tname, fname);
177 switch (type)
178 {
179 case 1: gLog << all << "Calibrated Data File detected..." << endl; break;
180 case 2: gLog << all << "Raw Data File detected..." << endl; break;
181 case -2: gLog << err << "File type unknown... abort." << endl; return;
182 case -1: gLog << err << "Tree " << tname << " not found... abort." << endl; return;
183 }
184
185 if (type==2)
186 {
187 TFile file(pname, "READ");
188 if (!file.IsZombie())
189 pcam = new MPedPhotCam;
190 if (pcam)
191 {
192 if (pcam->Read()<=0)
193 {
194 delete pcam;
195 pcam = NULL;
196 }
197
198 if (file.FindKey("MBadPixelsCam"))
199 {
200 MBadPixelsCam bad;
201 if (bad.Read()>0)
202 badpix->Merge(bad);
203 }
204 }
205 file.Close();
206 file.Open(cname, "READ");
207 if (!file.IsZombie())
208 ccam = new MCalibrationChargeCam;
209 if (ccam)
210 {
211 if (ccam->Read()<=0)
212 {
213 delete ccam;
214 ccam = NULL;
215 }
216
217 if (file.FindKey("MBadPixelsCam"))
218 {
219 MBadPixelsCam bad;
220 if (bad.Read()>0)
221 badpix->Merge(bad);
222 }
223 }
224 file.Close();
225 }
226
227 //
228 // Setup an empty job, with a reader task only.
229 // All tasks and parameter containers are deleted automatically
230 // (via SetOwner())
231 //
232 MTaskList *tlist = new MTaskList;
233 tlist->SetOwner();
234
235 MReadMarsFile *read = new MReadMarsFile(tname, fname);
236 read->DisableAutoScheme();
237 tlist->AddToList(read);
238
239 MGeomApply *apl = new MGeomApply;
240 tlist->AddToList(apl);
241
242 MParList *plist = new MParList;
243 plist->SetOwner();
244 plist->AddToList(tlist);
245 plist->AddToList(badpix);
246 if (pcam)
247 plist->AddToList(pcam);
248 if (ccam)
249 plist->AddToList(ccam);
250
251 MHEvent *evt01 = new MHEvent(MHEvent::kEvtSignalRaw);
252 MHEvent *evt02 = new MHEvent(MHEvent::kEvtSignalRaw);
253 MHEvent *evt03 = new MHEvent(MHEvent::kEvtPedestal);
254 MHEvent *evt04 = new MHEvent(MHEvent::kEvtPedestalRMS);
255 MHEvent *evt06a= new MHEvent(MHEvent::kEvtCleaningData);
256 MHEvent *evt06b= new MHEvent(MHEvent::kEvtCleaningLevels);
257 MHEvent *evt07 = new MHEvent(MHEvent::kEvtIdxMax);
258 MHEvent *evt08 = new MHEvent(MHEvent::kEvtArrTime);
259 MHEvent *evt09 = new MHEvent(MHEvent::kEvtTrigPix);
260 MHEvent *evt10 = new MHEvent(MHEvent::kEvtIslandIndex);
261
262 evt01->SetName("Signal");
263 evt02->SetName("Cleaned");
264 evt03->SetName("Pedestal");
265 evt04->SetName("PedRMS");
266 evt06a->SetName("CleanData");
267 evt06b->SetName("CleanLevels");
268 evt07->SetName("Max Slice Idx");
269 evt08->SetName("Arrival Time");
270 evt09->SetName("Trigger");
271 evt10->SetName("Islands");
272
273 // This makes sure, that the containers are deleted...
274 plist->AddToList(evt01);
275 plist->AddToList(evt02);
276 plist->AddToList(evt03);
277 plist->AddToList(evt04);
278 plist->AddToList(evt06a);
279 plist->AddToList(evt06b);
280 plist->AddToList(evt07);
281 plist->AddToList(evt08);
282 plist->AddToList(evt09);
283 plist->AddToList(evt10);
284
285 MCerPhotAnal2 *nanal = new MCerPhotAnal2;
286 MFillH *fill01 = new MFillH(evt01, "MCerPhotEvt", "MFillH01");
287 MImgCleanStd *clean = new MImgCleanStd;
288 MFillH *fill02 = new MFillH(evt02, "MCerPhotEvt", "MFillH02");
289 MFillH *fill03 = new MFillH(evt03, "MPedPhotCam", "MFillH03");
290 MFillH *fill04 = new MFillH(evt04, "MPedPhotCam", "MFillH04");
291 MFillH *fill06a= new MFillH(evt06a, "MCameraData", "MFillH06a");
292 MFillH *fill06b= new MFillH(evt06b, "MCameraData", "MFillH06b");
293 //MBlindPixelCalc *blind = new MBlindPixelCalc;
294 MHillasCalc *hcalc = new MHillasCalc;
295 MMcTriggerLvl2Calc *trcal = new MMcTriggerLvl2Calc;
296 MFillH *fill09 = new MFillH(evt09, "MMcTriggerLvl2", "MFillH09");
297 MFillH *fill10 = new MFillH(evt10, "MCerPhotEvt", "MFillH10");
298
299 // If no pedestal or no calibration file is availble
300 if (type==2)
301 {
302 if (!pcam || !ccam)
303 {
304 MFilter *f1=new MFDataMember("MRawRunHeader.fRunType", '>', 255.5);
305 MFilter *f2=new MFDataMember("MRawRunHeader.fRunType", '<', 255.5);
306 f1->SetName("MFMonteCarlo");
307 f2->SetName("MFRealData");
308
309 MMcPedestalCopy *pcopy = new MMcPedestalCopy;
310 MMcPedestalNSBAdd *pdnsb = new MMcPedestalNSBAdd;
311
312 MExtractSignal* extra = new MExtractSignal();
313 extra->SetRange(0, 14, 0, 14);
314 extra->SetSaturationLimit(240);
315
316 MMcCalibrationUpdate* mcupd = new MMcCalibrationUpdate;
317 mcupd->SetOuterPixelsGainScaling(kFALSE);
318
319 MCalibrate* mccal = new MCalibrate;
320 mccal->SetCalibrationMode(MCalibrate::kFfactor);
321
322 // MC
323 extra->SetFilter(f1);
324 mcupd->SetFilter(f1);
325 mccal->SetFilter(f1);
326 trcal->SetFilter(f1);
327 //fill09->SetFilter(f1);
328
329 // Data
330 nanal->SetFilter(f2);
331
332 // TaskList
333 tlist->AddToList(f1);
334 tlist->AddToList(f2);
335 tlist->AddToList(pcopy);
336 tlist->AddToList(pdnsb);
337
338 tlist->AddToList(extra);
339 tlist->AddToList(mcupd);
340 tlist->AddToList(mccal);
341
342
343 tlist->AddToList(nanal);
344 }
345 else
346 {
347 MCalibrate* calib = new MCalibrate;
348 tlist->AddToList(calib);
349 }
350 }
351
352 tlist->AddToList(fill01);
353 tlist->AddToList(clean);
354 tlist->AddToList(fill02);
355 tlist->AddToList(fill03);
356 tlist->AddToList(fill04);
357 tlist->AddToList(fill06a);
358 tlist->AddToList(fill06b);
359 //tlist->AddToList(blind);
360 tlist->AddToList(fill10);
361 tlist->AddToList(hcalc);
362 if ((!pcam || !ccam) && type==2)
363 {
364 MArrivalTimeCalc *tcalc = new MArrivalTimeCalc;
365 MFillH *fill07 = new MFillH(evt07, "MRawEvtData", "MFillH7");
366 MFillH *fill08 = new MFillH(evt08, "MArrivalTimeCam", "MFillH8");
367 tlist->AddToList(tcalc);
368 tlist->AddToList(fill07);
369 tlist->AddToList(fill08);
370
371 tlist->AddToList(trcal);
372 tlist->AddToList(fill09);
373 }
374 if (type==1)
375 {
376 MFillH *fill08 = new MFillH(evt08, "MArrivalTime", "MFillH8");
377 tlist->AddToList(fill08);
378 }
379
380 //
381 // Now distribute Display to all tasks
382 //
383 fEvtLoop = new MEvtLoop(Form("Mars %s", gSystem->BaseName(fname)));
384 fEvtLoop->SetOwner();
385 fEvtLoop->SetParList(plist);
386 fEvtLoop->SetDisplay(this);
387 fEvtLoop->ReadEnv("mars.rc");
388
389}
390
391// --------------------------------------------------------------------------
392//
393// Add the top part of the frame: This is filename and treename display
394//
395void MEventDisplay::AddTopFramePart1(TGCompositeFrame *frame,
396 const char *filename,
397 const char *treename)
398{
399 //
400 // --- the top1 part of the window ---
401 //
402 TGHorizontalFrame *top1 = new TGHorizontalFrame(frame, 1, 1);
403 fList->Add(top1);
404
405 //
406 // create gui elements
407 //
408 TGLabel *file = new TGLabel(top1, new TGString(Form("%s#%s", filename, treename)));
409 fList->Add(file);
410
411 //
412 // layout and add gui elements in/to frame
413 //
414 TGLayoutHints *laystd = new TGLayoutHints(kLHintsCenterX, 5, 5);
415 fList->Add(laystd);
416
417 top1->AddFrame(file, laystd);
418
419 //
420 // --- the top1 part of the window ---
421 //
422 TGHorizontalFrame *top2 = new TGHorizontalFrame(frame, 1, 1);
423 fList->Add(top2);
424
425 //
426 // layout and add frames
427 //
428 TGLayoutHints *laytop1 = new TGLayoutHints(kLHintsCenterX, 5, 5, 5);
429 fList->Add(laytop1);
430 frame->AddFrame(top1, laytop1);
431 frame->AddFrame(top2, laytop1);
432}
433
434// --------------------------------------------------------------------------
435//
436// Add the second part of the top frame: This are the event number controls
437//
438void MEventDisplay::AddTopFramePart2(TGCompositeFrame *frame)
439{
440 //
441 // --- the top2 part of the window ---
442 //
443 TGHorizontalFrame *top2 = new TGHorizontalFrame(frame, 1, 1);
444 fList->Add(top2);
445
446 //
447 // Create the gui elements
448 //
449 TGTextButton *prevevt = new TGTextButton(top2, " << ", kEvtPrev);
450 prevevt->Associate(this);
451
452 TGLabel *evtnr = new TGLabel(top2, new TGString("Event:"));
453
454 TGTextEntry *entry=new TGTextEntry(top2, new TGTextBuffer(100), kEvtNumber);
455 entry->Resize(60, entry->GetDefaultHeight());
456 entry->Associate(this);
457
458 fNumOfEvts = new TGLabel(top2, "of .");
459
460 TGTextButton *nextevt = new TGTextButton (top2, " >> ", kEvtNext);
461 nextevt->Associate(this);
462
463 //
464 // Add gui elements to 'atotodel'
465 //
466 fList->Add(prevevt);
467 fList->Add(evtnr);
468 fList->Add(entry);
469 fList->Add(fNumOfEvts);
470 fList->Add(nextevt);
471
472 //
473 // add the gui elements to the frame
474 //
475 TGLayoutHints *laystd = new TGLayoutHints(kLHintsLeft|kLHintsCenterY, 5, 5);
476 fList->Add(laystd);
477
478 top2->AddFrame(prevevt, laystd);
479 top2->AddFrame(evtnr, laystd);
480 top2->AddFrame(entry, laystd);
481 top2->AddFrame(fNumOfEvts, laystd);
482 top2->AddFrame(nextevt, laystd);
483
484 TGLayoutHints *laystd2 = new TGLayoutHints(kLHintsCenterX, 5, 5, 5, 5);
485 fList->Add(laystd2);
486 frame->AddFrame(top2, laystd2);
487
488 //
489 // Add trailing line...
490 //
491 TGHorizontal3DLine *line = new TGHorizontal3DLine(frame);
492 TGLayoutHints *layline = new TGLayoutHints(kLHintsExpandX);
493 fList->Add(line);
494 fList->Add(layline);
495 frame->AddFrame(line, layline);
496}
497
498// --------------------------------------------------------------------------
499//
500// Add the user frame part of the display
501//
502void MEventDisplay::AddUserFrame(const char* filename)
503{
504 fUserFrame->ChangeOptions(kHorizontalFrame);
505
506 TGCompositeFrame *vf1 = new TGVerticalFrame(fUserFrame, 1, 1);
507 TGCompositeFrame *vf2 = new TGVerticalFrame(fUserFrame, 1, 1);
508
509 AddTopFramePart1(vf1, filename, "Events");
510 AddTopFramePart2(vf1);
511
512 // create root embedded canvas and add it to the tab
513 TRootEmbeddedCanvas *ec = new TRootEmbeddedCanvas("Slices", vf2, vf1->GetDefaultHeight()*3/2, vf1->GetDefaultHeight(), 0);
514 vf2->AddFrame(ec);
515 fList->Add(ec);
516
517 // set background and border mode of the canvas
518 fCanvas = ec->GetCanvas();
519 fCanvas->SetBorderMode(0);
520 gROOT->GetListOfCanvases()->Add(fCanvas);
521 //fCanvas->SetBorderSize(1);
522 //fCanvas->SetBit(kNoContextMenu);
523 //fCanvas->SetFillColor(16);
524
525 TGLayoutHints *lay = new TGLayoutHints(kLHintsExpandX);
526 fUserFrame->AddFrame(vf1, lay);
527 fUserFrame->AddFrame(vf2);
528}
529
530// --------------------------------------------------------------------------
531//
532// Checks if the event number is valid, and if so reads the new event
533// and updates the display
534//
535void MEventDisplay::ReadinEvent(Int_t dir)
536{
537 MParList *plist = (MParList*) fEvtLoop->GetParList();
538 MTaskList *tlist = (MTaskList*) fEvtLoop->GetTaskList();
539 MGeomCam *geom = (MGeomCam*) plist->FindObject("MGeomCam");
540 MRawEvtData *raw = (MRawEvtData*)plist->FindObject("MRawEvtData");
541
542 //
543 // Read first event.
544 //
545 MReadTree *reader = (MReadTree*)fEvtLoop->FindTask("MRead");
546
547 const Int_t num = reader->GetNumEntry();
548
549 do
550 {
551 if (dir<0 && !reader->DecEventNum())
552 {
553 reader->SetEventNum(num);
554 return;
555 }
556 if (dir>0 && !reader->IncEventNum())
557 {
558 reader->SetEventNum(num);
559 return;
560 }
561
562 if (!tlist->Process())
563 return;
564
565 reader->DecEventNum();
566
567 } while (raw && raw->GetNumPixels()<1 && dir!=0);
568
569 //
570 // Cleare the 'FADC canvas'
571 //
572 fCanvas->Clear();
573 fCanvas->Modified();
574 fCanvas->Update();
575
576 //
577 // Print parameters
578 //
579 ((MHillas*) plist->FindObject("MHillas"))->Print(*geom);
580 ((MHillasExt*) plist->FindObject("MHillasExt"))->Print(*geom);
581 ((MHillasSrc*) plist->FindObject("MHillasSrc"))->Print(*geom);
582 plist->FindObject("MImagePar")->Print();
583 ((MNewImagePar*)plist->FindObject("MNewImagePar"))->Print(*geom);
584
585 //
586 // UpdateDisplay
587 //
588 gStyle->SetOptStat(1101);
589 Update();
590
591 TGTextEntry *entry = (TGTextEntry*)fList->FindWidget(kEvtNumber);
592 entry->SetText(Form("%d", reader->GetNumEntry()+1));
593}
594
595// --------------------------------------------------------------------------
596//
597// Read first event to get display booted
598//
599void MEventDisplay::ReadFirstEvent()
600{
601 if (!fEvtLoop->PreProcess())
602 return;
603
604 UnLock();
605
606 //
607 // Get parlist
608 //
609 MParList *plist = (MParList*)fEvtLoop->GetParList();
610
611 //
612 // Add Geometry tab
613 //
614 AddGeometryTabs();
615
616 //
617 // Now read event...
618 //
619 ReadinEvent();
620
621 MReadTree *reader = (MReadTree*)fEvtLoop->FindTask("MRead");
622 TGString *txt = new TGString(Form("of %d", reader->GetEntries()));
623 fNumOfEvts->SetText(txt);
624
625 //
626 // Draw ellipse on top of all pads
627 //
628 TObject *hillas1 = plist->FindObject("MHillas");
629 TObject *hillas2 = plist->FindObject("MHillasSrc");
630 for (int i=1; i<7;i++)
631 {
632 TCanvas *c = GetCanvas(i);
633 c->GetPad(1)->cd(1);
634 hillas1->Draw();
635 hillas2->Draw();
636 }
637}
638
639// --------------------------------------------------------------------------
640//
641// Adds the geometry tab
642//
643void MEventDisplay::AddGeometryTabs()
644{
645 MGeomCam *geom = (MGeomCam*)fEvtLoop->GetParList()->FindObject("MGeomCam");
646 if (!geom)
647 return;
648
649 TCanvas &c1=AddTab("Geometry");
650
651 MHCamera *cam = new MHCamera(*geom);
652 cam->SetBit(TH1::kNoStats|MHCamera::kNoLegend|kCanDelete);
653 cam->Draw("pixelindex");
654
655 c1.Modified();
656 c1.Update();
657
658 TCanvas &c2=AddTab("Sectors");
659
660 cam = new MHCamera(*geom);
661 cam->SetBit(TH1::kNoStats|MHCamera::kNoLegend|kCanDelete);
662 cam->Draw("sectorindex");
663
664 c2.Modified();
665 c2.Update();
666}
667
668// --------------------------------------------------------------------------
669//
670// ProcessMessage(Long_t msg, Long_t parm1, Long_t parm2)
671//
672// Processes information from all GUI items.
673// Selecting an item usually generates an event with 4 parameters.
674// The first two are packed into msg (first and second bytes).
675// The other two are parm1 and parm2.
676//
677Bool_t MEventDisplay::ProcessMessage(Long_t msg, Long_t mp1, Long_t mp2)
678{
679 switch (GET_MSG(msg))
680 {
681 case kC_TEXTENTRY:
682 switch(GET_SUBMSG(msg))
683 {
684 case kTE_ENTER:
685 switch(GET_SUBMSG(msg))
686 {
687 case kTE_ENTER:
688 {
689 TGTextEntry *entry = (TGTextEntry*)fList->FindWidget(kEvtNumber);
690 MReadTree *reader = (MReadTree*)fEvtLoop->FindTask("MRead");
691 if (reader->SetEventNum(atoi(entry->GetText())-1))
692 ReadinEvent();
693 }
694 return kTRUE;
695 }
696 return kTRUE;
697 }
698 break;
699
700 case kC_COMMAND:
701 switch (GET_SUBMSG(msg))
702 {
703 case kCM_TAB:
704 {
705 //
706 // Set name for 'FADC canvas'. The name is the anchor for MHCamera.
707 // and clear the canvas
708 TCanvas *c = GetCanvas(mp1);
709 if (!c)
710 break;
711 MHEvent *o = (MHEvent*)fEvtLoop->GetParList()->FindObject(c->GetName());
712 if (!o)
713 break;
714 fCanvas->SetName(Form("%p;%p;PixelContent", o->GetHist(), c->GetPad(1)->GetPad(1)));
715 }
716 break;
717
718 case kCM_BUTTON:
719 switch (mp1)
720 {
721 case kEvtPrev:
722 ReadinEvent(-1);
723 return kTRUE;
724
725 case kEvtNext:
726 ReadinEvent(+1);
727 return kTRUE;
728 }
729 return kTRUE;
730 }
731 break;
732 }
733
734 return MStatusDisplay::ProcessMessage(msg, mp1, mp2);
735}
Note: See TracBrowser for help on using the repository browser.