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

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