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 "AliMUONTriggerSegmentationV2.h"
83 #include "AliMpDEIterator.h"
84 #include "AliMpSegFactory.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 "AliMpDEManager.h"
97 #include "AliHeader.h"
104 #include <TPaveLabel.h>
105 #include <TPaveText.h>
106 #include <TDiamond.h>
110 #include <TVirtualX.h>
112 #include <TGeometry.h>
113 #include <TMarker3DBox.h>
114 #include <TParticle.h>
115 #include <TPolyLine3D.h>
119 ClassImp(AliMUONDisplay)
122 //_____________________________________________________________________________
123 AliMUONDisplay::AliMUONDisplay()
128 fDrawClusters(kTRUE),
140 /// Default constructor
143 //_____________________________________________________________________________
144 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
149 fDrawClusters(kTRUE),
157 fNextCathode(kFALSE),
162 /// Standard constructor to create an event display object.
166 gAlice->SetDisplay(this);
168 // Initialize display default parameters
171 // Set front view by default
180 // Create display canvas
182 if (ysize < 100) ysize = 750;
183 Int_t xsize = Int_t(size*830./ysize);
184 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
185 fCanvas->ToggleEventStatus();
187 // Create main display pad
188 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
191 fPad->SetFillColor(30);
192 fPad->SetBorderSize(2);
197 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
199 fColPad->SetFillColor(17);
200 fColPad->SetBorderSize(2);
205 // Create user interface control pad
209 // Create Range and mode pad
212 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
215 fTrigPad->SetFillColor(22);
216 fTrigPad->SetBorderSize(2);
217 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
218 fRangeSlider->SetObject(this);
219 char pickmode[] = "gAlice->Display()->SetPickMode()";
221 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
222 fPickButton->SetFillColor(38);
224 char zoommode[] = "gAlice->Display()->SetZoomMode()";
225 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
226 fZoomButton->SetFillColor(38);
228 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
229 fArcButton->SetFillColor(kGreen);
231 char butUnzoom[] = "gAlice->Display()->UnZoom()";
232 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
233 button->SetFillColor(38);
235 AppendPad(); // append display object as last object to force selection
238 fTrigPad->SetEditable(kFALSE);
239 fButtons->SetEditable(kFALSE);
242 // initialize container
244 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
249 //_____________________________________________________________________________
250 AliMUONDisplay::~AliMUONDisplay()
254 // Delete space point structure
255 if (fPoints) fPoints->Delete();
259 if (fPhits) fPhits->Delete();
263 if (fRpoints) fRpoints->Delete();
268 //_____________________________________________________________________________
269 void AliMUONDisplay::Clear(Option_t *)
271 /// Delete graphics temporary objects
274 //_____________________________________________________________________________
275 void AliMUONDisplay::DisplayButtons()
277 /// Create the user interface buttons
280 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
282 fButtons->SetFillColor(38);
283 fButtons->SetBorderSize(2);
286 // Int_t butcolor = 33;
287 Float_t dbutton = 0.08;
294 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
295 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
296 button->SetFillColor(38);
300 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
301 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
302 button->SetFillColor(38);
306 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
307 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
308 button->SetFillColor(38);
312 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
313 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
314 button->SetFillColor(38);
318 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
319 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
320 button->SetFillColor(38);
324 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
325 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
326 button->SetFillColor(38);
330 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
331 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
332 button->SetFillColor(38);
336 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
337 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
338 button->SetFillColor(38);
342 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
343 diamond->SetFillColor(50);
344 diamond->SetTextAlign(22);
345 diamond->SetTextColor(5);
346 diamond->SetTextSize(0.11);
348 diamond->AddText(".. ");
349 diamond->AddText("ROOT");
350 diamond->AddText("MUON");
351 diamond->AddText("... ");
352 diamond->AddText(" ");
355 //_____________________________________________________________________________
356 void AliMUONDisplay::CreateColors() const
358 /// Create the colors palette used to display clusters
373 new TColor(color,r,g,b);
383 new TColor(color,r,g,b);
393 new TColor(color,r,g,b);
403 new TColor(color,r,g,b);
413 new TColor(color,r,g,b);
420 //_____________________________________________________________________________
421 void AliMUONDisplay::DisplayColorScale()
423 /// Display pulse height color scale
427 Float_t xlow, ylow, xup, yup, hs;
428 Float_t x1, y1, x2, y2;
432 TText *text = new TText(0,0,"");
433 text->SetTextFont(61);
434 text->SetTextSize(0.2);
435 text->SetTextAlign(22);
438 Int_t adcmax=4096; // default 12 bits ADC
444 //*-* draw colortable boxes
445 hs = (y2-y1)/Float_t(22);
449 ylow = y1 + hs*(Float_t(i));
450 yup = y1 + hs*(Float_t(i+1));
452 Double_t logscale=Double_t(i+1)*(TMath::Log(adcmax)/22);
453 Int_t scale=(Int_t)TMath::Exp(logscale);
454 sprintf(label,"%d",scale);
455 box = new TBox(xlow, ylow, xup, yup);
457 box->SetFillColor(color);
458 text->DrawText(xlow+0.7, 0.5*(ylow+yup),label);
462 //______________________________________________________________________________
463 Int_t AliMUONDisplay::DistancetoPrimitive(Int_t px, Int_t)
465 /// Compute distance from point px,py to objects in event
467 gPad->SetCursor(kCross);
469 if (gPad == fTrigPad) return 9999;
471 const Int_t kBig = 9999;
473 Float_t xmin = gPad->GetX1();
474 Float_t xmax = gPad->GetX2();
475 Float_t dx = 0.02*(xmax - xmin);
476 Float_t x = gPad->AbsPixeltoX(px);
477 if (x < xmin+dx || x > xmax-dx) return dist;
479 if (fZoomMode) return 0;
483 //_____________________________________________________________________________
484 void AliMUONDisplay::Draw(Option_t *)
486 /// Display current event
494 //_____________________________________________________________________________
495 void AliMUONDisplay::DrawChamber()
497 /// Display current event
499 fDrawTracks = kFALSE;
501 DrawView(fTheta, fPhi, fPsi);
502 // Display the event number and title
507 //_____________________________________________________________________________
508 void AliMUONDisplay::DrawReco(Option_t *)
510 /// Display current event
513 // print kinematics of generated particles
515 // Draw global view of muon system
517 DrawGlobalView(135, -50, -140);
519 // Display the event number and title
524 //_____________________________________________________________________________
525 void AliMUONDisplay::PrintKinematics()
527 /// Print kinematic tree
529 AliRunLoader * runLoader;
530 TParticle *particle = new TParticle();
532 Float_t vertex[3], momentum[3];
535 runLoader = fLoader->GetRunLoader();
539 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
540 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
541 nPart = (Int_t)runLoader->TreeK()->GetEntries();
542 for(Int_t iPart = 0; iPart < nPart; iPart++) {
543 runLoader->TreeK()->GetEvent(iPart);
544 vertex[0] = particle->Vx();
545 vertex[1] = particle->Vy();
546 vertex[2] = particle->Vz();
547 momentum[0] = particle->Px();
548 momentum[1] = particle->Py();
549 momentum[2] = particle->Pz();
551 printf("===================================================\n");
552 printf(" Generated particle # %d \n",iPart);
553 printf(" name: %s \n",particle->GetName());
554 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
555 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
560 //_____________________________________________________________________________
561 void AliMUONDisplay::DrawSegmentation()
563 /// \todo to be re-written for new seg
564 /// Draw graphical representation of segmenatation
565 /// Attention: still experimental code
568 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
569 // AliMUONChamber* iChamber;
571 // AliSegmentation* seg;
572 // iChamber = &(pMUON->Chamber(fChamber));
573 // seg=iChamber->SegmentationModel(icat);
575 // Float_t zpos=iChamber->Z();
576 // Float_t r=iChamber->ROuter();
578 // TMarker3DBox *marker;
580 // for (Int_t j=0; j<seg->Npy(); j++) {
582 // y0=j*seg->Dpy()-seg->Dpy()/2.;
583 // for (seg->FirstPad(0.,y0,0,300,0.);
587 // if (seg->ISector()==0) continue;
589 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
590 // Float_t dpx=seg->Dpx(seg->ISector())/2;
591 // Float_t dpy=seg->Dpy(seg->ISector())/2;
592 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
593 // marker->SetLineColor(seg->ISector()+1);
594 // marker->SetFillStyle(1001);
595 // marker->SetFillColor(0);
600 // for (Int_t j=0; j<250; j++) {
601 // Float_t x0=j*seg->Dpx();
602 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
604 // for (seg->FirstPad(x0,0,0,0,y0);
608 // if (seg->ISector()==0) continue;
611 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
612 // Float_t dpx=seg->Dpx(seg->ISector())/2;
613 // Float_t dpy=seg->Dpy(seg->ISector())/2;
614 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
615 // marker->SetLineColor(seg->ISector()+1);
616 // marker->SetFillStyle(1001);
617 // marker->SetFillColor(0);
624 //_____________________________________________________________________________
625 void AliMUONDisplay::DrawClusters()
627 /// Draw clusters for MUON chambers
629 Int_t ndigits, digit;
637 ndigits = points->GetEntriesFast();
638 for (digit=0;digit<ndigits;digit++){
639 pm = (AliMUONPoints*)points->UncheckedAt(digit);
643 for (Int_t im=0;im<3;im++) {
644 TMarker3DBox *marker=pm->GetMarker(im);
649 fClustersCuts += pm->GetN();
653 //_____________________________________________________________________________
654 void AliMUONDisplay::DrawHits()
656 /// Draw hits for MUON chambers
660 Int_t ntracks, track;
667 ntracks = points->GetEntriesFast();
668 for (track=0;track<ntracks;track++) {
669 pm = (AliMUONPoints*)points->UncheckedAt(track);
672 fHitsCuts += pm->GetN();
677 //_____________________________________________________________________________
678 void AliMUONDisplay::DrawCoG()
680 /// Draw hits for MUON chambers
682 if (!fDrawCoG) return;
683 if (fChamber > 10) return;
684 LoadCoG(fChamber,fCathode);
692 ncog = points->GetEntriesFast();
693 for (icog=0;icog<ncog;icog++) {
694 pm = (AliMUONPoints*)points->UncheckedAt(icog);
699 //_____________________________________________________________________________
700 void AliMUONDisplay::DrawTracks()
704 if (!fDrawTracks) return;
707 Int_t nTrack, iTrack;
713 nTrack = points->GetEntriesFast();
714 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
715 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
720 //_____________________________________________________________________________
722 void AliMUONDisplay::DrawTitle(Option_t *option)
724 /// Draw the event title
726 Float_t xmin = gPad->GetX1();
727 Float_t xmax = gPad->GetX2();
728 Float_t ymin = gPad->GetY1();
729 Float_t ymax = gPad->GetY2();
730 Float_t dx = xmax-xmin;
731 Float_t dy = ymax-ymin;
733 AliRunLoader * runLoader;
735 runLoader = fLoader->GetRunLoader();
740 if (strlen(option) == 0) {
741 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
742 // title->SetTextSize(0.023932);
743 title->SetTextSize(0.02);
744 title->SetBit(kCanDelete);
745 title->SetFillColor(42);
748 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
749 runLoader->GetEventNumber(),
750 gAlice->GetHeader()->GetRun(),
753 title->AddText(ptitle);
754 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
755 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
756 nparticles, fHitsCuts,fClustersCuts);
757 title->AddText(ptitle);
759 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
760 label->SetBit(kCanDelete);
761 label->SetFillColor(42);
766 //_____________________________________________________________________________
767 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
769 /// Draw a view of MUON clusters
771 AliInfo(" Draw View");
773 gPad->SetCursor(kWatch);
774 // gPad->SetFillColor(39);
775 gPad->SetFillColor(1);
777 // gPad->SetFillColor(39);
778 gPad->SetFillColor(1);
781 TView *view = new TView(1);
783 Float_t range = fRrange*fRangeSlider->GetMaximum();
784 view->SetRange(-range,-range,-range,range, range, range);
785 // zoom back to full scale only if DrawView not called from NextCathode
794 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
796 // Recovering the chamber
797 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
799 const AliMUONGeometryTransformer* kGeomTransformer
800 = pMUON->GetGeometryTransformer();
802 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
803 AliMpSegFactory segFactory;
804 // Mapping segmentation factory will be used only to create mapping
805 // segmentations if not present in the DE segmentations
807 // Display MUON Chamber Geometry
809 sprintf(nodeName,"MUON%d",100+fChamber);
810 printf(">>>> chamber is %d\n",fChamber);
814 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() ) {
816 Int_t detElemId = it.CurrentDE();
817 AliMpSectorSegmentation * seg =
818 (AliMpSectorSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
819 const AliMpSector * sector = seg->GetSector();
821 // get sector measurements
822 TVector2 position = sector->Position();
823 TVector2 dimension = sector->Dimensions(); // half length
825 Float_t xlocal1 = position.Px(); // FIXME: not really needed as it's 0 ?
826 Float_t ylocal1 = position.Py(); // FIXME: not really needed as it's 0 ?
827 Float_t xlocal2 = dimension.Px() * 2.;
828 Float_t ylocal2 = dimension.Px() * 2.;
830 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
831 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
834 TPolyLine3D* poly = new TPolyLine3D();
837 poly->SetPoint(nPoint++, xg1, yg1, 0.);
838 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
839 Float_t x = xg1 + xg2 * TMath::Cos(d);
840 Float_t y = yg1 + yg2 * TMath::Sin(d);
841 poly->SetPoint(nPoint++, x, y, 0.);
843 poly->SetPoint(nPoint++, xg1, yg1, 0.);
845 poly->SetLineColor(2);
854 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() )
856 Int_t detElemId = it.CurrentDE();
857 AliMpStationType stationType = AliMpDEManager::GetStationType(detElemId);
859 if ( segmentation->HasDE(detElemId) )
861 const AliMpVSegmentation* seg = segmentation->GetMpSegmentation(detElemId, 0);
863 // Create mapping segmentation if old trigger segmentation
864 // (not using mapping)
865 seg = segFactory.CreateMpSegmentation(detElemId, 0);
869 Float_t deltax = seg->Dimensions().X();
870 Float_t deltay = seg->Dimensions().Y();
871 Float_t xlocal1 = -deltax;
872 Float_t ylocal1 = -deltay;
873 Float_t xlocal2 = +deltax;
874 Float_t ylocal2 = +deltay;
875 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
876 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
878 // drawing slat active volumes
879 Float_t xCenter = (xg1 + xg2)/2.;
880 Float_t yCenter = (yg1 + yg2)/2.;
882 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
884 box->SetFillStyle(0);
885 box->SetLineColor( stationType == kStationTrigger ? 4 : 2);
888 if ( stationType == kStation345 )
890 // drawing inner circle + disc
891 TPolyLine3D* poly = new TPolyLine3D();
892 TPolyLine3D* poly1 = new TPolyLine3D();
896 for (Float_t d = 0; d < 6.24; d+= 0.005)
898 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
899 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
900 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
901 poly->SetPoint(nPoint++, x, y, 0.);
902 poly1->SetPoint(nPoint1++, x, y, 0.);
904 poly->SetLineColor(1);
906 poly1->SetLineColor(2);
914 // // Delete mapping segmentations created with factory
915 // segFactory.DeleteSegmentations();
917 //add clusters to the pad
921 // DrawSegmentation();
922 // add itself to the list (must be last)
924 view->SetView(phi, theta, psi, iret);
927 //_____________________________________________________________________________
928 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
930 /// Draw a view of muons chambers with tracks
932 gPad->SetCursor(kWatch);
933 // gPad->SetFillColor(39);
934 gPad->SetFillColor(1);
936 // gPad->SetFillColor(39);
937 gPad->SetFillColor(1);
941 TView *view = new TView(1);
943 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
944 view->SetRange(-range,-range,-range,range,range,range);
946 // Display all MUON Chambers segmentation
949 sprintf(nodeName,"alice");
951 node1=gAlice->GetGeometry()->GetNode(nodeName);
952 if (node1) node1->Draw("same");
955 // Draw clusters for all chambers
956 Int_t chamberSave = fChamber;
957 for (fChamber = 1; fChamber <= 10; fChamber++){
960 fChamber = chamberSave;
961 // Draw reconstructed tracks
968 Float_t x0 = (-1+shift)/zoom;
969 Float_t y0 = (-1+shift)/zoom;
970 Float_t x1 = (1+shift)/zoom;
971 Float_t y1 = (1+shift)/zoom;
972 gPad->Range(x0,y0,x1,y1);
973 view->SetView(phi, theta, psi, iret);
977 //______________________________________________________________________________
978 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
980 /// Execute action corresponding to the mouse event
982 static Float_t x0, y0, x1, y1;
984 static Int_t pxold, pyold;
985 static Int_t px0, py0;
986 static Int_t linedrawn;
989 if (px == 0 && py == 0) { //when called by sliders
990 if (event == kButton1Up) {
995 if (!fZoomMode && gPad->GetView()) {
996 gPad->GetView()->ExecuteRotateView(event, px, py);
1000 // something to zoom ?
1001 gPad->SetCursor(kCross);
1006 gVirtualX->SetLineColor(-1);
1007 gPad->TAttLine::Modify(); //Change line attributes only if necessary
1008 x0 = gPad->AbsPixeltoX(px);
1009 y0 = gPad->AbsPixeltoY(py);
1011 pxold = px; pyold = py;
1015 case kButton1Motion:
1016 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1020 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1024 gPad->GetCanvas()->FeedbackMode(kFALSE);
1025 if (px == px0) return;
1026 if (py == py0) return;
1027 x1 = gPad->AbsPixeltoX(px);
1028 y1 = gPad->AbsPixeltoY(py);
1030 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1031 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1032 gPad->Range(x0,y0,x1,y1);
1033 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1035 fZoomX0[fZooms] = x0;
1036 fZoomY0[fZooms] = y0;
1037 fZoomX1[fZooms] = x1;
1038 fZoomY1[fZooms] = y1;
1040 gPad->Modified(kTRUE);
1045 //___________________________________________
1046 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1048 /// Read digits info and store x,y,z info in arrays fPoints.
1049 /// Loop on all detectors
1051 if (chamber > 14) return;
1057 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1059 AliMUONGeometrySegmentation* segmentation2 = 0x0;
1061 GetMUONData()->SetTreeAddress("D");
1063 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
1064 if (muonDigits == 0) return;
1066 gAlice->ResetDigits();
1069 if (GetLoader()->TreeD()) {
1070 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
1071 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1072 GetMUONData()->GetDigits();
1075 Int_t ndigits = muonDigits->GetEntriesFast();
1076 if (ndigits == 0) return;
1077 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1079 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1082 AliMUONPoints *points = 0;
1083 TMarker3DBox *marker = 0;
1086 Float_t adcmax = 1024; // default
1087 if (chamber<11) adcmax = 4096;
1089 // check if trigger is using new or old segmentation
1091 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
1092 if ( segmentation->GetMpSegmentation(1100, cathode-1) ) old = false;
1094 if ( old && chamber > 10) {
1095 if (chamber > 10) printf(">>> old segmentation for trigger \n");
1096 else printf(">>> old segmentation for tracking \n");
1099 = pMUON->GetSegmentation()->GetModuleSegmentation(chamber-1, cathode-1);
1100 for (Int_t digit = 0; digit < ndigits; digit++) {
1101 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1102 if (mdig->Cathode() != cathode-1) continue;
1105 // First get all needed parameters
1107 Int_t charge = mdig->Signal();
1108 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1109 Int_t color = 261+index;
1110 Int_t colorTrigger = 2;
1111 if (color > 282) color = 282;
1113 if (chamber > 10) { // trigger chamber
1115 Int_t sumCharge = 0;
1116 for (Int_t icharge = 0; icharge < 10; icharge++) {
1117 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1119 Int_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1120 if(sumCharge <= 10 || testCharge > 0) {
1121 colorTrigger = color;
1127 // get the center of the pad - add on x and y half of pad size
1128 Float_t xpad, ypad, zpad;
1132 Int_t detElemId = mdig->DetElemId();
1133 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1134 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1135 dpx = segmentation2->Dpx(detElemId, isec)/2;
1136 dpy = segmentation2->Dpy(detElemId, isec)/2;
1138 // Then set the objects
1139 points = new AliMUONPoints(npoints);
1140 fPoints->AddAt(points,digit);
1142 points->SetMarkerColor(colorTrigger);
1144 points->SetMarkerColor(color);
1146 points->SetMarkerStyle(21);
1147 points->SetMarkerSize(0.5);
1148 points->SetParticle(-1);
1149 points->SetHitIndex(-1);
1150 points->SetTrackIndex(-1);
1151 points->SetDigitIndex(digit);
1152 points->SetPoint(0,xpad,ypad,zpos);
1154 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1155 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1158 marker->SetLineColor(lineColor);
1159 marker->SetFillStyle(1001);
1160 marker->SetFillColor(color);
1161 marker->SetRefObject((TObject*)points);
1162 points->Set3DMarker(0, marker);
1163 } // end loop on digits
1166 if (chamber > 10) printf(">>> new segmentation for trigger \n");
1167 else printf(">>> new segmentation for tracking \n");
1169 const AliMUONGeometryTransformer* kGeomTransformer
1170 = pMUON->GetGeometryTransformer();
1172 //loop over all digits and store their position
1173 for (Int_t digit = 0; digit < ndigits; digit++) {
1174 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1175 if (mdig->Cathode() != cathode-1) continue;
1177 // get all needed parameters
1178 Int_t ix=mdig->PadX();
1179 Int_t iy=mdig->PadY();
1180 Int_t detElemId=mdig->DetElemId();
1181 Int_t charge = mdig->Signal();
1182 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1183 Int_t color = 261+index;
1184 Int_t colorTrigger = 2;
1185 if (color > 282) color = 282;
1187 const AliMpVSegmentation* seg =
1188 pMUON->GetSegmentation()->GetMpSegmentation(detElemId,cathode-1);
1190 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1192 if (chamber > 10) { // trigger chamber
1193 Int_t sumCharge = 0;
1194 Int_t n = mdig->Ntracks();
1195 for (Int_t icharge = 0; icharge < n; icharge++) {
1196 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1198 assert(sumCharge==mdig->Signal());
1199 Int_t testCharge = sumCharge-(Int_t(sumCharge/n))*n;
1200 if(sumCharge <= n || testCharge > 0) {
1201 colorTrigger = color;
1207 // get the pad position and dimensions
1208 Float_t xlocal1 = pad.Position().X();
1209 Float_t ylocal1 = pad.Position().Y();
1210 Float_t xlocal2 = pad.Dimensions().X();
1211 Float_t ylocal2 = pad.Dimensions().Y();
1213 Float_t xg1, xg2, yg1, yg2, zg1;
1215 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
1216 // (no transformation for pad dimensions)
1220 // Then set the objects
1221 points = new AliMUONPoints(npoints);
1222 fPoints->AddAt(points,digit);
1224 points->SetMarkerColor(colorTrigger);
1226 points->SetMarkerColor(color);
1228 points->SetMarkerStyle(21);
1229 points->SetMarkerSize(0.5);
1230 points->SetParticle(-1);
1231 points->SetHitIndex(-1);
1232 points->SetTrackIndex(-1);
1233 points->SetDigitIndex(digit);
1234 points->SetPoint(0,xg1,yg1,zpos);
1236 Int_t lineColor = (zg1-zpos > 0) ? 2:3;
1237 marker=new TMarker3DBox(xg1,yg1,zpos,xg2,yg2,0,0,0);
1239 marker->SetLineColor(lineColor);
1240 marker->SetFillStyle(1001);
1241 marker->SetFillColor(color);
1242 marker->SetRefObject((TObject*)points);
1243 points->Set3DMarker(0, marker);
1245 } // end loop on digits
1246 } // end of new segmentation
1248 //___________________________________________
1249 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1251 /// Read raw clusters info and store x,y,z info in arrays fRpoints.
1252 /// Loop on all detectors
1254 if (chamber > 10) return;
1258 GetMUONData()->SetTreeAddress("RC");
1259 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1261 if (muonRawClusters == 0) return;
1264 if (GetMUONData()->TreeR()) {
1265 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1266 GetMUONData()->TreeR()->GetEvent(0);
1269 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1270 if (nrawcl == 0) return;
1271 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1273 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1274 AliMUONRawCluster *mRaw;
1275 AliMUONPoints *points = 0;
1277 //loop over all raw clusters and store their position
1278 points = new AliMUONPoints(nrawcl);
1279 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1280 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1281 points->SetMarkerColor(51);
1282 points->SetMarkerStyle(2);
1283 points->SetMarkerSize(1.);
1284 points->SetParticle(-1);
1285 points->SetHitIndex(-1);
1286 points->SetTrackIndex(-1);
1287 points->SetDigitIndex(-1);
1288 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1289 fRpoints->AddAt(points,iraw);
1290 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1294 //___________________________________________
1295 void AliMUONDisplay::LoadTracks()
1299 AliMUONTrack* recTrack = 0;
1300 AliMUONTrackParam* trackParam = 0;
1301 TClonesArray * trackParamAtHit = 0;
1305 GetMUONData()->SetTreeAddress("RT");
1306 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1307 if (recTracksArray == NULL) return;
1308 GetMUONData()->GetRecTracks();
1310 Int_t nRecTracks = 0;
1312 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1315 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1317 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1318 // reading info from tracks
1319 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1321 Int_t nTrackHits = recTrack->GetNTrackHits();
1323 if (nTrackHits == 0) continue;
1326 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1327 points->SetLineColor(6);
1328 points->SetLineWidth(1);
1329 fRpoints->AddAt(points,iRecTracks);
1335 trackParam = recTrack->GetTrackParamAtVertex();
1336 xRec = trackParam->GetNonBendingCoor();
1337 yRec = trackParam->GetBendingCoor();
1338 zRec = trackParam->GetZ();
1339 points->SetPoint(iPoint,xRec,yRec,zRec);
1342 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1343 trackParamAtHit = recTrack->GetTrackParamAtHit();
1344 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1345 xRec = trackParam->GetNonBendingCoor();
1346 yRec = trackParam->GetBendingCoor();
1347 zRec = trackParam->GetZ();
1348 points->SetPoint(iPoint,xRec,yRec,zRec);
1350 } // end loop rec. hits
1351 PrintTrack(iRecTracks,recTrack);
1352 } // end loop tracks
1357 //___________________________________________
1358 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1360 /// Print reconstructed track
1362 AliMUONTrackParam *trackParam;
1363 Float_t vertex[3], momentum[3];
1364 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1367 trackParam = recTrack->GetTrackParamAtVertex();
1368 vertex[0] = trackParam->GetNonBendingCoor();
1369 vertex[1] = trackParam->GetBendingCoor();
1370 vertex[2] = trackParam->GetZ();
1371 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1372 bendingSlope = trackParam->GetBendingSlope();
1373 nonBendingSlope = trackParam->GetNonBendingSlope();
1374 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1375 momentum[0] = momentum[2] * nonBendingSlope;
1376 momentum[1] = momentum[2] * bendingSlope;
1377 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1378 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1380 printf("===================================================\n");
1381 printf(" Reconstructed track # %d \n",iRecTracks);
1382 printf(" charge: %d \n",charge);
1383 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1384 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1385 printf(" track chi2/dof: %f \n",chi2dof);
1389 //___________________________________________
1390 void AliMUONDisplay::LoadHits(Int_t chamber)
1392 /// Read hits info and store x,y,z info in arrays fPhits.
1393 /// Loop on all detectors
1395 if (chamber > 14) return;
1402 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1404 if (GetMUONData()->TreeH()) {
1405 GetMUONData()->SetTreeAddress("H");
1406 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
1408 for (track = 0; track < ntracks; track++) {
1409 GetMUONData()->ResetHits();
1410 GetMUONData()->GetTrack(track);//skowron
1411 TClonesArray *muonHits = GetMUONData()->Hits();
1412 if (muonHits == 0) return;
1413 nthits += muonHits->GetEntriesFast();
1415 if (fPhits == 0) fPhits = new TObjArray(nthits);
1417 for (track=0; track<ntracks;track++) {
1418 GetMUONData()->ResetHits();
1419 GetMUONData()->GetTrack(track);//skowron
1420 TClonesArray *muonHits = GetMUONData()->Hits();
1421 if (muonHits == 0) return;
1422 Int_t nhits = muonHits->GetEntriesFast();
1423 if (nhits == 0) continue;
1425 AliMUONPoints *points = 0;
1427 for (Int_t hit=0;hit<nhits;hit++) {
1428 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1429 Int_t nch = mHit->Chamber(); // chamber number
1430 if (nch != chamber) continue;
1432 // Retrieve info and set the objects
1434 points = new AliMUONPoints(npoints);
1435 fPhits->AddAt(points,nhold+hit);
1436 points->SetMarkerColor(kRed);
1437 points->SetMarkerStyle(5);
1438 points->SetMarkerSize(1.);
1439 points->SetParticle(mHit->Track());
1440 points->SetHitIndex(hit);
1441 points->SetTrackIndex(track);
1442 points->SetDigitIndex(-1);
1443 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1444 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1451 //_____________________________________________________________________________
1452 void AliMUONDisplay::Paint(Option_t *)
1454 /// Paint miscellaneous items
1457 //_____________________________________________________________________________
1458 void AliMUONDisplay::SetPickMode()
1460 /// Set parameters for pick mode.
1464 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1465 fTrigPad->Modified();
1468 //_____________________________________________________________________________
1469 void AliMUONDisplay::SetZoomMode()
1471 /// Set parameters for zoom mode
1475 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1476 fTrigPad->Modified();
1479 //_____________________________________________________________________________
1480 void AliMUONDisplay::NextChamber(Int_t delta)
1482 /// To go from chamber to next chamber if delta = 1
1483 /// or previous chamber otherwise
1485 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1487 if (fChamber > 1) fChamber--;
1491 LoadDigits(fChamber, fCathode);
1495 //_____________________________________________________________________________
1496 void AliMUONDisplay::NextCathode()
1498 /// To switch to other cathode plane
1502 if (fCathode == 1) {
1503 LoadDigits(fChamber, 2);
1505 LoadDigits(fChamber, 1);
1507 fNextCathode = kTRUE; // to keep the same zoom
1509 fNextCathode = kFALSE;
1510 TPad *pad = (TPad*)gPad->GetPadSave();
1511 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1512 fZoomX1[fZooms], fZoomY1[fZooms]);
1518 //_____________________________________________________________________________
1519 void AliMUONDisplay::Trigger()
1521 /// Print global trigger output
1523 AliMUONGlobalTrigger* globalTrig;
1525 GetMUONData()->SetTreeAddress("GLT");
1526 GetMUONData()->GetTriggerD();
1528 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1529 if (globalTrig == 0) return;
1531 globalTrig->Print("full");
1533 // // returns Trigger Decision for current event
1534 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1536 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1537 // AliMUONData* muonData = decision->GetMUONData();
1538 // muonData->SetTreeAddress("D");
1539 // decision->Trigger();
1541 //_____________________________________________________________________________
1542 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1544 /// Set chamber and cathode number
1551 LoadDigits(chamber,cathode);
1555 //_____________________________________________________________________________
1556 void AliMUONDisplay::SetEvent(Int_t newevent)
1560 gAlice->GetEvent(newevent);
1562 if (!gAlice->TreeD()) return;
1565 LoadDigits(fChamber,fCathode);
1569 //_____________________________________________________________________________
1570 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1572 /// Set view range along R and Z
1582 //_____________________________________________________________________________
1583 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1585 /// Change viewing angles for current event
1593 TView *view = gPad->GetView();
1594 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1599 //_____________________________________________________________________________
1600 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1602 /// Display (current event_number + delta)
1603 /// - delta = 1 shown next event
1604 /// - delta = -1 show previous event
1606 AliRunLoader * runLoader;
1608 runLoader = fLoader->GetRunLoader();
1613 //runLoader->CleanDetectors();
1614 //runLoader->CleanKinematics();
1615 Int_t currentEvent = runLoader->GetEventNumber();
1616 Int_t newEvent = currentEvent + delta;
1617 runLoader->GetEvent(newEvent);
1620 LoadDigits(fChamber, fCathode);
1625 //______________________________________________________________________________
1626 void AliMUONDisplay::UnZoom()
1630 if (fZooms <= 0) return;
1632 TPad *pad = (TPad*)gPad->GetPadSave();
1633 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1637 //_____________________________________________________________________________
1638 void AliMUONDisplay::ResetPoints()
1640 /// Reset array of points
1648 //_____________________________________________________________________________
1649 void AliMUONDisplay::ResetPhits()
1651 /// Reset array of points
1659 //_____________________________________________________________________________
1660 void AliMUONDisplay::ResetRpoints()
1662 /// Reset array of points