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