]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDDisplay.cxx
Misalignment-related bug fixed
[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
32 #include <TSystem.h>
33 #include <TApplication.h>
34 #include <TButton.h>
35 #include <TCanvas.h>
36 #include <TGeoManager.h>
37 #include <TH1D.h>
38 #include <TMarker3DBox.h>
39 #include <TMath.h>
40 #include <TSlider.h>
41 #include <TStyle.h>
42 #include <TView.h>
43 #include <TVirtualX.h>
44 // #include <TArrayF.h>
45 // #include <TParticle.h>
46
47 #include "AliFMDDisplay.h"      // ALIFMDDISPLAY_H
48 #include "AliFMDHit.h"          // ALIFMDHIT_H
49 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
50 #include "AliFMDRecPoint.h"     // ALIFMDRECPOINT_H
51 #include "AliFMDGeometry.h"     // ALIFMDGEOMETRY_H
52 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
53 #include <AliESDFMD.h>          // ALIESDFMD_H
54 // #include <AliLog.h>
55 #include "AliFMDDebug.h" // Better debug macros
56
57 //____________________________________________________________________
58 ClassImp(AliFMDDisplay)
59 #if 0
60   ; // This is here to keep Emacs for indenting the next line
61 #endif
62
63 //____________________________________________________________________
64 AliFMDDisplay* AliFMDDisplay::fgInstance = 0;
65
66 //____________________________________________________________________
67 AliFMDDisplay* 
68 AliFMDDisplay::Instance()
69 {
70   // Return static instance 
71   // If the instance does not exist
72   // it is not created!
73   return fgInstance;
74 }
75
76 //____________________________________________________________________
77 AliFMDDisplay::~AliFMDDisplay()
78 {
79   // Destructor. 
80   // Cleans 
81   // up 
82   if (fMarkers) {
83     fMarkers->Delete();
84     delete fMarkers;
85   }
86   if (fHits) {
87     fHits->Clear();
88     delete fHits;
89   }
90   if (fPad)     delete fPad;
91   fButtons.Delete();
92   if (fSlider)  delete fSlider;
93   if (fCanvas)  delete fCanvas;
94 }
95   
96 //____________________________________________________________________
97 AliFMDDisplay::AliFMDDisplay(Bool_t onlyFMD, const char* gAliceFile)
98   : AliFMDInput(gAliceFile),
99     fWait(kFALSE),
100     fMarkers(0),
101     fHits(0),
102     fCanvas(0), 
103     fPad(0), 
104     fSlider(0),
105     fZoomMode(kFALSE),
106     fX0(0),
107     fY0(0),
108     fX1(0),
109     fY1(0),
110     fMultCut(0),
111     fPedestalFactor(0),
112     fXPixel(0),
113     fYPixel(0),
114     fOldXPixel(0),
115     fOldYPixel(0),
116     fLineDrawn(0),
117     fOnlyFMD(onlyFMD)
118 {
119   // Constructor of an FMD display object. 
120   // Must be called 
121   // before Instance 
122   AddLoad(kGeometry);
123   if (fgInstance) delete fgInstance;
124   fgInstance = this;
125   SetMultiplicityCut();
126   SetPedestalFactor();
127 }
128
129 //____________________________________________________________________
130 void           
131 AliFMDDisplay::MakeCanvas(const char** which)
132 {
133   // Make a canvas 
134   // Parameters: 
135   //   which   Which button to put up. 
136   gStyle->SetPalette(1);
137   Double_t y1 = .10;
138   Int_t    w  = 700;
139   fCanvas = new TCanvas("display", "Display", w, Int_t(w / (1-y1)));
140   fCanvas->SetFillColor(1);
141   fCanvas->ToggleEventStatus();
142   fCanvas->cd();
143   fPad = new TPad("view", "3DView", 0.0, y1, 1.0, 1.0, 1, 0, 0);
144   fPad->Draw();
145
146   Double_t yb = 0;
147   fCanvas->cd();
148   if (TESTBIT(fTreeMask, kESD) || 
149       TESTBIT(fTreeMask, kDigits) || 
150       TESTBIT(fTreeMask, kRaw)) {
151     yb = .05;
152     fSlider = new TSlider("multCut", "Multiplicity cut", 0, 0, 1, yb);
153     fSlider->SetMethod("AliFMDDisplay::Instance()->ChangeCut()");
154     fSlider->Draw();
155     fSlider->SetMinimum(TESTBIT(fTreeMask, kESD) ? fMultCut * 10 :
156                         fPedestalFactor * 10);
157   }
158   const char** p = which;
159   const char*  m;
160   Int_t        n  = 0;
161   Int_t        j  = 0;
162   while (*(p++)) n++;
163   AliInfo(Form("Got %d buttons", n));
164   Float_t      x0 = 0;
165   Float_t      dx = 1. / n;
166   p               = which;
167   while ((m = *(p++))) {
168     fCanvas->cd();
169     AliInfo(Form("Adding button %s", m));
170     TButton* b = new TButton(m, Form("AliFMDDisplay::Instance()->%s()", m),
171                              x0, yb, x0 + dx, y1);
172     b->Draw();
173     fButtons.Add(b);
174     x0 += dx;
175     j++;
176   }
177 }
178
179 //____________________________________________________________________
180 void           
181 AliFMDDisplay::ShowOnlyFMD()
182 {
183   // Show only the FMD 
184   // Do not show 
185   // other volumes 
186   if (!fGeoManager) return;
187   static bool once = false;
188   if (once) return;
189   once = true;
190   AliInfo("Will only show the FMD");
191   TGeoVolume* top = gGeoManager->GetTopVolume();
192   top->InvisibleAll(kTRUE);
193   TGeoIterator next(top);
194   TGeoNode* node;
195   TGeoVolume* v = 0;
196   Bool_t hasFMD1 = kFALSE;
197   Bool_t hasFMD2 = kFALSE;
198   Bool_t hasFMD3 = kFALSE;
199   TObjArray toshow;
200   while ((node = static_cast<TGeoNode*>(next()))) {
201     const char* name = node->GetName();
202     if (!name) continue;
203     if (!(v = node->GetVolume())) continue;
204
205     if (name[0] == 'F') {
206       if (name[2] == 'M' && (name[3] == 'T' || name[3] == 'B')) {
207         // Virtual Master half-ring volume - top-level
208         Int_t det = node->GetNumber();
209         switch (det) {
210         case 1: hasFMD1 = true; break;
211         case 2: hasFMD2 = true; break;
212         case 3: hasFMD3 = true; break;
213         default: continue;
214         }
215         toshow.Add(v);
216       }
217       else if (name[3] == 'V' && (name[2] == 'T' || name[2] == 'B')) 
218         toshow.Add(v); // Virtual Half-ring, bare detectors
219       // else if (name[3] == 'H' && (name[2] == 'F' || name[2] == 'B')) 
220       //  toshow.Add(v); // Virtual Hybrid container 
221     }
222     v->SetVisibility(kFALSE);
223     v->SetVisDaughters(kFALSE);
224     v->InvisibleAll(kTRUE);
225   }
226   TIter i(&toshow);
227   while ((v = static_cast<TGeoVolume*>(i()))) {
228     v->SetVisibility(kTRUE);
229     v->SetVisDaughters(kTRUE);
230     v->InvisibleAll(kFALSE);
231   }  
232 }
233
234     
235 //____________________________________________________________________
236 void           
237 AliFMDDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py) 
238 {
239   // Execute an event on canvas 
240   // Parameters: 
241   //   event   What happened
242   //   px, py  Pixel coordinates 
243   if (px == 0 && py == 0) return;
244   if (!fZoomMode && fPad->GetView()) {
245     fPad->GetView()->ExecuteRotateView(event, px, py);
246     return;
247   }
248   fPad->SetCursor(kCross);
249   switch (event) {
250   case kButton1Down: 
251     fPad->TAttLine::Modify();
252     fX0        = fPad->AbsPixeltoX(px);
253     fY0        = fPad->AbsPixeltoY(py);
254     fXPixel    = fOldXPixel = px;
255     fYPixel    = fOldYPixel = py;
256     fLineDrawn = kFALSE;
257     return;
258   case kButton1Motion:
259     if (fLineDrawn) 
260       gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
261                          TVirtualX::kHollow);
262     fOldXPixel = px;
263     fOldYPixel = py;
264     fLineDrawn = kTRUE;
265     gVirtualX->DrawBox(fXPixel, fYPixel, fOldXPixel, fOldYPixel, 
266                        TVirtualX::kHollow);
267     return;
268   case kButton1Up:
269     fPad->GetCanvas()->FeedbackMode(kFALSE);
270     if (px == fXPixel || py == fYPixel) return;
271     fX1 = fPad->AbsPixeltoX(px);
272     fY1 = fPad->AbsPixeltoY(py);
273     if (fX1 < fX0) std::swap(fX0, fX1); 
274     if (fY1 < fY0) std::swap(fY0, fY1); 
275     fPad->Range(fX0, fY0, fX1, fY1);
276     fPad->Modified();
277     return;
278   }
279 }
280
281 //____________________________________________________________________
282 Int_t          
283 AliFMDDisplay::DistancetoPrimitive(Int_t px, Int_t) 
284 {
285   // Calculate the distance from point to 
286   // something in the canvas. 
287   // Depends on the zoom mode. 
288   fPad->SetCursor(kCross);
289   Float_t xmin = fPad->GetX1();
290   Float_t xmax = fPad->GetX2();
291   Float_t dx   = .02 * (xmax - xmin);
292   Float_t x    = fPad->AbsPixeltoX(px);
293   if (x < xmin + dx || x > xmax - dx) return 9999;
294   return (fZoomMode ? 0 : 7);
295 }
296 //____________________________________________________________________
297 Bool_t 
298 AliFMDDisplay::Init()
299 {
300   // Initialize.  GEt transforms and such,
301   // so that we can draw thins properly 
302   // Returns true on success 
303   if (!AliFMDInput::Init()) return kFALSE;
304   AliFMDGeometry* geom = AliFMDGeometry::Instance();
305   geom->Init();
306   geom->InitTransformations();
307   if (TESTBIT(fTreeMask, kDigits) || TESTBIT(fTreeMask, kRaw)) 
308     AliFMDParameters::Instance()->Init();
309
310   fMarkers = new TObjArray;
311   fHits    = new TObjArray;
312   fMarkers->SetOwner(kTRUE);
313   fHits->SetOwner(kFALSE);
314   return kTRUE;
315 }
316
317 //____________________________________________________________________
318 void
319 AliFMDDisplay::MakeAux()
320 {
321   // MAke the aux canvas 
322   // This is used to display spectra
323   // etc, 
324   if ((TESTBIT(fTreeMask, kESD) || 
325        TESTBIT(fTreeMask, kDigits) || 
326        TESTBIT(fTreeMask, kRaw))) {
327     if (!fAux) {
328       fAux = new TCanvas("aux", "Aux");
329       fAux->SetLogy();
330       if (TESTBIT(fTreeMask, kESD)) 
331         fSpec = new TH1D("spec", "Mult spectra", 500, 0, 10);
332       else 
333         fSpec = new TH1D("spec", "Adc spectra", 1024, -.5, 1023.5);
334       fSpecCut = static_cast<TH1*>(fSpec->Clone("specCut"));
335       fSpec->SetFillColor(2);
336       fSpec->SetFillStyle(3001);
337       fSpecCut->SetFillColor(4);
338       fSpecCut->SetFillStyle(3001);
339     }
340     else {
341       fSpec->Reset();
342       fSpecCut->Reset();
343     }
344   }
345 }
346
347 //____________________________________________________________________
348 void
349 AliFMDDisplay::DrawAux()
350 {
351   // Draw in the Aux the canvas 
352   // For example draw the spectra 
353   // or such stuff 
354   if (!fAux) return;
355   fAux->cd();
356   fAux->Clear();
357   fSpec->Draw();
358   fSpecCut->Draw("same");
359   fAux->Modified();
360   fAux->Update();
361   fAux->cd();
362 }
363
364 //____________________________________________________________________
365 Bool_t 
366 AliFMDDisplay::Begin(Int_t event) 
367 {
368   // Begin of event.  Make canvas is not already done 
369   // Parameters: 
370   //   event   The event number 
371   if (!fCanvas) {
372     const char* m[] = { "Continue", "Zoom", "Pick", "Redisplay", 0 }; 
373     MakeCanvas(m);
374   }
375   MakeAux();
376
377   // AliInfo("Clearing canvas");
378   // fCanvas->Clear();
379   if (!fGeoManager) {
380     Warning("End", "No geometry manager");
381     return kFALSE;
382   }
383   AliInfo("Drawing geometry");
384   fPad->cd();
385   fGeoManager->GetTopVolume()->Draw();
386   if (fOnlyFMD) ShowOnlyFMD();
387   AliInfo("Adjusting view");
388   Int_t irep;
389   if (fPad->GetView()) {
390     fPad->GetView()->SetView(-200, -40, 80, irep);
391     fPad->GetView()->Zoom();
392     fPad->Modified();
393     fPad->cd();
394   }
395   return AliFMDInput::Begin(event);
396 }
397
398 //____________________________________________________________________
399 void
400 AliFMDDisplay::AtEnd() 
401 {
402   // Called at of the event. 
403   // Draw stuff. 
404   // Draw spectrum. 
405   fPad->cd();
406   fMarkers->Draw();
407   fPad->cd();
408   AppendPad();
409   fPad->cd();
410   DrawAux();
411 }
412
413 //____________________________________________________________________
414 void
415 AliFMDDisplay::Idle() 
416 {
417   // Idle loop. 
418   // Sends the ROOT loop into the idle loop,
419   // so that we can go on. 
420   fWait = kTRUE;
421   while (fWait) {
422     gApplication->StartIdleing();
423     gSystem->InnerLoop();
424     gApplication->StopIdleing();
425   }
426   AliInfo("After idle loop");
427   if (fMarkers) fMarkers->Delete();
428   if (fHits)    fHits->Clear();
429   AliInfo("After clearing caches");
430 }
431
432 //____________________________________________________________________
433 Bool_t 
434 AliFMDDisplay::End()
435 {
436   // End of event.  Draw everything 
437   AtEnd();
438   Idle();
439   return AliFMDInput::End();
440 }
441
442 //____________________________________________________________________
443 Int_t
444 AliFMDDisplay::LookupColor(Float_t x, Float_t max) const
445 {
446   // Look-up color.  
447   // Get a colour from the  current palette depending 
448   // on the ratio x/max 
449   Int_t idx = Int_t(x / max * gStyle->GetNumberOfColors());
450   return gStyle->GetColorPalette(idx);
451 }
452
453 //____________________________________________________________________
454 void
455 AliFMDDisplay::ChangeCut() 
456 {
457   // Change the cut on the slider. 
458   // The factor depends on what is 
459   // drawn in the AUX canvas
460   fMultCut        = fSlider->GetMinimum() * 10;
461   fPedestalFactor = fSlider->GetMinimum() * 10;
462   AliInfo(Form("Multiplicity cut: %7.5f, Pedestal factor: %7.4f (%6.5f)", 
463                fMultCut, fPedestalFactor, fSlider->GetMinimum()));
464   Redisplay();
465 }
466
467 //____________________________________________________________________
468 void
469 AliFMDDisplay::Redisplay()
470 {
471   // Redisplay stuff. 
472   // Redraw markers, hits, 
473   // spectra 
474   if (fMarkers) fMarkers->Delete();
475   if (fHits)    fHits->Clear();
476   if (fSpec)    fSpec->Reset();
477   if (fSpecCut) fSpecCut->Reset();
478   Event();
479   AtEnd();
480 }
481
482 //____________________________________________________________________
483 void
484 AliFMDDisplay::AddMarker(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
485                          TObject* o, Float_t s, Float_t max)
486 {
487   // Add a marker to the display
488   //
489   //    det     Detector
490   //    rng     Ring
491   //    sec     Sector 
492   //    str     Strip
493   //    o       Object to refer to
494   //    s       Signal 
495   //    max     Maximum of signal 
496   //
497   AliFMDGeometry*   geom = AliFMDGeometry::Instance();
498   Double_t x, y, z;
499   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
500   Float_t  size  = .1;
501   Float_t  zsize = s / max * 10;
502   Float_t  r     = TMath::Sqrt(x * x + y * y);
503   Float_t  theta = TMath::ATan2(r, z);
504   Float_t  phi   = TMath::ATan2(y, x);
505   Float_t  rz    = z + (z < 0 ? 1 : -1) * zsize;
506   TMarker3DBox* marker = new  TMarker3DBox(x,y,rz,size,size,zsize,theta,phi);
507   if (o) marker->SetRefObject(o);
508   marker->SetLineColor(LookupColor(s, max));
509   fMarkers->Add(marker);
510 }
511   
512
513 //____________________________________________________________________
514 Bool_t 
515 AliFMDDisplay::ProcessHit(AliFMDHit* hit, TParticle* /* p */) 
516 {
517   // Process a hit. 
518   // Parameters: 
519   //   hit   Hit data
520   if (!hit) { AliError("No hit");   return kFALSE; }
521   // if (!p)   { AliError("No track"); return kFALSE; }
522
523   if (fHits) fHits->Add(hit);
524   Float_t  size  = .1;
525   Float_t  zsize = TMath::Sqrt(hit->Edep() * 20);
526   Float_t  z     = hit->Z() + (hit->Z() < 0 ? 1 : -1) * zsize; 
527   Float_t  pt    = TMath::Sqrt(hit->Py()*hit->Py()+hit->Px()*hit->Px());
528   Float_t  theta = TMath::ATan2(pt, hit->Pz());
529   Float_t  phi   = TMath::ATan2(hit->Py(), hit->Px());
530   TMarker3DBox* marker = new  TMarker3DBox(hit->X(), hit->Y(), z,
531                                            size, size, zsize, theta, phi);
532   marker->SetLineColor(LookupColor(hit->Edep(), 1));
533   marker->SetRefObject(hit);
534   fMarkers->Add(marker);
535   return kTRUE;
536 }
537
538 //____________________________________________________________________
539 Bool_t 
540 AliFMDDisplay::ProcessDigit(AliFMDDigit* digit)
541 {
542   // Process a digit 
543   // Parameters: 
544   //   digit Digit information 
545   if (!digit) { AliError("No digit");   return kFALSE; }
546
547   AliFMDParameters* parm = AliFMDParameters::Instance();
548   UShort_t det           =  digit->Detector();
549   Char_t   ring          =  digit->Ring();
550   UShort_t sec           =  digit->Sector();
551   UShort_t str           =  digit->Strip();
552   Double_t ped           =  parm->GetPedestal(det,ring, sec, str);
553   Double_t pedW          =  parm->GetPedestalWidth(det,ring, sec, str);
554   Double_t threshold     =  ped + fPedestalFactor * pedW;
555   Float_t  counts        =  digit->Counts();
556   AliFMDDebug(10, ("FMD%d%c[%2d,%3d] ADC: %d > %d (=%4.2f+%4.2f*%4.2f)", 
557                     digit->Detector(), digit->Ring(), digit->Sector(), 
558                     digit->Strip(), Int_t(counts), Int_t(threshold), 
559                     ped, fPedestalFactor, pedW));
560   if (fSpec) fSpec->Fill(counts);
561   if (counts < threshold) return kTRUE;
562   if (fHits) fHits->Add(digit);
563   if (fSpecCut) fSpecCut->Fill(counts);
564
565   AddMarker(det, ring, sec, str, digit, counts, 1024);
566   return kTRUE;
567 }
568
569 //____________________________________________________________________
570 Bool_t 
571 AliFMDDisplay::ProcessRaw(AliFMDDigit* digit)
572 {
573   // PRocess raw data 
574   // Parameters: 
575   //   digit Digit information 
576   return ProcessDigit(digit);
577 }
578
579 //____________________________________________________________________
580 Bool_t 
581 AliFMDDisplay::ProcessRecPoint(AliFMDRecPoint* recpoint)
582 {
583   // Process reconstructed point 
584   // Parameters: 
585   //  recpoint  Reconstructed multiplicity/energy
586   if (!recpoint) { AliError("No recpoint");   return kFALSE; }
587   if (recpoint->Particles() < fMultCut) return kTRUE;
588   if (fHits) fHits->Add(recpoint);
589   AddMarker(recpoint->Detector(), recpoint->Ring(), recpoint->Sector(),  
590             recpoint->Strip(), recpoint, recpoint->Particles(), 20);
591   return kTRUE;
592 }
593
594 //____________________________________________________________________
595 Bool_t 
596 AliFMDDisplay::ProcessESD(UShort_t det, Char_t rng, UShort_t sec, UShort_t str,
597                           Float_t, Float_t mult)
598 {
599   // Process data from ESD 
600   // Parameters 
601   //   det,rng,sec,str   Detector coordinates. 
602   //   mult              Multiplicity. 
603   Double_t cmult = mult;
604   if (fSpec) fSpec->Fill(cmult);
605   if (cmult < fMultCut || cmult == AliESDFMD::kInvalidMult) return kTRUE;
606   AddMarker(det,rng,sec,str, 0, cmult, 20);
607   if (fSpecCut) fSpecCut->Fill(cmult);
608   return kTRUE;
609 }
610
611 //____________________________________________________________________
612 //
613 // EOF
614 //