]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDSpectraDisplay.cxx
Some changes to make it possible to run the DA's
[u/mrichter/AliRoot.git] / FMD / AliFMDSpectraDisplay.cxx
1 // 
2 //____________________________________________________________________
3 //
4 // $Id: DrawHits.C 22496 2007-11-26 13:50:44Z cholm $
5 //
6 // Script that contains a class to draw hits, using the
7 // AliFMDInputHits class in the util library. 
8 //
9 // Use the script `Compile.C' to compile this class using ACLic. 
10 //
11 #include <AliCDBManager.h>
12 #include <AliFMDParameters.h>
13 #include <AliFMDHit.h>
14 #include <AliFMDDigit.h>
15 #include <AliFMDSDigit.h>
16 #include <AliFMDRecPoint.h>
17 #include <AliESDFMD.h>
18 #include <AliFMDPattern.h>
19 #include <AliFMDSpectraDisplay.h>
20 #include <TBrowser.h>
21 #include <TDirectory.h>
22 #include <TObjArray.h>
23 #include <TH1.h>
24 #include <TH2.h>
25 #include <TH3.h>
26 #include <TGraph.h>
27 #include <TF1.h>
28 #include <iostream>
29 #include <TStyle.h>
30 #include <TEnv.h>
31 #include <TCanvas.h>
32 #include <TGFrame.h>
33 #include <TGCanvas.h>
34 #include <TGListTree.h>
35 #include <TGClient.h>
36 #include <TSystem.h>
37 #include <KeySymbols.h>
38 #include <TClass.h>
39 #include <RQ_OBJECT.h>
40 #include <TSlider.h>
41 #define NESTED(X) AliFMDSpectraDisplay::AliFMDSpectraDisplay # X
42
43 //==================================================================
44 void 
45 AliFMDSpectraDisplay::AliFMDSpectraDisplayElement::MakeHistograms(const 
46                                                                   TAxis* 
47                                                                   axis) 
48 {
49   // Create the 
50   // needed histograms
51   // for this element
52   if (fFull) return;
53   if (axis->IsVariableBinSize()) {
54     fFull = new TH1F(Form("f_%s", GetName()), Form("%s spectra", GetName()),
55                      axis->GetXbins()->fN, axis->GetXbins()->fArray);
56     fCut = new TH1F(Form("c_%s", GetName()), 
57                     Form("%s restricted spectra", GetName()),
58                     axis->GetXbins()->fN, axis->GetXbins()->fArray);
59   }
60   else { 
61     fFull = new TH1F(Form("f_%s", GetName()), Form("%s spectra", GetName()),
62                      axis->GetNbins(), axis->GetXmin(), axis->GetXmax());
63     fCut = new TH1F(Form("c_%s", GetName()), 
64                     Form("%s restricted spectra", GetName()),
65                     axis->GetNbins(), axis->GetXmin(), axis->GetXmax());
66   }
67   fFull->SetFillColor(kRed);
68   fFull->SetFillStyle(3001);
69   fCut->SetFillColor(kBlue);
70   fCut->SetFillStyle(3001);
71 }
72 //__________________________________________________________________
73 void AliFMDSpectraDisplay::AliFMDSpectraDisplayElement::DoFill(Double_t v) 
74 {
75   // Fill into histograms
76   if (fFull) fFull->Fill(v);
77 }
78 //__________________________________________________________________
79 void AliFMDSpectraDisplay::AliFMDSpectraDisplayElement::Show(Option_t* option, 
80                                        Double_t l, Double_t h) 
81 {
82   // Show this element
83   // 
84   if (!fFull) return;
85   gPad->SetLogy(fFull->GetMaximum() > 10);
86   fFull->Draw(option);
87   Double_t lx = fFull->GetXaxis()->GetXmin();
88   Double_t hx = fFull->GetXaxis()->GetXmax();
89   Double_t rr = (hx-lx);
90   Double_t ll = rr * l + lx;
91   Double_t hh = rr * h + lx;
92   for (Int_t i = 1; i <= fFull->GetNbinsX(); i++) { 
93     if (fFull->GetBinCenter(i) <= ll || 
94         fFull->GetBinCenter(i) >  hh) { 
95       fCut->SetBinContent(i, 0); 
96       continue;
97     }
98     fCut->SetBinContent(i, fFull->GetBinContent(i));
99   }
100   fCut->Draw(Form("%s same", option));
101 }
102
103 //__________________________________________________________________
104 Int_t
105 AliFMDSpectraDisplay::AliFMDSpectraDisplayElement::Compare(const TObject*) const
106 {
107   return -1;
108 }
109
110
111 //==================================================================
112 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::AliFMDSpectraDisplayTop(TGCompositeFrame& frame, 
113                                                  TCanvas* canvas) 
114   : AliFMDSpectraDisplayElement("All", "Everything"), 
115     fHints(kLHintsExpandX|kLHintsExpandY,3, 3, 3, 3),
116     fContainer(&frame, 200, 350), 
117     fList(&fContainer, kHorizontalFrame),
118     fChildren(0),
119     fCurrentEntry(0),
120     fCanvas(canvas),
121     fkHist1DIcon(gClient->GetPicture("h1_t.xpm")),
122     fkHist2DIcon(gClient->GetPicture("h2_t.xpm")),
123     fkHist3DIcon(gClient->GetPicture("h3_t.xpm")),
124     fkGraphIcon(gClient->GetPicture("graph.xpm")), 
125     fAxis(0),
126     fEntry(*(fList.AddItem(0, GetName(), this)))
127 {
128   // Constructor 
129   // Parameters: 
130   //    frame    PArent frame 
131   //    canvas   Canvas to draw in
132   fContainer.AddFrame(&fList, &fHints);
133   frame.AddFrame(&fContainer, &fHints);
134
135   fList.Connect("Clicked(TGListTreeItem*,Int_t)", 
136                 "AliFMDSpectraDisplay::AliFMDSpectraDisplayTop", this, 
137                 "HandleEntry(TGListTreeItem*,Int_t)");
138   fList.Connect("KeyPressed(TGListTreeItem*,ULong_t,ULong_t)", 
139                 "AliFMDSpectraDisplay::AliFMDSpectraDisplayTop", this, 
140                 "HandleKey(TGListTreeItem*,UInt_t,UInt_t)");
141   fList.Connect("ReturnPressed(TGListTreeItem*)", 
142                 "AliFMDSpectraDisplay::AliFMDSpectraDisplayTop", this, 
143                 "HandleReturn(TGListTreeItem*)");
144 }
145 //____________________________________________________________________
146 void
147 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::SetAxis(TAxis* axis) 
148 {
149   // Set the axis of histograms
150   fAxis = axis;
151   MakeHistograms(axis);
152 }
153 //____________________________________________________________________
154 void
155 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::ClearCanvas()
156 {
157   // Clear the canvas
158   if (!fCanvas) return;
159   fCanvas->Clear();
160 }
161   
162 //____________________________________________________________________
163 void
164 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::ClearList()
165 {
166   // Clear the lsit
167   fList.DeleteItem(fList.GetFirstItem());
168   UpdateList();
169 }
170   
171   
172 //____________________________________________________________________
173 void
174 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::HandleReturn(TGListTreeItem * f)
175 {
176   // HAndle when return is pressed
177   if (!f) { 
178     fList.UnselectAll(kFALSE);
179     fList.SetSelected(0);
180     return;
181   }
182   fList.ToggleItem(f);
183   UpdateList();
184 }
185   
186   
187 //____________________________________________________________________
188 void
189 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::HandleKey(TGListTreeItem * f, UInt_t keysym, UInt_t)
190 {
191   // Handle a key stroke
192   if (!f) { 
193     fList.UnselectAll(kFALSE);
194     fList.SetSelected(0);
195     return;
196   }
197   TGListTreeItem* next = 0;
198   switch (keysym) {
199   case kKey_Up:
200     next = f->GetPrevSibling();
201     if (!next) { 
202       next = f->GetParent();
203       if (next) fList.CloseItem(next);
204     }
205     break;
206   case kKey_Down:
207     next = f->GetNextSibling();
208     if (!next && f->GetParent()) {
209       next = f->GetParent()->GetNextSibling();
210       fList.CloseItem(f->GetParent());
211     }
212     break;
213   case kKey_Left:
214     next = f->GetParent();
215     if (next) fList.CloseItem(next);
216     break;
217   case kKey_Right:
218     next = f->GetFirstChild();
219     if (next) fList.OpenItem(f);
220     break;
221   case kKey_PageUp:
222     fList.PageUp(kTRUE);
223     next = fList.GetSelected();
224     break;
225   case kKey_PageDown:
226     fList.PageDown(kTRUE);
227     next = fList.GetSelected();
228     break;
229   }
230   // if (next) gClient->NeedRedraw(&fList);
231   if (next && next != f) {
232     fList.ClearHighlighted();
233     fList.SetSelected(next);
234     HandleEntry(next,0);
235   }
236   if (next) gClient->NeedRedraw(&fList);
237 }
238
239 //____________________________________________________________________
240 void
241 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::HandleEntry(TGListTreeItem* entry, Int_t /*id*/) 
242 {
243   // Handle selection of entries
244   TGListTreeItem* old = fCurrentEntry;
245   if (entry) {
246     if (!entry->GetUserData()) return;
247     fCurrentEntry = entry;
248   }
249   else {
250     fCurrentEntry = 0;
251     ClearCanvas();
252   }
253   if (old != fCurrentEntry && fCanvas) fCanvas->cd();
254   SelectionChanged();
255 }
256
257 //____________________________________________________________________
258 void
259 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::UpdateList() 
260 {
261   // Update list
262   gClient->NeedRedraw(&fList);
263 }
264
265 //____________________________________________________________________
266 void
267 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::UpdateCanvas() 
268 {
269   // update canvas
270   if (!fCanvas) return;
271   fCanvas->Modified();
272   fCanvas->Update();
273   fCanvas->cd();
274 }
275 //____________________________________________________________________
276 TObject* AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::Current() const
277 {
278   // Get currently selected entry if any
279   if (!fCurrentEntry) return 0;
280   if (!fCurrentEntry->GetUserData()) return 0;
281   return static_cast<TObject*>(fCurrentEntry->GetUserData());
282 }
283   
284 //__________________________________________________________________
285 AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector& 
286 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::GetOrAdd(UShort_t id) 
287
288   // Get or add a sub-element
289   Int_t     idx = id - 1;
290   AliFMDSpectraDisplayDetector* d   = 0;
291   if (fChildren.GetEntriesFast() <= idx ||
292       !(d = static_cast<AliFMDSpectraDisplayDetector*>(fChildren.At(idx))))  { 
293     d = new AliFMDSpectraDisplayDetector(id, *this);
294     fChildren.AddAtAndExpand(d, idx);
295     // GetTop().GetList().SortChildren(&fEntry);
296     // GetTop().GetList().Sort(&(d->GetEntry()));
297   }
298   return *d;
299 }
300 //__________________________________________________________________
301 void 
302 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::Fill(UShort_t det, Char_t ring, 
303                                                     UShort_t sec, UShort_t str,
304                                                     Double_t v) 
305
306   AliFMDSpectraDisplayDetector& d = GetOrAdd(det);
307   d.Fill(ring, sec, str, v);
308   DoFill(v);
309 }
310 //__________________________________________________________________
311 Int_t
312 AliFMDSpectraDisplay::AliFMDSpectraDisplayTop::Compare(const TObject*) const
313 {
314   // Compare to another object
315   return -1;
316 }
317 //==================================================================
318 AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector::AliFMDSpectraDisplayDetector(UShort_t det, 
319                                            AliFMDSpectraDisplayTop& tree) 
320   : AliFMDSpectraDisplayElement(Form("FMD%d", det), "FMD Sub-detector"), 
321     fId(det), 
322     fParent(tree),
323     fChildren(0),
324     fEntry(*(tree.GetList().AddItem(&(tree.GetEntry()), GetName())))
325 {
326   // Constructor
327   fEntry.SetUserData(this);
328   fEntry.SetText(GetName());
329   if (GetTop().GetAxis()) MakeHistograms(GetTop().GetAxis());
330 }
331 //__________________________________________________________________
332 AliFMDSpectraDisplay::AliFMDSpectraDisplayRing& 
333 AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector::GetOrAdd(Char_t id) 
334
335   // Get or add an element
336   Int_t idx = (id == 'I' || id == 'i') ? 0 : 1;
337   AliFMDSpectraDisplayRing* r   = 0;;
338   if (fChildren.GetEntriesFast() <= idx ||
339       !(r = static_cast<AliFMDSpectraDisplayRing*>(fChildren.At(idx))))  { 
340     r = new AliFMDSpectraDisplayRing(id, *this);
341     fChildren.AddAtAndExpand(r, idx);
342     // GetTop().GetList().SortChildren(&fEntry);
343     // GetTop().GetList().Sort(&(r->GetEntry()));
344   }
345   return *r;
346 }
347 //__________________________________________________________________
348 void 
349 AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector::Fill(Char_t ring, 
350                                                          UShort_t sec, 
351                                                          UShort_t str, 
352                                                          Double_t v) 
353
354   // Fill values
355   AliFMDSpectraDisplayRing& r = GetOrAdd(ring);
356   r.Fill(sec, str, v);
357   DoFill(v);
358 }  
359 //__________________________________________________________________
360 Int_t
361 AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector::Compare(const TObject* o) const
362 {
363   // Compare to other element
364   std::cout << "Comparing detector to a " << o->ClassName() << std::endl;
365   if (o->IsA() == AliFMDSpectraDisplay::AliFMDSpectraDisplayDetector::Class()) { 
366     const AliFMDSpectraDisplayDetector* ro = 
367       static_cast<const AliFMDSpectraDisplayDetector*>(o);
368     return (Id() <  ro->Id() ? -1 : 
369             Id() == ro->Id() ? 0 : 1);
370   }
371   return -1;
372 }
373 //==================================================================
374 AliFMDSpectraDisplay::AliFMDSpectraDisplayRing::AliFMDSpectraDisplayRing(Char_t id,
375                             AliFMDSpectraDisplayDetector& d) 
376   : AliFMDSpectraDisplayElement(Form("FMD%d%c", d.Id(), id), "FMD Ring"), 
377     fParent(d),
378     fId(id),
379     fChildren(0),
380     fEntry(*(GetTop().GetList().AddItem(&(d.GetEntry()), GetName(), this)))
381 {
382   // Constructor
383   fEntry.SetText(GetName());
384   if (GetTop().GetAxis()) MakeHistograms(GetTop().GetAxis());
385 }
386 //__________________________________________________________________
387 AliFMDSpectraDisplay::AliFMDSpectraDisplaySector& 
388 AliFMDSpectraDisplay::AliFMDSpectraDisplayRing::GetOrAdd(UShort_t id) 
389
390   // Get or add another element
391   AliFMDSpectraDisplaySector* s = 0;
392   if (fChildren.GetEntriesFast() <= id ||
393       !(s = static_cast<AliFMDSpectraDisplaySector*>(fChildren.At(id)))) {
394     s = new AliFMDSpectraDisplaySector(id, *this);
395     fChildren.AddAtAndExpand(s, id);
396     // GetTop().GetList().SortChildren(&fEntry);
397     // GetTop().GetList().Sort(&(s->GetEntry()));
398   }
399   return *s;
400 }
401 //__________________________________________________________________
402 void 
403 AliFMDSpectraDisplay::AliFMDSpectraDisplayRing::Fill(UShort_t sec, 
404                                                      UShort_t str, 
405                                                      Double_t v) 
406
407   // Fill values
408   AliFMDSpectraDisplaySector& s = GetOrAdd(sec);
409   s.Fill(str, v);
410   DoFill(v);
411 }
412 //__________________________________________________________________
413 Int_t
414 AliFMDSpectraDisplay::AliFMDSpectraDisplayRing::Compare(const TObject* o) const
415 {
416   // Compare to other element
417   std::cout << "Comparing ring to a " << o->ClassName() << std::endl;
418   if (o->IsA() == AliFMDSpectraDisplay::AliFMDSpectraDisplayRing::Class()) { 
419     const AliFMDSpectraDisplayRing* ro = 
420       static_cast<const AliFMDSpectraDisplayRing*>(o);
421     return (Id() <  ro->Id() ? -1 : 
422             Id() == ro->Id() ? 0 : 1);
423   }
424   return -1;
425 }
426 //==================================================================
427 AliFMDSpectraDisplay::AliFMDSpectraDisplaySector::AliFMDSpectraDisplaySector(UShort_t id, 
428                                       AliFMDSpectraDisplayRing& r) 
429   : AliFMDSpectraDisplayElement(Form("FMD%d%c_%02d",r.DetectorId(),r.Id(),id), 
430                                 "FMD Sector"), 
431     fParent(r),
432     fId(id),
433     fChildren(0),
434     fEntry(*(GetTop().GetList().AddItem(&(r.GetEntry()), GetName(), this)))
435 {
436   // Constructor 
437   fEntry.SetText(GetName());
438   if (GetTop().GetAxis()) MakeHistograms(GetTop().GetAxis());
439 }
440 //__________________________________________________________________
441 AliFMDSpectraDisplay::AliFMDSpectraDisplayStrip& 
442 AliFMDSpectraDisplay::AliFMDSpectraDisplaySector::GetOrAdd(UShort_t id) 
443
444   // Get or add another element
445   AliFMDSpectraDisplayStrip* s = 0;
446   if (fChildren.GetEntriesFast() <= id || 
447       !(s = static_cast<AliFMDSpectraDisplayStrip*>(fChildren.At(id)))) {
448     s = new AliFMDSpectraDisplayStrip(id, *this);
449     fChildren.AddAtAndExpand(s, id);
450     // GetTop().GetList().SortChildren(&fEntry);
451     // GetTop().GetList().Sort(&(s->GetEntry()));
452   }
453   return *s;
454 }
455 //__________________________________________________________________
456 void 
457 AliFMDSpectraDisplay::AliFMDSpectraDisplaySector::Fill(UShort_t str, 
458                                                        Double_t v) 
459
460   // Fill values
461   AliFMDSpectraDisplayStrip& s = GetOrAdd(str);
462   s.Fill(v);
463   DoFill(v);
464 }
465 //__________________________________________________________________
466 Int_t
467 AliFMDSpectraDisplay::AliFMDSpectraDisplaySector::Compare(const TObject* o) const
468 {
469   // Compare to another elemnt
470   std::cout << "Comparing sector to a " << o->ClassName() << std::endl;
471   if (o->IsA() == AliFMDSpectraDisplay::AliFMDSpectraDisplaySector::Class()) { 
472     const AliFMDSpectraDisplaySector* ro = 
473       static_cast<const AliFMDSpectraDisplaySector*>(o);
474     return (Id() <  ro->Id() ? -1 : 
475             Id() == ro->Id() ? 0 : 1);
476   }
477   return -1;
478 }
479 //==================================================================
480 AliFMDSpectraDisplay::AliFMDSpectraDisplayStrip::AliFMDSpectraDisplayStrip(UShort_t id, 
481                                   AliFMDSpectraDisplaySector& s) 
482   : AliFMDSpectraDisplayElement(Form("FMD%d%c_%02d_%03d", 
483                                      s.DetectorId(), s.RingId(), 
484                                      s.Id(), id), "FMD Strip"), 
485     fParent(s),
486     fId(id),
487     fEntry(*(GetTop().GetList().AddItem(&(s.GetEntry()), GetName(), this)))
488 {
489   // Constructor
490   fEntry.SetText(GetName());
491   fEntry.SetPictures(GetTop().GetH1Pic(), GetTop().GetH1Pic());
492   if (GetTop().GetAxis()) MakeHistograms(GetTop().GetAxis());
493 }
494 //__________________________________________________________________
495 void 
496 AliFMDSpectraDisplay::AliFMDSpectraDisplayStrip::Fill(Double_t v) 
497
498   // Fill values
499   DoFill(v);
500 }
501 //__________________________________________________________________
502 Int_t
503 AliFMDSpectraDisplay::AliFMDSpectraDisplayStrip::Compare(const TObject* o) const
504 {
505   // Compare to another element
506   std::cout << "Comparing strip to a " << o->ClassName() << std::endl;
507   if (o->IsA() == AliFMDSpectraDisplay::AliFMDSpectraDisplayStrip::Class()) { 
508     const AliFMDSpectraDisplayStrip* ro = 
509       static_cast<const AliFMDSpectraDisplayStrip*>(o);
510     return (Id() <  ro->Id() ? -1 : 
511             Id() == ro->Id() ? 0 : 1);
512   }
513   return -1;
514 }
515
516 //==================================================================
517 AliFMDSpectraDisplay::AliFMDSpectraDisplay()
518   : AliFMDPattern(), 
519     fSelector(gClient->GetRoot(), 100, 100), 
520     fTop(fSelector, fAux)
521 {
522   // Constructor
523   // AddLoad(AliFMDInput::kRaw);
524   SetName("RAW");
525   SetTitle("RAW");
526   
527   fTop.Connect("SelectionChanged()", 
528                "AliFMDSpectraDisplay", this, "HandleDraw()");
529   
530   fSelector.MapSubwindows();
531   fSelector.Resize(fSelector.GetDefaultSize());
532   fSelector.MapWindow();
533 }
534
535 //__________________________________________________________________
536 Bool_t 
537 AliFMDSpectraDisplay::HandleDraw()
538 {
539   // Handle draw request
540   TObject* user = fTop.Current();
541   if (!user) return kFALSE;
542   if (!user->InheritsFrom(AliFMDSpectraDisplay::AliFMDSpectraDisplayElement::Class())) { 
543     Warning("HandleDraw", "%s does not inherit from Spectra::Element", 
544             user->GetName());
545     return kFALSE;
546   }
547   fAux->cd();
548   AliFMDSpectraDisplayElement* e 
549     = static_cast<AliFMDSpectraDisplayElement*>(user);
550   e->Show("hist", fSlider->GetMinimum(), fSlider->GetMaximum());
551   fAux->Modified();
552   fAux->Update();
553   fAux->cd();
554   return kTRUE;
555 }
556 //__________________________________________________________________
557 void 
558 AliFMDSpectraDisplay::MakeAux()
559 {
560   // MAke auxilary canvas
561   AliFMDPattern::MakeAux();
562   if (!fAux) return;
563   fTop.SetAxis(fSpec->GetXaxis());
564 }
565 //__________________________________________________________________
566 void 
567 AliFMDSpectraDisplay::DrawAux()
568 {
569   // Draw in the Aux the canvas 
570   // For example draw the spectra 
571   // or such stuff 
572   if (fTop.Current() != &fTop && HandleDraw()) return;
573   AliFMDPattern::DrawAux();
574 }
575 //__________________________________________________________________
576 Bool_t 
577 AliFMDSpectraDisplay::ProcessHit(AliFMDHit* hit, TParticle* p) 
578 {
579   // Process a hit
580   fTop.Fill(hit->Detector(), 
581             hit->Ring(), 
582             hit->Sector(),
583             hit->Strip(), 
584             hit->Edep());
585   return AliFMDPattern::ProcessHit(hit, p);
586 }
587 //__________________________________________________________________
588 Bool_t 
589 AliFMDSpectraDisplay::ProcessDigit(AliFMDDigit* digit)
590 {
591   // Process a digit
592   fTop.Fill(digit->Detector(), 
593             digit->Ring(), 
594             digit->Sector(),
595             digit->Strip(), 
596             digit->Counts());
597   return AliFMDDisplay::ProcessDigit(digit);
598 }
599 //__________________________________________________________________
600 Bool_t 
601 AliFMDSpectraDisplay::ProcessSDigit(AliFMDSDigit* sdigit)
602 {
603   // Process a summable digit
604   fTop.Fill(sdigit->Detector(), 
605             sdigit->Ring(), 
606             sdigit->Sector(),
607             sdigit->Strip(), 
608             sdigit->Counts());
609   return AliFMDDisplay::ProcessSDigit(sdigit);
610 }
611 //__________________________________________________________________
612 Bool_t 
613 AliFMDSpectraDisplay::ProcessRawDigit(AliFMDDigit* digit)
614 {
615   // Process a raw digit
616   return ProcessDigit(digit);
617 }
618 //__________________________________________________________________
619 Bool_t 
620 AliFMDSpectraDisplay::ProcessRecPoint(AliFMDRecPoint* recpoint)
621 {
622   // Process a rec-point
623   fTop.Fill(recpoint->Detector(), 
624             recpoint->Ring(), 
625             recpoint->Sector(),
626             recpoint->Strip(), 
627             recpoint->Particles());
628   return AliFMDDisplay::ProcessRecPoint(recpoint);
629 }
630 //__________________________________________________________________
631 Bool_t 
632 AliFMDSpectraDisplay::ProcessESD(UShort_t det, Char_t rng, UShort_t sec, 
633                                  UShort_t str, Float_t x, Float_t mult)
634 {
635   // Process ESD entry
636   fTop.Fill(det, rng, sec, str, mult);
637   return AliFMDDisplay::ProcessESD(det, rng, sec, str, x, mult);
638 }
639 //__________________________________________________________________
640 //
641 // EOF
642 //