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 //////////////////////////////////////////////////////////////////////////
22 // Utility class to display ALICE outline, tracks, hits,.. //
24 //////////////////////////////////////////////////////////////////////////
31 #include <TPaveLabel.h>
32 #include <TPaveText.h>
37 #include <TVirtualX.h>
39 #include <TGeometry.h>
40 #include <TMarker3DBox.h>
41 #include <TParticle.h>
42 #include <TPolyLine3D.h>
45 #include "AliMUONDisplay.h"
48 #include "AliMUONPoints.h"
49 #include "AliMUONGlobalTrigger.h"
50 #include "AliHeader.h"
52 #include "AliMUONHit.h"
53 #include "AliMUONDigit.h"
54 #include "AliMUONRawCluster.h"
55 #include "AliMUONTrack.h"
56 #include "AliMUONTrackParam.h"
58 #include "AliMUONGeometryTransformer.h"
59 #include "AliMUONGeometryModule.h"
60 #include "AliMpSlatSegmentation.h"
61 #include "AliMpSlat.h"
62 #include "AliMpSectorSegmentation.h"
63 #include "AliMpSector.h"
65 #include "AliMpTriggerSegmentation.h"
66 #include "AliMpTrigger.h"
68 #include "AliMUONSegmentation.h"
69 #include "AliMUONGeometrySegmentation.h"
70 #include "AliMUONConstants.h"
73 // to manage the same zoom on both cathodes
77 ClassImp(AliMUONDisplay)
80 //_____________________________________________________________________________
81 AliMUONDisplay::AliMUONDisplay()
89 fNextCathode = kFALSE;
93 //_____________________________________________________________________________
94 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
97 // Create an event display object.
98 // A canvas named "edisplay" is created with a vertical size in pixels
100 // A QUICK Overview of the Event Display functions
101 // ===============================================
103 // The event display can ve invoked by executing the macro "display.C"
104 // A canvas like in the picture below will appear.
106 // On the left side of the canvas, the following buttons appear:
107 // *Next* to move to the next event
108 // *Previous* to move to the previous event
110 // *Pick* Select this option to be able to point on a track with the
111 // mouse. Once on the track, use the right button to select
112 // an action. For example, select SetMarkerAttributes to
113 // change the marker type/color/size for the track.
114 // *Zoom* Select this option (default) if you want to zoom.
115 // To zoom, simply select the selected area with the left button.
116 // *UnZoom* To revert to the previous picture size.
118 // slider R On the left side, the vertical slider can be used to
119 // set the default picture size.
121 // When you are in Zoom mode, you can click on the black part of the canvas
122 // to select special options with the right mouse button.
125 // When you are in pick mode, you can "Inspect" the object pointed by the mouse.
126 // When you are on a track, select the menu item "InspectParticle"
127 // to display the current particle attributes.
129 // You can activate the Root browser by selecting the Inspect menu
130 // in the canvas tool bar menu. Then select "Start Browser"
131 // This will open a new canvas with the browser. At this point, you may want
132 // to display some histograms (from the Trees). Go to the "File" menu
133 // of the browser and click on "New canvas".
134 // In the browser, click on item "ROOT files" in the left pane.
135 // Click on galice.root.
137 // Click on TPC for example
138 // Click on any variable (eg TPC.fX) to histogram the variable.
140 // If you are lost, you can click on HELP in any Root canvas or browser.
143 <img src="gif/AliMUONDisplay.gif">
150 gAlice->SetDisplay(this);
152 // Initialize display default parameters
154 // Set front view by default
161 fDrawClusters = kTRUE;
163 fDrawTracks = kFALSE;
172 // Create display canvas
174 if (ysize < 100) ysize = 750;
175 Int_t xsize = Int_t(size*830./ysize);
176 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
177 fCanvas->ToggleEventStatus();
179 // Create main display pad
180 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
183 fPad->SetFillColor(30);
184 fPad->SetBorderSize(2);
189 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
191 fColPad->SetFillColor(17);
192 fColPad->SetBorderSize(2);
197 // Create user interface control pad
201 // Create Range and mode pad
204 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
207 fTrigPad->SetFillColor(22);
208 fTrigPad->SetBorderSize(2);
209 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
210 fRangeSlider->SetObject(this);
211 char pickmode[] = "gAlice->Display()->SetPickMode()";
213 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
214 fPickButton->SetFillColor(38);
216 char zoommode[] = "gAlice->Display()->SetZoomMode()";
217 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
218 fZoomButton->SetFillColor(38);
220 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
221 fArcButton->SetFillColor(kGreen);
223 char butUnzoom[] = "gAlice->Display()->UnZoom()";
224 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
225 button->SetFillColor(38);
227 AppendPad(); // append display object as last object to force selection
230 fTrigPad->SetEditable(kFALSE);
231 fButtons->SetEditable(kFALSE);
233 fNextCathode = kFALSE;
235 // initialize container
237 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
242 AliMUONDisplay::AliMUONDisplay(const AliMUONDisplay & display)
243 : AliDisplay(display)
245 // Protected copy constructor
247 AliFatal("Not implemented.");
252 //_____________________________________________________________________________
253 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
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()
498 fDrawTracks = kFALSE;
500 DrawView(fTheta, fPhi, fPsi);
501 // Display the event number and title
506 //_____________________________________________________________________________
507 void AliMUONDisplay::DrawReco(Option_t *)
509 // Display current event
512 // print kinematics of generated particles
514 // Draw global view of muon system
516 DrawGlobalView(135, -50, -140);
518 // Display the event number and title
523 //_____________________________________________________________________________
524 void AliMUONDisplay::PrintKinematics()
526 AliRunLoader * runLoader;
527 TParticle *particle = new TParticle();
529 Float_t vertex[3], momentum[3];
532 runLoader = fLoader->GetRunLoader();
536 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
537 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
538 nPart = (Int_t)runLoader->TreeK()->GetEntries();
539 for(Int_t iPart = 0; iPart < nPart; iPart++) {
540 runLoader->TreeK()->GetEvent(iPart);
541 vertex[0] = particle->Vx();
542 vertex[1] = particle->Vy();
543 vertex[2] = particle->Vz();
544 momentum[0] = particle->Px();
545 momentum[1] = particle->Py();
546 momentum[2] = particle->Pz();
548 printf("===================================================\n");
549 printf(" Generated particle # %d \n",iPart);
550 printf(" name: %s \n",particle->GetName());
551 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
552 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
557 void AliMUONDisplay::DrawSegmentation()
559 // to be re-written for new seg
560 // Draw graphical representation of segmenatation
561 // Attention: still experimental code
564 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
565 // AliMUONChamber* iChamber;
567 // AliSegmentation* seg;
568 // iChamber = &(pMUON->Chamber(fChamber));
569 // seg=iChamber->SegmentationModel(icat);
571 // Float_t zpos=iChamber->Z();
572 // Float_t r=iChamber->ROuter();
574 // TMarker3DBox *marker;
576 // for (Int_t j=0; j<seg->Npy(); j++) {
578 // y0=j*seg->Dpy()-seg->Dpy()/2.;
579 // for (seg->FirstPad(0.,y0,0,300,0.);
583 // if (seg->ISector()==0) continue;
585 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
586 // Float_t dpx=seg->Dpx(seg->ISector())/2;
587 // Float_t dpy=seg->Dpy(seg->ISector())/2;
588 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
589 // marker->SetLineColor(seg->ISector()+1);
590 // marker->SetFillStyle(1001);
591 // marker->SetFillColor(0);
596 // for (Int_t j=0; j<250; j++) {
597 // Float_t x0=j*seg->Dpx();
598 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
600 // for (seg->FirstPad(x0,0,0,0,y0);
604 // if (seg->ISector()==0) continue;
607 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
608 // Float_t dpx=seg->Dpx(seg->ISector())/2;
609 // Float_t dpy=seg->Dpy(seg->ISector())/2;
610 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
611 // marker->SetLineColor(seg->ISector()+1);
612 // marker->SetFillStyle(1001);
613 // marker->SetFillColor(0);
620 //_____________________________________________________________________________
621 void AliMUONDisplay::DrawClusters()
623 // Draw clusters for MUON chambers
625 Int_t ndigits, digit;
633 ndigits = points->GetEntriesFast();
634 for (digit=0;digit<ndigits;digit++){
635 pm = (AliMUONPoints*)points->UncheckedAt(digit);
639 for (Int_t im=0;im<3;im++) {
640 TMarker3DBox *marker=pm->GetMarker(im);
645 fClustersCuts += pm->GetN();
649 //_____________________________________________________________________________
650 void AliMUONDisplay::DrawHits()
652 // Draw hits for MUON chambers
656 Int_t ntracks, track;
663 ntracks = points->GetEntriesFast();
664 for (track=0;track<ntracks;track++) {
665 pm = (AliMUONPoints*)points->UncheckedAt(track);
668 fHitsCuts += pm->GetN();
673 //_____________________________________________________________________________
674 void AliMUONDisplay::DrawCoG()
676 // Draw hits for MUON chambers
677 if (!fDrawCoG) return;
678 if (fChamber > 10) return;
679 LoadCoG(fChamber,fCathode);
687 ncog = points->GetEntriesFast();
688 for (icog=0;icog<ncog;icog++) {
689 pm = (AliMUONPoints*)points->UncheckedAt(icog);
694 //_____________________________________________________________________________
695 void AliMUONDisplay::DrawTracks()
698 if (!fDrawTracks) return;
701 Int_t nTrack, iTrack;
707 nTrack = points->GetEntriesFast();
708 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
709 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
714 //_____________________________________________________________________________
716 void AliMUONDisplay::DrawTitle(Option_t *option)
718 // Draw the event title
720 Float_t xmin = gPad->GetX1();
721 Float_t xmax = gPad->GetX2();
722 Float_t ymin = gPad->GetY1();
723 Float_t ymax = gPad->GetY2();
724 Float_t dx = xmax-xmin;
725 Float_t dy = ymax-ymin;
727 AliRunLoader * runLoader;
729 runLoader = fLoader->GetRunLoader();
734 if (strlen(option) == 0) {
735 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
736 // title->SetTextSize(0.023932);
737 title->SetTextSize(0.02);
738 title->SetBit(kCanDelete);
739 title->SetFillColor(42);
742 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
743 runLoader->GetEventNumber(),
744 gAlice->GetHeader()->GetRun(),
747 title->AddText(ptitle);
748 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
749 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
750 nparticles, fHitsCuts,fClustersCuts);
751 title->AddText(ptitle);
753 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
754 label->SetBit(kCanDelete);
755 label->SetFillColor(42);
760 //_____________________________________________________________________________
761 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
763 // Draw a view of MUON clusters
764 AliInfo(" Draw View");
766 gPad->SetCursor(kWatch);
767 // gPad->SetFillColor(39);
768 gPad->SetFillColor(1);
770 // gPad->SetFillColor(39);
771 gPad->SetFillColor(1);
775 TView *view = new TView(1);
777 Float_t range = fRrange*fRangeSlider->GetMaximum();
778 view->SetRange(-range,-range,-range,range, range, range);
779 // zoom back to full scale only if DrawView not called from NextCathode
788 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
790 // Recovering the chamber
791 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
793 const AliMUONGeometryTransformer* kGeomTransformer
794 = pMUON->GetGeometryTransformer();
796 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
798 // Display MUON Chamber Geometry
800 sprintf(nodeName,"MUON%d",100+fChamber);
801 printf(">>>> chamber is %d\n",fChamber);
804 for(Int_t id = 0; id < 4; id++) {
806 Int_t detElemId = fChamber*100+id;
807 AliMpSectorSegmentation * seg =
808 (AliMpSectorSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
809 const AliMpSector * sector = seg->GetSector();
811 // get sector measurements
812 TVector2 position = sector->Position();
813 TVector2 dimension = sector->Dimensions(); // half length
815 Float_t xlocal1 = position.Px(); // mm to cm
816 Float_t ylocal1 = position.Py();
817 Float_t xlocal2 = dimension.Px() * 2.;
818 Float_t ylocal2 = dimension.Px() * 2.;
820 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
821 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
824 TPolyLine3D* poly = new TPolyLine3D();
827 poly->SetPoint(nPoint++, xg1, yg1, 0.);
828 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
829 Float_t x = xg1 + xg2 * TMath::Cos(d);
830 Float_t y = yg1 + yg2 * TMath::Sin(d);
831 poly->SetPoint(nPoint++, x, y, 0.);
833 poly->SetPoint(nPoint++, xg1, yg1, 0.);
835 poly->SetLineColor(2);
843 if(fChamber >4 && fChamber <11) {
845 for(id=0; id<26; id++) {
846 Int_t detElemId = fChamber*100+id;
847 if ( segmentation->HasDE(detElemId) ) {
848 AliMpSlatSegmentation * seg =
849 (AliMpSlatSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
851 const AliMpSlat* slat = seg->Slat();
852 Float_t deltax = slat->DX();
853 Float_t deltay = slat->DY();
854 Float_t xlocal1 = -deltax;
855 Float_t ylocal1 = -deltay;
856 Float_t xlocal2 = +deltax;
857 Float_t ylocal2 = +deltay;
858 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
859 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
861 // drawing slat active volumes
862 Float_t xCenter = (xg1 + xg2)/2.;
863 Float_t yCenter = (yg1 + yg2)/2.;
865 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
867 box->SetFillStyle(0);
868 box->SetLineColor(2);
871 // drawing inner circle + disc
872 TPolyLine3D* poly = new TPolyLine3D();
873 TPolyLine3D* poly1 = new TPolyLine3D();
877 for (Float_t d = 0; d < 6.24; d+= 0.005) {
878 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
879 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
880 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
881 poly->SetPoint(nPoint++, x, y, 0.);
882 poly1->SetPoint(nPoint1++, x, y, 0.);
885 poly->SetLineColor(1);
887 poly1->SetLineColor(2);
894 /*-- Trigger Chambers ---------------------------------------*/
895 if(fChamber >10 && fChamber <15) {
897 for(id=0; id<18; id++) {
898 Int_t detElemId = fChamber*100+id;
899 AliMpTriggerSegmentation * seg
900 = (AliMpTriggerSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
902 const AliMpTrigger* slat = seg->Slat();
903 Float_t deltax = slat->DX();
904 Float_t deltay = slat->DY();
905 Float_t xlocal1 = -deltax;
906 Float_t ylocal1 = -deltay;
907 Float_t xlocal2 = +deltax;
908 Float_t ylocal2 = +deltay;
910 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
911 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
913 // drawing slat active volumes
914 Float_t xCenter = (xg1 + xg2)/2.;
915 Float_t yCenter = (yg1 + yg2)/2.;
917 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
919 box->SetFillStyle(0);
920 box->SetLineColor(4);
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 TView *view = new TView(1);
951 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
952 view->SetRange(-range,-range,-range,range,range,range);
954 // Display all MUON Chambers segmentation
957 sprintf(nodeName,"alice");
959 node1=gAlice->GetGeometry()->GetNode(nodeName);
960 if (node1) node1->Draw("same");
963 // Draw clusters for all chambers
964 Int_t chamberSave = fChamber;
965 for (fChamber = 1; fChamber <= 10; fChamber++){
968 fChamber = chamberSave;
969 // Draw reconstructed tracks
976 Float_t x0 = (-1+shift)/zoom;
977 Float_t y0 = (-1+shift)/zoom;
978 Float_t x1 = (1+shift)/zoom;
979 Float_t y1 = (1+shift)/zoom;
980 gPad->Range(x0,y0,x1,y1);
981 view->SetView(phi, theta, psi, iret);
985 //______________________________________________________________________________
986 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
988 // Execute action corresponding to the mouse event
990 static Float_t x0, y0, x1, y1;
992 static Int_t pxold, pyold;
993 static Int_t px0, py0;
994 static Int_t linedrawn;
997 if (px == 0 && py == 0) { //when called by sliders
998 if (event == kButton1Up) {
1003 if (!fZoomMode && gPad->GetView()) {
1004 gPad->GetView()->ExecuteRotateView(event, px, py);
1008 // something to zoom ?
1009 gPad->SetCursor(kCross);
1014 gVirtualX->SetLineColor(-1);
1015 gPad->TAttLine::Modify(); //Change line attributes only if necessary
1016 x0 = gPad->AbsPixeltoX(px);
1017 y0 = gPad->AbsPixeltoY(py);
1019 pxold = px; pyold = py;
1023 case kButton1Motion:
1024 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1028 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1032 gPad->GetCanvas()->FeedbackMode(kFALSE);
1033 if (px == px0) return;
1034 if (py == py0) return;
1035 x1 = gPad->AbsPixeltoX(px);
1036 y1 = gPad->AbsPixeltoY(py);
1038 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1039 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1040 gPad->Range(x0,y0,x1,y1);
1041 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1043 fZoomX0[fZooms] = x0;
1044 fZoomY0[fZooms] = y0;
1045 fZoomX1[fZooms] = x1;
1046 fZoomY1[fZooms] = y1;
1048 gPad->Modified(kTRUE);
1053 //___________________________________________
1054 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1056 // Read digits info and store x,y,z info in arrays fPoints
1057 // Loop on all detectors
1059 if (chamber > 14) return;
1065 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1066 AliMUONGeometrySegmentation* segmentation2 = 0x0;
1068 GetMUONData()->SetTreeAddress("D");
1070 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
1071 if (muonDigits == 0) return;
1073 gAlice->ResetDigits();
1076 if (GetLoader()->TreeD()) {
1077 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
1078 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1079 GetMUONData()->GetDigits();
1082 Int_t ndigits = muonDigits->GetEntriesFast();
1083 if (ndigits == 0) return;
1084 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1086 //segmentation2 = iChamber->SegmentationModel2(cathode);
1088 = pMUON->GetSegmentation()->GetModuleSegmentation(chamber-1, cathode-1);
1090 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1093 AliMUONPoints *points = 0;
1094 TMarker3DBox *marker = 0;
1096 //loop over all digits and store their position
1099 Float_t adcmax = 1024; // default
1100 if (chamber<11) adcmax = 4096;
1102 for (Int_t digit = 0; digit < ndigits; digit++) {
1103 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1104 if (mdig->Cathode() != cathode-1) continue;
1107 // First get all needed parameters
1109 Int_t charge = mdig->Signal();
1110 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1111 Int_t color = 261+index;
1112 Int_t colorTrigger = 2;
1113 if (color > 282) color = 282;
1115 if (chamber > 10) { // trigger chamber
1117 Int_t sumCharge = 0;
1118 for (Int_t icharge = 0; icharge < 10; icharge++) {
1119 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1121 Int_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1122 if(sumCharge <= 10 || testCharge > 0) {
1123 colorTrigger = color;
1129 // get the center of the pad - add on x and y half of pad size
1130 Float_t xpad, ypad, zpad;
1134 Int_t detElemId = mdig->DetElemId();
1135 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1136 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1137 dpx = segmentation2->Dpx(detElemId, isec)/2;
1138 dpy = segmentation2->Dpy(detElemId, isec)/2;
1140 // segmentation->Dump();
1143 // Then set the objects
1145 points = new AliMUONPoints(npoints);
1146 fPoints->AddAt(points,digit);
1148 points->SetMarkerColor(colorTrigger);
1150 points->SetMarkerColor(color);
1152 points->SetMarkerStyle(21);
1153 points->SetMarkerSize(0.5);
1154 points->SetParticle(-1);
1155 points->SetHitIndex(-1);
1156 points->SetTrackIndex(-1);
1157 points->SetDigitIndex(digit);
1158 points->SetPoint(0,xpad,ypad,zpos);
1160 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1161 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1164 marker->SetLineColor(lineColor);
1165 marker->SetFillStyle(1001);
1166 marker->SetFillColor(color);
1167 marker->SetRefObject((TObject*)points);
1168 points->Set3DMarker(0, marker);
1171 //___________________________________________
1172 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1174 // Read raw clusters info and store x,y,z info in arrays fRpoints
1175 // Loop on all detectors
1177 if (chamber > 10) return;
1181 GetMUONData()->SetTreeAddress("RC");
1182 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1184 if (muonRawClusters == 0) return;
1187 if (GetMUONData()->TreeR()) {
1188 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1189 GetMUONData()->TreeR()->GetEvent(0);
1192 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1193 if (nrawcl == 0) return;
1194 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1196 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1197 AliMUONRawCluster *mRaw;
1198 AliMUONPoints *points = 0;
1200 //loop over all raw clusters and store their position
1201 points = new AliMUONPoints(nrawcl);
1202 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1203 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1204 points->SetMarkerColor(51);
1205 points->SetMarkerStyle(2);
1206 points->SetMarkerSize(1.);
1207 points->SetParticle(-1);
1208 points->SetHitIndex(-1);
1209 points->SetTrackIndex(-1);
1210 points->SetDigitIndex(-1);
1211 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1212 fRpoints->AddAt(points,iraw);
1213 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1217 //___________________________________________
1218 void AliMUONDisplay::LoadTracks()
1221 AliMUONTrack* recTrack = 0;
1222 AliMUONTrackParam* trackParam = 0;
1223 TClonesArray * trackParamAtHit = 0;
1227 GetMUONData()->SetTreeAddress("RT");
1228 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1229 if (recTracksArray == NULL) return;
1230 GetMUONData()->GetRecTracks();
1232 Int_t nRecTracks = 0;
1234 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1237 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1239 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1240 // reading info from tracks
1241 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1243 Int_t nTrackHits = recTrack->GetNTrackHits();
1245 if (nTrackHits == 0) continue;
1248 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1249 points->SetLineColor(6);
1250 points->SetLineWidth(1);
1251 fRpoints->AddAt(points,iRecTracks);
1257 trackParam = recTrack->GetTrackParamAtVertex();
1258 xRec = trackParam->GetNonBendingCoor();
1259 yRec = trackParam->GetBendingCoor();
1260 zRec = trackParam->GetZ();
1261 points->SetPoint(iPoint,xRec,yRec,zRec);
1264 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1265 trackParamAtHit = recTrack->GetTrackParamAtHit();
1266 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1267 xRec = trackParam->GetNonBendingCoor();
1268 yRec = trackParam->GetBendingCoor();
1269 zRec = trackParam->GetZ();
1270 points->SetPoint(iPoint,xRec,yRec,zRec);
1272 } // end loop rec. hits
1273 PrintTrack(iRecTracks,recTrack);
1274 } // end loop tracks
1279 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1281 AliMUONTrackParam *trackParam;
1282 Float_t vertex[3], momentum[3];
1283 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1286 trackParam = recTrack->GetTrackParamAtVertex();
1287 vertex[0] = trackParam->GetNonBendingCoor();
1288 vertex[1] = trackParam->GetBendingCoor();
1289 vertex[2] = trackParam->GetZ();
1290 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1291 bendingSlope = trackParam->GetBendingSlope();
1292 nonBendingSlope = trackParam->GetNonBendingSlope();
1293 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1294 momentum[0] = momentum[2] * nonBendingSlope;
1295 momentum[1] = momentum[2] * bendingSlope;
1296 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1297 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1299 printf("===================================================\n");
1300 printf(" Reconstructed track # %d \n",iRecTracks);
1301 printf(" charge: %d \n",charge);
1302 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1303 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1304 printf(" track chi2/dof: %f \n",chi2dof);
1308 //___________________________________________
1309 void AliMUONDisplay::LoadHits(Int_t chamber)
1311 // Read hits info and store x,y,z info in arrays fPhits
1312 // Loop on all detectors
1314 if (chamber > 14) return;
1321 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1323 if (GetMUONData()->TreeH()) {
1324 GetMUONData()->SetTreeAddress("H");
1325 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
1327 for (track = 0; track < ntracks; track++) {
1328 GetMUONData()->ResetHits();
1329 GetMUONData()->GetTrack(track);//skowron
1330 TClonesArray *muonHits = GetMUONData()->Hits();
1331 if (muonHits == 0) return;
1332 nthits += muonHits->GetEntriesFast();
1334 if (fPhits == 0) fPhits = new TObjArray(nthits);
1336 for (track=0; track<ntracks;track++) {
1337 GetMUONData()->ResetHits();
1338 GetMUONData()->GetTrack(track);//skowron
1339 TClonesArray *muonHits = GetMUONData()->Hits();
1340 if (muonHits == 0) return;
1341 Int_t nhits = muonHits->GetEntriesFast();
1342 if (nhits == 0) continue;
1344 AliMUONPoints *points = 0;
1346 for (Int_t hit=0;hit<nhits;hit++) {
1347 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1348 Int_t nch = mHit->Chamber(); // chamber number
1349 if (nch != chamber) continue;
1351 // Retrieve info and set the objects
1353 points = new AliMUONPoints(npoints);
1354 fPhits->AddAt(points,nhold+hit);
1355 points->SetMarkerColor(kRed);
1356 points->SetMarkerStyle(5);
1357 points->SetMarkerSize(1.);
1358 points->SetParticle(mHit->Track());
1359 points->SetHitIndex(hit);
1360 points->SetTrackIndex(track);
1361 points->SetDigitIndex(-1);
1362 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1363 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1370 //_____________________________________________________________________________
1371 void AliMUONDisplay::Paint(Option_t *)
1373 // Paint miscellaneous items
1376 //_____________________________________________________________________________
1377 void AliMUONDisplay::SetPickMode()
1379 // Set parameters for pick mode.
1383 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1384 fTrigPad->Modified();
1387 //_____________________________________________________________________________
1388 void AliMUONDisplay::SetZoomMode()
1390 // Set parameters for zoom mode
1393 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1394 fTrigPad->Modified();
1397 //_____________________________________________________________________________
1398 void AliMUONDisplay::NextChamber(Int_t delta)
1400 // to go from chamber to next chamber if delta = 1
1401 // or previous chamber otherwise
1403 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1405 if (fChamber > 1) fChamber--;
1409 LoadDigits(fChamber, fCathode);
1413 //_____________________________________________________________________________
1414 void AliMUONDisplay::NextCathode()
1416 // to switch to other cathode plane
1419 if (fCathode == 1) {
1420 LoadDigits(fChamber, 2);
1422 LoadDigits(fChamber, 1);
1424 fNextCathode = kTRUE; // to keep the same zoom
1426 fNextCathode = kFALSE;
1427 TPad *pad = (TPad*)gPad->GetPadSave();
1428 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1429 fZoomX1[fZooms], fZoomY1[fZooms]);
1435 //_____________________________________________________________________________
1436 void AliMUONDisplay::Trigger()
1439 AliMUONGlobalTrigger* globalTrig;
1441 GetMUONData()->SetTreeAddress("GLT");
1442 GetMUONData()->GetTriggerD();
1444 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1445 if (globalTrig == 0) return;
1447 printf("===================================================\n");
1448 printf(" Global Trigger output Low pt High pt All\n");
1450 printf(" number of Single Plus :\t");
1451 printf("%i\t",globalTrig->SinglePlusLpt());
1452 printf("%i\t",globalTrig->SinglePlusHpt());
1453 printf("%i\t",globalTrig->SinglePlusApt());
1456 printf(" number of Single Minus :\t");
1457 printf("%i\t",globalTrig->SingleMinusLpt());
1458 printf("%i\t",globalTrig->SingleMinusHpt());
1459 printf("%i\t",globalTrig->SingleMinusApt());
1462 printf(" number of Single Undefined :\t");
1463 printf("%i\t",globalTrig->SingleUndefLpt());
1464 printf("%i\t",globalTrig->SingleUndefHpt());
1465 printf("%i\t",globalTrig->SingleUndefApt());
1468 printf(" number of UnlikeSign pair :\t");
1469 printf("%i\t",globalTrig->PairUnlikeLpt());
1470 printf("%i\t",globalTrig->PairUnlikeHpt());
1471 printf("%i\t",globalTrig->PairUnlikeApt());
1474 printf(" number of LikeSign pair :\t");
1475 printf("%i\t",globalTrig->PairLikeLpt());
1476 printf("%i\t",globalTrig->PairLikeHpt());
1477 printf("%i\t",globalTrig->PairLikeApt());
1479 printf("===================================================\n");
1483 // // returns Trigger Decision for current event
1484 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1486 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1487 // AliMUONData* muonData = decision->GetMUONData();
1488 // muonData->SetTreeAddress("D");
1489 // decision->Trigger();
1491 //_____________________________________________________________________________
1492 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1494 // Set chamber and cathode number
1500 LoadDigits(chamber,cathode);
1504 void AliMUONDisplay::SetEvent(Int_t newevent)
1507 gAlice->GetEvent(newevent);
1509 if (!gAlice->TreeD()) return;
1512 LoadDigits(fChamber,fCathode);
1516 //_____________________________________________________________________________
1517 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1519 // Set view range along R and Z
1528 //_____________________________________________________________________________
1529 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1531 // change viewing angles for current event
1539 TView *view = gPad->GetView();
1540 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1545 //_____________________________________________________________________________
1546 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1548 AliRunLoader * runLoader;
1550 runLoader = fLoader->GetRunLoader();
1554 // Display (current event_number + delta)
1555 // delta = 1 shown next event
1556 // delta = -1 show previous event
1558 //runLoader->CleanDetectors();
1559 //runLoader->CleanKinematics();
1560 Int_t currentEvent = runLoader->GetEventNumber();
1561 Int_t newEvent = currentEvent + delta;
1562 runLoader->GetEvent(newEvent);
1565 LoadDigits(fChamber, fCathode);
1570 //______________________________________________________________________________
1571 void AliMUONDisplay::UnZoom()
1574 if (fZooms <= 0) return;
1576 TPad *pad = (TPad*)gPad->GetPadSave();
1577 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1581 //_____________________________________________________________________________
1582 void AliMUONDisplay::ResetPoints()
1585 // Reset array of points
1593 //_____________________________________________________________________________
1594 void AliMUONDisplay::ResetPhits()
1597 // Reset array of points
1605 //_____________________________________________________________________________
1606 void AliMUONDisplay::ResetRpoints()
1609 // Reset array of points
1618 AliMUONDisplay & AliMUONDisplay::operator = (const AliMUONDisplay & rhs)
1620 // Protected assignement operator
1622 if (this == &rhs) return *this;
1624 AliFatal("Not implemented.");