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