]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDDisplay.cxx
Use Alignable volume names
[u/mrichter/AliRoot.git] / FMD / AliFMDDisplay.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id$ */
16 /** @file    AliFMDDisplay.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:39:09 2006
19     @brief   FMD Event display 
20 */
21 //___________________________________________________________________
22 //
23 // The classes defined here, are utility classes for reading in data
24 // for the FMD.  They are  put in a seperate library to not polute the
25 // normal libraries.  The classes are intended to be used as base
26 // classes for customized class that do some sort of analysis on the
27 // various types of data produced by the FMD. 
28 //
29 // Latest changes by Christian Holm Christensen
30 //
31 #include "AliFMDDisplay.h"      // ALIFMDDISPLAY_H
32 #include "AliFMDHit.h"          // ALIFMDHIT_H
33 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
34 #include "AliFMDRecPoint.h"     // ALIFMDRECPOINT_H
35 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
36 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
37 #include <AliESDFMD.h>          // ALIESDFMD_H
38 #include <AliLog.h>
39 #include <TStyle.h>
40 // #include <TArrayF.h>
41 #include <TMarker3DBox.h>
42 #include <TGeoManager.h>
43 // #include <TMath.h>
44 #include <TApplication.h>
45 #include <TButton.h>
46 // #include <TParticle.h>
47 #include <TCanvas.h>
48 #include <TView.h>
49 #include <TVirtualX.h>
50 #include <TSlider.h>
51 #include <TH1D.h>
52
53 //____________________________________________________________________
54 ClassImp(AliFMDDisplay)
55 #if 0
56   ; // This is here to keep Emacs for indenting the next line
57 #endif
58
59 //____________________________________________________________________
60 AliFMDDisplay* AliFMDDisplay::fgInstance = 0;
61
62 //____________________________________________________________________
63 AliFMDDisplay* 
64 AliFMDDisplay::Instance()
65 {
66   // Return static instance 
67   return fgInstance;
68 }
69
70 //____________________________________________________________________
71 AliFMDDisplay::~AliFMDDisplay()
72 {
73   if (fMarkers) {
74     fMarkers->Delete();
75     delete fMarkers;
76   }
77   if (fHits) {
78     fHits->Clear();
79     delete fHits;
80   }
81   if (fPad)     delete fPad;
82   fButtons.Delete();
83   if (fSlider)  delete fSlider;
84   if (fCanvas)  delete fCanvas;
85 }
86   
87 //____________________________________________________________________
88 AliFMDDisplay::AliFMDDisplay(Bool_t onlyFMD, const char* gAliceFile)
89   : AliFMDInput(gAliceFile),
90     fWait(kFALSE),
91     fMarkers(0),
92     fHits(0),
93     fCanvas(0), 
94     fPad(0), 
95     fSlider(0),
96     fZoomMode(kFALSE),
97     fX0(0),
98     fY0(0),
99     fX1(0),
100     fY1(0),
101     fMultCut(0),
102     fPedestalFactor(0),
103     fXPixel(0),
104     fYPixel(0),
105     fOldXPixel(0),
106     fOldYPixel(0),
107     fLineDrawn(0),
108     fOnlyFMD(onlyFMD)
109 {
110   // Constructor of an FMD display object. 
111   AddLoad(kGeometry);
112   if (fgInstance) delete fgInstance;
113   fgInstance = this;
114   SetMultiplicityCut();
115   SetPedestalFactor();
116 }
117
118 //____________________________________________________________________
119 void           
120 AliFMDDisplay::MakeCanvas(const char** which)
121 {
122   gStyle->SetPalette(1);
123   Double_t y1 = .10;
124   Int_t    w  = 700;
125   fCanvas = new TCanvas("display", "Display", w, Int_t(w / (1-y1)));
126   fCanvas->SetFillColor(1);
127   fCanvas->ToggleEventStatus();
128   fCanvas->cd();
129   fPad = new TPad("view", "3DView", 0.0, y1, 1.0, 1.0, 1, 0, 0);
130   fPad->Draw();
131
132   Double_t yb = 0;
133   fCanvas->cd();
134   if (TESTBIT(fTreeMask, kESD) || 
135       TESTBIT(fTreeMask, kDigits) || 
136       TESTBIT(fTreeMask, kRaw)) {
137     yb = .05;
138     fSlider = new TSlider("multCut", "Multiplicity cut", 0, 0, 1, yb);
139     fSlider->SetMethod("AliFMDDisplay::Instance()->ChangeCut()");
140     fSlider->SetMinimum(TESTBIT(fTreeMask, kESD) ? fMultCut :
141                         fPedestalFactor * 10);
142     fSlider->Draw();
143   }
144   const char** p = which;
145   const char*  m;
146   Int_t        n  = 0;
147   Int_t        j  = 0;
148   while (*(p++)) n++;
149   AliInfo(Form("Got %d buttons", n));
150   Float_t      x0 = 0;
151   Float_t      dx = 1. / n;
152   p               = which;
153   while ((m = *(p++))) {
154     fCanvas->cd();
155     AliInfo(Form("Adding button %s", m));
156     TButton* b = new TButton(m, Form("AliFMDDisplay::Instance()->%s()", m),
157                              x0, yb, x0 + dx, y1);
158     b->Draw();
159     fButtons.Add(b);
160     x0 += dx;
161     j++;
162   }
163 }
164
165 //____________________________________________________________________
166 void           
167 AliFMDDisplay::ShowOnlyFMD()
168 {
169   if (!fGeoManager) return;
170   static bool once = false;
171   if (once) return;
172   once = true;
173   AliInfo("Will only show the FMD");
174   TGeoVolume* top = gGeoManager->GetTopVolume();
175   top->InvisibleAll(kTRUE);
176   TGeoIterator next(top);
177   TGeoNode* node;
178   TGeoVolume* v = 0;
179   Bool_t hasFMD1 = kFALSE;
180   Bool_t hasFMD2 = kFALSE;
181   Bool_t hasFMD3 = kFALSE;
182   TObjArray toshow;
183   while ((node = static_cast<TGeoNode*>(next()))) {
184     const char* name = node->GetName();
185     if (!name) continue;
186     if (!(v = node->GetVolume())) continue;
187
188     if (name[0] == 'F') {
189       if (name[2] == 'M' && (name[3] == 'T' || name[3] == 'B')) {
190         // Virtual Master half-ring volume - top-level
191         Int_t det = node->GetNumber();
192         switch (det) {
193         case 1: hasFMD1 = true; break;
194         case 2: hasFMD2 = true; break;
195         case 3: hasFMD3 = true; break;
196         default: continue;
197         }
198         toshow.Add(v);
199       }
200       else if (name[3] == 'V' && (name[2] == 'T' || name[2] == 'B')) 
201         toshow.Add(v); // Virtual Half-ring, bare detectors
202       // else if (name[3] == 'H' && (name[2] == 'F' || name[2] == 'B')) 
203       //  toshow.Add(v); // Virtual Hybrid container 
204     }
205     v->SetVisibility(kFALSE);
206     v->SetVisDaughters(kFALSE);
207     v->InvisibleAll(kTRUE);
208   }
209   TIter i(&toshow);
210   while ((v = static_cast<TGeoVolume*>(i()))) {
211     v->SetVisibility(kTRUE);
212     v->SetVisDaughters(kTRUE);
213     v->InvisibleAll(kFALSE);
214   }  
215 }
216
217     
218 //____________________________________________________________________
219 void           
220 AliFMDDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py) 
221 {
222   // AliInfo(Form("Event %d, at (%d,%d)", px, py));
223   if (px == 0 && py == 0) return;
224   if (!fZoomMode && fPad->GetView()) {
225     fPad->GetView()->ExecuteRotateView(event, px, py);
226     return;
227   }
228   fPad->SetCursor(kCross);
229   switch (event) {
230   case kButton1Down: 
231     fPad->TAttLine::Modify();
232     fX0        = fPad->AbsPixeltoX(px);
233     fY0        = fPad->AbsPixeltoY(py);
234     fXPixel    = fOldXPixel = px;
235     fYPixel    = fOldYPixel = py;
236     fLineDrawn = kFALSE;
237     return;
238   case kButton1Motion:
239     if (fLineDrawn) 
240       gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
241                          TVirtualX::kHollow);
242     fOldXPixel = px;
243     fOldYPixel = py;
244     fLineDrawn = kTRUE;
245     gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
246                        TVirtualX::kHollow);
247     return;
248   case kButton1Up:
249     fPad->GetCanvas()->FeedbackMode(kFALSE);
250     if (px == fXPixel || py == fYPixel) return;
251     fX1 = fPad->AbsPixeltoX(px);
252     fY1 = fPad->AbsPixeltoY(py);
253     if (fX1 < fX0) std::swap(fX0, fX1); 
254     if (fY1 < fY0) std::swap(fY0, fY1); 
255     fPad->Range(fX0, fY0, fX1, fY1);
256     fPad->Modified();
257     return;
258   }
259 }
260
261 //____________________________________________________________________
262 Int_t          
263 AliFMDDisplay::DistancetoPrimitive(Int_t px, Int_t) 
264 {
265   // AliInfo(Form("@ (%d,%d)", px, py));
266   fPad->SetCursor(kCross);
267   Float_t xmin = fPad->GetX1();
268   Float_t xmax = fPad->GetX2();
269   Float_t dx   = .02 * (xmax - xmin);
270   Float_t x    = fPad->AbsPixeltoX(px);
271   if (x < xmin + dx || x > xmax - dx) return 9999;
272   return (fZoomMode ? 0 : 7);
273 }
274 //____________________________________________________________________
275 Bool_t 
276 AliFMDDisplay::Init()
277 {
278   // Initialize.  GEt transforms and such, 
279   if (!AliFMDInput::Init()) return kFALSE;
280   AliFMDGeometry* geom = AliFMDGeometry::Instance();
281   geom->Init();
282   geom->InitTransformations();
283   fMarkers = new TObjArray;
284   fHits    = new TObjArray;
285   fMarkers->SetOwner(kTRUE);
286   fHits->SetOwner(kFALSE);
287   return kTRUE;
288 }
289
290 //____________________________________________________________________
291 void
292 AliFMDDisplay::MakeAux()
293 {
294   if ((TESTBIT(fTreeMask, kESD) || 
295        TESTBIT(fTreeMask, kDigits) || 
296        TESTBIT(fTreeMask, kRaw))) {
297     if (!fAux) {
298       fAux = new TCanvas("aux", "Aux");
299       fAux->SetLogy();
300       if (TESTBIT(fTreeMask, kESD)) 
301         fSpec = new TH1D("spec", "Mult spectra", 150, 0, 3);
302       else 
303         fSpec = new TH1D("spec", "Adc spectra", 1024, -.5, 1023.5);
304       fSpecCut = static_cast<TH1*>(fSpec->Clone("specCut"));
305       fSpec->SetFillColor(2);
306       fSpec->SetFillStyle(3001);
307       fSpecCut->SetFillColor(4);
308       fSpecCut->SetFillStyle(3001);
309     }
310     else {
311       fSpec->Reset();
312       fSpecCut->Reset();
313     }
314   }
315 }
316
317 //____________________________________________________________________
318 void
319 AliFMDDisplay::DrawAux()
320 {
321   if (!fAux) return;
322   fAux->cd();
323   fAux->Clear();
324   fSpec->Draw();
325   fSpecCut->Draw("same");
326   fAux->Modified();
327   fAux->Update();
328   fAux->cd();
329 }
330
331 //____________________________________________________________________
332 Bool_t 
333 AliFMDDisplay::Begin(Int_t event) 
334 {
335   // Begin of event.  Make canvas is not already done 
336   if (!fCanvas) {
337     const char* m[] = { "Continue", "Zoom", "Pick", "Redisplay", 0 }; 
338     MakeCanvas(m);
339     MakeAux();
340   }
341   // AliInfo("Clearing canvas");
342   // fCanvas->Clear();
343   if (!fGeoManager) {
344     Warning("End", "No geometry manager");
345     return kFALSE;
346   }
347   AliInfo("Drawing geometry");
348   fPad->cd();
349   fGeoManager->GetTopVolume()->Draw();
350   if (fOnlyFMD) ShowOnlyFMD();
351   AliInfo("Adjusting view");
352   Int_t irep;
353   if (fPad->GetView()) {
354     fPad->GetView()->SetView(-200, -40, 80, irep);
355     fPad->GetView()->Zoom();
356     fPad->Modified();
357     fPad->cd();
358   }
359   return AliFMDInput::Begin(event);
360 }
361
362 //____________________________________________________________________
363 void
364 AliFMDDisplay::AtEnd() 
365 {
366   fPad->cd();
367   fMarkers->Draw();
368   fPad->cd();
369   AppendPad();
370   fPad->cd();
371   DrawAux();
372 }
373
374 //____________________________________________________________________
375 void
376 AliFMDDisplay::Idle() 
377 {
378   fWait = kTRUE;
379   while (fWait) {
380     gApplication->StartIdleing();
381     gSystem->InnerLoop();
382     gApplication->StopIdleing();
383   }
384   AliInfo("After idle loop");
385   if (fMarkers) fMarkers->Delete();
386   if (fHits)    fHits->Clear();
387   AliInfo("After clearing caches");
388 }
389
390 //____________________________________________________________________
391 Bool_t 
392 AliFMDDisplay::End()
393 {
394   // End of event.  Draw everything 
395   AtEnd();
396   Idle();
397   return AliFMDInput::End();
398 }
399
400 //____________________________________________________________________
401 Int_t
402 AliFMDDisplay::LookupColor(Float_t x, Float_t max) const
403 {
404   // Look-up color 
405   Int_t idx = Int_t(x / max * gStyle->GetNumberOfColors());
406   return gStyle->GetColorPalette(idx);
407 }
408
409 //____________________________________________________________________
410 void
411 AliFMDDisplay::ChangeCut() 
412 {
413   fMultCut        = fSlider->GetMinimum();
414   fPedestalFactor = fSlider->GetMinimum() * 10;
415   AliInfo(Form("Multiplicity cut: %7.5f, Pedestal factor: %7.4f (%6.5f)", 
416                fMultCut, fPedestalFactor, fSlider->GetMinimum()));
417   Redisplay();
418 }
419
420 //____________________________________________________________________
421 void
422 AliFMDDisplay::Redisplay()
423 {
424   if (fMarkers) fMarkers->Delete();
425   if (fHits)    fHits->Clear();
426   if (fSpec)    fSpec->Reset();
427   if (fSpecCut) fSpecCut->Reset();
428   Event();
429   AtEnd();
430 }
431
432 //____________________________________________________________________
433 void
434 AliFMDDisplay::AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
435                          TObject* o, Float_t s, Float_t max)
436 {
437   // Add a marker to the display
438   //
439   //    det     Detector
440   //    rng     Ring
441   //    sec     Sector 
442   //    str     Strip
443   //    o       Object to refer to
444   //    s       Signal 
445   //    max     Maximum of signal 
446   //
447   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
448   Double_t x, y, z;
449   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
450   Float_t  size  = .1;
451   Float_t  zsize = s / max * 10;
452   Float_t  r     = TMath::Sqrt(x * x + y * y);
453   Float_t  theta = TMath::ATan2(r, z);
454   Float_t  phi   = TMath::ATan2(y, x);
455   Float_t  rz    = z + (z < 0 ? 1 : -1) * zsize;
456   TMarker3DBox* marker = new  TMarker3DBox(x,y,rz,size,size,zsize,theta,phi);
457   if (o) marker->SetRefObject(o);
458   marker->SetLineColor(LookupColor(s, max));
459   fMarkers->Add(marker);
460 }
461   
462
463 //____________________________________________________________________
464 Bool_t 
465 AliFMDDisplay::ProcessHit(AliFMDHit* hit, TParticle* /* p */) 
466 {
467   // Process a hit 
468   if (!hit) { AliError("No hit");   return kFALSE; }
469   // if (!p)   { AliError("No track"); return kFALSE; }
470
471   if (fHits) fHits->Add(hit);
472   Float_t  size  = .1;
473   Float_t  zsize = TMath::Sqrt(hit->Edep() * 20);
474   Float_t  z     = hit->Z() + (hit->Z() < 0 ? 1 : -1) * zsize; 
475   Float_t  pt    = TMath::Sqrt(hit->Py()*hit->Py()+hit->Px()*hit->Px());
476   Float_t  theta = TMath::ATan2(pt, hit->Pz());
477   Float_t  phi   = TMath::ATan2(hit->Py(), hit->Px());
478   TMarker3DBox* marker = new  TMarker3DBox(hit->X(), hit->Y(), z,
479                                            size, size, zsize, theta, phi);
480   marker->SetLineColor(LookupColor(hit->Edep(), 1));
481   marker->SetRefObject(hit);
482   fMarkers->Add(marker);
483   return kTRUE;
484 }
485
486 //____________________________________________________________________
487 Bool_t 
488 AliFMDDisplay::ProcessDigit(AliFMDDigit* digit)
489 {
490   // Process a digit 
491   if (!digit) { AliError("No digit");   return kFALSE; }
492
493   AliFMDParameters* parm = AliFMDParameters::Instance();
494   UShort_t det           =  digit->Detector();
495   Char_t   ring          =  digit->Ring();
496   UShort_t sec           =  digit->Sector();
497   UShort_t str           =  digit->Strip();
498   Double_t ped           =  parm->GetPedestal(det,ring, sec, str);
499   Double_t pedW          =  parm->GetPedestalWidth(det,ring, sec, str);
500   Double_t threshold     =  ped * fPedestalFactor * pedW;
501   Float_t  counts        = digit->Counts();
502   if (fSpec) fSpec->Fill(counts);
503   if (counts < threshold) return kTRUE;
504   if (fHits) fHits->Add(digit);
505   if (fSpecCut) fSpecCut->Fill(counts);
506
507   AddMarker(det, ring, sec, str, digit, counts, 1024);
508   return kTRUE;
509 }
510
511 //____________________________________________________________________
512 Bool_t 
513 AliFMDDisplay::ProcessRaw(AliFMDDigit* digit)
514 {
515   // PRocess raw data 
516   return ProcessDigit(digit);
517 }
518
519 //____________________________________________________________________
520 Bool_t 
521 AliFMDDisplay::ProcessRecPoint(AliFMDRecPoint* recpoint)
522 {
523   // Process reconstructed point 
524   if (!recpoint) { AliError("No recpoint");   return kFALSE; }
525   if (recpoint->Particles() < fMultCut) return kTRUE;
526   if (fHits) fHits->Add(recpoint);
527   AddMarker(recpoint->Detector(), recpoint->Ring(), recpoint->Sector(),  
528             recpoint->Strip(), recpoint, recpoint->Particles(), 20);
529   return kTRUE;
530 }
531
532 //____________________________________________________________________
533 Bool_t 
534 AliFMDDisplay::ProcessESD(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
535                           Float_t eta, Float_t mult)
536 {
537   Double_t cmult = (mult * 
538                     TMath::Abs(TMath::Cos(2.*TMath::ATan(TMath::Exp(-eta)))));
539   if (fSpec) fSpec->Fill(cmult);
540   if (cmult < fMultCut || cmult == AliESDFMD::kInvalidMult) return kTRUE;
541   AddMarker(det,rng,sec,str, 0, cmult, 20);
542   if (fSpecCut) fSpecCut->Fill(cmult);
543   return kTRUE;
544 }
545
546 //____________________________________________________________________
547 //
548 // EOF
549 //