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 <TList.h> // TList::Add
|
---|
35 | #include <TStyle.h> // gStyle->SetOptStat
|
---|
36 | #include <TCanvas.h> // TCanvas::cd
|
---|
37 |
|
---|
38 | //
|
---|
39 | // root GUI
|
---|
40 | //
|
---|
41 | #include <TGLabel.h> // TGLabel
|
---|
42 | #include <TGButton.h> // TGPictureButton
|
---|
43 | #include <TG3DLine.h> // TGHorizontal3DLine
|
---|
44 | #include <TGTextEntry.h> // TGTextEntry
|
---|
45 | #include <TGButtonGroup.h> // TGVButtonGroup
|
---|
46 | #include <TRootEmbeddedCanvas.h> // TRootEmbeddedCanvas
|
---|
47 |
|
---|
48 | //
|
---|
49 | // General
|
---|
50 | //
|
---|
51 | #include "MGList.h" // MGList
|
---|
52 |
|
---|
53 | #include "MParList.h" // MParList::AddToList
|
---|
54 | #include "MEvtLoop.h" // MEvtLoop::GetParList
|
---|
55 | #include "MTaskList.h" // MTaskList::AddToList
|
---|
56 |
|
---|
57 | //
|
---|
58 | // Tasks
|
---|
59 | //
|
---|
60 | #include "MReadMarsFile.h" // MReadMarsFile
|
---|
61 | #include "MGeomApply.h" // MGeomApply
|
---|
62 | #include "MFDataMember.h" // MFDataMember
|
---|
63 | #include "MMcPedestalCopy.h" // MMcPedestalCopy
|
---|
64 | #include "MMcPedestalNSBAdd.h" // MMcPedestalNSBAdd
|
---|
65 | #include "MCerPhotCalc.h" // MCerPhotCalc
|
---|
66 | #include "MCerPhotAnal2.h" // MCerPhotAnal2
|
---|
67 | #include "MImgCleanStd.h" // MImgCleanStd
|
---|
68 | #include "MHillasCalc.h" // MHillasCalc
|
---|
69 | #include "MHillasSrcCalc.h" // MHillasSrcCalc
|
---|
70 | #include "MBlindPixelCalc.h" // MBlindPixelCalc
|
---|
71 | #include "MFillH.h" // MFillH
|
---|
72 |
|
---|
73 | //
|
---|
74 | // Container
|
---|
75 | //
|
---|
76 | #include "MHillas.h" // MHillas::Print(const MGeomCam&)
|
---|
77 | #include "MHillasExt.h" // MHillasExt::Print(const MGeomCam&)
|
---|
78 | #include "MHillasSrc.h" // MHillasSrc::Print(const MGeomCam&)
|
---|
79 | #include "MHEvent.h" // MHEvent
|
---|
80 | #include "MHCamera.h" // MHCamera
|
---|
81 | #include "MRawEvtData.h" // MRawEvtData
|
---|
82 |
|
---|
83 | ClassImp(MEventDisplay);
|
---|
84 |
|
---|
85 | // --------------------------------------------------------------------------
|
---|
86 | //
|
---|
87 | // Constructor.
|
---|
88 | //
|
---|
89 | MEventDisplay::MEventDisplay(const char *filename) : MStatusDisplay()
|
---|
90 | {
|
---|
91 | //
|
---|
92 | // Setup Task list for hillas calculation
|
---|
93 | //
|
---|
94 | SetupTaskList("Events", filename);
|
---|
95 |
|
---|
96 | //
|
---|
97 | // Add MEventDisplay GUI elements to the display
|
---|
98 | //
|
---|
99 | AddUserFrame(filename);
|
---|
100 |
|
---|
101 | //
|
---|
102 | // Show new part of the window, resize to correct aspect ratio
|
---|
103 | //
|
---|
104 | // FIXME: This should be done by MStatusDisplay automatically
|
---|
105 | Resize(GetWidth(), GetHeight() + fUserFrame->GetDefaultHeight());
|
---|
106 | SetWindowName("Event Display");
|
---|
107 | MapSubwindows();
|
---|
108 |
|
---|
109 | //
|
---|
110 | // Readin first event and display it
|
---|
111 | //
|
---|
112 | ReadFirstEvent();
|
---|
113 | }
|
---|
114 |
|
---|
115 | // --------------------------------------------------------------------------
|
---|
116 | //
|
---|
117 | // Destructor: PostProcess eventloop, delete eventloop with all containers
|
---|
118 | //
|
---|
119 | MEventDisplay::~MEventDisplay()
|
---|
120 | {
|
---|
121 | fEvtLoop->PostProcess();
|
---|
122 | delete fEvtLoop;
|
---|
123 | }
|
---|
124 |
|
---|
125 | // --------------------------------------------------------------------------
|
---|
126 | //
|
---|
127 | // Return reading task
|
---|
128 | //
|
---|
129 | MReadTree *MEventDisplay::GetReader() const
|
---|
130 | {
|
---|
131 | MParList *plist = (MParList*)fEvtLoop->GetParList();
|
---|
132 | MTaskList *tlist = (MTaskList*)plist->FindObject("MTaskList");
|
---|
133 | MReadTree *reader = (MReadTree*)tlist->FindObject("MRead");
|
---|
134 | return reader;
|
---|
135 | }
|
---|
136 |
|
---|
137 | // --------------------------------------------------------------------------
|
---|
138 | //
|
---|
139 | // Setup Task and parameter list for hillas calculation,
|
---|
140 | // preprocess tasks and read in first event (process)
|
---|
141 | //
|
---|
142 | void MEventDisplay::SetupTaskList(const char *tname, const char *fname)
|
---|
143 | {
|
---|
144 | //
|
---|
145 | // Setup an empty job, with a reader task only.
|
---|
146 | // All tasks and parameter containers are deleted automatically
|
---|
147 | // (via SetOwner())
|
---|
148 | //
|
---|
149 | MTaskList *tlist = new MTaskList;
|
---|
150 | tlist->SetOwner();
|
---|
151 |
|
---|
152 | MReadMarsFile *read = new MReadMarsFile(tname, fname);
|
---|
153 | read->DisableAutoScheme();
|
---|
154 | tlist->AddToList(read);
|
---|
155 |
|
---|
156 | MGeomApply *apl = new MGeomApply;
|
---|
157 | tlist->AddToList(apl);
|
---|
158 |
|
---|
159 | MParList *plist = new MParList;
|
---|
160 | plist->SetOwner();
|
---|
161 | plist->AddToList(tlist);
|
---|
162 |
|
---|
163 | fEvtLoop = new MEvtLoop;
|
---|
164 | fEvtLoop->SetOwner();
|
---|
165 | fEvtLoop->SetParList(plist);
|
---|
166 |
|
---|
167 | MHEvent *evt1 = new MHEvent(MHEvent::kEvtSignal);
|
---|
168 | MHEvent *evt2 = new MHEvent(MHEvent::kEvtSignal);
|
---|
169 | MHEvent *evt3 = new MHEvent(MHEvent::kEvtPedestal);
|
---|
170 | MHEvent *evt4 = new MHEvent(MHEvent::kEvtPedestalRMS);
|
---|
171 | MHEvent *evt5 = new MHEvent(MHEvent::kEvtRelativeSignal);
|
---|
172 | MHEvent *evt6 = new MHEvent(MHEvent::kEvtCleaningLevels);
|
---|
173 | evt1->SetName("Signal");
|
---|
174 | evt2->SetName("Cleaned");
|
---|
175 | evt3->SetName("Pedestal");
|
---|
176 | evt4->SetName("PedRMS");
|
---|
177 | evt5->SetName("Signal/PedRMS");
|
---|
178 | evt6->SetName("CleanLevels");
|
---|
179 | plist->AddToList(evt1);
|
---|
180 | plist->AddToList(evt2);
|
---|
181 | plist->AddToList(evt3);
|
---|
182 | plist->AddToList(evt4);
|
---|
183 | plist->AddToList(evt5);
|
---|
184 | plist->AddToList(evt6);
|
---|
185 |
|
---|
186 | MMcPedestalCopy *pcopy = new MMcPedestalCopy;
|
---|
187 | MMcPedestalNSBAdd *pdnsb = new MMcPedestalNSBAdd;
|
---|
188 | MCerPhotCalc *ncalc = new MCerPhotCalc;
|
---|
189 | MCerPhotAnal2 *nanal = new MCerPhotAnal2;
|
---|
190 | MFillH *fill1 = new MFillH(evt1, "MCerPhotEvt", "MFillH1");
|
---|
191 | MImgCleanStd *clean = new MImgCleanStd;
|
---|
192 | MFillH *fill2 = new MFillH(evt2, "MCerPhotEvt", "MFillH2");
|
---|
193 | MFillH *fill3 = new MFillH(evt3, "MPedestalCam", "MFillH3");
|
---|
194 | MFillH *fill4 = new MFillH(evt4, "MPedestalCam", "MFillH4");
|
---|
195 | MFillH *fill5 = new MFillH(evt5, "MCameraData", "MFillH5");
|
---|
196 | MFillH *fill6 = new MFillH(evt6, "MCameraData", "MFillH6");
|
---|
197 | MBlindPixelCalc *blind = new MBlindPixelCalc;
|
---|
198 | MHillasCalc *hcalc = new MHillasCalc;
|
---|
199 | MHillasSrcCalc *scalc = new MHillasSrcCalc;
|
---|
200 |
|
---|
201 | MFilter *f1=new MFDataMember("MRawRunHeader.fRunType", '>', 255.5);
|
---|
202 | MFilter *f2=new MFDataMember("MRawRunHeader.fRunType", '<', 255.5);
|
---|
203 |
|
---|
204 | ncalc->SetFilter(f1);
|
---|
205 | nanal->SetFilter(f2);
|
---|
206 |
|
---|
207 | tlist->AddToList(f1);
|
---|
208 | tlist->AddToList(f2);
|
---|
209 | tlist->AddToList(pcopy);
|
---|
210 | tlist->AddToList(pdnsb);
|
---|
211 | tlist->AddToList(ncalc);
|
---|
212 | tlist->AddToList(nanal);
|
---|
213 | tlist->AddToList(fill1);
|
---|
214 | tlist->AddToList(clean);
|
---|
215 | tlist->AddToList(fill2);
|
---|
216 | tlist->AddToList(fill3);
|
---|
217 | tlist->AddToList(fill4);
|
---|
218 | tlist->AddToList(fill5);
|
---|
219 | tlist->AddToList(fill6);
|
---|
220 | tlist->AddToList(blind);
|
---|
221 | tlist->AddToList(hcalc);
|
---|
222 | tlist->AddToList(scalc);
|
---|
223 |
|
---|
224 | //
|
---|
225 | // Now distribute Display to all tasks
|
---|
226 | //
|
---|
227 | tlist->SetDisplay(this);
|
---|
228 | }
|
---|
229 |
|
---|
230 | // --------------------------------------------------------------------------
|
---|
231 | //
|
---|
232 | // Add the top part of the frame: This is filename and treename display
|
---|
233 | //
|
---|
234 | void MEventDisplay::AddTopFramePart1(TGCompositeFrame *frame,
|
---|
235 | const char *filename,
|
---|
236 | const char *treename)
|
---|
237 | {
|
---|
238 | //
|
---|
239 | // --- the top1 part of the window ---
|
---|
240 | //
|
---|
241 | TGHorizontalFrame *top1 = new TGHorizontalFrame(frame, 1, 1);
|
---|
242 | fList->Add(top1);
|
---|
243 |
|
---|
244 | //
|
---|
245 | // create gui elements
|
---|
246 | //
|
---|
247 | TGLabel *file = new TGLabel(top1, new TGString(Form("%s#%s", filename, treename)));
|
---|
248 | fList->Add(file);
|
---|
249 |
|
---|
250 | //
|
---|
251 | // layout and add gui elements in/to frame
|
---|
252 | //
|
---|
253 | TGLayoutHints *laystd = new TGLayoutHints(kLHintsCenterX, 5, 5);
|
---|
254 | fList->Add(laystd);
|
---|
255 |
|
---|
256 | top1->AddFrame(file, laystd);
|
---|
257 |
|
---|
258 | //
|
---|
259 | // --- the top1 part of the window ---
|
---|
260 | //
|
---|
261 | TGHorizontalFrame *top2 = new TGHorizontalFrame(frame, 1, 1);
|
---|
262 | fList->Add(top2);
|
---|
263 |
|
---|
264 | //
|
---|
265 | // layout and add frames
|
---|
266 | //
|
---|
267 | TGLayoutHints *laytop1 = new TGLayoutHints(kLHintsCenterX, 5, 5, 5);
|
---|
268 | fList->Add(laytop1);
|
---|
269 | frame->AddFrame(top1, laytop1);
|
---|
270 | frame->AddFrame(top2, laytop1);
|
---|
271 | }
|
---|
272 |
|
---|
273 | // --------------------------------------------------------------------------
|
---|
274 | //
|
---|
275 | // Add the second part of the top frame: This are the event number controls
|
---|
276 | //
|
---|
277 | void MEventDisplay::AddTopFramePart2(TGCompositeFrame *frame)
|
---|
278 | {
|
---|
279 | //
|
---|
280 | // --- the top2 part of the window ---
|
---|
281 | //
|
---|
282 | TGHorizontalFrame *top2 = new TGHorizontalFrame(frame, 1, 1);
|
---|
283 | fList->Add(top2);
|
---|
284 |
|
---|
285 | //
|
---|
286 | // Create the gui elements
|
---|
287 | //
|
---|
288 | TGTextButton *prevevt = new TGTextButton(top2, " << ", kEvtPrev);
|
---|
289 | prevevt->Associate(this);
|
---|
290 |
|
---|
291 | TGLabel *evtnr = new TGLabel(top2, new TGString("Event:"));
|
---|
292 |
|
---|
293 | TGTextEntry *entry=new TGTextEntry(top2, new TGTextBuffer(100), kEvtNumber);
|
---|
294 | entry->Resize(60, entry->GetDefaultHeight());
|
---|
295 | entry->Associate(this);
|
---|
296 |
|
---|
297 | fNumOfEvts = new TGLabel(top2, "of .");
|
---|
298 |
|
---|
299 | TGTextButton *nextevt = new TGTextButton (top2, " >> ", kEvtNext);
|
---|
300 | nextevt->Associate(this);
|
---|
301 |
|
---|
302 | //
|
---|
303 | // Add gui elements to 'atotodel'
|
---|
304 | //
|
---|
305 | fList->Add(prevevt);
|
---|
306 | fList->Add(evtnr);
|
---|
307 | fList->Add(entry);
|
---|
308 | fList->Add(fNumOfEvts);
|
---|
309 | fList->Add(nextevt);
|
---|
310 |
|
---|
311 | //
|
---|
312 | // add the gui elements to the frame
|
---|
313 | //
|
---|
314 | TGLayoutHints *laystd = new TGLayoutHints(kLHintsLeft|kLHintsCenterY, 5, 5);
|
---|
315 | fList->Add(laystd);
|
---|
316 |
|
---|
317 | top2->AddFrame(prevevt, laystd);
|
---|
318 | top2->AddFrame(evtnr, laystd);
|
---|
319 | top2->AddFrame(entry, laystd);
|
---|
320 | top2->AddFrame(fNumOfEvts, laystd);
|
---|
321 | top2->AddFrame(nextevt, laystd);
|
---|
322 |
|
---|
323 | TGLayoutHints *laystd2 = new TGLayoutHints(kLHintsCenterX, 5, 5, 5, 5);
|
---|
324 | fList->Add(laystd2);
|
---|
325 | frame->AddFrame(top2, laystd2);
|
---|
326 |
|
---|
327 | //
|
---|
328 | // Add trailing line...
|
---|
329 | //
|
---|
330 | TGHorizontal3DLine *line = new TGHorizontal3DLine(frame);
|
---|
331 | TGLayoutHints *layline = new TGLayoutHints(kLHintsExpandX);
|
---|
332 | fList->Add(line);
|
---|
333 | fList->Add(layline);
|
---|
334 | frame->AddFrame(line, layline);
|
---|
335 | }
|
---|
336 |
|
---|
337 | // --------------------------------------------------------------------------
|
---|
338 | //
|
---|
339 | // Add the user frame part of the display
|
---|
340 | //
|
---|
341 | void MEventDisplay::AddUserFrame(const char* filename)
|
---|
342 | {
|
---|
343 | fUserFrame->ChangeOptions(kHorizontalFrame);
|
---|
344 |
|
---|
345 | TGCompositeFrame *vf1 = new TGVerticalFrame(fUserFrame, 1, 1);
|
---|
346 | TGCompositeFrame *vf2 = new TGVerticalFrame(fUserFrame, 1, 1);
|
---|
347 |
|
---|
348 | AddTopFramePart1(vf1, filename, "Events");
|
---|
349 | AddTopFramePart2(vf1);
|
---|
350 |
|
---|
351 | // create root embedded canvas and add it to the tab
|
---|
352 | TRootEmbeddedCanvas *ec = new TRootEmbeddedCanvas("Slices", vf2, vf1->GetDefaultHeight()*3/2, vf1->GetDefaultHeight(), 0);
|
---|
353 | vf2->AddFrame(ec);
|
---|
354 | fList->Add(ec);
|
---|
355 |
|
---|
356 | // set background and border mode of the canvas
|
---|
357 | fCanvas = ec->GetCanvas();
|
---|
358 | fCanvas->SetBorderMode(0);
|
---|
359 | gROOT->GetListOfCanvases()->Add(fCanvas);
|
---|
360 | //fCanvas->SetBorderSize(1);
|
---|
361 | //fCanvas->SetBit(kNoContextMenu);
|
---|
362 | //fCanvas->SetFillColor(16);
|
---|
363 |
|
---|
364 | TGLayoutHints *lay = new TGLayoutHints(kLHintsExpandX);
|
---|
365 | fUserFrame->AddFrame(vf1, lay);
|
---|
366 | fUserFrame->AddFrame(vf2);
|
---|
367 | }
|
---|
368 |
|
---|
369 | // --------------------------------------------------------------------------
|
---|
370 | //
|
---|
371 | // Checks if the event number is valid, and if so reads the new event
|
---|
372 | // and updates the display
|
---|
373 | //
|
---|
374 | void MEventDisplay::ReadinEvent(Int_t dir)
|
---|
375 | {
|
---|
376 | MParList *plist = (MParList*)fEvtLoop->GetParList();
|
---|
377 | MGeomCam *geom = (MGeomCam*) plist->FindObject("MGeomCam");
|
---|
378 | MTaskList *tlist = (MTaskList*)plist->FindObject("MTaskList");
|
---|
379 | MRawEvtData *raw = (MRawEvtData*)plist->FindObject("MRawEvtData");
|
---|
380 | if (!raw)
|
---|
381 | return;
|
---|
382 |
|
---|
383 | //
|
---|
384 | // Read first event.
|
---|
385 | //
|
---|
386 | MReadTree *reader = GetReader();
|
---|
387 |
|
---|
388 | const Int_t num = reader->GetNumEntry();
|
---|
389 |
|
---|
390 | do
|
---|
391 | {
|
---|
392 | if (dir<0 && !reader->DecEventNum())
|
---|
393 | {
|
---|
394 | reader->SetEventNum(num);
|
---|
395 | return;
|
---|
396 | }
|
---|
397 | if (dir>0 && !reader->IncEventNum())
|
---|
398 | {
|
---|
399 | reader->SetEventNum(num);
|
---|
400 | return;
|
---|
401 | }
|
---|
402 |
|
---|
403 | if (!tlist->Process())
|
---|
404 | return;
|
---|
405 |
|
---|
406 | reader->DecEventNum();
|
---|
407 |
|
---|
408 | } while (raw->GetNumPixels()<1 && dir!=0);
|
---|
409 |
|
---|
410 | //
|
---|
411 | // Cleare the 'FADC canvas'
|
---|
412 | //
|
---|
413 | fCanvas->Clear();
|
---|
414 | fCanvas->Modified();
|
---|
415 | fCanvas->Update();
|
---|
416 |
|
---|
417 | //
|
---|
418 | // Print parameters
|
---|
419 | //
|
---|
420 | ((MHillas*) plist->FindObject("MHillas"))->Print(*geom);
|
---|
421 | ((MHillasExt*)plist->FindObject("MHillasExt"))->Print(*geom);
|
---|
422 | ((MHillasSrc*)plist->FindObject("MHillasSrc"))->Print(*geom);
|
---|
423 | plist->FindObject("MNewImagePar")->Print();
|
---|
424 |
|
---|
425 | //
|
---|
426 | // UpdateDisplay
|
---|
427 | //
|
---|
428 | gStyle->SetOptStat(1101);
|
---|
429 | Update();
|
---|
430 |
|
---|
431 | TGTextEntry *entry = (TGTextEntry*)fList->FindWidget(kEvtNumber);
|
---|
432 | entry->SetText(Form("%d", reader->GetNumEntry()+1));
|
---|
433 | }
|
---|
434 |
|
---|
435 | // --------------------------------------------------------------------------
|
---|
436 | //
|
---|
437 | // Read first event to get display booted
|
---|
438 | //
|
---|
439 | void MEventDisplay::ReadFirstEvent()
|
---|
440 | {
|
---|
441 | if (!fEvtLoop->PreProcess())
|
---|
442 | return;
|
---|
443 |
|
---|
444 | //
|
---|
445 | // Get parlist
|
---|
446 | //
|
---|
447 | MParList *plist = (MParList*)fEvtLoop->GetParList();
|
---|
448 |
|
---|
449 | //
|
---|
450 | // Add Geometry tab
|
---|
451 | //
|
---|
452 | AddGeometryTabs();
|
---|
453 |
|
---|
454 | //
|
---|
455 | // Now read event...
|
---|
456 | //
|
---|
457 | ReadinEvent();
|
---|
458 |
|
---|
459 | TGString *txt = new TGString(Form("of %d", GetReader()->GetEntries()));
|
---|
460 | fNumOfEvts->SetText(txt);
|
---|
461 |
|
---|
462 | //
|
---|
463 | // Draw ellipse on top of all pads
|
---|
464 | //
|
---|
465 | TObject *hillas = plist->FindObject("MHillas");
|
---|
466 | for (int i=1; i<7;i++)
|
---|
467 | {
|
---|
468 | TCanvas *c = GetCanvas(i);
|
---|
469 | c->GetPad(1)->cd(1);
|
---|
470 | hillas->Draw();
|
---|
471 | }
|
---|
472 | }
|
---|
473 |
|
---|
474 | // --------------------------------------------------------------------------
|
---|
475 | //
|
---|
476 | // Adds the geometry tab
|
---|
477 | //
|
---|
478 | void MEventDisplay::AddGeometryTabs()
|
---|
479 | {
|
---|
480 | MGeomCam *geom = (MGeomCam*)fEvtLoop->GetParList()->FindObject("MGeomCam");
|
---|
481 | if (!geom)
|
---|
482 | return;
|
---|
483 |
|
---|
484 | TCanvas &c1=AddTab("Geometry");
|
---|
485 |
|
---|
486 | MHCamera *cam = new MHCamera(*geom);
|
---|
487 | cam->SetBit(TH1::kNoStats|MHCamera::kNoLegend|kCanDelete);
|
---|
488 | cam->Draw("pixelindex");
|
---|
489 | fList->Add(cam);
|
---|
490 |
|
---|
491 | c1.Modified();
|
---|
492 | c1.Update();
|
---|
493 |
|
---|
494 | TCanvas &c2=AddTab("Sectors");
|
---|
495 |
|
---|
496 | cam = new MHCamera(*geom);
|
---|
497 | cam->SetBit(TH1::kNoStats|MHCamera::kNoLegend|kCanDelete);
|
---|
498 | cam->Draw("sectorindex");
|
---|
499 | fList->Add(cam);
|
---|
500 |
|
---|
501 | c2.Modified();
|
---|
502 | c2.Update();
|
---|
503 | }
|
---|
504 |
|
---|
505 | // --------------------------------------------------------------------------
|
---|
506 | //
|
---|
507 | // ProcessMessage(Long_t msg, Long_t parm1, Long_t parm2)
|
---|
508 | //
|
---|
509 | // Processes information from all GUI items.
|
---|
510 | // Selecting an item usually generates an event with 4 parameters.
|
---|
511 | // The first two are packed into msg (first and second bytes).
|
---|
512 | // The other two are parm1 and parm2.
|
---|
513 | //
|
---|
514 | Bool_t MEventDisplay::ProcessMessage(Long_t msg, Long_t mp1, Long_t mp2)
|
---|
515 | {
|
---|
516 | switch (GET_MSG(msg))
|
---|
517 | {
|
---|
518 | case kC_TEXTENTRY:
|
---|
519 | switch(GET_SUBMSG(msg))
|
---|
520 | {
|
---|
521 | case kTE_ENTER:
|
---|
522 | switch(GET_SUBMSG(msg))
|
---|
523 | {
|
---|
524 | case kTE_ENTER:
|
---|
525 | {
|
---|
526 | TGTextEntry *entry = (TGTextEntry*)fList->FindWidget(kEvtNumber);
|
---|
527 | if (GetReader()->SetEventNum(atoi(entry->GetText())-1))
|
---|
528 | ReadinEvent();
|
---|
529 | }
|
---|
530 | return kTRUE;
|
---|
531 | }
|
---|
532 | return kTRUE;
|
---|
533 | }
|
---|
534 | break;
|
---|
535 |
|
---|
536 | case kC_COMMAND:
|
---|
537 | switch (GET_SUBMSG(msg))
|
---|
538 | {
|
---|
539 | case kCM_TAB:
|
---|
540 | {
|
---|
541 | //
|
---|
542 | // Set name for 'FADC canvas'. The name is the anchor for MHCamera.
|
---|
543 | // and clear the canvas
|
---|
544 | TCanvas *c = GetCanvas(mp1);
|
---|
545 | if (!c)
|
---|
546 | break;
|
---|
547 | MHEvent *o = (MHEvent*)fEvtLoop->GetParList()->FindObject(c->GetName());
|
---|
548 | if (!o)
|
---|
549 | break;
|
---|
550 | fCanvas->SetName(Form("%p;%p;PixelContent", o->GetHist(), c->GetPad(1)));
|
---|
551 | }
|
---|
552 | break;
|
---|
553 |
|
---|
554 | case kCM_BUTTON:
|
---|
555 | switch (mp1)
|
---|
556 | {
|
---|
557 | case kEvtPrev:
|
---|
558 | ReadinEvent(-1);
|
---|
559 | return kTRUE;
|
---|
560 |
|
---|
561 | case kEvtNext:
|
---|
562 | ReadinEvent(+1);
|
---|
563 | return kTRUE;
|
---|
564 | }
|
---|
565 | return kTRUE;
|
---|
566 | }
|
---|
567 | break;
|
---|
568 | }
|
---|
569 |
|
---|
570 | return MStatusDisplay::ProcessMessage(msg, mp1, mp2);
|
---|
571 | }
|
---|