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