1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 /// \class AliMUONDisplay
19 /// Create an event display object.
20 /// A canvas named "edisplay" is created with a vertical size in pixels \n
22 /// A QUICK Overview of the Event Display functions \n
23 /// =============================================== \n
25 /// The event display can ve invoked by executing the macro "display.C" \n
26 /// A canvas like in the picture below will appear.
28 /// On the left side of the canvas, the following buttons appear:
29 /// - *Next* to move to the next event
30 /// - *Previous* to move to the previous event
31 /// - *Pick* Select this option to be able to point on a track with the
32 /// mouse. Once on the track, use the right button to select
33 /// an action. For example, select SetMarkerAttributes to
34 /// change the marker type/color/size for the track.
35 /// - *Zoom* Select this option (default) if you want to zoom.
36 /// To zoom, simply select the selected area with the left button.
37 /// - *UnZoom* To revert to the previous picture size.
39 /// - slider R On the left side, the vertical slider can be used to
40 /// set the default picture size.
42 /// When you are in Zoom mode, you can click on the black part of the canvas
43 /// to select special options with the right mouse button.
45 /// When you are in pick mode, you can "Inspect" the object pointed by the mouse.
46 /// When you are on a track, select the menu item "InspectParticle"
47 /// to display the current particle attributes.
49 /// You can activate the Root browser by selecting the Inspect menu
50 /// in the canvas tool bar menu. Then select "Start Browser"
51 /// This will open a new canvas with the browser. At this point, you may want
52 /// to display some histograms (from the Trees). Go to the "File" menu
53 /// of the browser and click on "New canvas".
54 /// In the browser, click on item "ROOT files" in the left pane.
55 /// Click on galice.root.
57 /// Click on TPC for example
58 /// Click on any variable (eg TPC.fX) to histogram the variable.
60 /// If you are lost, you can click on HELP in any Root canvas or browser.
64 <img src="gif/AliMUONDisplay.gif">
68 #include "AliMUONDisplay.h"
70 #include "AliMUONPoints.h"
71 #include "AliMUONGlobalTrigger.h"
72 #include "AliMUONHit.h"
73 #include "AliMUONDigit.h"
74 #include "AliMUONRawCluster.h"
75 #include "AliMUONTrack.h"
76 #include "AliMUONTrackParam.h"
77 #include "AliMUONGeometryTransformer.h"
78 #include "AliMUONSegmentation.h"
79 #include "AliMUONGeometrySegmentation.h"
80 #include "AliMUONConstants.h"
81 #include "AliMUONTriggerSegmentation.h"
83 #include "AliMpDEIterator.h"
84 #include "AliMpSegmentation.h"
85 #include "AliMpSlatSegmentation.h"
86 #include "AliMpSlat.h"
87 #include "AliMpSectorSegmentation.h"
88 #include "AliMpSector.h"
89 #include "AliMpTriggerSegmentation.h"
90 #include "AliMpTrigger.h"
91 #include "AliMpStationType.h"
92 #include "AliMpCathodType.h"
93 #include "AliMpDEManager.h"
98 #include "AliHeader.h"
103 #if ROOT_VERSION_CODE>= 331523
109 #include <TPaveLabel.h>
110 #include <TPaveText.h>
111 #include <TDiamond.h>
115 #include <TVirtualX.h>
117 #include <TGeometry.h>
118 #include <TMarker3DBox.h>
119 #include <TParticle.h>
120 #include <TPolyLine3D.h>
124 ClassImp(AliMUONDisplay)
127 //_____________________________________________________________________________
128 AliMUONDisplay::AliMUONDisplay()
133 fDrawClusters(kTRUE),
145 /// Default constructor
148 //_____________________________________________________________________________
149 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
154 fDrawClusters(kTRUE),
162 fNextCathode(kFALSE),
167 /// Standard constructor to create an event display object.
171 gAlice->SetDisplay(this);
173 // Initialize display default parameters
176 // Set front view by default
185 // Create display canvas
187 if (ysize < 100) ysize = 750;
188 Int_t xsize = Int_t(size*830./ysize);
189 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
190 fCanvas->ToggleEventStatus();
192 // Create main display pad
193 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
196 fPad->SetFillColor(30);
197 fPad->SetBorderSize(2);
202 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
204 fColPad->SetFillColor(17);
205 fColPad->SetBorderSize(2);
210 // Create user interface control pad
214 // Create Range and mode pad
217 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
220 fTrigPad->SetFillColor(22);
221 fTrigPad->SetBorderSize(2);
222 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
223 fRangeSlider->SetObject(this);
224 char pickmode[] = "gAlice->Display()->SetPickMode()";
226 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
227 fPickButton->SetFillColor(38);
229 char zoommode[] = "gAlice->Display()->SetZoomMode()";
230 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
231 fZoomButton->SetFillColor(38);
233 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
234 fArcButton->SetFillColor(kGreen);
236 char butUnzoom[] = "gAlice->Display()->UnZoom()";
237 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
238 button->SetFillColor(38);
240 AppendPad(); // append display object as last object to force selection
243 fTrigPad->SetEditable(kFALSE);
244 fButtons->SetEditable(kFALSE);
247 // initialize container
249 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
254 //_____________________________________________________________________________
255 AliMUONDisplay::~AliMUONDisplay()
259 // Delete space point structure
260 if (fPoints) fPoints->Delete();
264 if (fPhits) fPhits->Delete();
268 if (fRpoints) fRpoints->Delete();
273 //_____________________________________________________________________________
274 void AliMUONDisplay::Clear(Option_t *)
276 /// Delete graphics temporary objects
279 //_____________________________________________________________________________
280 void AliMUONDisplay::DisplayButtons()
282 /// Create the user interface buttons
285 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
287 fButtons->SetFillColor(38);
288 fButtons->SetBorderSize(2);
291 // Int_t butcolor = 33;
292 Float_t dbutton = 0.08;
299 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
300 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
301 button->SetFillColor(38);
305 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
306 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
307 button->SetFillColor(38);
311 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
312 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
313 button->SetFillColor(38);
317 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
318 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
319 button->SetFillColor(38);
323 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
324 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
325 button->SetFillColor(38);
329 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
330 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
331 button->SetFillColor(38);
335 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
336 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
337 button->SetFillColor(38);
341 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
342 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
343 button->SetFillColor(38);
347 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
348 diamond->SetFillColor(50);
349 diamond->SetTextAlign(22);
350 diamond->SetTextColor(5);
351 diamond->SetTextSize(0.11);
353 diamond->AddText(".. ");
354 diamond->AddText("ROOT");
355 diamond->AddText("MUON");
356 diamond->AddText("... ");
357 diamond->AddText(" ");
360 //_____________________________________________________________________________
361 void AliMUONDisplay::CreateColors() const
363 /// Create the colors palette used to display clusters
378 new TColor(color,r,g,b);
388 new TColor(color,r,g,b);
398 new TColor(color,r,g,b);
408 new TColor(color,r,g,b);
418 new TColor(color,r,g,b);
425 //_____________________________________________________________________________
426 void AliMUONDisplay::DisplayColorScale()
428 /// Display pulse height color scale
432 Float_t xlow, ylow, xup, yup, hs;
433 Float_t x1, y1, x2, y2;
437 TText *text = new TText(0,0,"");
438 text->SetTextFont(61);
439 text->SetTextSize(0.2);
440 text->SetTextAlign(22);
443 Int_t adcmax=4096; // default 12 bits ADC
449 //*-* draw colortable boxes
450 hs = (y2-y1)/Float_t(22);
454 ylow = y1 + hs*(Float_t(i));
455 yup = y1 + hs*(Float_t(i+1));
457 Double_t logscale=Double_t(i+1)*(TMath::Log(adcmax)/22);
458 Int_t scale=(Int_t)TMath::Exp(logscale);
459 sprintf(label,"%d",scale);
460 box = new TBox(xlow, ylow, xup, yup);
462 box->SetFillColor(color);
463 text->DrawText(xlow+0.7, 0.5*(ylow+yup),label);
467 //______________________________________________________________________________
468 Int_t AliMUONDisplay::DistancetoPrimitive(Int_t px, Int_t)
470 /// Compute distance from point px,py to objects in event
472 gPad->SetCursor(kCross);
474 if (gPad == fTrigPad) return 9999;
476 const Int_t kBig = 9999;
478 Float_t xmin = gPad->GetX1();
479 Float_t xmax = gPad->GetX2();
480 Float_t dx = 0.02*(xmax - xmin);
481 Float_t x = gPad->AbsPixeltoX(px);
482 if (x < xmin+dx || x > xmax-dx) return dist;
484 if (fZoomMode) return 0;
488 //_____________________________________________________________________________
489 void AliMUONDisplay::Draw(Option_t *)
491 /// Display current event
499 //_____________________________________________________________________________
500 void AliMUONDisplay::DrawChamber()
502 /// Display current event
504 fDrawTracks = kFALSE;
506 DrawView(fTheta, fPhi, fPsi);
507 // Display the event number and title
512 //_____________________________________________________________________________
513 void AliMUONDisplay::DrawReco(Option_t *)
515 /// Display current event
518 // print kinematics of generated particles
520 // Draw global view of muon system
522 DrawGlobalView(135, -50, -140);
524 // Display the event number and title
529 //_____________________________________________________________________________
530 void AliMUONDisplay::PrintKinematics()
532 /// Print kinematic tree
534 AliRunLoader * runLoader;
535 TParticle *particle = new TParticle();
537 Float_t vertex[3], momentum[3];
540 runLoader = fLoader->GetRunLoader();
544 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
545 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
546 nPart = (Int_t)runLoader->TreeK()->GetEntries();
547 for(Int_t iPart = 0; iPart < nPart; iPart++) {
548 runLoader->TreeK()->GetEvent(iPart);
549 vertex[0] = particle->Vx();
550 vertex[1] = particle->Vy();
551 vertex[2] = particle->Vz();
552 momentum[0] = particle->Px();
553 momentum[1] = particle->Py();
554 momentum[2] = particle->Pz();
556 printf("===================================================\n");
557 printf(" Generated particle # %d \n",iPart);
558 printf(" name: %s \n",particle->GetName());
559 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
560 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
565 //_____________________________________________________________________________
566 void AliMUONDisplay::DrawSegmentation()
568 /// \todo to be re-written for new seg
569 /// Draw graphical representation of segmenatation
570 /// Attention: still experimental code
573 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
574 // AliMUONChamber* iChamber;
576 // AliSegmentation* seg;
577 // iChamber = &(pMUON->Chamber(fChamber));
578 // seg=iChamber->SegmentationModel(icat);
580 // Float_t zpos=iChamber->Z();
581 // Float_t r=iChamber->ROuter();
583 // TMarker3DBox *marker;
585 // for (Int_t j=0; j<seg->Npy(); j++) {
587 // y0=j*seg->Dpy()-seg->Dpy()/2.;
588 // for (seg->FirstPad(0.,y0,0,300,0.);
592 // if (seg->ISector()==0) continue;
594 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
595 // Float_t dpx=seg->Dpx(seg->ISector())/2;
596 // Float_t dpy=seg->Dpy(seg->ISector())/2;
597 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
598 // marker->SetLineColor(seg->ISector()+1);
599 // marker->SetFillStyle(1001);
600 // marker->SetFillColor(0);
605 // for (Int_t j=0; j<250; j++) {
606 // Float_t x0=j*seg->Dpx();
607 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
609 // for (seg->FirstPad(x0,0,0,0,y0);
613 // if (seg->ISector()==0) continue;
616 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
617 // Float_t dpx=seg->Dpx(seg->ISector())/2;
618 // Float_t dpy=seg->Dpy(seg->ISector())/2;
619 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
620 // marker->SetLineColor(seg->ISector()+1);
621 // marker->SetFillStyle(1001);
622 // marker->SetFillColor(0);
629 //_____________________________________________________________________________
630 void AliMUONDisplay::DrawClusters()
632 /// Draw clusters for MUON chambers
634 Int_t ndigits, digit;
642 ndigits = points->GetEntriesFast();
643 for (digit=0;digit<ndigits;digit++){
644 pm = (AliMUONPoints*)points->UncheckedAt(digit);
648 for (Int_t im=0;im<3;im++) {
649 TMarker3DBox *marker=pm->GetMarker(im);
654 fClustersCuts += pm->GetN();
658 //_____________________________________________________________________________
659 void AliMUONDisplay::DrawHits()
661 /// Draw hits for MUON chambers
665 Int_t ntracks, track;
672 ntracks = points->GetEntriesFast();
673 for (track=0;track<ntracks;track++) {
674 pm = (AliMUONPoints*)points->UncheckedAt(track);
677 fHitsCuts += pm->GetN();
682 //_____________________________________________________________________________
683 void AliMUONDisplay::DrawCoG()
685 /// Draw hits for MUON chambers
687 if (!fDrawCoG) return;
688 if (fChamber > 10) return;
689 LoadCoG(fChamber,fCathode);
697 ncog = points->GetEntriesFast();
698 for (icog=0;icog<ncog;icog++) {
699 pm = (AliMUONPoints*)points->UncheckedAt(icog);
704 //_____________________________________________________________________________
705 void AliMUONDisplay::DrawTracks()
709 if (!fDrawTracks) return;
712 Int_t nTrack, iTrack;
718 nTrack = points->GetEntriesFast();
719 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
720 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
725 //_____________________________________________________________________________
727 void AliMUONDisplay::DrawTitle(Option_t *option)
729 /// Draw the event title
731 Float_t xmin = gPad->GetX1();
732 Float_t xmax = gPad->GetX2();
733 Float_t ymin = gPad->GetY1();
734 Float_t ymax = gPad->GetY2();
735 Float_t dx = xmax-xmin;
736 Float_t dy = ymax-ymin;
738 AliRunLoader * runLoader;
740 runLoader = fLoader->GetRunLoader();
745 if (strlen(option) == 0) {
746 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
747 // title->SetTextSize(0.023932);
748 title->SetTextSize(0.02);
749 title->SetBit(kCanDelete);
750 title->SetFillColor(42);
753 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
754 runLoader->GetEventNumber(),
755 gAlice->GetHeader()->GetRun(),
758 title->AddText(ptitle);
759 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
760 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
761 nparticles, fHitsCuts,fClustersCuts);
762 title->AddText(ptitle);
764 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
765 label->SetBit(kCanDelete);
766 label->SetFillColor(42);
771 //_____________________________________________________________________________
772 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
774 /// Draw a view of MUON clusters
776 AliInfo(" Draw View");
778 gPad->SetCursor(kWatch);
779 // gPad->SetFillColor(39);
780 gPad->SetFillColor(1);
782 // gPad->SetFillColor(39);
783 gPad->SetFillColor(1);
786 #if ROOT_VERSION_CODE>= 331523
787 Double_t rmin[]={-1,-1,-1};
788 Double_t rmax[]={ 1, 1, 1};
789 TView *view = new TView3D(1,rmin,rmax);
791 TView *view = new TView(1);
794 Float_t range = fRrange*fRangeSlider->GetMaximum();
795 view->SetRange(-range,-range,-range,range, range, range);
796 // zoom back to full scale only if DrawView not called from NextCathode
805 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
807 // Recovering the chamber
808 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
810 const AliMUONGeometryTransformer* kGeomTransformer
811 = pMUON->GetGeometryTransformer();
813 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
815 // Display MUON Chamber Geometry
817 sprintf(nodeName,"MUON%d",100+fChamber);
818 printf(">>>> chamber is %d\n",fChamber);
822 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() ) {
824 Int_t detElemId = it.CurrentDEId();
825 AliMpSectorSegmentation * seg =
826 (AliMpSectorSegmentation *) AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
827 const AliMpSector * sector = seg->GetSector();
829 // get sector measurements
830 TVector2 position = sector->Position();
831 TVector2 dimension = sector->Dimensions(); // half length
833 Float_t xlocal1 = position.Px(); // FIXME: not really needed as it's 0 ?
834 Float_t ylocal1 = position.Py(); // FIXME: not really needed as it's 0 ?
835 Float_t xlocal2 = dimension.Px() * 2.;
836 Float_t ylocal2 = dimension.Px() * 2.;
838 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
839 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
842 TPolyLine3D* poly = new TPolyLine3D();
845 poly->SetPoint(nPoint++, xg1, yg1, 0.);
846 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
847 Float_t x = xg1 + xg2 * TMath::Cos(d);
848 Float_t y = yg1 + yg2 * TMath::Sin(d);
849 poly->SetPoint(nPoint++, x, y, 0.);
851 poly->SetPoint(nPoint++, xg1, yg1, 0.);
853 poly->SetLineColor(2);
862 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() )
864 Int_t detElemId = it.CurrentDEId();
865 AliMp::StationType stationType = AliMpDEManager::GetStationType(detElemId);
867 if ( segmentation->HasDE(detElemId) )
869 const AliMpVSegmentation* seg
870 = AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
873 Float_t deltax = seg->Dimensions().X();
874 Float_t deltay = seg->Dimensions().Y();
875 Float_t xlocal1 = -deltax;
876 Float_t ylocal1 = -deltay;
877 Float_t xlocal2 = +deltax;
878 Float_t ylocal2 = +deltay;
879 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
880 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
882 // drawing slat active volumes
883 Float_t xCenter = (xg1 + xg2)/2.;
884 Float_t yCenter = (yg1 + yg2)/2.;
886 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
888 box->SetFillStyle(0);
889 box->SetLineColor( stationType == AliMp::kStationTrigger ? 4 : 2);
892 if ( stationType == AliMp::kStation345 )
894 // drawing inner circle + disc
895 TPolyLine3D* poly = new TPolyLine3D();
896 TPolyLine3D* poly1 = new TPolyLine3D();
900 for (Float_t d = 0; d < 6.24; d+= 0.005)
902 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
903 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
904 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
905 poly->SetPoint(nPoint++, x, y, 0.);
906 poly1->SetPoint(nPoint1++, x, y, 0.);
908 poly->SetLineColor(1);
910 poly1->SetLineColor(2);
918 //add clusters to the pad
922 // DrawSegmentation();
923 // add itself to the list (must be last)
925 view->SetView(phi, theta, psi, iret);
928 //_____________________________________________________________________________
929 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
931 /// Draw a view of muons chambers with tracks
933 gPad->SetCursor(kWatch);
934 // gPad->SetFillColor(39);
935 gPad->SetFillColor(1);
937 // gPad->SetFillColor(39);
938 gPad->SetFillColor(1);
942 #if ROOT_VERSION_CODE>= 331523
943 Double_t rmin[]={-1,-1,-1};
944 Double_t rmax[]={ 1, 1, 1};
945 TView *view = new TView3D(1,rmin,rmax);
947 TView *view = new TView(1);
950 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
951 view->SetRange(-range,-range,-range,range,range,range);
953 // Display all MUON Chambers segmentation
956 sprintf(nodeName,"alice");
958 node1=gAlice->GetGeometry()->GetNode(nodeName);
959 if (node1) node1->Draw("same");
962 // Draw clusters for all chambers
963 Int_t chamberSave = fChamber;
964 for (fChamber = 1; fChamber <= 10; fChamber++){
967 fChamber = chamberSave;
968 // Draw reconstructed tracks
975 Float_t x0 = (-1+shift)/zoom;
976 Float_t y0 = (-1+shift)/zoom;
977 Float_t x1 = (1+shift)/zoom;
978 Float_t y1 = (1+shift)/zoom;
979 gPad->Range(x0,y0,x1,y1);
980 view->SetView(phi, theta, psi, iret);
984 //______________________________________________________________________________
985 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
987 /// Execute action corresponding to the mouse event
989 static Float_t x0, y0, x1, y1;
991 static Int_t pxold, pyold;
992 static Int_t px0, py0;
993 static Int_t linedrawn;
996 if (px == 0 && py == 0) { //when called by sliders
997 if (event == kButton1Up) {
1002 if (!fZoomMode && gPad->GetView()) {
1003 gPad->GetView()->ExecuteRotateView(event, px, py);
1007 // something to zoom ?
1008 gPad->SetCursor(kCross);
1013 gVirtualX->SetLineColor(-1);
1014 gPad->TAttLine::Modify(); //Change line attributes only if necessary
1015 x0 = gPad->AbsPixeltoX(px);
1016 y0 = gPad->AbsPixeltoY(py);
1018 pxold = px; pyold = py;
1022 case kButton1Motion:
1023 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1027 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1031 gPad->GetCanvas()->FeedbackMode(kFALSE);
1032 if (px == px0) return;
1033 if (py == py0) return;
1034 x1 = gPad->AbsPixeltoX(px);
1035 y1 = gPad->AbsPixeltoY(py);
1037 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1038 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1039 gPad->Range(x0,y0,x1,y1);
1040 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1042 fZoomX0[fZooms] = x0;
1043 fZoomY0[fZooms] = y0;
1044 fZoomX1[fZooms] = x1;
1045 fZoomY1[fZooms] = y1;
1047 gPad->Modified(kTRUE);
1052 //___________________________________________
1053 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1055 /// Read digits info and store x,y,z info in arrays fPoints.
1056 /// Loop on all detectors
1058 if (chamber > 14) return;
1064 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1066 GetMUONData()->SetTreeAddress("D");
1068 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
1069 if (muonDigits == 0) return;
1071 gAlice->ResetDigits();
1074 if (GetLoader()->TreeD()) {
1075 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
1076 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1077 GetMUONData()->GetDigits();
1080 Int_t ndigits = muonDigits->GetEntriesFast();
1081 if (ndigits == 0) return;
1082 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1084 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1087 AliMUONPoints *points = 0;
1088 TMarker3DBox *marker = 0;
1091 Float_t adcmax = 1024; // default
1092 if (chamber<11) adcmax = 4096;
1094 // check if trigger is using new or old segmentation
1096 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
1097 const AliMUONVGeometryDESegmentation* kdeSegmentation
1098 = segmentation->GetDESegmentation(1100, cathode-1);
1099 if ( dynamic_cast<const AliMUONTriggerSegmentation*>(kdeSegmentation) ) old = false;
1101 if ( old && chamber > 10) {
1102 if (chamber > 10) printf(">>> old segmentation for trigger \n");
1103 else printf(">>> old segmentation for tracking \n");
1105 for (Int_t digit = 0; digit < ndigits; digit++) {
1106 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1107 if (mdig->Cathode() != cathode-1) continue;
1110 // First get all needed parameters
1112 Float_t charge = mdig->Signal();
1113 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1114 Int_t color = 261+index;
1115 Int_t colorTrigger = 2;
1116 if (color > 282) color = 282;
1118 if (chamber > 10) { // trigger chamber
1120 Float_t sumCharge = 0;
1121 for (Int_t icharge = 0; icharge < 10; icharge++) {
1122 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1124 Float_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1125 if(sumCharge <= 10 || testCharge > 0) {
1126 colorTrigger = color;
1132 // get the center of the pad - add on x and y half of pad size
1133 Float_t xpad, ypad, zpad;
1137 Int_t detElemId = mdig->DetElemId();
1138 AliMUONGeometrySegmentation* segmentation2
1139 = pMUON->GetSegmentation()->GetModuleSegmentationByDEId(detElemId, cathode-1);
1140 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1141 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1142 dpx = segmentation2->Dpx(detElemId, isec)/2;
1143 dpy = segmentation2->Dpy(detElemId, isec)/2;
1145 // Then set the objects
1146 points = new AliMUONPoints(npoints);
1147 fPoints->AddAt(points,digit);
1149 points->SetMarkerColor(colorTrigger);
1151 points->SetMarkerColor(color);
1153 points->SetMarkerStyle(21);
1154 points->SetMarkerSize(0.5);
1155 points->SetParticle(-1);
1156 points->SetHitIndex(-1);
1157 points->SetTrackIndex(-1);
1158 points->SetDigitIndex(digit);
1159 points->SetPoint(0,xpad,ypad,zpos);
1161 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1162 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1165 marker->SetLineColor(lineColor);
1166 marker->SetFillStyle(1001);
1167 marker->SetFillColor(color);
1168 marker->SetRefObject((TObject*)points);
1169 points->Set3DMarker(0, marker);
1170 } // end loop on digits
1173 if (chamber > 10) printf(">>> new segmentation for trigger \n");
1174 else printf(">>> new segmentation for tracking \n");
1176 const AliMUONGeometryTransformer* kGeomTransformer
1177 = pMUON->GetGeometryTransformer();
1179 //loop over all digits and store their position
1180 for (Int_t digit = 0; digit < ndigits; digit++) {
1181 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1182 if (mdig->Cathode() != cathode-1) continue;
1184 // get all needed parameters
1185 Int_t ix=mdig->PadX();
1186 Int_t iy=mdig->PadY();
1187 Int_t detElemId=mdig->DetElemId();
1188 Float_t charge = mdig->Signal();
1189 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1190 Int_t color = 261+index;
1191 Int_t colorTrigger = 2;
1192 if (color > 282) color = 282;
1194 const AliMpVSegmentation* seg =
1195 AliMpSegmentation::Instance()
1196 ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode-1));
1198 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1200 if (chamber > 10) { // trigger chamber
1201 Float_t sumCharge = 0;
1202 Int_t n = mdig->Ntracks();
1203 for (Int_t icharge = 0; icharge < n; icharge++) {
1204 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1206 Float_t testCharge = sumCharge-(Int_t(sumCharge/n))*n;
1207 if(sumCharge <= n || testCharge > 0) {
1208 colorTrigger = color;
1214 // get the pad position and dimensions
1215 Float_t xlocal1 = pad.Position().X();
1216 Float_t ylocal1 = pad.Position().Y();
1217 Float_t xlocal2 = pad.Dimensions().X();
1218 Float_t ylocal2 = pad.Dimensions().Y();
1220 Float_t xg1, xg2, yg1, yg2, zg1;
1222 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
1223 // (no transformation for pad dimensions)
1227 // Then set the objects
1228 points = new AliMUONPoints(npoints);
1229 fPoints->AddAt(points,digit);
1231 points->SetMarkerColor(colorTrigger);
1233 points->SetMarkerColor(color);
1235 points->SetMarkerStyle(21);
1236 points->SetMarkerSize(0.5);
1237 points->SetParticle(-1);
1238 points->SetHitIndex(-1);
1239 points->SetTrackIndex(-1);
1240 points->SetDigitIndex(digit);
1241 points->SetPoint(0,xg1,yg1,zpos);
1243 Int_t lineColor = (zg1-zpos > 0) ? 2:3;
1244 marker=new TMarker3DBox(xg1,yg1,zpos,xg2,yg2,0,0,0);
1246 marker->SetLineColor(lineColor);
1247 marker->SetFillStyle(1001);
1248 marker->SetFillColor(color);
1249 marker->SetRefObject((TObject*)points);
1250 points->Set3DMarker(0, marker);
1252 } // end loop on digits
1253 } // end of new segmentation
1255 //___________________________________________
1256 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1258 /// Read raw clusters info and store x,y,z info in arrays fRpoints.
1259 /// Loop on all detectors
1261 if (chamber > 10) return;
1265 GetMUONData()->SetTreeAddress("RC");
1266 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1268 if (muonRawClusters == 0) return;
1271 if (GetMUONData()->TreeR()) {
1272 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1273 GetMUONData()->TreeR()->GetEvent(0);
1276 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1277 if (nrawcl == 0) return;
1278 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1280 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1281 AliMUONRawCluster *mRaw;
1282 AliMUONPoints *points = 0;
1284 //loop over all raw clusters and store their position
1285 points = new AliMUONPoints(nrawcl);
1286 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1287 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1288 points->SetMarkerColor(51);
1289 points->SetMarkerStyle(2);
1290 points->SetMarkerSize(1.);
1291 points->SetParticle(-1);
1292 points->SetHitIndex(-1);
1293 points->SetTrackIndex(-1);
1294 points->SetDigitIndex(-1);
1295 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1296 fRpoints->AddAt(points,iraw);
1297 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1301 //___________________________________________
1302 void AliMUONDisplay::LoadTracks()
1306 AliMUONTrack* recTrack = 0;
1307 AliMUONTrackParam* trackParam = 0;
1308 TClonesArray * trackParamAtHit = 0;
1312 GetMUONData()->SetTreeAddress("RT");
1313 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1314 if (recTracksArray == NULL) return;
1315 GetMUONData()->GetRecTracks();
1317 Int_t nRecTracks = 0;
1319 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1322 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1324 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1325 // reading info from tracks
1326 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1328 Int_t nTrackHits = recTrack->GetNTrackHits();
1330 if (nTrackHits == 0) continue;
1333 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1334 points->SetLineColor(6);
1335 points->SetLineWidth(1);
1336 fRpoints->AddAt(points,iRecTracks);
1342 // vertex unknown at the tracking level -> put it at (0,0,0)
1343 points->SetPoint(iPoint,0.,0.,0.);
1346 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1347 trackParamAtHit = recTrack->GetTrackParamAtHit();
1348 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1349 xRec = trackParam->GetNonBendingCoor();
1350 yRec = trackParam->GetBendingCoor();
1351 zRec = trackParam->GetZ();
1352 points->SetPoint(iPoint,xRec,yRec,zRec);
1354 } // end loop rec. hits
1355 PrintTrack(iRecTracks,recTrack);
1356 } // end loop tracks
1361 //___________________________________________
1362 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1364 /// Print reconstructed track
1366 AliMUONTrackParam *trackParam;
1367 Float_t vertex[3], momentum[3];
1368 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1371 trackParam = recTrack->GetTrackParamAtVertex(); // meaningless since the vertex is not known at the tracking level
1372 vertex[0] = trackParam->GetNonBendingCoor();
1373 vertex[1] = trackParam->GetBendingCoor();
1374 vertex[2] = trackParam->GetZ();
1375 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1376 bendingSlope = trackParam->GetBendingSlope();
1377 nonBendingSlope = trackParam->GetNonBendingSlope();
1378 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1379 momentum[0] = momentum[2] * nonBendingSlope;
1380 momentum[1] = momentum[2] * bendingSlope;
1381 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1382 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1384 printf("===================================================\n");
1385 printf("//*****************************************************************//\n");
1386 printf("// meaningless since the vertex is not known at the tracking level //\n");
1387 printf("//*****************************************************************//\n");
1388 printf(" Reconstructed track # %d \n",iRecTracks);
1389 printf(" charge: %d \n",charge);
1390 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1391 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1392 printf(" track chi2/dof: %f \n",chi2dof);
1396 //___________________________________________
1397 void AliMUONDisplay::LoadHits(Int_t chamber)
1399 /// Read hits info and store x,y,z info in arrays fPhits.
1400 /// Loop on all detectors
1402 if (chamber > 14) return;
1409 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1411 if (GetMUONData()->TreeH()) {
1412 GetMUONData()->SetTreeAddress("H");
1413 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
1415 for (track = 0; track < ntracks; track++) {
1416 GetMUONData()->ResetHits();
1417 GetMUONData()->GetTrack(track);//skowron
1418 TClonesArray *muonHits = GetMUONData()->Hits();
1419 if (muonHits == 0) return;
1420 nthits += muonHits->GetEntriesFast();
1422 if (fPhits == 0) fPhits = new TObjArray(nthits);
1424 for (track=0; track<ntracks;track++) {
1425 GetMUONData()->ResetHits();
1426 GetMUONData()->GetTrack(track);//skowron
1427 TClonesArray *muonHits = GetMUONData()->Hits();
1428 if (muonHits == 0) return;
1429 Int_t nhits = muonHits->GetEntriesFast();
1430 if (nhits == 0) continue;
1432 AliMUONPoints *points = 0;
1434 for (Int_t hit=0;hit<nhits;hit++) {
1435 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1436 Int_t nch = mHit->Chamber(); // chamber number
1437 if (nch != chamber) continue;
1439 // Retrieve info and set the objects
1441 points = new AliMUONPoints(npoints);
1442 fPhits->AddAt(points,nhold+hit);
1443 points->SetMarkerColor(kRed);
1444 points->SetMarkerStyle(5);
1445 points->SetMarkerSize(1.);
1446 points->SetParticle(mHit->Track());
1447 points->SetHitIndex(hit);
1448 points->SetTrackIndex(track);
1449 points->SetDigitIndex(-1);
1450 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1451 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1458 //_____________________________________________________________________________
1459 void AliMUONDisplay::Paint(Option_t *)
1461 /// Paint miscellaneous items
1464 //_____________________________________________________________________________
1465 void AliMUONDisplay::SetPickMode()
1467 /// Set parameters for pick mode.
1471 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1472 fTrigPad->Modified();
1475 //_____________________________________________________________________________
1476 void AliMUONDisplay::SetZoomMode()
1478 /// Set parameters for zoom mode
1482 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1483 fTrigPad->Modified();
1486 //_____________________________________________________________________________
1487 void AliMUONDisplay::NextChamber(Int_t delta)
1489 /// To go from chamber to next chamber if delta = 1
1490 /// or previous chamber otherwise
1492 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1494 if (fChamber > 1) fChamber--;
1498 LoadDigits(fChamber, fCathode);
1502 //_____________________________________________________________________________
1503 void AliMUONDisplay::NextCathode()
1505 /// To switch to other cathode plane
1509 if (fCathode == 1) {
1510 LoadDigits(fChamber, 2);
1512 LoadDigits(fChamber, 1);
1514 fNextCathode = kTRUE; // to keep the same zoom
1516 fNextCathode = kFALSE;
1517 TPad *pad = (TPad*)gPad->GetPadSave();
1518 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1519 fZoomX1[fZooms], fZoomY1[fZooms]);
1525 //_____________________________________________________________________________
1526 void AliMUONDisplay::Trigger()
1528 /// Print global trigger output
1530 AliMUONGlobalTrigger* globalTrig;
1532 GetMUONData()->SetTreeAddress("GLT");
1533 GetMUONData()->GetTriggerD();
1535 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1536 if (globalTrig == 0) return;
1538 globalTrig->Print("full");
1540 // // returns Trigger Decision for current event
1541 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1543 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1544 // AliMUONData* muonData = decision->GetMUONData();
1545 // muonData->SetTreeAddress("D");
1546 // decision->Trigger();
1548 //_____________________________________________________________________________
1549 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1551 /// Set chamber and cathode number
1558 LoadDigits(chamber,cathode);
1562 //_____________________________________________________________________________
1563 void AliMUONDisplay::SetEvent(Int_t newevent)
1567 gAlice->GetEvent(newevent);
1569 if (!gAlice->TreeD()) return;
1572 LoadDigits(fChamber,fCathode);
1576 //_____________________________________________________________________________
1577 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1579 /// Set view range along R and Z
1589 //_____________________________________________________________________________
1590 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1592 /// Change viewing angles for current event
1600 TView *view = gPad->GetView();
1601 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1606 //_____________________________________________________________________________
1607 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1609 /// Display (current event_number + delta)
1610 /// - delta = 1 shown next event
1611 /// - delta = -1 show previous event
1613 AliRunLoader * runLoader;
1615 runLoader = fLoader->GetRunLoader();
1620 //runLoader->CleanDetectors();
1621 //runLoader->CleanKinematics();
1622 Int_t currentEvent = runLoader->GetEventNumber();
1623 Int_t newEvent = currentEvent + delta;
1624 runLoader->GetEvent(newEvent);
1627 LoadDigits(fChamber, fCathode);
1632 //______________________________________________________________________________
1633 void AliMUONDisplay::UnZoom()
1637 if (fZooms <= 0) return;
1639 TPad *pad = (TPad*)gPad->GetPadSave();
1640 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1644 //_____________________________________________________________________________
1645 void AliMUONDisplay::ResetPoints()
1647 /// Reset array of points
1655 //_____________________________________________________________________________
1656 void AliMUONDisplay::ResetPhits()
1658 /// Reset array of points
1666 //_____________________________________________________________________________
1667 void AliMUONDisplay::ResetRpoints()
1669 /// Reset array of points