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"
82 #include "AliMUONSimData.h"
83 #include "AliMUONRecData.h"
85 #include "AliMpDEIterator.h"
86 #include "AliMpSegmentation.h"
87 #include "AliMpSlatSegmentation.h"
88 #include "AliMpSlat.h"
89 #include "AliMpSectorSegmentation.h"
90 #include "AliMpSector.h"
91 #include "AliMpTriggerSegmentation.h"
92 #include "AliMpTrigger.h"
93 #include "AliMpStationType.h"
94 #include "AliMpCathodType.h"
95 #include "AliMpDEManager.h"
100 #include "AliHeader.h"
105 #if ROOT_VERSION_CODE>= 331523
111 #include <TPaveLabel.h>
112 #include <TPaveText.h>
113 #include <TDiamond.h>
117 #include <TVirtualX.h>
119 #include <TGeometry.h>
120 #include <TMarker3DBox.h>
121 #include <TParticle.h>
122 #include <TPolyLine3D.h>
126 ClassImp(AliMUONDisplay)
129 //_____________________________________________________________________________
130 AliMUONDisplay::AliMUONDisplay()
135 fDrawClusters(kTRUE),
149 /// Default constructor
152 //_____________________________________________________________________________
153 AliMUONDisplay::AliMUONDisplay(Int_t size,
154 AliLoader * simLoader, AliLoader * recLoader)
159 fDrawClusters(kTRUE),
167 fNextCathode(kFALSE),
168 fSimLoader(simLoader),
169 fRecLoader(recLoader),
174 /// Standard constructor to create an event display object.
178 gAlice->SetDisplay(this);
180 // Initialize display default parameters
183 // Set front view by default
192 // Create display canvas
194 if (ysize < 100) ysize = 750;
195 Int_t xsize = Int_t(size*830./ysize);
196 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
197 fCanvas->ToggleEventStatus();
199 // Create main display pad
200 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
203 fPad->SetFillColor(30);
204 fPad->SetBorderSize(2);
209 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
211 fColPad->SetFillColor(17);
212 fColPad->SetBorderSize(2);
217 // Create user interface control pad
221 // Create Range and mode pad
224 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
227 fTrigPad->SetFillColor(22);
228 fTrigPad->SetBorderSize(2);
229 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
230 fRangeSlider->SetObject(this);
231 char pickmode[] = "gAlice->Display()->SetPickMode()";
233 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
234 fPickButton->SetFillColor(38);
236 char zoommode[] = "gAlice->Display()->SetZoomMode()";
237 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
238 fZoomButton->SetFillColor(38);
240 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
241 fArcButton->SetFillColor(kGreen);
243 char butUnzoom[] = "gAlice->Display()->UnZoom()";
244 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
245 button->SetFillColor(38);
247 AppendPad(); // append display object as last object to force selection
250 fTrigPad->SetEditable(kFALSE);
251 fButtons->SetEditable(kFALSE);
254 // initialize container
256 fMUONSimData = new AliMUONSimData(fSimLoader,"MUON","MUON");
258 fMUONRecData = new AliMUONRecData(fRecLoader,"MUON","MUON");
261 //_____________________________________________________________________________
262 AliMUONDisplay::~AliMUONDisplay()
266 // Delete space point structure
267 if (fPoints) fPoints->Delete();
271 if (fPhits) fPhits->Delete();
275 if (fRpoints) fRpoints->Delete();
280 //_____________________________________________________________________________
281 void AliMUONDisplay::Clear(Option_t *)
283 /// Delete graphics temporary objects
286 //_____________________________________________________________________________
287 void AliMUONDisplay::DisplayButtons()
289 /// Create the user interface buttons
292 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
294 fButtons->SetFillColor(38);
295 fButtons->SetBorderSize(2);
298 // Int_t butcolor = 33;
299 Float_t dbutton = 0.08;
306 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
307 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
308 button->SetFillColor(38);
312 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
313 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
314 button->SetFillColor(38);
318 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
319 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
320 button->SetFillColor(38);
324 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
325 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
326 button->SetFillColor(38);
330 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
331 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
332 button->SetFillColor(38);
336 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
337 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
338 button->SetFillColor(38);
342 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
343 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
344 button->SetFillColor(38);
348 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
349 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
350 button->SetFillColor(38);
354 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
355 diamond->SetFillColor(50);
356 diamond->SetTextAlign(22);
357 diamond->SetTextColor(5);
358 diamond->SetTextSize(0.11);
360 diamond->AddText(".. ");
361 diamond->AddText("ROOT");
362 diamond->AddText("MUON");
363 diamond->AddText("... ");
364 diamond->AddText(" ");
367 //_____________________________________________________________________________
368 void AliMUONDisplay::CreateColors() const
370 /// Create the colors palette used to display clusters
385 new TColor(color,r,g,b);
395 new TColor(color,r,g,b);
405 new TColor(color,r,g,b);
415 new TColor(color,r,g,b);
425 new TColor(color,r,g,b);
432 //_____________________________________________________________________________
433 void AliMUONDisplay::DisplayColorScale()
435 /// Display pulse height color scale
439 Float_t xlow, ylow, xup, yup, hs;
440 Float_t x1, y1, x2, y2;
444 TText *text = new TText(0,0,"");
445 text->SetTextFont(61);
446 text->SetTextSize(0.2);
447 text->SetTextAlign(22);
450 Int_t adcmax=4096; // default 12 bits ADC
456 //*-* draw colortable boxes
457 hs = (y2-y1)/Float_t(22);
461 ylow = y1 + hs*(Float_t(i));
462 yup = y1 + hs*(Float_t(i+1));
464 Double_t logscale=Double_t(i+1)*(TMath::Log(adcmax)/22);
465 Int_t scale=(Int_t)TMath::Exp(logscale);
466 sprintf(label,"%d",scale);
467 box = new TBox(xlow, ylow, xup, yup);
469 box->SetFillColor(color);
470 text->DrawText(xlow+0.7, 0.5*(ylow+yup),label);
474 //______________________________________________________________________________
475 Int_t AliMUONDisplay::DistancetoPrimitive(Int_t px, Int_t)
477 /// Compute distance from point px,py to objects in event
479 gPad->SetCursor(kCross);
481 if (gPad == fTrigPad) return 9999;
483 const Int_t kBig = 9999;
485 Float_t xmin = gPad->GetX1();
486 Float_t xmax = gPad->GetX2();
487 Float_t dx = 0.02*(xmax - xmin);
488 Float_t x = gPad->AbsPixeltoX(px);
489 if (x < xmin+dx || x > xmax-dx) return dist;
491 if (fZoomMode) return 0;
495 //_____________________________________________________________________________
496 void AliMUONDisplay::Draw(Option_t *)
498 /// Display current event
506 //_____________________________________________________________________________
507 void AliMUONDisplay::DrawChamber()
509 /// Display current event
511 fDrawTracks = kFALSE;
513 DrawView(fTheta, fPhi, fPsi);
514 // Display the event number and title
519 //_____________________________________________________________________________
520 void AliMUONDisplay::DrawReco(Option_t *)
522 /// Display current event
525 // print kinematics of generated particles
527 // Draw global view of muon system
529 DrawGlobalView(135, -50, -140);
531 // Display the event number and title
536 //_____________________________________________________________________________
537 void AliMUONDisplay::PrintKinematics()
539 /// Print kinematic tree
541 if ( !fSimLoader || !fRecLoader ) {
542 AliErrorStream() << "Detector loaders are not defined." << endl;
546 AliRunLoader* runLoader = fSimLoader->GetRunLoader();
547 TParticle *particle = new TParticle();
549 Float_t vertex[3], momentum[3];
551 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
552 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
553 nPart = (Int_t)runLoader->TreeK()->GetEntries();
554 for(Int_t iPart = 0; iPart < nPart; iPart++) {
555 runLoader->TreeK()->GetEvent(iPart);
556 vertex[0] = particle->Vx();
557 vertex[1] = particle->Vy();
558 vertex[2] = particle->Vz();
559 momentum[0] = particle->Px();
560 momentum[1] = particle->Py();
561 momentum[2] = particle->Pz();
563 printf("===================================================\n");
564 printf(" Generated particle # %d \n",iPart);
565 printf(" name: %s \n",particle->GetName());
566 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
567 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
572 //_____________________________________________________________________________
573 void AliMUONDisplay::DrawSegmentation()
575 /// \todo to be re-written for new seg
576 /// Draw graphical representation of segmenatation
577 /// Attention: still experimental code
580 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
581 // AliMUONChamber* iChamber;
583 // AliSegmentation* seg;
584 // iChamber = &(pMUON->Chamber(fChamber));
585 // seg=iChamber->SegmentationModel(icat);
587 // Float_t zpos=iChamber->Z();
588 // Float_t r=iChamber->ROuter();
590 // TMarker3DBox *marker;
592 // for (Int_t j=0; j<seg->Npy(); j++) {
594 // y0=j*seg->Dpy()-seg->Dpy()/2.;
595 // for (seg->FirstPad(0.,y0,0,300,0.);
599 // if (seg->ISector()==0) continue;
601 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
602 // Float_t dpx=seg->Dpx(seg->ISector())/2;
603 // Float_t dpy=seg->Dpy(seg->ISector())/2;
604 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
605 // marker->SetLineColor(seg->ISector()+1);
606 // marker->SetFillStyle(1001);
607 // marker->SetFillColor(0);
612 // for (Int_t j=0; j<250; j++) {
613 // Float_t x0=j*seg->Dpx();
614 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
616 // for (seg->FirstPad(x0,0,0,0,y0);
620 // if (seg->ISector()==0) continue;
623 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
624 // Float_t dpx=seg->Dpx(seg->ISector())/2;
625 // Float_t dpy=seg->Dpy(seg->ISector())/2;
626 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
627 // marker->SetLineColor(seg->ISector()+1);
628 // marker->SetFillStyle(1001);
629 // marker->SetFillColor(0);
636 //_____________________________________________________________________________
637 void AliMUONDisplay::DrawClusters()
639 /// Draw clusters for MUON chambers
641 Int_t ndigits, digit;
649 ndigits = points->GetEntriesFast();
650 for (digit=0;digit<ndigits;digit++){
651 pm = (AliMUONPoints*)points->UncheckedAt(digit);
655 for (Int_t im=0;im<3;im++) {
656 TMarker3DBox *marker=pm->GetMarker(im);
661 fClustersCuts += pm->GetN();
665 //_____________________________________________________________________________
666 void AliMUONDisplay::DrawHits()
668 /// Draw hits for MUON chambers
672 Int_t ntracks, track;
679 ntracks = points->GetEntriesFast();
680 for (track=0;track<ntracks;track++) {
681 pm = (AliMUONPoints*)points->UncheckedAt(track);
684 fHitsCuts += pm->GetN();
689 //_____________________________________________________________________________
690 void AliMUONDisplay::DrawCoG()
692 /// Draw hits for MUON chambers
694 if (!fDrawCoG) return;
695 if (fChamber > 10) return;
696 LoadCoG(fChamber,fCathode);
704 ncog = points->GetEntriesFast();
705 for (icog=0;icog<ncog;icog++) {
706 pm = (AliMUONPoints*)points->UncheckedAt(icog);
711 //_____________________________________________________________________________
712 void AliMUONDisplay::DrawTracks()
716 if (!fDrawTracks) return;
719 Int_t nTrack, iTrack;
725 nTrack = points->GetEntriesFast();
726 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
727 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
732 //_____________________________________________________________________________
734 void AliMUONDisplay::DrawTitle(Option_t *option)
736 /// Draw the event title
738 if ( !fSimLoader || !fRecLoader ) {
739 AliErrorStream() << "Detector loaders are not defined." << endl;
743 Float_t xmin = gPad->GetX1();
744 Float_t xmax = gPad->GetX2();
745 Float_t ymin = gPad->GetY1();
746 Float_t ymax = gPad->GetY2();
747 Float_t dx = xmax-xmin;
748 Float_t dy = ymax-ymin;
750 AliRunLoader * runLoader = fSimLoader->GetRunLoader();
752 if (strlen(option) == 0) {
753 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
754 // title->SetTextSize(0.023932);
755 title->SetTextSize(0.02);
756 title->SetBit(kCanDelete);
757 title->SetFillColor(42);
760 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
761 runLoader->GetEventNumber(),
762 gAlice->GetHeader()->GetRun(),
765 title->AddText(ptitle);
766 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
767 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
768 nparticles, fHitsCuts,fClustersCuts);
769 title->AddText(ptitle);
771 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
772 label->SetBit(kCanDelete);
773 label->SetFillColor(42);
778 //_____________________________________________________________________________
779 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
781 /// Draw a view of MUON clusters
783 AliInfo(" Draw View");
785 gPad->SetCursor(kWatch);
786 // gPad->SetFillColor(39);
787 gPad->SetFillColor(1);
789 // gPad->SetFillColor(39);
790 gPad->SetFillColor(1);
793 #if ROOT_VERSION_CODE>= 331523
794 Double_t rmin[]={-1,-1,-1};
795 Double_t rmax[]={ 1, 1, 1};
796 TView *view = new TView3D(1,rmin,rmax);
798 TView *view = new TView(1);
801 Float_t range = fRrange*fRangeSlider->GetMaximum();
802 view->SetRange(-range,-range,-range,range, range, range);
803 // zoom back to full scale only if DrawView not called from NextCathode
812 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
814 // Recovering the chamber
815 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
817 const AliMUONGeometryTransformer* kGeomTransformer
818 = pMUON->GetGeometryTransformer();
820 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
822 // Display MUON Chamber Geometry
824 sprintf(nodeName,"MUON%d",100+fChamber);
825 printf(">>>> chamber is %d\n",fChamber);
829 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() ) {
831 Int_t detElemId = it.CurrentDEId();
832 AliMpSectorSegmentation * seg =
833 (AliMpSectorSegmentation *) AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
834 const AliMpSector * sector = seg->GetSector();
836 // get sector measurements
837 TVector2 position = sector->Position();
838 TVector2 dimension = sector->Dimensions(); // half length
840 Float_t xlocal1 = position.Px(); // FIXME: not really needed as it's 0 ?
841 Float_t ylocal1 = position.Py(); // FIXME: not really needed as it's 0 ?
842 Float_t xlocal2 = dimension.Px() * 2.;
843 Float_t ylocal2 = dimension.Px() * 2.;
845 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
846 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
849 TPolyLine3D* poly = new TPolyLine3D();
852 poly->SetPoint(nPoint++, xg1, yg1, 0.);
853 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
854 Float_t x = xg1 + xg2 * TMath::Cos(d);
855 Float_t y = yg1 + yg2 * TMath::Sin(d);
856 poly->SetPoint(nPoint++, x, y, 0.);
858 poly->SetPoint(nPoint++, xg1, yg1, 0.);
860 poly->SetLineColor(2);
869 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() )
871 Int_t detElemId = it.CurrentDEId();
872 AliMp::StationType stationType = AliMpDEManager::GetStationType(detElemId);
874 if ( segmentation->HasDE(detElemId) )
876 const AliMpVSegmentation* seg
877 = AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
880 Float_t deltax = seg->Dimensions().X();
881 Float_t deltay = seg->Dimensions().Y();
882 Float_t xlocal1 = -deltax;
883 Float_t ylocal1 = -deltay;
884 Float_t xlocal2 = +deltax;
885 Float_t ylocal2 = +deltay;
886 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
887 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
889 // drawing slat active volumes
890 Float_t xCenter = (xg1 + xg2)/2.;
891 Float_t yCenter = (yg1 + yg2)/2.;
893 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
895 box->SetFillStyle(0);
896 box->SetLineColor( stationType == AliMp::kStationTrigger ? 4 : 2);
899 if ( stationType == AliMp::kStation345 )
901 // drawing inner circle + disc
902 TPolyLine3D* poly = new TPolyLine3D();
903 TPolyLine3D* poly1 = new TPolyLine3D();
907 for (Float_t d = 0; d < 6.24; d+= 0.005)
909 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
910 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
911 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
912 poly->SetPoint(nPoint++, x, y, 0.);
913 poly1->SetPoint(nPoint1++, x, y, 0.);
915 poly->SetLineColor(1);
917 poly1->SetLineColor(2);
925 //add clusters to the pad
929 // DrawSegmentation();
930 // add itself to the list (must be last)
932 view->SetView(phi, theta, psi, iret);
935 //_____________________________________________________________________________
936 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
938 /// Draw a view of muons chambers with tracks
940 gPad->SetCursor(kWatch);
941 // gPad->SetFillColor(39);
942 gPad->SetFillColor(1);
944 // gPad->SetFillColor(39);
945 gPad->SetFillColor(1);
949 #if ROOT_VERSION_CODE>= 331523
950 Double_t rmin[]={-1,-1,-1};
951 Double_t rmax[]={ 1, 1, 1};
952 TView *view = new TView3D(1,rmin,rmax);
954 TView *view = new TView(1);
957 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
958 view->SetRange(-range,-range,-range,range,range,range);
960 // Display all MUON Chambers segmentation
963 sprintf(nodeName,"alice");
965 node1=gAlice->GetGeometry()->GetNode(nodeName);
966 if (node1) node1->Draw("same");
969 // Draw clusters for all chambers
970 Int_t chamberSave = fChamber;
971 for (fChamber = 1; fChamber <= 10; fChamber++){
974 fChamber = chamberSave;
975 // Draw reconstructed tracks
982 Float_t x0 = (-1+shift)/zoom;
983 Float_t y0 = (-1+shift)/zoom;
984 Float_t x1 = (1+shift)/zoom;
985 Float_t y1 = (1+shift)/zoom;
986 gPad->Range(x0,y0,x1,y1);
987 view->SetView(phi, theta, psi, iret);
991 //______________________________________________________________________________
992 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
994 /// Execute action corresponding to the mouse event
996 static Float_t x0, y0, x1, y1;
998 static Int_t pxold, pyold;
999 static Int_t px0, py0;
1000 static Int_t linedrawn;
1003 if (px == 0 && py == 0) { //when called by sliders
1004 if (event == kButton1Up) {
1009 if (!fZoomMode && gPad->GetView()) {
1010 gPad->GetView()->ExecuteRotateView(event, px, py);
1014 // something to zoom ?
1015 gPad->SetCursor(kCross);
1020 gVirtualX->SetLineColor(-1);
1021 gPad->TAttLine::Modify(); //Change line attributes only if necessary
1022 x0 = gPad->AbsPixeltoX(px);
1023 y0 = gPad->AbsPixeltoY(py);
1025 pxold = px; pyold = py;
1029 case kButton1Motion:
1030 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1034 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1038 gPad->GetCanvas()->FeedbackMode(kFALSE);
1039 if (px == px0) return;
1040 if (py == py0) return;
1041 x1 = gPad->AbsPixeltoX(px);
1042 y1 = gPad->AbsPixeltoY(py);
1044 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1045 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1046 gPad->Range(x0,y0,x1,y1);
1047 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1049 fZoomX0[fZooms] = x0;
1050 fZoomY0[fZooms] = y0;
1051 fZoomX1[fZooms] = x1;
1052 fZoomY1[fZooms] = y1;
1054 gPad->Modified(kTRUE);
1059 //___________________________________________
1060 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1062 /// Read digits info and store x,y,z info in arrays fPoints.
1063 /// Loop on all detectors
1065 if (chamber > 14) return;
1071 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1073 GetMUONSimData()->SetTreeAddress("D");
1075 TClonesArray *muonDigits = GetMUONSimData()->Digits(chamber-1);
1076 if (muonDigits == 0) return;
1078 gAlice->ResetDigits();
1081 if (GetSimLoader()->TreeD()) {
1082 nent = (Int_t) GetSimLoader()->TreeD()->GetEntries();
1083 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1084 GetMUONSimData()->GetDigits();
1087 Int_t ndigits = muonDigits->GetEntriesFast();
1088 if (ndigits == 0) return;
1089 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1091 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1094 AliMUONPoints *points = 0;
1095 TMarker3DBox *marker = 0;
1098 Float_t adcmax = 1024; // default
1099 if (chamber<11) adcmax = 4096;
1101 // check if trigger is using new or old segmentation
1103 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
1104 const AliMUONVGeometryDESegmentation* kdeSegmentation
1105 = segmentation->GetDESegmentation(1100, cathode-1);
1106 if ( dynamic_cast<const AliMUONTriggerSegmentation*>(kdeSegmentation) ) old = false;
1108 if ( old && chamber > 10) {
1109 if (chamber > 10) printf(">>> old segmentation for trigger \n");
1110 else printf(">>> old segmentation for tracking \n");
1112 for (Int_t digit = 0; digit < ndigits; digit++) {
1113 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1114 if (mdig->Cathode() != cathode-1) continue;
1117 // First get all needed parameters
1119 Float_t charge = mdig->Signal();
1120 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1121 Int_t color = 261+index;
1122 Int_t colorTrigger = 2;
1123 if (color > 282) color = 282;
1125 if (chamber > 10) { // trigger chamber
1127 Float_t sumCharge = 0;
1128 for (Int_t icharge = 0; icharge < 10; icharge++) {
1129 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1131 Float_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1132 if(sumCharge <= 10 || testCharge > 0) {
1133 colorTrigger = color;
1139 // get the center of the pad - add on x and y half of pad size
1140 Float_t xpad, ypad, zpad;
1144 Int_t detElemId = mdig->DetElemId();
1145 AliMUONGeometrySegmentation* segmentation2
1146 = pMUON->GetSegmentation()->GetModuleSegmentationByDEId(detElemId, cathode-1);
1147 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1148 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1149 dpx = segmentation2->Dpx(detElemId, isec)/2;
1150 dpy = segmentation2->Dpy(detElemId, isec)/2;
1152 // Then set the objects
1153 points = new AliMUONPoints(npoints);
1154 fPoints->AddAt(points,digit);
1156 points->SetMarkerColor(colorTrigger);
1158 points->SetMarkerColor(color);
1160 points->SetMarkerStyle(21);
1161 points->SetMarkerSize(0.5);
1162 points->SetParticle(-1);
1163 points->SetHitIndex(-1);
1164 points->SetTrackIndex(-1);
1165 points->SetDigitIndex(digit);
1166 points->SetPoint(0,xpad,ypad,zpos);
1168 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1169 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1172 marker->SetLineColor(lineColor);
1173 marker->SetFillStyle(1001);
1174 marker->SetFillColor(color);
1175 marker->SetRefObject((TObject*)points);
1176 points->Set3DMarker(0, marker);
1177 } // end loop on digits
1180 if (chamber > 10) printf(">>> new segmentation for trigger \n");
1181 else printf(">>> new segmentation for tracking \n");
1183 const AliMUONGeometryTransformer* kGeomTransformer
1184 = pMUON->GetGeometryTransformer();
1186 //loop over all digits and store their position
1187 for (Int_t digit = 0; digit < ndigits; digit++) {
1188 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1189 if (mdig->Cathode() != cathode-1) continue;
1191 // get all needed parameters
1192 Int_t ix=mdig->PadX();
1193 Int_t iy=mdig->PadY();
1194 Int_t detElemId=mdig->DetElemId();
1195 Float_t charge = mdig->Signal();
1196 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1197 Int_t color = 261+index;
1198 Int_t colorTrigger = 2;
1199 if (color > 282) color = 282;
1201 const AliMpVSegmentation* seg =
1202 AliMpSegmentation::Instance()
1203 ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode-1));
1205 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1207 if (chamber > 10) { // trigger chamber
1208 Float_t sumCharge = 0;
1209 Int_t n = mdig->Ntracks();
1210 for (Int_t icharge = 0; icharge < n; icharge++) {
1211 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1213 Float_t testCharge = sumCharge-(Int_t(sumCharge/n))*n;
1214 if(sumCharge <= n || testCharge > 0) {
1215 colorTrigger = color;
1221 // get the pad position and dimensions
1222 Float_t xlocal1 = pad.Position().X();
1223 Float_t ylocal1 = pad.Position().Y();
1224 Float_t xlocal2 = pad.Dimensions().X();
1225 Float_t ylocal2 = pad.Dimensions().Y();
1227 Float_t xg1, xg2, yg1, yg2, zg1;
1229 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
1230 // (no transformation for pad dimensions)
1234 // Then set the objects
1235 points = new AliMUONPoints(npoints);
1236 fPoints->AddAt(points,digit);
1238 points->SetMarkerColor(colorTrigger);
1240 points->SetMarkerColor(color);
1242 points->SetMarkerStyle(21);
1243 points->SetMarkerSize(0.5);
1244 points->SetParticle(-1);
1245 points->SetHitIndex(-1);
1246 points->SetTrackIndex(-1);
1247 points->SetDigitIndex(digit);
1248 points->SetPoint(0,xg1,yg1,zpos);
1250 Int_t lineColor = (zg1-zpos > 0) ? 2:3;
1251 marker=new TMarker3DBox(xg1,yg1,zpos,xg2,yg2,0,0,0);
1253 marker->SetLineColor(lineColor);
1254 marker->SetFillStyle(1001);
1255 marker->SetFillColor(color);
1256 marker->SetRefObject((TObject*)points);
1257 points->Set3DMarker(0, marker);
1259 } // end loop on digits
1260 } // end of new segmentation
1262 //___________________________________________
1263 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1265 /// Read raw clusters info and store x,y,z info in arrays fRpoints.
1266 /// Loop on all detectors
1268 if (chamber > 10) return;
1272 GetMUONRecData()->SetTreeAddress("RC");
1273 TClonesArray *muonRawClusters = GetMUONRecData()->RawClusters(chamber-1);
1275 if (muonRawClusters == 0) return;
1278 if (GetMUONRecData()->TreeR()) {
1279 nent=(Int_t) GetMUONRecData()->TreeR()->GetEntries();
1280 GetMUONRecData()->TreeR()->GetEvent(0);
1283 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1284 if (nrawcl == 0) return;
1285 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1287 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1288 AliMUONRawCluster *mRaw;
1289 AliMUONPoints *points = 0;
1291 //loop over all raw clusters and store their position
1292 points = new AliMUONPoints(nrawcl);
1293 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1294 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1295 points->SetMarkerColor(51);
1296 points->SetMarkerStyle(2);
1297 points->SetMarkerSize(1.);
1298 points->SetParticle(-1);
1299 points->SetHitIndex(-1);
1300 points->SetTrackIndex(-1);
1301 points->SetDigitIndex(-1);
1302 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1303 fRpoints->AddAt(points,iraw);
1304 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1308 //___________________________________________
1309 void AliMUONDisplay::LoadTracks()
1313 AliMUONTrack* recTrack = 0;
1314 AliMUONTrackParam* trackParam = 0;
1315 TClonesArray * trackParamAtHit = 0;
1319 GetMUONRecData()->SetTreeAddress("RT");
1320 TClonesArray* recTracksArray = GetMUONRecData()->RecTracks();
1321 if (recTracksArray == NULL) return;
1322 GetMUONRecData()->GetRecTracks();
1324 Int_t nRecTracks = 0;
1326 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1329 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1331 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1332 // reading info from tracks
1333 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1335 Int_t nTrackHits = recTrack->GetNTrackHits();
1337 if (nTrackHits == 0) continue;
1340 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1341 points->SetLineColor(6);
1342 points->SetLineWidth(1);
1343 fRpoints->AddAt(points,iRecTracks);
1349 // vertex unknown at the tracking level -> put it at (0,0,0)
1350 points->SetPoint(iPoint,0.,0.,0.);
1353 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1354 trackParamAtHit = recTrack->GetTrackParamAtHit();
1355 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1356 xRec = trackParam->GetNonBendingCoor();
1357 yRec = trackParam->GetBendingCoor();
1358 zRec = trackParam->GetZ();
1359 points->SetPoint(iPoint,xRec,yRec,zRec);
1361 } // end loop rec. hits
1362 PrintTrack(iRecTracks,recTrack);
1363 } // end loop tracks
1368 //___________________________________________
1369 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1371 /// Print reconstructed track
1373 AliMUONTrackParam *trackParam;
1374 Float_t vertex[3], momentum[3];
1375 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1378 trackParam = recTrack->GetTrackParamAtVertex(); // meaningless since the vertex is not known at the tracking level
1379 vertex[0] = trackParam->GetNonBendingCoor();
1380 vertex[1] = trackParam->GetBendingCoor();
1381 vertex[2] = trackParam->GetZ();
1382 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
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]);
1392 if ( trackParam->GetInverseBendingMomentum() != 0. ) {
1393 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1394 bendingSlope = trackParam->GetBendingSlope();
1395 nonBendingSlope = trackParam->GetNonBendingSlope();
1396 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1397 momentum[0] = momentum[2] * nonBendingSlope;
1398 momentum[1] = momentum[2] * bendingSlope;
1399 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1402 AliErrorStream() << "Cannot calculate momentum: pYZ = inf" << endl;
1405 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1406 printf(" track chi2/dof: %f \n",chi2dof);
1409 //___________________________________________
1410 void AliMUONDisplay::LoadHits(Int_t chamber)
1412 /// Read hits info and store x,y,z info in arrays fPhits.
1413 /// Loop on all detectors
1415 if (chamber > 14) return;
1422 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1424 if (GetMUONSimData()->TreeH()) {
1425 GetMUONSimData()->SetTreeAddress("H");
1426 Int_t ntracks = (Int_t)GetMUONSimData()->TreeH()->GetEntries(); //skowron
1428 for (track = 0; track < ntracks; track++) {
1429 GetMUONSimData()->ResetHits();
1430 GetMUONSimData()->GetTrack(track);//skowron
1431 TClonesArray *muonHits = GetMUONSimData()->Hits();
1432 if (muonHits == 0) return;
1433 nthits += muonHits->GetEntriesFast();
1435 if (fPhits == 0) fPhits = new TObjArray(nthits);
1437 for (track=0; track<ntracks;track++) {
1438 GetMUONSimData()->ResetHits();
1439 GetMUONSimData()->GetTrack(track);//skowron
1440 TClonesArray *muonHits = GetMUONSimData()->Hits();
1441 if (muonHits == 0) return;
1442 Int_t nhits = muonHits->GetEntriesFast();
1443 if (nhits == 0) continue;
1445 AliMUONPoints *points = 0;
1447 for (Int_t hit=0;hit<nhits;hit++) {
1448 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1449 Int_t nch = mHit->Chamber(); // chamber number
1450 if (nch != chamber) continue;
1452 // Retrieve info and set the objects
1454 points = new AliMUONPoints(npoints);
1455 fPhits->AddAt(points,nhold+hit);
1456 points->SetMarkerColor(kRed);
1457 points->SetMarkerStyle(5);
1458 points->SetMarkerSize(1.);
1459 points->SetParticle(mHit->Track());
1460 points->SetHitIndex(hit);
1461 points->SetTrackIndex(track);
1462 points->SetDigitIndex(-1);
1463 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1464 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1471 //_____________________________________________________________________________
1472 void AliMUONDisplay::Paint(Option_t *)
1474 /// Paint miscellaneous items
1477 //_____________________________________________________________________________
1478 void AliMUONDisplay::SetPickMode()
1480 /// Set parameters for pick mode.
1484 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1485 fTrigPad->Modified();
1488 //_____________________________________________________________________________
1489 void AliMUONDisplay::SetZoomMode()
1491 /// Set parameters for zoom mode
1495 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1496 fTrigPad->Modified();
1499 //_____________________________________________________________________________
1500 void AliMUONDisplay::NextChamber(Int_t delta)
1502 /// To go from chamber to next chamber if delta = 1
1503 /// or previous chamber otherwise
1505 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1507 if (fChamber > 1) fChamber--;
1511 LoadDigits(fChamber, fCathode);
1515 //_____________________________________________________________________________
1516 void AliMUONDisplay::NextCathode()
1518 /// To switch to other cathode plane
1522 if (fCathode == 1) {
1523 LoadDigits(fChamber, 2);
1525 LoadDigits(fChamber, 1);
1527 fNextCathode = kTRUE; // to keep the same zoom
1529 fNextCathode = kFALSE;
1530 TPad *pad = (TPad*)gPad->GetPadSave();
1531 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1532 fZoomX1[fZooms], fZoomY1[fZooms]);
1538 //_____________________________________________________________________________
1539 void AliMUONDisplay::Trigger()
1541 /// Print global trigger output
1543 AliMUONGlobalTrigger* globalTrig;
1545 GetMUONRecData()->SetTreeAddress("GLT");
1546 GetMUONRecData()->GetTriggerD();
1548 globalTrig = (AliMUONGlobalTrigger*)GetMUONRecData()->GlobalTrigger()->UncheckedAt(0);
1549 if (globalTrig == 0) return;
1551 globalTrig->Print("full");
1553 // // returns Trigger Decision for current event
1554 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1556 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1557 // AliMUONData* muonData = decision->GetMUONRecData();
1558 // muonData->SetTreeAddress("D");
1559 // decision->Trigger();
1561 //_____________________________________________________________________________
1562 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1564 /// Set chamber and cathode number
1571 LoadDigits(chamber,cathode);
1575 //_____________________________________________________________________________
1576 void AliMUONDisplay::SetEvent(Int_t newevent)
1580 gAlice->GetEvent(newevent);
1582 if (!gAlice->TreeD()) return;
1585 LoadDigits(fChamber,fCathode);
1589 //_____________________________________________________________________________
1590 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1592 /// Set view range along R and Z
1602 //_____________________________________________________________________________
1603 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1605 /// Change viewing angles for current event
1613 TView *view = gPad->GetView();
1614 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1619 //_____________________________________________________________________________
1620 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1622 /// Display (current event_number + delta)
1623 /// - delta = 1 shown next event
1624 /// - delta = -1 show previous event
1626 if ( !fSimLoader || !fRecLoader ) {
1627 AliErrorStream() << "Detector loaders are not defined." << endl;
1631 AliRunLoader* runSimLoader = fSimLoader->GetRunLoader();
1632 AliRunLoader* runRecLoader = fRecLoader->GetRunLoader();
1636 //runLoader->CleanDetectors();
1637 //runLoader->CleanKinematics();
1638 Int_t currentEvent = runSimLoader->GetEventNumber();
1639 Int_t newEvent = currentEvent + delta;
1640 runSimLoader->GetEvent(newEvent);
1641 runRecLoader->GetEvent(newEvent);
1644 LoadDigits(fChamber, fCathode);
1649 //______________________________________________________________________________
1650 void AliMUONDisplay::UnZoom()
1654 if (fZooms <= 0) return;
1656 TPad *pad = (TPad*)gPad->GetPadSave();
1657 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1661 //_____________________________________________________________________________
1662 void AliMUONDisplay::ResetPoints()
1664 /// Reset array of points
1672 //_____________________________________________________________________________
1673 void AliMUONDisplay::ResetPhits()
1675 /// Reset array of points
1683 //_____________________________________________________________________________
1684 void AliMUONDisplay::ResetRpoints()
1686 /// Reset array of points