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"
105 #include <TPaveLabel.h>
106 #include <TPaveText.h>
107 #include <TDiamond.h>
111 #include <TVirtualX.h>
113 #include <TGeometry.h>
114 #include <TMarker3DBox.h>
115 #include <TParticle.h>
116 #include <TPolyLine3D.h>
120 ClassImp(AliMUONDisplay)
123 //_____________________________________________________________________________
124 AliMUONDisplay::AliMUONDisplay()
129 fDrawClusters(kTRUE),
141 /// Default constructor
144 //_____________________________________________________________________________
145 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
150 fDrawClusters(kTRUE),
158 fNextCathode(kFALSE),
163 /// Standard constructor to create an event display object.
167 gAlice->SetDisplay(this);
169 // Initialize display default parameters
172 // Set front view by default
181 // Create display canvas
183 if (ysize < 100) ysize = 750;
184 Int_t xsize = Int_t(size*830./ysize);
185 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
186 fCanvas->ToggleEventStatus();
188 // Create main display pad
189 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
192 fPad->SetFillColor(30);
193 fPad->SetBorderSize(2);
198 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
200 fColPad->SetFillColor(17);
201 fColPad->SetBorderSize(2);
206 // Create user interface control pad
210 // Create Range and mode pad
213 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
216 fTrigPad->SetFillColor(22);
217 fTrigPad->SetBorderSize(2);
218 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
219 fRangeSlider->SetObject(this);
220 char pickmode[] = "gAlice->Display()->SetPickMode()";
222 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
223 fPickButton->SetFillColor(38);
225 char zoommode[] = "gAlice->Display()->SetZoomMode()";
226 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
227 fZoomButton->SetFillColor(38);
229 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
230 fArcButton->SetFillColor(kGreen);
232 char butUnzoom[] = "gAlice->Display()->UnZoom()";
233 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
234 button->SetFillColor(38);
236 AppendPad(); // append display object as last object to force selection
239 fTrigPad->SetEditable(kFALSE);
240 fButtons->SetEditable(kFALSE);
243 // initialize container
245 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
250 //_____________________________________________________________________________
251 AliMUONDisplay::~AliMUONDisplay()
255 // Delete space point structure
256 if (fPoints) fPoints->Delete();
260 if (fPhits) fPhits->Delete();
264 if (fRpoints) fRpoints->Delete();
269 //_____________________________________________________________________________
270 void AliMUONDisplay::Clear(Option_t *)
272 /// Delete graphics temporary objects
275 //_____________________________________________________________________________
276 void AliMUONDisplay::DisplayButtons()
278 /// Create the user interface buttons
281 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
283 fButtons->SetFillColor(38);
284 fButtons->SetBorderSize(2);
287 // Int_t butcolor = 33;
288 Float_t dbutton = 0.08;
295 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
296 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
297 button->SetFillColor(38);
301 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
302 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
303 button->SetFillColor(38);
307 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
308 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
309 button->SetFillColor(38);
313 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
314 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
315 button->SetFillColor(38);
319 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
320 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
321 button->SetFillColor(38);
325 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
326 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
327 button->SetFillColor(38);
331 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
332 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
333 button->SetFillColor(38);
337 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
338 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
339 button->SetFillColor(38);
343 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
344 diamond->SetFillColor(50);
345 diamond->SetTextAlign(22);
346 diamond->SetTextColor(5);
347 diamond->SetTextSize(0.11);
349 diamond->AddText(".. ");
350 diamond->AddText("ROOT");
351 diamond->AddText("MUON");
352 diamond->AddText("... ");
353 diamond->AddText(" ");
356 //_____________________________________________________________________________
357 void AliMUONDisplay::CreateColors() const
359 /// Create the colors palette used to display clusters
374 new TColor(color,r,g,b);
384 new TColor(color,r,g,b);
394 new TColor(color,r,g,b);
404 new TColor(color,r,g,b);
414 new TColor(color,r,g,b);
421 //_____________________________________________________________________________
422 void AliMUONDisplay::DisplayColorScale()
424 /// Display pulse height color scale
428 Float_t xlow, ylow, xup, yup, hs;
429 Float_t x1, y1, x2, y2;
433 TText *text = new TText(0,0,"");
434 text->SetTextFont(61);
435 text->SetTextSize(0.2);
436 text->SetTextAlign(22);
439 Int_t adcmax=4096; // default 12 bits ADC
445 //*-* draw colortable boxes
446 hs = (y2-y1)/Float_t(22);
450 ylow = y1 + hs*(Float_t(i));
451 yup = y1 + hs*(Float_t(i+1));
453 Double_t logscale=Double_t(i+1)*(TMath::Log(adcmax)/22);
454 Int_t scale=(Int_t)TMath::Exp(logscale);
455 sprintf(label,"%d",scale);
456 box = new TBox(xlow, ylow, xup, yup);
458 box->SetFillColor(color);
459 text->DrawText(xlow+0.7, 0.5*(ylow+yup),label);
463 //______________________________________________________________________________
464 Int_t AliMUONDisplay::DistancetoPrimitive(Int_t px, Int_t)
466 /// Compute distance from point px,py to objects in event
468 gPad->SetCursor(kCross);
470 if (gPad == fTrigPad) return 9999;
472 const Int_t kBig = 9999;
474 Float_t xmin = gPad->GetX1();
475 Float_t xmax = gPad->GetX2();
476 Float_t dx = 0.02*(xmax - xmin);
477 Float_t x = gPad->AbsPixeltoX(px);
478 if (x < xmin+dx || x > xmax-dx) return dist;
480 if (fZoomMode) return 0;
484 //_____________________________________________________________________________
485 void AliMUONDisplay::Draw(Option_t *)
487 /// Display current event
495 //_____________________________________________________________________________
496 void AliMUONDisplay::DrawChamber()
498 /// Display current event
500 fDrawTracks = kFALSE;
502 DrawView(fTheta, fPhi, fPsi);
503 // Display the event number and title
508 //_____________________________________________________________________________
509 void AliMUONDisplay::DrawReco(Option_t *)
511 /// Display current event
514 // print kinematics of generated particles
516 // Draw global view of muon system
518 DrawGlobalView(135, -50, -140);
520 // Display the event number and title
525 //_____________________________________________________________________________
526 void AliMUONDisplay::PrintKinematics()
528 /// Print kinematic tree
530 AliRunLoader * runLoader;
531 TParticle *particle = new TParticle();
533 Float_t vertex[3], momentum[3];
536 runLoader = fLoader->GetRunLoader();
540 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
541 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
542 nPart = (Int_t)runLoader->TreeK()->GetEntries();
543 for(Int_t iPart = 0; iPart < nPart; iPart++) {
544 runLoader->TreeK()->GetEvent(iPart);
545 vertex[0] = particle->Vx();
546 vertex[1] = particle->Vy();
547 vertex[2] = particle->Vz();
548 momentum[0] = particle->Px();
549 momentum[1] = particle->Py();
550 momentum[2] = particle->Pz();
552 printf("===================================================\n");
553 printf(" Generated particle # %d \n",iPart);
554 printf(" name: %s \n",particle->GetName());
555 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
556 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
561 //_____________________________________________________________________________
562 void AliMUONDisplay::DrawSegmentation()
564 /// \todo to be re-written for new seg
565 /// Draw graphical representation of segmenatation
566 /// Attention: still experimental code
569 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
570 // AliMUONChamber* iChamber;
572 // AliSegmentation* seg;
573 // iChamber = &(pMUON->Chamber(fChamber));
574 // seg=iChamber->SegmentationModel(icat);
576 // Float_t zpos=iChamber->Z();
577 // Float_t r=iChamber->ROuter();
579 // TMarker3DBox *marker;
581 // for (Int_t j=0; j<seg->Npy(); j++) {
583 // y0=j*seg->Dpy()-seg->Dpy()/2.;
584 // for (seg->FirstPad(0.,y0,0,300,0.);
588 // if (seg->ISector()==0) continue;
590 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
591 // Float_t dpx=seg->Dpx(seg->ISector())/2;
592 // Float_t dpy=seg->Dpy(seg->ISector())/2;
593 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
594 // marker->SetLineColor(seg->ISector()+1);
595 // marker->SetFillStyle(1001);
596 // marker->SetFillColor(0);
601 // for (Int_t j=0; j<250; j++) {
602 // Float_t x0=j*seg->Dpx();
603 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
605 // for (seg->FirstPad(x0,0,0,0,y0);
609 // if (seg->ISector()==0) continue;
612 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
613 // Float_t dpx=seg->Dpx(seg->ISector())/2;
614 // Float_t dpy=seg->Dpy(seg->ISector())/2;
615 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
616 // marker->SetLineColor(seg->ISector()+1);
617 // marker->SetFillStyle(1001);
618 // marker->SetFillColor(0);
625 //_____________________________________________________________________________
626 void AliMUONDisplay::DrawClusters()
628 /// Draw clusters for MUON chambers
630 Int_t ndigits, digit;
638 ndigits = points->GetEntriesFast();
639 for (digit=0;digit<ndigits;digit++){
640 pm = (AliMUONPoints*)points->UncheckedAt(digit);
644 for (Int_t im=0;im<3;im++) {
645 TMarker3DBox *marker=pm->GetMarker(im);
650 fClustersCuts += pm->GetN();
654 //_____________________________________________________________________________
655 void AliMUONDisplay::DrawHits()
657 /// Draw hits for MUON chambers
661 Int_t ntracks, track;
668 ntracks = points->GetEntriesFast();
669 for (track=0;track<ntracks;track++) {
670 pm = (AliMUONPoints*)points->UncheckedAt(track);
673 fHitsCuts += pm->GetN();
678 //_____________________________________________________________________________
679 void AliMUONDisplay::DrawCoG()
681 /// Draw hits for MUON chambers
683 if (!fDrawCoG) return;
684 if (fChamber > 10) return;
685 LoadCoG(fChamber,fCathode);
693 ncog = points->GetEntriesFast();
694 for (icog=0;icog<ncog;icog++) {
695 pm = (AliMUONPoints*)points->UncheckedAt(icog);
700 //_____________________________________________________________________________
701 void AliMUONDisplay::DrawTracks()
705 if (!fDrawTracks) return;
708 Int_t nTrack, iTrack;
714 nTrack = points->GetEntriesFast();
715 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
716 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
721 //_____________________________________________________________________________
723 void AliMUONDisplay::DrawTitle(Option_t *option)
725 /// Draw the event title
727 Float_t xmin = gPad->GetX1();
728 Float_t xmax = gPad->GetX2();
729 Float_t ymin = gPad->GetY1();
730 Float_t ymax = gPad->GetY2();
731 Float_t dx = xmax-xmin;
732 Float_t dy = ymax-ymin;
734 AliRunLoader * runLoader;
736 runLoader = fLoader->GetRunLoader();
741 if (strlen(option) == 0) {
742 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
743 // title->SetTextSize(0.023932);
744 title->SetTextSize(0.02);
745 title->SetBit(kCanDelete);
746 title->SetFillColor(42);
749 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
750 runLoader->GetEventNumber(),
751 gAlice->GetHeader()->GetRun(),
754 title->AddText(ptitle);
755 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
756 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
757 nparticles, fHitsCuts,fClustersCuts);
758 title->AddText(ptitle);
760 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
761 label->SetBit(kCanDelete);
762 label->SetFillColor(42);
767 //_____________________________________________________________________________
768 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
770 /// Draw a view of MUON clusters
772 AliInfo(" Draw View");
774 gPad->SetCursor(kWatch);
775 // gPad->SetFillColor(39);
776 gPad->SetFillColor(1);
778 // gPad->SetFillColor(39);
779 gPad->SetFillColor(1);
782 TView *view = new TView(1);
784 Float_t range = fRrange*fRangeSlider->GetMaximum();
785 view->SetRange(-range,-range,-range,range, range, range);
786 // zoom back to full scale only if DrawView not called from NextCathode
795 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
797 // Recovering the chamber
798 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
800 const AliMUONGeometryTransformer* kGeomTransformer
801 = pMUON->GetGeometryTransformer();
803 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
805 // Display MUON Chamber Geometry
807 sprintf(nodeName,"MUON%d",100+fChamber);
808 printf(">>>> chamber is %d\n",fChamber);
812 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() ) {
814 Int_t detElemId = it.CurrentDEId();
815 AliMpSectorSegmentation * seg =
816 (AliMpSectorSegmentation *) AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
817 const AliMpSector * sector = seg->GetSector();
819 // get sector measurements
820 TVector2 position = sector->Position();
821 TVector2 dimension = sector->Dimensions(); // half length
823 Float_t xlocal1 = position.Px(); // FIXME: not really needed as it's 0 ?
824 Float_t ylocal1 = position.Py(); // FIXME: not really needed as it's 0 ?
825 Float_t xlocal2 = dimension.Px() * 2.;
826 Float_t ylocal2 = dimension.Px() * 2.;
828 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
829 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
832 TPolyLine3D* poly = new TPolyLine3D();
835 poly->SetPoint(nPoint++, xg1, yg1, 0.);
836 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
837 Float_t x = xg1 + xg2 * TMath::Cos(d);
838 Float_t y = yg1 + yg2 * TMath::Sin(d);
839 poly->SetPoint(nPoint++, x, y, 0.);
841 poly->SetPoint(nPoint++, xg1, yg1, 0.);
843 poly->SetLineColor(2);
852 for ( it.First(fChamber-1); ! it.IsDone(); it.Next() )
854 Int_t detElemId = it.CurrentDEId();
855 AliMp::StationType stationType = AliMpDEManager::GetStationType(detElemId);
857 if ( segmentation->HasDE(detElemId) )
859 const AliMpVSegmentation* seg
860 = AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::kCath0);
863 Float_t deltax = seg->Dimensions().X();
864 Float_t deltay = seg->Dimensions().Y();
865 Float_t xlocal1 = -deltax;
866 Float_t ylocal1 = -deltay;
867 Float_t xlocal2 = +deltax;
868 Float_t ylocal2 = +deltay;
869 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
870 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
872 // drawing slat active volumes
873 Float_t xCenter = (xg1 + xg2)/2.;
874 Float_t yCenter = (yg1 + yg2)/2.;
876 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
878 box->SetFillStyle(0);
879 box->SetLineColor( stationType == AliMp::kStationTrigger ? 4 : 2);
882 if ( stationType == AliMp::kStation345 )
884 // drawing inner circle + disc
885 TPolyLine3D* poly = new TPolyLine3D();
886 TPolyLine3D* poly1 = new TPolyLine3D();
890 for (Float_t d = 0; d < 6.24; d+= 0.005)
892 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
893 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
894 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
895 poly->SetPoint(nPoint++, x, y, 0.);
896 poly1->SetPoint(nPoint1++, x, y, 0.);
898 poly->SetLineColor(1);
900 poly1->SetLineColor(2);
908 //add clusters to the pad
912 // DrawSegmentation();
913 // add itself to the list (must be last)
915 view->SetView(phi, theta, psi, iret);
918 //_____________________________________________________________________________
919 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
921 /// Draw a view of muons chambers with tracks
923 gPad->SetCursor(kWatch);
924 // gPad->SetFillColor(39);
925 gPad->SetFillColor(1);
927 // gPad->SetFillColor(39);
928 gPad->SetFillColor(1);
932 TView *view = new TView(1);
934 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
935 view->SetRange(-range,-range,-range,range,range,range);
937 // Display all MUON Chambers segmentation
940 sprintf(nodeName,"alice");
942 node1=gAlice->GetGeometry()->GetNode(nodeName);
943 if (node1) node1->Draw("same");
946 // Draw clusters for all chambers
947 Int_t chamberSave = fChamber;
948 for (fChamber = 1; fChamber <= 10; fChamber++){
951 fChamber = chamberSave;
952 // Draw reconstructed tracks
959 Float_t x0 = (-1+shift)/zoom;
960 Float_t y0 = (-1+shift)/zoom;
961 Float_t x1 = (1+shift)/zoom;
962 Float_t y1 = (1+shift)/zoom;
963 gPad->Range(x0,y0,x1,y1);
964 view->SetView(phi, theta, psi, iret);
968 //______________________________________________________________________________
969 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
971 /// Execute action corresponding to the mouse event
973 static Float_t x0, y0, x1, y1;
975 static Int_t pxold, pyold;
976 static Int_t px0, py0;
977 static Int_t linedrawn;
980 if (px == 0 && py == 0) { //when called by sliders
981 if (event == kButton1Up) {
986 if (!fZoomMode && gPad->GetView()) {
987 gPad->GetView()->ExecuteRotateView(event, px, py);
991 // something to zoom ?
992 gPad->SetCursor(kCross);
997 gVirtualX->SetLineColor(-1);
998 gPad->TAttLine::Modify(); //Change line attributes only if necessary
999 x0 = gPad->AbsPixeltoX(px);
1000 y0 = gPad->AbsPixeltoY(py);
1002 pxold = px; pyold = py;
1006 case kButton1Motion:
1007 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1011 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1015 gPad->GetCanvas()->FeedbackMode(kFALSE);
1016 if (px == px0) return;
1017 if (py == py0) return;
1018 x1 = gPad->AbsPixeltoX(px);
1019 y1 = gPad->AbsPixeltoY(py);
1021 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1022 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1023 gPad->Range(x0,y0,x1,y1);
1024 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1026 fZoomX0[fZooms] = x0;
1027 fZoomY0[fZooms] = y0;
1028 fZoomX1[fZooms] = x1;
1029 fZoomY1[fZooms] = y1;
1031 gPad->Modified(kTRUE);
1036 //___________________________________________
1037 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1039 /// Read digits info and store x,y,z info in arrays fPoints.
1040 /// Loop on all detectors
1042 if (chamber > 14) return;
1048 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1050 GetMUONData()->SetTreeAddress("D");
1052 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
1053 if (muonDigits == 0) return;
1055 gAlice->ResetDigits();
1058 if (GetLoader()->TreeD()) {
1059 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
1060 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1061 GetMUONData()->GetDigits();
1064 Int_t ndigits = muonDigits->GetEntriesFast();
1065 if (ndigits == 0) return;
1066 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1068 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1071 AliMUONPoints *points = 0;
1072 TMarker3DBox *marker = 0;
1075 Float_t adcmax = 1024; // default
1076 if (chamber<11) adcmax = 4096;
1078 // check if trigger is using new or old segmentation
1080 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
1081 const AliMUONVGeometryDESegmentation* kdeSegmentation
1082 = segmentation->GetDESegmentation(1100, cathode-1);
1083 if ( dynamic_cast<const AliMUONTriggerSegmentation*>(kdeSegmentation) ) old = false;
1085 if ( old && chamber > 10) {
1086 if (chamber > 10) printf(">>> old segmentation for trigger \n");
1087 else printf(">>> old segmentation for tracking \n");
1089 for (Int_t digit = 0; digit < ndigits; digit++) {
1090 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1091 if (mdig->Cathode() != cathode-1) continue;
1094 // First get all needed parameters
1096 Float_t charge = mdig->Signal();
1097 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1098 Int_t color = 261+index;
1099 Int_t colorTrigger = 2;
1100 if (color > 282) color = 282;
1102 if (chamber > 10) { // trigger chamber
1104 Float_t sumCharge = 0;
1105 for (Int_t icharge = 0; icharge < 10; icharge++) {
1106 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1108 Float_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1109 if(sumCharge <= 10 || testCharge > 0) {
1110 colorTrigger = color;
1116 // get the center of the pad - add on x and y half of pad size
1117 Float_t xpad, ypad, zpad;
1121 Int_t detElemId = mdig->DetElemId();
1122 AliMUONGeometrySegmentation* segmentation2
1123 = pMUON->GetSegmentation()->GetModuleSegmentationByDEId(detElemId, cathode-1);
1124 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1125 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1126 dpx = segmentation2->Dpx(detElemId, isec)/2;
1127 dpy = segmentation2->Dpy(detElemId, isec)/2;
1129 // Then set the objects
1130 points = new AliMUONPoints(npoints);
1131 fPoints->AddAt(points,digit);
1133 points->SetMarkerColor(colorTrigger);
1135 points->SetMarkerColor(color);
1137 points->SetMarkerStyle(21);
1138 points->SetMarkerSize(0.5);
1139 points->SetParticle(-1);
1140 points->SetHitIndex(-1);
1141 points->SetTrackIndex(-1);
1142 points->SetDigitIndex(digit);
1143 points->SetPoint(0,xpad,ypad,zpos);
1145 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1146 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1149 marker->SetLineColor(lineColor);
1150 marker->SetFillStyle(1001);
1151 marker->SetFillColor(color);
1152 marker->SetRefObject((TObject*)points);
1153 points->Set3DMarker(0, marker);
1154 } // end loop on digits
1157 if (chamber > 10) printf(">>> new segmentation for trigger \n");
1158 else printf(">>> new segmentation for tracking \n");
1160 const AliMUONGeometryTransformer* kGeomTransformer
1161 = pMUON->GetGeometryTransformer();
1163 //loop over all digits and store their position
1164 for (Int_t digit = 0; digit < ndigits; digit++) {
1165 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1166 if (mdig->Cathode() != cathode-1) continue;
1168 // get all needed parameters
1169 Int_t ix=mdig->PadX();
1170 Int_t iy=mdig->PadY();
1171 Int_t detElemId=mdig->DetElemId();
1172 Float_t charge = mdig->Signal();
1173 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1174 Int_t color = 261+index;
1175 Int_t colorTrigger = 2;
1176 if (color > 282) color = 282;
1178 const AliMpVSegmentation* seg =
1179 AliMpSegmentation::Instance()
1180 ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode-1));
1182 AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
1184 if (chamber > 10) { // trigger chamber
1185 Float_t sumCharge = 0;
1186 Int_t n = mdig->Ntracks();
1187 for (Int_t icharge = 0; icharge < n; icharge++) {
1188 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1190 Float_t testCharge = sumCharge-(Int_t(sumCharge/n))*n;
1191 if(sumCharge <= n || testCharge > 0) {
1192 colorTrigger = color;
1198 // get the pad position and dimensions
1199 Float_t xlocal1 = pad.Position().X();
1200 Float_t ylocal1 = pad.Position().Y();
1201 Float_t xlocal2 = pad.Dimensions().X();
1202 Float_t ylocal2 = pad.Dimensions().Y();
1204 Float_t xg1, xg2, yg1, yg2, zg1;
1206 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
1207 // (no transformation for pad dimensions)
1211 // Then set the objects
1212 points = new AliMUONPoints(npoints);
1213 fPoints->AddAt(points,digit);
1215 points->SetMarkerColor(colorTrigger);
1217 points->SetMarkerColor(color);
1219 points->SetMarkerStyle(21);
1220 points->SetMarkerSize(0.5);
1221 points->SetParticle(-1);
1222 points->SetHitIndex(-1);
1223 points->SetTrackIndex(-1);
1224 points->SetDigitIndex(digit);
1225 points->SetPoint(0,xg1,yg1,zpos);
1227 Int_t lineColor = (zg1-zpos > 0) ? 2:3;
1228 marker=new TMarker3DBox(xg1,yg1,zpos,xg2,yg2,0,0,0);
1230 marker->SetLineColor(lineColor);
1231 marker->SetFillStyle(1001);
1232 marker->SetFillColor(color);
1233 marker->SetRefObject((TObject*)points);
1234 points->Set3DMarker(0, marker);
1236 } // end loop on digits
1237 } // end of new segmentation
1239 //___________________________________________
1240 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1242 /// Read raw clusters info and store x,y,z info in arrays fRpoints.
1243 /// Loop on all detectors
1245 if (chamber > 10) return;
1249 GetMUONData()->SetTreeAddress("RC");
1250 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1252 if (muonRawClusters == 0) return;
1255 if (GetMUONData()->TreeR()) {
1256 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1257 GetMUONData()->TreeR()->GetEvent(0);
1260 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1261 if (nrawcl == 0) return;
1262 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1264 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1265 AliMUONRawCluster *mRaw;
1266 AliMUONPoints *points = 0;
1268 //loop over all raw clusters and store their position
1269 points = new AliMUONPoints(nrawcl);
1270 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1271 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1272 points->SetMarkerColor(51);
1273 points->SetMarkerStyle(2);
1274 points->SetMarkerSize(1.);
1275 points->SetParticle(-1);
1276 points->SetHitIndex(-1);
1277 points->SetTrackIndex(-1);
1278 points->SetDigitIndex(-1);
1279 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1280 fRpoints->AddAt(points,iraw);
1281 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1285 //___________________________________________
1286 void AliMUONDisplay::LoadTracks()
1290 AliMUONTrack* recTrack = 0;
1291 AliMUONTrackParam* trackParam = 0;
1292 TClonesArray * trackParamAtHit = 0;
1296 GetMUONData()->SetTreeAddress("RT");
1297 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1298 if (recTracksArray == NULL) return;
1299 GetMUONData()->GetRecTracks();
1301 Int_t nRecTracks = 0;
1303 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1306 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1308 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1309 // reading info from tracks
1310 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1312 Int_t nTrackHits = recTrack->GetNTrackHits();
1314 if (nTrackHits == 0) continue;
1317 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1318 points->SetLineColor(6);
1319 points->SetLineWidth(1);
1320 fRpoints->AddAt(points,iRecTracks);
1326 trackParam = recTrack->GetTrackParamAtVertex();
1327 xRec = trackParam->GetNonBendingCoor();
1328 yRec = trackParam->GetBendingCoor();
1329 zRec = trackParam->GetZ();
1330 points->SetPoint(iPoint,xRec,yRec,zRec);
1333 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1334 trackParamAtHit = recTrack->GetTrackParamAtHit();
1335 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1336 xRec = trackParam->GetNonBendingCoor();
1337 yRec = trackParam->GetBendingCoor();
1338 zRec = trackParam->GetZ();
1339 points->SetPoint(iPoint,xRec,yRec,zRec);
1341 } // end loop rec. hits
1342 PrintTrack(iRecTracks,recTrack);
1343 } // end loop tracks
1348 //___________________________________________
1349 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1351 /// Print reconstructed track
1353 AliMUONTrackParam *trackParam;
1354 Float_t vertex[3], momentum[3];
1355 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1358 trackParam = recTrack->GetTrackParamAtVertex();
1359 vertex[0] = trackParam->GetNonBendingCoor();
1360 vertex[1] = trackParam->GetBendingCoor();
1361 vertex[2] = trackParam->GetZ();
1362 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1363 bendingSlope = trackParam->GetBendingSlope();
1364 nonBendingSlope = trackParam->GetNonBendingSlope();
1365 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1366 momentum[0] = momentum[2] * nonBendingSlope;
1367 momentum[1] = momentum[2] * bendingSlope;
1368 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1369 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1371 printf("===================================================\n");
1372 printf(" Reconstructed track # %d \n",iRecTracks);
1373 printf(" charge: %d \n",charge);
1374 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1375 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1376 printf(" track chi2/dof: %f \n",chi2dof);
1380 //___________________________________________
1381 void AliMUONDisplay::LoadHits(Int_t chamber)
1383 /// Read hits info and store x,y,z info in arrays fPhits.
1384 /// Loop on all detectors
1386 if (chamber > 14) return;
1393 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1395 if (GetMUONData()->TreeH()) {
1396 GetMUONData()->SetTreeAddress("H");
1397 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
1399 for (track = 0; track < ntracks; track++) {
1400 GetMUONData()->ResetHits();
1401 GetMUONData()->GetTrack(track);//skowron
1402 TClonesArray *muonHits = GetMUONData()->Hits();
1403 if (muonHits == 0) return;
1404 nthits += muonHits->GetEntriesFast();
1406 if (fPhits == 0) fPhits = new TObjArray(nthits);
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 Int_t nhits = muonHits->GetEntriesFast();
1414 if (nhits == 0) continue;
1416 AliMUONPoints *points = 0;
1418 for (Int_t hit=0;hit<nhits;hit++) {
1419 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1420 Int_t nch = mHit->Chamber(); // chamber number
1421 if (nch != chamber) continue;
1423 // Retrieve info and set the objects
1425 points = new AliMUONPoints(npoints);
1426 fPhits->AddAt(points,nhold+hit);
1427 points->SetMarkerColor(kRed);
1428 points->SetMarkerStyle(5);
1429 points->SetMarkerSize(1.);
1430 points->SetParticle(mHit->Track());
1431 points->SetHitIndex(hit);
1432 points->SetTrackIndex(track);
1433 points->SetDigitIndex(-1);
1434 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1435 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1442 //_____________________________________________________________________________
1443 void AliMUONDisplay::Paint(Option_t *)
1445 /// Paint miscellaneous items
1448 //_____________________________________________________________________________
1449 void AliMUONDisplay::SetPickMode()
1451 /// Set parameters for pick mode.
1455 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1456 fTrigPad->Modified();
1459 //_____________________________________________________________________________
1460 void AliMUONDisplay::SetZoomMode()
1462 /// Set parameters for zoom mode
1466 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1467 fTrigPad->Modified();
1470 //_____________________________________________________________________________
1471 void AliMUONDisplay::NextChamber(Int_t delta)
1473 /// To go from chamber to next chamber if delta = 1
1474 /// or previous chamber otherwise
1476 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1478 if (fChamber > 1) fChamber--;
1482 LoadDigits(fChamber, fCathode);
1486 //_____________________________________________________________________________
1487 void AliMUONDisplay::NextCathode()
1489 /// To switch to other cathode plane
1493 if (fCathode == 1) {
1494 LoadDigits(fChamber, 2);
1496 LoadDigits(fChamber, 1);
1498 fNextCathode = kTRUE; // to keep the same zoom
1500 fNextCathode = kFALSE;
1501 TPad *pad = (TPad*)gPad->GetPadSave();
1502 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1503 fZoomX1[fZooms], fZoomY1[fZooms]);
1509 //_____________________________________________________________________________
1510 void AliMUONDisplay::Trigger()
1512 /// Print global trigger output
1514 AliMUONGlobalTrigger* globalTrig;
1516 GetMUONData()->SetTreeAddress("GLT");
1517 GetMUONData()->GetTriggerD();
1519 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1520 if (globalTrig == 0) return;
1522 globalTrig->Print("full");
1524 // // returns Trigger Decision for current event
1525 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1527 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1528 // AliMUONData* muonData = decision->GetMUONData();
1529 // muonData->SetTreeAddress("D");
1530 // decision->Trigger();
1532 //_____________________________________________________________________________
1533 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1535 /// Set chamber and cathode number
1542 LoadDigits(chamber,cathode);
1546 //_____________________________________________________________________________
1547 void AliMUONDisplay::SetEvent(Int_t newevent)
1551 gAlice->GetEvent(newevent);
1553 if (!gAlice->TreeD()) return;
1556 LoadDigits(fChamber,fCathode);
1560 //_____________________________________________________________________________
1561 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1563 /// Set view range along R and Z
1573 //_____________________________________________________________________________
1574 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1576 /// Change viewing angles for current event
1584 TView *view = gPad->GetView();
1585 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1590 //_____________________________________________________________________________
1591 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1593 /// Display (current event_number + delta)
1594 /// - delta = 1 shown next event
1595 /// - delta = -1 show previous event
1597 AliRunLoader * runLoader;
1599 runLoader = fLoader->GetRunLoader();
1604 //runLoader->CleanDetectors();
1605 //runLoader->CleanKinematics();
1606 Int_t currentEvent = runLoader->GetEventNumber();
1607 Int_t newEvent = currentEvent + delta;
1608 runLoader->GetEvent(newEvent);
1611 LoadDigits(fChamber, fCathode);
1616 //______________________________________________________________________________
1617 void AliMUONDisplay::UnZoom()
1621 if (fZooms <= 0) return;
1623 TPad *pad = (TPad*)gPad->GetPadSave();
1624 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1628 //_____________________________________________________________________________
1629 void AliMUONDisplay::ResetPoints()
1631 /// Reset array of points
1639 //_____________________________________________________________________________
1640 void AliMUONDisplay::ResetPhits()
1642 /// Reset array of points
1650 //_____________________________________________________________________________
1651 void AliMUONDisplay::ResetRpoints()
1653 /// Reset array of points