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

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