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