]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDisplay.cxx
Dummy default constructor
[u/mrichter/AliRoot.git] / RICH / AliRICHDisplay.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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
16
17 //////////////////////////////////////////////////////////////////////////
18 //                                                                      //
19 // AliDisplay                                                           //
20 //                                                                      //
21 // Utility class to display ALICE outline, tracks, hits,..              //
22 //                                                                      //
23 //////////////////////////////////////////////////////////////////////////
24
25 #include <TROOT.h>
26 #include <TTree.h>
27 #include <TButton.h>
28 #include <TColor.h>
29 #include <TCanvas.h>
30 #include <TView.h>
31 #include <TText.h>
32 #include <TPolyMarker3D.h>
33 #include <TAtt3D.h>
34 #include <TAttLine.h>
35 #include <TPolyMarker.h>
36 #include <TPaveLabel.h>
37 #include <TPaveText.h>
38 #include <TList.h>
39 #include <TDiamond.h>
40 #include <TNode.h>
41 #include <TArc.h>
42 #include <TTUBE.h>
43 #include <TSlider.h>
44 #include <TGeometry.h>
45 #include <TSliderBox.h>
46 #include <TGaxis.h>
47 #include <TVirtualX.h>
48 #include <TMath.h>
49 #include <TRandom.h>
50 #include <X3DBuffer.h>
51 #include <TParticle.h>
52
53 #include "AliRun.h"
54 #include "AliPDG.h"
55 #include "AliDetector.h"
56 #include "AliRICH.h"
57 #include "AliRICHConst.h"
58 #include "AliRICHDisplay.h"
59 #include "AliRICHPoints.h"
60 #include "AliHeader.h"
61
62 #include "AliRICHSDigit.h"
63 #include "AliMC.h"
64
65 ClassImp(AliRICHDisplay)
66     
67
68 //____________________________________________________________________________
69 AliRICHDisplay::AliRICHDisplay()
70
71
72 // default constructor
73     fColPad = 0;
74     fPoints = 0;
75     fPhits = 0;
76     fPCerenkovs = 0;
77     fCanvas = 0;
78     fRpoints = 0;
79     fRecpoints = 0; 
80 }
81
82 //_____________________________________________________________________________
83 AliRICHDisplay::AliRICHDisplay(Int_t size)
84 {
85 // Create an event display object.
86 // A canvas named "edisplay" is created with a vertical size in pixels
87 //
88 //    A QUICK Overview of the Event Display functions
89 //    ===============================================
90 //
91 //  The event display can ve invoked by executing the macro "display.C"
92 // A canvas like in the picture below will appear.
93 //
94 //  On the left side of the canvas, the following buttons appear:
95 //   *Next*       to move to the next event
96 //   *Previous*   to move to the previous event
97
98 //   *Pick*       Select this option to be able to point on a track with the
99 //                mouse. Once on the track, use the right button to select
100 //                an action. For example, select SetMarkerAttributes to
101 //                change the marker type/color/size for the track.
102 //   *Zoom*       Select this option (default) if you want to zoom.
103 //                To zoom, simply select the selected area with the left button.
104 //   *UnZoom*     To revert to the previous picture size.
105 //
106 //   slider R     On the left side, the vertical slider can be used to
107 //                set the default picture size.
108 //
109 //    When you are in Zoom mode, you can click on the black part of the canvas
110 //  to select special options with the right mouse button.
111
112 //
113 //  When you are in pick mode, you can "Inspect" the object pointed by the mouse.
114 //  When you are on a track, select the menu item "InspectParticle"
115 //  to display the current particle attributes.
116 //
117 //  You can activate the Root browser by selecting the Inspect menu
118 //  in the canvas tool bar menu. Then select "Start Browser"
119 //  This will open a new canvas with the browser. At this point, you may want
120 //  to display some histograms (from the Trees). Go to the "File" menu
121 //  of the browser and click on "New canvas".
122 //  In the browser, click on item "ROOT files" in the left pane.
123 //  Click on galice.root.
124 //  Click on TH
125 //  Click on TPC for example
126 //  Click on any variable (eg TPC.fX) to histogram the variable.
127 //
128 //   If you are lost, you can click on HELP in any Root canvas or browser.
129 //Begin_Html
130 /*
131   <img src="gif/AliRICHDisplay.gif">
132 */
133 //End_Html
134     
135     
136     fPad = 0;
137     
138     gAlice->SetDisplay(this);
139     
140     // Initialize display default parameters
141     SetRange();
142     
143     // Set front view by default
144     fTheta = 90;              //inclined HMPID
145     fPhi   = 30;              //inclined HMPID
146     //fTheta = 90;               //normal HMPID
147     //fPhi   = 90;                //normal HMPID
148     fPsi   = 0;
149     fChamber = 1;
150     fCathode = 1;
151     //   fRzone   = 1.e10;
152     fDrawClusters  = kTRUE;
153     fDrawCoG       = kTRUE;
154     fDrawRecHits   = kTRUE;
155     fZoomMode      = 1;
156     fZooms         = 0;
157     fClustersCuts  = 0;
158     fPoints        = 0;
159     fPCerenkovs    = 0;
160     fPhits         = 0;
161     fRpoints       = 0;   
162     fRecpoints       = 0;
163     // Create colors
164     CreateColors();
165     // Create display canvas
166     Int_t ysize = size;
167     if (ysize < 100) ysize = 750;
168     Int_t xsize = Int_t(size*830./ysize);
169     fCanvas = new TCanvas("Canvas", "RICH Clusters Display",14,47,xsize,ysize);
170     fCanvas->ToggleEventStatus();
171     
172     // Create main display pad
173     fPad = new TPad("viewpad", "RICH display",0.15,0,0.9,1);
174     fPad->Draw();
175     fPad->Modified();
176     fPad->SetFillColor(1);
177     fPad->SetBorderSize(2);
178     
179     fCanvas->cd();
180     
181     // Create colors pad
182     fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
183     fColPad->Draw();
184     fColPad->Modified();
185     fColPad->SetFillColor(19);
186     fColPad->SetBorderSize(2);
187     fColPad->cd();
188     DisplayColorScale();
189     
190     fCanvas->cd();
191     
192     // Create user interface control pad
193     DisplayButtons();
194     fCanvas->cd();
195     
196     // Create Range and mode pad
197     Float_t dxtr     = 0.15;
198     Float_t dytr     = 0.45;
199     fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
200     fTrigPad->Draw();
201     fTrigPad->cd();
202     fTrigPad->SetFillColor(22);
203     fTrigPad->SetBorderSize(2);
204     fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
205     fRangeSlider->SetObject(this);
206     char pickmode[] = "gAlice->Display()->SetPickMode()";
207     Float_t db = 0.09;
208     fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
209     fPickButton->SetFillColor(38);
210     fPickButton->Draw();
211     char zoommode[] = "gAlice->Display()->SetZoomMode()";
212     fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
213     fZoomButton->SetFillColor(38);
214     fZoomButton->Draw();
215     fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
216     fArcButton->SetFillColor(kGreen);
217     fArcButton->Draw();
218     char butUnzoom[] = "gAlice->Display()->UnZoom()";
219     TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
220     button->SetFillColor(38);
221     button->Draw();
222     AppendPad(); // append display object as last object to force selection
223     
224     fTrigPad->SetEditable(kFALSE);
225     fButtons->SetEditable(kFALSE);
226     fCanvas->cd();
227     fCanvas->Update();
228 }
229
230
231 //_____________________________________________________________________________
232 AliRICHDisplay::~AliRICHDisplay()
233 {
234     // Delete space point structure
235     if (fPoints) fPoints->Delete();
236     delete fPoints;
237     fPoints     = 0;
238     //
239     if (fPhits) fPhits->Delete();
240     delete fPhits;
241     fPhits     = 0;
242     //
243     if (fRpoints) fRpoints->Delete();
244     delete fRpoints;
245     fRpoints     = 0;
246     //
247     if (fRecpoints) fRecpoints->Delete();
248     delete fRecpoints;
249     fRecpoints     = 0;
250     //
251     if (fPCerenkovs) fPCerenkovs->Delete();
252     delete fPCerenkovs;
253     fPCerenkovs     = 0;
254 }
255
256 //_____________________________________________________________________________
257 void AliRICHDisplay::Clear(Option_t *)
258 {
259 //    Delete graphics temporary objects
260 }
261
262 //_____________________________________________________________________________
263 void AliRICHDisplay::DisplayButtons()
264 {
265 //    Create the user interface buttons
266     
267     
268     fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
269     fButtons->Draw();
270     fButtons->SetFillColor(38);
271     fButtons->SetBorderSize(2);
272     fButtons->cd();
273     
274     //   Int_t butcolor = 33;
275     Float_t dbutton = 0.08;
276     Float_t y  = 0.96;
277     Float_t dy = 0.014;
278     Float_t x0 = 0.05;
279     Float_t x1 = 0.95;
280     
281     TButton *button;
282     char but1[] = "gAlice->Display()->ShowNextEvent(1)";
283     button = new TButton("Next",but1,x0,y-dbutton,x1,y);
284     button->SetFillColor(38);
285     button->Draw();
286     
287     y -= dbutton +dy;
288     char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
289     button = new TButton("Previous",but2,x0,y-dbutton,x1,y);
290     button->SetFillColor(38);
291     button->Draw();
292
293     y -= dbutton +dy;
294     char but7[] = "gAlice->Display()->DrawViewGL()";
295     button = new TButton("OpenGL",but7,x0,y-dbutton,x1,y);
296     button->SetFillColor(38);
297     button->Draw();
298     
299     y -= dbutton +dy;
300     char but8[] = "gAlice->Display()->DrawViewX3D()";
301     button = new TButton("X3D",but8,x0,y-dbutton,x1,y);
302     button->SetFillColor(38);
303     button->Draw();
304     
305     // display logo
306     TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
307     diamond->SetFillColor(50);
308     diamond->SetTextAlign(22);
309     diamond->SetTextColor(5);
310     diamond->SetTextSize(0.11);
311     diamond->Draw();
312     diamond->AddText(".. ");
313     diamond->AddText("ROOT");
314     diamond->AddText("RICH");
315     diamond->AddText("... ");
316     diamond->AddText(" ");
317 }
318
319 //_____________________________________________________________________________
320 void AliRICHDisplay::CreateColors()
321 {
322 //    Create the colors palette used to display clusters
323     
324     Int_t k,i;
325     Int_t color;
326     Float_t r,g,b;
327     
328     for (k=1;k<=5;k++) {
329         switch(k) {
330         case 1:
331             for (i=1;i<=5;i++) {
332                 r=1.;
333                 g=i*0.2;  
334                 b=0.;
335                 color=i;
336                 color=700+23-color;
337                 new TColor(color,r,g,b);
338             } 
339             break;
340         case 2:
341             for (i=1;i<=4;i++) {
342                 r=1.1-i*0.2;
343                 g=1.;  
344                 b=0.;
345                 color=i+5;
346                 color=700+23-color;
347                 new TColor(color,r,g,b);
348             } 
349             break;
350     case 3:
351         for (i=1;i<=4;i++) {
352             r=0.;
353             g=1.;  
354             b=i*0.2+0.2;
355             color=i+9;
356             color=700+23-color;
357             new TColor(color,r,g,b);
358         } 
359         break;
360         case 4:
361             for (i=1;i<=4;i++) {
362                 r=0.;
363                 g=1.1-i*0.2;  
364                 b=1.;
365                 color=i+13;
366                 color=700+23-color;
367                 new TColor(color,r,g,b);
368             } 
369             break;
370         case 5:
371             for (i=1;i<=5;i++) {
372                 r=i*0.2;
373                 g=0.;  
374                 b=1.;
375                 color=i+17;
376                 color=700+23-color;
377                 new TColor(color,r,g,b);
378             } 
379             break;
380         }
381     }
382 }
383
384 //_____________________________________________________________________________
385 void AliRICHDisplay::DisplayColorScale()
386 {
387
388 // Draw the color scale in the RICH display canvas
389     
390     Int_t i;
391     Int_t color;
392     Float_t xlow, ylow, xup, yup, hs;
393     Float_t x1, y1, x2, y2;
394     x1 = y1 = 0;
395     x2 = y2 = 20;
396     
397     gPad->SetFillColor(0);
398     gPad->Clear();
399     gPad->Range(x1,y1,x2,y2);
400     TText *text = new TText(0,0,"");
401     text->SetTextFont(61);
402     text->SetTextSize(0.2);
403     text->SetTextAlign(22);
404     
405     TBox *box;
406     char label[8];
407 //*-* draw colortable boxes
408     hs = (y2-y1)/Float_t(22);
409     xlow=x1+1;
410     xup=x2-9;
411     for (i=0;i<22;i++) {
412         ylow = y1 + hs*(Float_t(i));
413         yup  = y1 + hs*(Float_t(i+1));
414         color = 701+i;
415         Double_t logscale=Double_t(i+1)*(TMath::Log(kadc_satm)/22);
416         Int_t scale=(Int_t)TMath::Exp(logscale);
417         sprintf(label,"%d",scale);
418         box = new TBox(xlow, ylow, xup, yup);
419         box->SetFillColor(color);
420         box->Draw();
421         text->DrawText(xup+4, 0.5*(ylow+yup),label);
422     }
423 }
424
425 //______________________________________________________________________________
426 Int_t AliRICHDisplay::DistancetoPrimitive(Int_t px, Int_t)
427 {
428 // Compute distance from point px,py to objects in event
429     
430     gPad->SetCursor(kCross);
431     
432     if (gPad == fTrigPad) return 9999;
433     
434     const Int_t kBig = 9999;
435     Int_t dist   = kBig;
436     Float_t xmin = gPad->GetX1();
437     Float_t xmax = gPad->GetX2();
438     Float_t dx   = 0.02*(xmax - xmin);
439     Float_t x    = gPad->AbsPixeltoX(px);
440     if (x < xmin+dx || x > xmax-dx) return dist;
441     if (fZoomMode) return 0;
442     else           return 7;
443 }
444
445 //_____________________________________________________________________________
446 void AliRICHDisplay::Draw(Option_t *)
447 {
448 //    Display current event
449
450     fPad->cd();
451     
452     DrawView(fTheta, fPhi, fPsi);      // see how to draw PGON+inner frames
453     
454     // Display the event number and title
455     fPad->cd();
456     DrawTitle();
457 }
458
459 //_____________________________________________________________________________
460 void AliRICHDisplay::DrawCoG()
461 {
462 //    Draw hits for RICH chambers
463
464     if (!fDrawCoG) return;
465     ResetRpoints();
466     for (Int_t chamber=0;chamber<kNCH;chamber++) {
467         LoadCoG(chamber,1);
468     }
469     
470     Int_t ncog, icog;
471     TObjArray *points;
472     AliRICHPoints *pm;
473     points = fRpoints;
474     if (!points) return;
475     ncog = points->GetEntriesFast();
476     for (icog=0; icog < ncog; icog++) {
477         pm = (AliRICHPoints*)points->UncheckedAt(icog);
478         if (!pm) continue;
479         pm->Draw();
480     }
481 }
482
483 void AliRICHDisplay::DrawRecHits()
484 {
485 //    Draw rec hits for RICH chambers
486
487     if (!fDrawRecHits) return;
488     //ResetRecpoints();
489     for (Int_t chamber=0;chamber<kNCH;chamber++) {
490         LoadRecHits(chamber,1);
491     }
492     
493     Int_t nrec, irec;
494     TObjArray *points;
495     AliRICHPoints *pm;
496     points = fRecpoints;
497     if (!points) return;
498     nrec = points->GetEntriesFast();
499     //printf("Nrec %d\n",nrec);
500     for (irec=0; irec < nrec; irec++) {
501         pm = (AliRICHPoints*)points->UncheckedAt(irec);
502         if (!pm) continue;
503         pm->Draw();
504     }
505 }
506
507 //_____________________________________________________________________________
508
509 void AliRICHDisplay::DrawCerenkovs()
510 {
511 //    Draw cerenkovs hits for RICH chambers
512     
513     LoadCerenkovs(fChamber);
514     //printf("\nDrawCerenkovs\n");
515     
516     Int_t ntracks, track;
517     TObjArray *cpoints;
518     AliRICHPoints *pm;
519     
520     fHitsCuts = 0;
521     cpoints = fPCerenkovs;
522     //printf ("Cpoints: %p",cpoints);
523     if (!cpoints) return;
524     ntracks = cpoints->GetEntriesFast();
525     //printf("DrawCerenkovs - ntracks %d \n",ntracks);
526     for (track=0;track<ntracks;track++) {
527         pm = (AliRICHPoints*)cpoints->UncheckedAt(track);
528         if (!pm) continue;
529         pm->Draw();
530         fHitsCuts += pm->GetN();
531     }
532 }
533
534 //_____________________________________________________________________________
535
536 void AliRICHDisplay::DrawClusters()
537 {
538 //    Draw clusterss for RICH chambers
539     
540     Int_t ndigits, digit;
541     TObjArray *points;
542     AliRICHPoints *pm;
543     
544     fClustersCuts = 0;
545     points = fPoints;
546     if (!points) return;
547     ndigits = points->GetEntriesFast();
548     for (digit=0;digit<ndigits;digit++){
549         pm = (AliRICHPoints*)points->UncheckedAt(digit);
550         if (!pm) continue;
551         pm->Draw();
552         Float_t *pxyz;
553         pxyz=pm->GetP();
554         for (Int_t im=0;im<3;im++) {
555             TMarker3DBox *marker=pm->GetMarker(im);
556             if (marker)
557                 marker->Draw();
558         }
559         fClustersCuts +=pm->GetN();
560     }
561 }
562
563 //_____________________________________________________________________________
564 void AliRICHDisplay::DrawHits()
565 {
566 //    Draw hits for RICH chambers
567     
568     LoadHits(fChamber);
569     
570     Int_t ntracks, track;
571     TObjArray *points;
572     AliRICHPoints *pm;
573     
574     fHitsCuts = 0;
575     points = Phits();
576     if (!fPhits) return;
577 //    ntracks = points->GetEntriesFast();
578     ntracks = fPhits->GetEntriesFast();
579     
580     //printf("DrawHits - ntracks %d \n",ntracks);
581     for (track=0;track<ntracks;track++) {
582         pm = (AliRICHPoints*)fPhits->UncheckedAt(track);
583         if (!pm) continue;
584         pm->Draw();
585         fHitsCuts += pm->GetN();
586     }
587 }
588
589
590 //_____________________________________________________________________________
591 void AliRICHDisplay::DrawTitle(Option_t *option)
592 {
593 //    Draw the event title
594     
595     Float_t xmin = gPad->GetX1();
596     Float_t xmax = gPad->GetX2();
597     Float_t ymin = gPad->GetY1();
598     Float_t ymax = gPad->GetY2();
599     Float_t dx   = xmax-xmin;
600     Float_t dy   = ymax-ymin;
601     
602     if (strlen(option) == 0) {
603         TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
604         title->SetBit(kCanDelete);
605         title->SetFillColor(42);
606         title->Draw();
607         char ptitle[100];
608         sprintf(ptitle,"Alice event: %d, Run:%d",
609                 gAlice->GetHeader()->GetEvent(), gAlice->GetHeader()->GetRun());
610         title->AddText(ptitle);
611         Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
612         sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
613                 nparticles, fHitsCuts,fClustersCuts);
614         title->AddText(ptitle);
615     } else {
616         TPaveLabel *label = 
617             new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
618         label->SetBit(kCanDelete);
619         label->SetFillColor(42);
620         label->Draw();
621     }
622 }
623
624 //_____________________________________________________________________________
625 void AliRICHDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
626 {
627 //    Draw a view of RICH clusters
628     
629    gPad->SetCursor(kWatch);
630    gPad->SetFillColor(1);
631    gPad->Clear();
632    
633    Int_t iret;
634    TView *view = new TView(1);
635    Float_t range = fRrange*fRangeSlider->GetMaximum();
636    view->SetRange(-range,-range,-range,range, range, range);
637    fZoomX0[0] = -1;
638    fZoomY0[0] = -1;
639    fZoomX1[0] =  1;
640    fZoomY1[0] =  1;
641    fZooms = 0;
642    
643    //Display RICH Chamber Geometry
644    gAlice->GetGeometry()->Draw("same");
645    
646    //add clusters to the pad
647    DrawClusters();
648    DrawHits();
649    //DrawCerenkovs();
650    if (gAlice->TreeR())
651      {
652        //printf("Calling DrawCoG\n");
653          DrawCoG();
654        //printf("Calling DrawRecHits\n");
655          DrawRecHits();
656      }
657    /*for (Int_t i=0;i<7;i++)
658      LoadRecHits(i,1);*/
659    
660    // add itself to the list (must be last)
661    AppendPad();
662    
663    view->SetView(phi, theta, psi, iret);
664 }
665
666 //______________________________________________________________________________
667 void AliRICHDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
668 {
669 //  Execute action corresponding to the mouse event
670     
671     static Float_t x0, y0, x1, y1;
672     
673    static Int_t pxold, pyold;
674    static Int_t px0, py0;
675    static Int_t linedrawn;
676    Float_t temp;
677    
678    if (px == 0 && py == 0) { //when called by sliders
679        if (event == kButton1Up) {
680          printf("Drawing event %d\n",event);
681          Draw();
682        }
683        return;
684    }
685    if (!fZoomMode && gPad->GetView()) {
686        gPad->GetView()->ExecuteRotateView(event, px, py);
687        return;
688    }
689
690    // something to zoom ?
691    gPad->SetCursor(kCross);
692    
693    switch (event) {
694        
695    case kButton1Down:
696        gVirtualX->SetLineColor(-1);
697        gPad->TAttLine::Modify();  //Change line attributes only if necessary
698        x0 = gPad->AbsPixeltoX(px);
699        y0 = gPad->AbsPixeltoY(py);
700        px0   = px; py0   = py;
701        pxold = px; pyold = py;
702        linedrawn = 0;
703        return;
704        
705    case kButton1Motion:
706        if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
707        pxold = px;
708        pyold = py;
709        linedrawn = 1;
710        gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
711        return;
712        
713    case kButton1Up:
714        gPad->GetCanvas()->FeedbackMode(kFALSE);
715        if (px == px0) return;
716        if (py == py0) return;
717        x1 = gPad->AbsPixeltoX(px);
718        y1 = gPad->AbsPixeltoY(py);
719        
720        if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
721        if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
722        gPad->Range(x0,y0,x1,y1);
723        if (fZooms < kMAXZOOM-1) {
724            fZooms++;
725            fZoomX0[fZooms] = x0;
726            fZoomY0[fZooms] = y0;
727            fZoomX1[fZooms] = x1;
728            fZoomY1[fZooms] = y1;
729        }
730        gPad->Modified(kTRUE);
731        return;
732    }
733    
734 }
735 //___________________________________________
736 void AliRICHDisplay::LoadCoG(Int_t chamber, Int_t cathode)
737 {
738 // Read raw clusters info and store x,y,z info in arrays fRpoints
739 // Loop on all detectors
740
741    if (chamber > 7) return;
742
743    AliRICH *pRICH  = (AliRICH*)gAlice->GetModule("RICH");
744    AliRICHChamber*  iChamber;
745
746    TClonesArray *pRICHrawclust  = pRICH->Clusters(chamber);
747    if (pRICHrawclust == 0) return;
748
749    pRICH->ResetClusters();
750
751
752    Int_t nent=(Int_t)gAlice->TreeR()->GetEntries();
753    gAlice->TreeR()->GetEvent(nent-1+cathode-1);
754    Int_t nrawcl = pRICHrawclust->GetEntriesFast();
755    //printf ("nrawcl:%d\n",nrawcl);
756    if (nrawcl == 0) return;
757    if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
758    
759    iChamber = &(pRICH->Chamber(chamber));
760    AliRICHcluster  *mRaw;
761    AliRICHPoints *points = 0;
762    //
763    //loop over all raw clusters and store their position
764    points = new AliRICHPoints(nrawcl);
765    for (Int_t iraw=0;iraw<nrawcl;iraw++) {
766        mRaw   = (AliRICHcluster*)pRICHrawclust->UncheckedAt(iraw);
767        fRpoints->AddAt(points,iraw);
768        points->SetMarkerColor(3);
769        points->SetMarkerStyle(3);
770        points->SetMarkerSize(1.);
771        points->SetParticle(-1);
772        points->SetHitIndex(-1);
773        points->SetTrackIndex(-1);
774        points->SetDigitIndex(-1);
775        Float_t  vectorLoc[3]={mRaw->X(),5,mRaw->Y()};
776        Float_t  vectorGlob[3];
777        iChamber->LocaltoGlobal(vectorLoc,vectorGlob);
778        points->SetPoint(iraw,vectorGlob[0],vectorGlob[1],vectorGlob[2]);
779    }
780 }//LoadCoG()
781 //___________________________________________
782 void AliRICHDisplay::LoadRecHits(Int_t chamber, Int_t cathode)
783 {
784   chamber++;cathode++;
785 }
786 //___________________________________________
787 void AliRICHDisplay::LoadDigits()
788 {
789 // Read digits info and store x,y,z info in arrays fPoints
790 // Loop on all detectors
791     
792    ResetPoints();
793    AliRICH *pRICH  = (AliRICH*)gAlice->GetDetector("RICH");
794    AliRICHChamber*       iChamber;
795    AliSegmentation*      segmentation;
796    Int_t nAllDigits=0;
797    Int_t ich;
798
799    //printf("Entering LoadDigits\n");
800
801    if (gAlice->TreeD())
802      {
803    
804        for (ich=1; ich<=kNCH; ich++) {
805          TClonesArray *pRICHdigits  = pRICH->Digits(ich);
806          if (pRICHdigits == 0) continue;
807          gAlice->ResetDigits();
808          gAlice->TreeD()->GetEvent(0);
809          Int_t ndigits = pRICHdigits->GetEntriesFast();
810          nAllDigits+=ndigits;
811        }
812        
813        if (fPoints == 0) fPoints = new TObjArray(nAllDigits);   
814        Int_t counter=0;
815        for (ich=1; ich<=kNCH; ich++) {
816          TClonesArray *pRICHdigits  = pRICH->Digits(ich);
817          if (pRICHdigits == 0) continue;
818          gAlice->ResetDigits();
819          gAlice->TreeD()->GetEvent(0);
820          Int_t ndigits = pRICHdigits->GetEntriesFast();
821          if (ndigits == 0) continue;
822          iChamber = &(pRICH->Chamber(ich));
823          segmentation=iChamber->GetSegmentationModel();
824          Float_t dpx  = segmentation->Dpx();
825          Float_t dpy  = segmentation->Dpy();
826          
827          //printf("Dpx:%d, Dpy:%d\n",dpx,dpy);
828          
829          AliRICHdigit  *mdig;
830          AliRICHPoints *points = 0;
831          TMarker3DBox  *marker = 0;
832          //
833          //loop over all digits and store their position
834          Int_t npoints=1;
835          
836          for (Int_t digit=0;digit<ndigits;digit++) {
837            mdig    = (AliRICHdigit*)pRICHdigits->UncheckedAt(digit);
838            points = new AliRICHPoints(npoints);
839            fPoints->AddAt(points,counter);
840            counter++;
841            Int_t charge=(Int_t)mdig->Q();
842            Int_t index=Int_t(TMath::Log(charge)/(TMath::Log(kadc_satm)/22));
843            Int_t color=701+index;
844            if (color>722) color=722;
845            points->SetMarkerColor(color);
846            points->SetMarkerStyle(21);
847            points->SetMarkerSize(0.5);
848            Float_t xpad, ypad, zpad;
849            segmentation->GetPadC(mdig->X(), mdig->Y(),xpad, ypad, zpad);
850            Float_t vectorLoc[3]={xpad,5,ypad};
851            Float_t  vectorGlob[3];
852            iChamber->LocaltoGlobal(vectorLoc,vectorGlob);
853            points->SetParticle(-1);
854            points->SetHitIndex(-1);
855            points->SetTrackIndex(-1);
856            points->SetDigitIndex(digit);
857            points->SetPoint(0,vectorGlob[0],vectorGlob[1],vectorGlob[2]);
858            //printf("Y position (digit): %f\n", vectorGlob[1]);
859            
860            segmentation->GetPadC(mdig->X(), mdig->Y(), xpad, ypad, zpad);
861            Float_t theta = iChamber->GetRotMatrix()->GetTheta();
862            Float_t phi   = iChamber->GetRotMatrix()->GetPhi();     
863            marker=new TMarker3DBox(vectorGlob[0],vectorGlob[1],vectorGlob[2],
864                                    dpy/2,0,dpx/2,theta,phi);
865            marker->SetLineColor(2);
866            marker->SetFillStyle(1001);
867            marker->SetFillColor(color);
868            marker->SetRefObject((TObject*)points);
869            points->Set3DMarker(0, marker);
870          } // loop over digits
871        } // loop over chambers 
872      } //if TreeD
873 }//LoadDigits();
874 //__________________________________________________________________________________________________
875 void AliRICHDisplay::LoadHits(Int_t chamber)
876 {
877 // Read hits info and store x,y,z info in arrays fPhits
878 // Loop on all detectors
879     
880     
881     fChamber=chamber; 
882     ResetPhits();
883     
884     AliRICH *pRICH  = (AliRICH*)gAlice->GetDetector("RICH");
885     AliRICHChamber*  iChamber;
886     
887     iChamber = &(pRICH->Chamber(chamber-1));
888     Int_t ntracks = (Int_t)pRICH->TreeH()->GetEntries();
889     Int_t track;
890     
891     if (fPhits == 0) fPhits = new TObjArray(ntracks);
892     //TVector *xp = new TVector(1000);
893     //TVector *yp = new TVector(1000);
894     //TVector *zp = new TVector(1000);
895     //TVector *ptrk = new TVector(1000);
896     //TVector *phit = new TVector(1000);
897     Int_t nAllHits=0;
898     for (track=0; track<ntracks;track++) {
899         gAlice->ResetHits();
900         pRICH->TreeH()->GetEvent(track);
901         TClonesArray *pRICHhits  = pRICH->Hits();
902         if (pRICHhits == 0) return;
903         Int_t nhits = pRICHhits->GetEntriesFast();
904         nAllHits+=nhits;
905     }
906
907     fPhits = new TObjArray(nAllHits);
908
909     Int_t npoints=0;
910     for (track=0; track<ntracks;track++) {
911         gAlice->ResetHits();
912         pRICH->TreeH()->GetEvent(track);
913         TClonesArray *pRICHhits  = pRICH->Hits();
914         if (pRICHhits == 0) return;
915         Int_t nhits = pRICHhits->GetEntriesFast();
916         if (nhits == 0) continue;
917         AliRICHhit *mHit;
918         AliRICHPoints *points = 0;
919         for (Int_t hit=0;hit<nhits;hit++) {
920             points = new AliRICHPoints(1);
921             fPhits->AddAt(points,npoints);
922             mHit = (AliRICHhit*)pRICHhits->UncheckedAt(hit);
923             TParticle *current = (TParticle*)gAlice->GetMCApp()->Particle(mHit->Track());
924             if (current->GetPdgCode() == 50000050) {
925                 points->SetMarkerColor(kBlue);
926             } else if (current->GetPdgCode() == 50000051) {
927                 points->SetMarkerColor(kYellow);
928             } else {
929                 points->SetMarkerColor(kRed);
930             }
931             points->SetMarkerStyle(5);
932             points->SetMarkerSize(1.);
933             points->SetParticle(mHit->Track());
934             points->SetHitIndex(hit);
935             points->SetTrackIndex(track);
936             points->SetDigitIndex(-1);
937             points->SetPoint(0,mHit->X(), mHit->Y(), mHit->Z());
938             npoints++;
939         }
940     }
941 }//LoadHits()
942 //__________________________________________________________________________________________________
943 void AliRICHDisplay::LoadCerenkovs(Int_t chamber)
944 {
945 // Read cerenkov hits info and store x,y,z info in array fPCerenkovs
946 // Loop on all detectors
947     
948     fChamber=chamber; 
949     ResetPCerenkovs();
950     
951     AliRICH *pRICH  = (AliRICH*)gAlice->GetDetector("RICH");
952     AliRICHChamber*  iChamber;
953     
954     iChamber = &(pRICH->Chamber(chamber-1));
955     
956     pRICH->SetTreeAddress();
957     Int_t ntracks = (Int_t)pRICH->TreeH()->GetEntries();
958     
959     if (fPCerenkovs == 0) fPCerenkovs = new TObjArray(ntracks);
960     TVector *xp = new TVector(1000);
961     TVector *yp = new TVector(1000);
962     TVector *zp = new TVector(1000);
963     TVector *ptrk = new TVector(1000);
964     TVector *phit = new TVector(1000);
965     for (Int_t track=0; track<ntracks;track++) {
966         gAlice->ResetHits();
967         pRICH->TreeH()->GetEvent(track);
968         TClonesArray *pRICHCerenkovs  = pRICH->Cerenkovs();
969         if (pRICHCerenkovs == 0) return;
970         Int_t nhits = pRICHCerenkovs->GetEntriesFast();
971         if (nhits == 0) continue;
972         AliRICHCerenkov *mCerenkov;
973         AliRICHPoints *cpoints = 0;
974         Int_t npoints=0;
975         
976         
977 //Display Cerenkov hits in blue
978         
979         for (Int_t hit=0;hit<nhits;hit++) {
980             mCerenkov = (AliRICHCerenkov*)pRICHCerenkovs->UncheckedAt(hit);
981             (*xp)(npoints)=mCerenkov->X();
982             (*yp)(npoints)=mCerenkov->Y();
983             (*zp)(npoints)=mCerenkov->Z();
984             (*ptrk)(npoints)=Float_t(mCerenkov->GetTrack());
985             (*phit)(npoints)=Float_t(hit);
986             npoints++;
987         }
988         if (npoints == 0) continue;
989         cpoints = new AliRICHPoints(npoints);
990         for (Int_t p=0;p<npoints;p++) {
991             cpoints->SetMarkerColor(kBlue);
992             cpoints->SetMarkerStyle(5);
993             cpoints->SetMarkerSize(1.);
994             cpoints->SetParticle(Int_t((*ptrk)(p)));
995             cpoints->SetHitIndex(Int_t((*phit)(p)));
996             cpoints->SetTrackIndex(track);
997             cpoints->SetDigitIndex(-1);
998             cpoints->SetPoint(p,(*xp)(p),(*yp)(p),(*zp)(p));
999         }
1000         xp->Zero();
1001         yp->Zero();
1002         ptrk->Zero();
1003         phit->Zero();
1004         fPCerenkovs->AddAt(cpoints,track);
1005     }
1006 }
1007
1008 //_____________________________________________________________________________
1009 void AliRICHDisplay::Paint(Option_t *)
1010 {
1011 //    Paint miscellaneous items
1012
1013 }
1014
1015 //_____________________________________________________________________________
1016 void AliRICHDisplay::SetPickMode()
1017 {
1018
1019 // Toggle pick mode
1020
1021     fZoomMode = 0;
1022     
1023     fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1024     fTrigPad->Modified();
1025 }
1026
1027 //_____________________________________________________________________________
1028 void AliRICHDisplay::SetZoomMode()
1029 {
1030
1031 // Toggle Zoom mode
1032
1033     fZoomMode = 1;
1034     
1035     fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1036     fTrigPad->Modified();
1037 }
1038
1039 //_____________________________________________________________________________
1040 void AliRICHDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1041 {
1042 // Set chamber and cathode number
1043     fChamber = chamber;
1044     fCathode = cathode;
1045     
1046     //printf("SetChamberAndCathode - fChamber fCathode %d %d\n",fChamber,fCathode);
1047     if (!fPad) return;
1048     fPad->Clear();
1049     LoadDigits();
1050     Draw();
1051 }
1052
1053 //_____________________________________________________________________________
1054 void AliRICHDisplay::SetRange(Float_t rrange, Float_t zrange)
1055 {
1056 // Set view range along R and Z
1057     fRrange = rrange;
1058     fZrange = zrange;
1059     
1060    if (!fPad) return;
1061    fPad->Clear();
1062    Draw();
1063 }
1064
1065 //_____________________________________________________________________________
1066 void AliRICHDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1067 {
1068 //  change viewing angles for current event
1069     
1070     fPad->cd();
1071     fPhi   = phi;
1072     fTheta = theta;
1073     fPsi   = psi;
1074     Int_t iret = 0;
1075     
1076     TView *view = gPad->GetView();
1077     if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1078     else      Draw();
1079     
1080     gPad->Modified();
1081 }
1082
1083 //_____________________________________________________________________________
1084 void AliRICHDisplay::ShowNextEvent(Int_t delta)
1085 {
1086 //  Display (current event_number+delta)
1087 //    delta =  1  shown next event
1088 //    delta = -1 show previous event
1089     
1090     if (delta) {
1091         gAlice->Clear();
1092         Int_t currentEvent = gAlice->GetHeader()->GetEvent();
1093         Int_t newEvent     = currentEvent + delta;
1094         gAlice->GetEvent(newEvent);
1095         if (!gAlice->TreeD()) return; 
1096     }
1097     LoadDigits();
1098     DrawClusters();
1099     fPad->cd(); 
1100     Draw();
1101
1102   
1103 }
1104
1105 //______________________________________________________________________________
1106 void AliRICHDisplay::UnZoom()
1107 {
1108
1109 // Return to previous zoom factor
1110
1111     if (fZooms <= 0) return;
1112     fZooms--;
1113     TPad *pad = (TPad*)gPad->GetPadSave();
1114     pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1115     pad->Modified();
1116 }
1117
1118 //_____________________________________________________________________________
1119 void AliRICHDisplay::ResetPoints()
1120 {
1121     //
1122     // Reset array of points
1123     //
1124     if (fPoints) {
1125         fPoints->Delete();
1126         delete fPoints;
1127         fPoints = 0;
1128   }
1129 }
1130 //_____________________________________________________________________________
1131 void AliRICHDisplay::ResetRpoints()
1132 {
1133   //
1134   // Reset array of points
1135   //
1136   if (fRpoints) {
1137     fRpoints->Delete();
1138     delete fRpoints;
1139     fRpoints = 0;
1140   }
1141 }
1142 //_____________________________________________________________________________
1143 void AliRICHDisplay::ResetRecpoints()
1144 {
1145   //
1146   // Reset array of points
1147   //
1148   if (fRecpoints) {
1149     fRecpoints->Delete();
1150     delete fRecpoints;
1151     fRecpoints = 0;
1152   }
1153 }
1154 //_____________________________________________________________________________
1155 void AliRICHDisplay::ResetPhits()
1156 {
1157     //
1158     // Reset array of points
1159     //
1160     if (fPhits) {
1161         fPhits->Delete();
1162         delete fPhits;
1163         fPhits = 0;
1164     }
1165 }
1166 //_____________________________________________________________________________
1167 void AliRICHDisplay::ResetPCerenkovs()
1168 {
1169     //
1170     // Reset array of points
1171     //
1172     if (fPCerenkovs) {
1173         fPCerenkovs->Delete();
1174         delete fPCerenkovs;
1175         fPCerenkovs = 0;
1176     }
1177 }
1178
1179 //_____________________________________________________________________________
1180
1181 void AliRICHDisplay::DrawViewGL()
1182 {
1183 //    Draw current view using OPENGL
1184
1185    TPad *pad = (TPad*)gPad->GetPadSave();
1186    pad->cd();
1187    TView *view = pad->GetView();
1188    if (!view) return;
1189    pad->x3d("OPENGL");
1190 }
1191
1192 //_____________________________________________________________________________
1193 void AliRICHDisplay::DrawViewX3D()
1194 {
1195 //    Draw current view using X3D
1196
1197    TPad *pad = (TPad*)gPad->GetPadSave();
1198    pad->cd();
1199    TView *view = pad->GetView();
1200    if (!view) return;
1201    pad->x3d();
1202 }