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 "AliMUONSegmentationManager.h"
71 #include "AliMUONConstants.h"
74 // to manage the same zoom on both cathodes
78 ClassImp(AliMUONDisplay)
81 //_____________________________________________________________________________
82 AliMUONDisplay::AliMUONDisplay()
90 fNextCathode = kFALSE;
94 //_____________________________________________________________________________
95 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
98 // Create an event display object.
99 // A canvas named "edisplay" is created with a vertical size in pixels
101 // A QUICK Overview of the Event Display functions
102 // ===============================================
104 // The event display can ve invoked by executing the macro "display.C"
105 // A canvas like in the picture below will appear.
107 // On the left side of the canvas, the following buttons appear:
108 // *Next* to move to the next event
109 // *Previous* to move to the previous event
111 // *Pick* Select this option to be able to point on a track with the
112 // mouse. Once on the track, use the right button to select
113 // an action. For example, select SetMarkerAttributes to
114 // change the marker type/color/size for the track.
115 // *Zoom* Select this option (default) if you want to zoom.
116 // To zoom, simply select the selected area with the left button.
117 // *UnZoom* To revert to the previous picture size.
119 // slider R On the left side, the vertical slider can be used to
120 // set the default picture size.
122 // When you are in Zoom mode, you can click on the black part of the canvas
123 // to select special options with the right mouse button.
126 // When you are in pick mode, you can "Inspect" the object pointed by the mouse.
127 // When you are on a track, select the menu item "InspectParticle"
128 // to display the current particle attributes.
130 // You can activate the Root browser by selecting the Inspect menu
131 // in the canvas tool bar menu. Then select "Start Browser"
132 // This will open a new canvas with the browser. At this point, you may want
133 // to display some histograms (from the Trees). Go to the "File" menu
134 // of the browser and click on "New canvas".
135 // In the browser, click on item "ROOT files" in the left pane.
136 // Click on galice.root.
138 // Click on TPC for example
139 // Click on any variable (eg TPC.fX) to histogram the variable.
141 // If you are lost, you can click on HELP in any Root canvas or browser.
144 <img src="gif/AliMUONDisplay.gif">
151 gAlice->SetDisplay(this);
153 // Initialize display default parameters
155 // Set front view by default
162 fDrawClusters = kTRUE;
164 fDrawTracks = kFALSE;
173 // Create display canvas
175 if (ysize < 100) ysize = 750;
176 Int_t xsize = Int_t(size*830./ysize);
177 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
178 fCanvas->ToggleEventStatus();
180 // Create main display pad
181 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
184 fPad->SetFillColor(30);
185 fPad->SetBorderSize(2);
190 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
192 fColPad->SetFillColor(17);
193 fColPad->SetBorderSize(2);
198 // Create user interface control pad
202 // Create Range and mode pad
205 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
208 fTrigPad->SetFillColor(22);
209 fTrigPad->SetBorderSize(2);
210 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
211 fRangeSlider->SetObject(this);
212 char pickmode[] = "gAlice->Display()->SetPickMode()";
214 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
215 fPickButton->SetFillColor(38);
217 char zoommode[] = "gAlice->Display()->SetZoomMode()";
218 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
219 fZoomButton->SetFillColor(38);
221 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
222 fArcButton->SetFillColor(kGreen);
224 char butUnzoom[] = "gAlice->Display()->UnZoom()";
225 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
226 button->SetFillColor(38);
228 AppendPad(); // append display object as last object to force selection
231 fTrigPad->SetEditable(kFALSE);
232 fButtons->SetEditable(kFALSE);
234 fNextCathode = kFALSE;
236 // initialize container
238 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
243 AliMUONDisplay::AliMUONDisplay(const AliMUONDisplay & display)
244 : AliDisplay(display)
246 // Protected copy constructor
248 AliFatal("Not implemented.");
253 //_____________________________________________________________________________
254 AliMUONDisplay::~AliMUONDisplay()
256 // Delete space point structure
257 if (fPoints) fPoints->Delete();
261 if (fPhits) fPhits->Delete();
265 if (fRpoints) fRpoints->Delete();
270 //_____________________________________________________________________________
271 void AliMUONDisplay::Clear(Option_t *)
273 // Delete graphics temporary objects
276 //_____________________________________________________________________________
277 void AliMUONDisplay::DisplayButtons()
279 // Create the user interface buttons
282 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
284 fButtons->SetFillColor(38);
285 fButtons->SetBorderSize(2);
288 // Int_t butcolor = 33;
289 Float_t dbutton = 0.08;
296 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
297 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
298 button->SetFillColor(38);
302 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
303 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
304 button->SetFillColor(38);
308 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
309 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
310 button->SetFillColor(38);
314 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
315 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
316 button->SetFillColor(38);
320 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
321 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
322 button->SetFillColor(38);
326 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
327 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
328 button->SetFillColor(38);
332 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
333 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
334 button->SetFillColor(38);
338 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
339 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
340 button->SetFillColor(38);
344 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
345 diamond->SetFillColor(50);
346 diamond->SetTextAlign(22);
347 diamond->SetTextColor(5);
348 diamond->SetTextSize(0.11);
350 diamond->AddText(".. ");
351 diamond->AddText("ROOT");
352 diamond->AddText("MUON");
353 diamond->AddText("... ");
354 diamond->AddText(" ");
357 //_____________________________________________________________________________
358 void AliMUONDisplay::CreateColors() const
360 // Create the colors palette used to display clusters
375 new TColor(color,r,g,b);
385 new TColor(color,r,g,b);
395 new TColor(color,r,g,b);
405 new TColor(color,r,g,b);
415 new TColor(color,r,g,b);
422 //_____________________________________________________________________________
423 void AliMUONDisplay::DisplayColorScale()
425 // 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()
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 AliRunLoader * runLoader;
528 TParticle *particle = new TParticle();
530 Float_t vertex[3], momentum[3];
533 runLoader = fLoader->GetRunLoader();
537 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
538 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
539 nPart = (Int_t)runLoader->TreeK()->GetEntries();
540 for(Int_t iPart = 0; iPart < nPart; iPart++) {
541 runLoader->TreeK()->GetEvent(iPart);
542 vertex[0] = particle->Vx();
543 vertex[1] = particle->Vy();
544 vertex[2] = particle->Vz();
545 momentum[0] = particle->Px();
546 momentum[1] = particle->Py();
547 momentum[2] = particle->Pz();
549 printf("===================================================\n");
550 printf(" Generated particle # %d \n",iPart);
551 printf(" name: %s \n",particle->GetName());
552 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
553 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
558 void AliMUONDisplay::DrawSegmentation()
560 // to be re-written for new seg
561 // Draw graphical representation of segmenatation
562 // Attention: still experimental code
565 // AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
566 // AliMUONChamber* iChamber;
568 // AliSegmentation* seg;
569 // iChamber = &(pMUON->Chamber(fChamber));
570 // seg=iChamber->SegmentationModel(icat);
572 // Float_t zpos=iChamber->Z();
573 // Float_t r=iChamber->ROuter();
575 // TMarker3DBox *marker;
577 // for (Int_t j=0; j<seg->Npy(); j++) {
579 // y0=j*seg->Dpy()-seg->Dpy()/2.;
580 // for (seg->FirstPad(0.,y0,0,300,0.);
584 // if (seg->ISector()==0) continue;
586 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
587 // Float_t dpx=seg->Dpx(seg->ISector())/2;
588 // Float_t dpy=seg->Dpy(seg->ISector())/2;
589 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
590 // marker->SetLineColor(seg->ISector()+1);
591 // marker->SetFillStyle(1001);
592 // marker->SetFillColor(0);
597 // for (Int_t j=0; j<250; j++) {
598 // Float_t x0=j*seg->Dpx();
599 // Float_t y0=TMath::Sqrt(r*r-x0*x0);
601 // for (seg->FirstPad(x0,0,0,0,y0);
605 // if (seg->ISector()==0) continue;
608 // seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
609 // Float_t dpx=seg->Dpx(seg->ISector())/2;
610 // Float_t dpy=seg->Dpy(seg->ISector())/2;
611 // marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
612 // marker->SetLineColor(seg->ISector()+1);
613 // marker->SetFillStyle(1001);
614 // marker->SetFillColor(0);
621 //_____________________________________________________________________________
622 void AliMUONDisplay::DrawClusters()
624 // Draw clusters for MUON chambers
626 Int_t ndigits, digit;
634 ndigits = points->GetEntriesFast();
635 for (digit=0;digit<ndigits;digit++){
636 pm = (AliMUONPoints*)points->UncheckedAt(digit);
640 for (Int_t im=0;im<3;im++) {
641 TMarker3DBox *marker=pm->GetMarker(im);
646 fClustersCuts += pm->GetN();
650 //_____________________________________________________________________________
651 void AliMUONDisplay::DrawHits()
653 // Draw hits for MUON chambers
657 Int_t ntracks, track;
664 ntracks = points->GetEntriesFast();
665 for (track=0;track<ntracks;track++) {
666 pm = (AliMUONPoints*)points->UncheckedAt(track);
669 fHitsCuts += pm->GetN();
674 //_____________________________________________________________________________
675 void AliMUONDisplay::DrawCoG()
677 // Draw hits for MUON chambers
678 if (!fDrawCoG) return;
679 if (fChamber > 10) return;
680 LoadCoG(fChamber,fCathode);
688 ncog = points->GetEntriesFast();
689 for (icog=0;icog<ncog;icog++) {
690 pm = (AliMUONPoints*)points->UncheckedAt(icog);
695 //_____________________________________________________________________________
696 void AliMUONDisplay::DrawTracks()
699 if (!fDrawTracks) return;
702 Int_t nTrack, iTrack;
708 nTrack = points->GetEntriesFast();
709 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
710 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
715 //_____________________________________________________________________________
717 void AliMUONDisplay::DrawTitle(Option_t *option)
719 // Draw the event title
721 Float_t xmin = gPad->GetX1();
722 Float_t xmax = gPad->GetX2();
723 Float_t ymin = gPad->GetY1();
724 Float_t ymax = gPad->GetY2();
725 Float_t dx = xmax-xmin;
726 Float_t dy = ymax-ymin;
728 AliRunLoader * runLoader;
730 runLoader = fLoader->GetRunLoader();
735 if (strlen(option) == 0) {
736 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
737 // title->SetTextSize(0.023932);
738 title->SetTextSize(0.02);
739 title->SetBit(kCanDelete);
740 title->SetFillColor(42);
743 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
744 runLoader->GetEventNumber(),
745 gAlice->GetHeader()->GetRun(),
748 title->AddText(ptitle);
749 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
750 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
751 nparticles, fHitsCuts,fClustersCuts);
752 title->AddText(ptitle);
754 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
755 label->SetBit(kCanDelete);
756 label->SetFillColor(42);
761 //_____________________________________________________________________________
762 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
764 // Draw a view of MUON clusters
765 AliInfo(" Draw View");
767 gPad->SetCursor(kWatch);
768 // gPad->SetFillColor(39);
769 gPad->SetFillColor(1);
771 // gPad->SetFillColor(39);
772 gPad->SetFillColor(1);
776 TView *view = new TView(1);
778 Float_t range = fRrange*fRangeSlider->GetMaximum();
779 view->SetRange(-range,-range,-range,range, range, range);
780 // zoom back to full scale only if DrawView not called from NextCathode
789 Float_t xg1, xg2, yg1, yg2, zg1, zg2;
791 // Recovering the chamber
792 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
794 const AliMUONGeometryTransformer* kGeomTransformer
795 = pMUON->GetGeometryTransformer();
797 AliMUONSegmentation* segmentation = pMUON->GetSegmentation();
799 // Display MUON Chamber Geometry
801 sprintf(nodeName,"MUON%d",100+fChamber);
802 printf(">>>> chamber is %d\n",fChamber);
805 for(Int_t id = 0; id < 4; id++) {
807 Int_t detElemId = fChamber*100+id;
808 AliMpSectorSegmentation * seg =
809 (AliMpSectorSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
810 const AliMpSector * sector = seg->GetSector();
812 // get sector measurements
813 TVector2 position = sector->Position();
814 TVector2 dimension = sector->Dimensions(); // half length
816 Float_t xlocal1 = position.Px(); // mm to cm
817 Float_t ylocal1 = position.Py();
818 Float_t xlocal2 = dimension.Px() * 2.;
819 Float_t ylocal2 = dimension.Px() * 2.;
821 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
822 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
825 TPolyLine3D* poly = new TPolyLine3D();
828 poly->SetPoint(nPoint++, xg1, yg1, 0.);
829 for (Float_t d = 0; d < TMath::Pi()/2.; d+= 0.01) {
830 Float_t x = xg1 + xg2 * TMath::Cos(d);
831 Float_t y = yg1 + yg2 * TMath::Sin(d);
832 poly->SetPoint(nPoint++, x, y, 0.);
834 poly->SetPoint(nPoint++, xg1, yg1, 0.);
836 poly->SetLineColor(2);
844 if(fChamber >4 && fChamber <11) {
846 for(id=0; id<26; id++) {
847 Int_t detElemId = fChamber*100+id;
848 if ( segmentation->HasDE(detElemId) ) {
849 AliMpSlatSegmentation * seg =
850 (AliMpSlatSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
852 seg = (AliMpSlatSegmentation *)
853 AliMUONSegmentationManager::Segmentation(detElemId, kBendingPlane);
856 const AliMpSlat* slat = seg->Slat();
857 Float_t deltax = slat->DX();
858 Float_t deltay = slat->DY();
859 Float_t xlocal1 = -deltax;
860 Float_t ylocal1 = -deltay;
861 Float_t xlocal2 = +deltax;
862 Float_t ylocal2 = +deltay;
863 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
864 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
866 // drawing slat active volumes
867 Float_t xCenter = (xg1 + xg2)/2.;
868 Float_t yCenter = (yg1 + yg2)/2.;
870 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
872 box->SetFillStyle(0);
873 box->SetLineColor(2);
876 // drawing inner circle + disc
877 TPolyLine3D* poly = new TPolyLine3D();
878 TPolyLine3D* poly1 = new TPolyLine3D();
882 for (Float_t d = 0; d < 6.24; d+= 0.005) {
883 Double_t x = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Cos(d)/2.;
884 Double_t y = AliMUONConstants::Dmin((fChamber-1)/2) * TMath::Sin(d)/2.;
885 if (nPoint % 2 == 0) poly->SetPoint(nPoint++, 0., 0., 0.);
886 poly->SetPoint(nPoint++, x, y, 0.);
887 poly1->SetPoint(nPoint1++, x, y, 0.);
890 poly->SetLineColor(1);
892 poly1->SetLineColor(2);
899 /*-- Trigger Chambers ---------------------------------------*/
900 if(fChamber >10 && fChamber <15) {
902 for(id=0; id<18; id++) {
903 Int_t detElemId = fChamber*100+id;
904 AliMpTriggerSegmentation * seg
905 = (AliMpTriggerSegmentation *) segmentation->GetMpSegmentation(detElemId, 0);
907 seg = (AliMpTriggerSegmentation *)
908 AliMUONSegmentationManager::Segmentation(detElemId, kBendingPlane);
911 const AliMpTrigger* slat = seg->Slat();
912 Float_t deltax = slat->DX();
913 Float_t deltay = slat->DY();
914 Float_t xlocal1 = -deltax;
915 Float_t ylocal1 = -deltay;
916 Float_t xlocal2 = +deltax;
917 Float_t ylocal2 = +deltay;
919 kGeomTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg1, yg1, zg1);
920 kGeomTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg2, yg2, zg2);
922 // drawing slat active volumes
923 Float_t xCenter = (xg1 + xg2)/2.;
924 Float_t yCenter = (yg1 + yg2)/2.;
926 TMarker3DBox* box = new TMarker3DBox(xCenter,yCenter,0,xlocal1,ylocal2,0,0,0);
928 box->SetFillStyle(0);
929 box->SetLineColor(4);
934 //add clusters to the pad
938 // DrawSegmentation();
939 // add itself to the list (must be last)
941 view->SetView(phi, theta, psi, iret);
944 //_____________________________________________________________________________
945 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
947 // Draw a view of muons chambers with tracks
949 gPad->SetCursor(kWatch);
950 // gPad->SetFillColor(39);
951 gPad->SetFillColor(1);
953 // gPad->SetFillColor(39);
954 gPad->SetFillColor(1);
958 TView *view = new TView(1);
960 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
961 view->SetRange(-range,-range,-range,range,range,range);
963 // Display all MUON Chambers segmentation
966 sprintf(nodeName,"alice");
968 node1=gAlice->GetGeometry()->GetNode(nodeName);
969 if (node1) node1->Draw("same");
972 // Draw clusters for all chambers
973 Int_t chamberSave = fChamber;
974 for (fChamber = 1; fChamber <= 10; fChamber++){
977 fChamber = chamberSave;
978 // Draw reconstructed tracks
985 Float_t x0 = (-1+shift)/zoom;
986 Float_t y0 = (-1+shift)/zoom;
987 Float_t x1 = (1+shift)/zoom;
988 Float_t y1 = (1+shift)/zoom;
989 gPad->Range(x0,y0,x1,y1);
990 view->SetView(phi, theta, psi, iret);
994 //______________________________________________________________________________
995 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
997 // Execute action corresponding to the mouse event
999 static Float_t x0, y0, x1, y1;
1001 static Int_t pxold, pyold;
1002 static Int_t px0, py0;
1003 static Int_t linedrawn;
1006 if (px == 0 && py == 0) { //when called by sliders
1007 if (event == kButton1Up) {
1012 if (!fZoomMode && gPad->GetView()) {
1013 gPad->GetView()->ExecuteRotateView(event, px, py);
1017 // something to zoom ?
1018 gPad->SetCursor(kCross);
1023 gVirtualX->SetLineColor(-1);
1024 gPad->TAttLine::Modify(); //Change line attributes only if necessary
1025 x0 = gPad->AbsPixeltoX(px);
1026 y0 = gPad->AbsPixeltoY(py);
1028 pxold = px; pyold = py;
1032 case kButton1Motion:
1033 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1037 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
1041 gPad->GetCanvas()->FeedbackMode(kFALSE);
1042 if (px == px0) return;
1043 if (py == py0) return;
1044 x1 = gPad->AbsPixeltoX(px);
1045 y1 = gPad->AbsPixeltoY(py);
1047 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
1048 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
1049 gPad->Range(x0,y0,x1,y1);
1050 if (fZooms < AliMUONConstants::MaxZoom()-1) {
1052 fZoomX0[fZooms] = x0;
1053 fZoomY0[fZooms] = y0;
1054 fZoomX1[fZooms] = x1;
1055 fZoomY1[fZooms] = y1;
1057 gPad->Modified(kTRUE);
1062 //___________________________________________
1063 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
1065 // Read digits info and store x,y,z info in arrays fPoints
1066 // Loop on all detectors
1068 if (chamber > 14) return;
1074 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1075 AliMUONGeometrySegmentation* segmentation2 = 0x0;
1077 GetMUONData()->SetTreeAddress("D");
1079 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
1080 if (muonDigits == 0) return;
1082 gAlice->ResetDigits();
1085 if (GetLoader()->TreeD()) {
1086 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
1087 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
1088 GetMUONData()->GetDigits();
1091 Int_t ndigits = muonDigits->GetEntriesFast();
1092 if (ndigits == 0) return;
1093 if (fPoints == 0) fPoints = new TObjArray(ndigits);
1095 //segmentation2 = iChamber->SegmentationModel2(cathode);
1097 = pMUON->GetSegmentation()->GetModuleSegmentation(chamber-1, cathode-1);
1099 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1102 AliMUONPoints *points = 0;
1103 TMarker3DBox *marker = 0;
1105 //loop over all digits and store their position
1108 Float_t adcmax = 1024; // default
1109 if (chamber<11) adcmax = 4096;
1111 for (Int_t digit = 0; digit < ndigits; digit++) {
1112 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
1113 if (mdig->Cathode() != cathode-1) continue;
1116 // First get all needed parameters
1118 Int_t charge = mdig->Signal();
1119 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
1120 Int_t color = 261+index;
1121 Int_t colorTrigger = 2;
1122 if (color > 282) color = 282;
1124 if (chamber > 10) { // trigger chamber
1126 Int_t sumCharge = 0;
1127 for (Int_t icharge = 0; icharge < 10; icharge++) {
1128 sumCharge = sumCharge+mdig->TrackCharge(icharge);
1130 Int_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
1131 if(sumCharge <= 10 || testCharge > 0) {
1132 colorTrigger = color;
1138 // get the center of the pad - add on x and y half of pad size
1139 Float_t xpad, ypad, zpad;
1143 Int_t detElemId = mdig->DetElemId();
1144 segmentation2->GetPadC(detElemId, mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
1145 isec = segmentation2->Sector(detElemId, mdig->PadX(), mdig->PadY());
1146 dpx = segmentation2->Dpx(detElemId, isec)/2;
1147 dpy = segmentation2->Dpy(detElemId, isec)/2;
1149 // segmentation->Dump();
1152 // Then set the objects
1154 points = new AliMUONPoints(npoints);
1155 fPoints->AddAt(points,digit);
1157 points->SetMarkerColor(colorTrigger);
1159 points->SetMarkerColor(color);
1161 points->SetMarkerStyle(21);
1162 points->SetMarkerSize(0.5);
1163 points->SetParticle(-1);
1164 points->SetHitIndex(-1);
1165 points->SetTrackIndex(-1);
1166 points->SetDigitIndex(digit);
1167 points->SetPoint(0,xpad,ypad,zpos);
1169 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1170 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1173 marker->SetLineColor(lineColor);
1174 marker->SetFillStyle(1001);
1175 marker->SetFillColor(color);
1176 marker->SetRefObject((TObject*)points);
1177 points->Set3DMarker(0, marker);
1180 //___________________________________________
1181 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1183 // Read raw clusters info and store x,y,z info in arrays fRpoints
1184 // Loop on all detectors
1186 if (chamber > 10) return;
1190 GetMUONData()->SetTreeAddress("RC");
1191 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1193 if (muonRawClusters == 0) return;
1196 if (GetMUONData()->TreeR()) {
1197 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1198 GetMUONData()->TreeR()->GetEvent(0);
1201 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1202 if (nrawcl == 0) return;
1203 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1205 Float_t zpos = AliMUONConstants::DefaultChamberZ(chamber-1);
1206 AliMUONRawCluster *mRaw;
1207 AliMUONPoints *points = 0;
1209 //loop over all raw clusters and store their position
1210 points = new AliMUONPoints(nrawcl);
1211 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1212 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1213 points->SetMarkerColor(51);
1214 points->SetMarkerStyle(2);
1215 points->SetMarkerSize(1.);
1216 points->SetParticle(-1);
1217 points->SetHitIndex(-1);
1218 points->SetTrackIndex(-1);
1219 points->SetDigitIndex(-1);
1220 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1221 fRpoints->AddAt(points,iraw);
1222 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1226 //___________________________________________
1227 void AliMUONDisplay::LoadTracks()
1230 AliMUONTrack* recTrack = 0;
1231 AliMUONTrackParam* trackParam = 0;
1232 TClonesArray * trackParamAtHit = 0;
1236 GetMUONData()->SetTreeAddress("RT");
1237 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1238 if (recTracksArray == NULL) return;
1239 GetMUONData()->GetRecTracks();
1241 Int_t nRecTracks = 0;
1243 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1246 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1248 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1249 // reading info from tracks
1250 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1252 Int_t nTrackHits = recTrack->GetNTrackHits();
1254 if (nTrackHits == 0) continue;
1257 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1258 points->SetLineColor(6);
1259 points->SetLineWidth(1);
1260 fRpoints->AddAt(points,iRecTracks);
1266 trackParam = recTrack->GetTrackParamAtVertex();
1267 xRec = trackParam->GetNonBendingCoor();
1268 yRec = trackParam->GetBendingCoor();
1269 zRec = trackParam->GetZ();
1270 points->SetPoint(iPoint,xRec,yRec,zRec);
1273 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1274 trackParamAtHit = recTrack->GetTrackParamAtHit();
1275 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1276 xRec = trackParam->GetNonBendingCoor();
1277 yRec = trackParam->GetBendingCoor();
1278 zRec = trackParam->GetZ();
1279 points->SetPoint(iPoint,xRec,yRec,zRec);
1281 } // end loop rec. hits
1282 PrintTrack(iRecTracks,recTrack);
1283 } // end loop tracks
1288 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1290 AliMUONTrackParam *trackParam;
1291 Float_t vertex[3], momentum[3];
1292 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1295 trackParam = recTrack->GetTrackParamAtVertex();
1296 vertex[0] = trackParam->GetNonBendingCoor();
1297 vertex[1] = trackParam->GetBendingCoor();
1298 vertex[2] = trackParam->GetZ();
1299 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1300 bendingSlope = trackParam->GetBendingSlope();
1301 nonBendingSlope = trackParam->GetNonBendingSlope();
1302 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1303 momentum[0] = momentum[2] * nonBendingSlope;
1304 momentum[1] = momentum[2] * bendingSlope;
1305 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1306 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1308 printf("===================================================\n");
1309 printf(" Reconstructed track # %d \n",iRecTracks);
1310 printf(" charge: %d \n",charge);
1311 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1312 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1313 printf(" track chi2/dof: %f \n",chi2dof);
1317 //___________________________________________
1318 void AliMUONDisplay::LoadHits(Int_t chamber)
1320 // Read hits info and store x,y,z info in arrays fPhits
1321 // Loop on all detectors
1323 if (chamber > 14) return;
1330 Float_t zpos=AliMUONConstants::DefaultChamberZ(chamber-1);
1332 if (GetMUONData()->TreeH()) {
1333 GetMUONData()->SetTreeAddress("H");
1334 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
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 nthits += muonHits->GetEntriesFast();
1343 if (fPhits == 0) fPhits = new TObjArray(nthits);
1345 for (track=0; track<ntracks;track++) {
1346 GetMUONData()->ResetHits();
1347 GetMUONData()->GetTrack(track);//skowron
1348 TClonesArray *muonHits = GetMUONData()->Hits();
1349 if (muonHits == 0) return;
1350 Int_t nhits = muonHits->GetEntriesFast();
1351 if (nhits == 0) continue;
1353 AliMUONPoints *points = 0;
1355 for (Int_t hit=0;hit<nhits;hit++) {
1356 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1357 Int_t nch = mHit->Chamber(); // chamber number
1358 if (nch != chamber) continue;
1360 // Retrieve info and set the objects
1362 points = new AliMUONPoints(npoints);
1363 fPhits->AddAt(points,nhold+hit);
1364 points->SetMarkerColor(kRed);
1365 points->SetMarkerStyle(5);
1366 points->SetMarkerSize(1.);
1367 points->SetParticle(mHit->Track());
1368 points->SetHitIndex(hit);
1369 points->SetTrackIndex(track);
1370 points->SetDigitIndex(-1);
1371 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1372 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1379 //_____________________________________________________________________________
1380 void AliMUONDisplay::Paint(Option_t *)
1382 // Paint miscellaneous items
1385 //_____________________________________________________________________________
1386 void AliMUONDisplay::SetPickMode()
1388 // Set parameters for pick mode.
1392 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1393 fTrigPad->Modified();
1396 //_____________________________________________________________________________
1397 void AliMUONDisplay::SetZoomMode()
1399 // Set parameters for zoom mode
1402 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1403 fTrigPad->Modified();
1406 //_____________________________________________________________________________
1407 void AliMUONDisplay::NextChamber(Int_t delta)
1409 // to go from chamber to next chamber if delta = 1
1410 // or previous chamber otherwise
1412 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1414 if (fChamber > 1) fChamber--;
1418 LoadDigits(fChamber, fCathode);
1422 //_____________________________________________________________________________
1423 void AliMUONDisplay::NextCathode()
1425 // to switch to other cathode plane
1428 if (fCathode == 1) {
1429 LoadDigits(fChamber, 2);
1431 LoadDigits(fChamber, 1);
1433 fNextCathode = kTRUE; // to keep the same zoom
1435 fNextCathode = kFALSE;
1436 TPad *pad = (TPad*)gPad->GetPadSave();
1437 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1438 fZoomX1[fZooms], fZoomY1[fZooms]);
1444 //_____________________________________________________________________________
1445 void AliMUONDisplay::Trigger()
1448 AliMUONGlobalTrigger* globalTrig;
1450 GetMUONData()->SetTreeAddress("GLT");
1451 GetMUONData()->GetTriggerD();
1453 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1454 if (globalTrig == 0) return;
1456 printf("===================================================\n");
1457 printf(" Global Trigger output Low pt High pt All\n");
1459 printf(" number of Single Plus :\t");
1460 printf("%i\t",globalTrig->SinglePlusLpt());
1461 printf("%i\t",globalTrig->SinglePlusHpt());
1462 printf("%i\t",globalTrig->SinglePlusApt());
1465 printf(" number of Single Minus :\t");
1466 printf("%i\t",globalTrig->SingleMinusLpt());
1467 printf("%i\t",globalTrig->SingleMinusHpt());
1468 printf("%i\t",globalTrig->SingleMinusApt());
1471 printf(" number of Single Undefined :\t");
1472 printf("%i\t",globalTrig->SingleUndefLpt());
1473 printf("%i\t",globalTrig->SingleUndefHpt());
1474 printf("%i\t",globalTrig->SingleUndefApt());
1477 printf(" number of UnlikeSign pair :\t");
1478 printf("%i\t",globalTrig->PairUnlikeLpt());
1479 printf("%i\t",globalTrig->PairUnlikeHpt());
1480 printf("%i\t",globalTrig->PairUnlikeApt());
1483 printf(" number of LikeSign pair :\t");
1484 printf("%i\t",globalTrig->PairLikeLpt());
1485 printf("%i\t",globalTrig->PairLikeHpt());
1486 printf("%i\t",globalTrig->PairLikeApt());
1488 printf("===================================================\n");
1492 // // returns Trigger Decision for current event
1493 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1495 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1496 // AliMUONData* muonData = decision->GetMUONData();
1497 // muonData->SetTreeAddress("D");
1498 // decision->Trigger();
1500 //_____________________________________________________________________________
1501 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1503 // Set chamber and cathode number
1509 LoadDigits(chamber,cathode);
1513 void AliMUONDisplay::SetEvent(Int_t newevent)
1516 gAlice->GetEvent(newevent);
1518 if (!gAlice->TreeD()) return;
1521 LoadDigits(fChamber,fCathode);
1525 //_____________________________________________________________________________
1526 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1528 // Set view range along R and Z
1537 //_____________________________________________________________________________
1538 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1540 // change viewing angles for current event
1548 TView *view = gPad->GetView();
1549 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1554 //_____________________________________________________________________________
1555 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1557 AliRunLoader * runLoader;
1559 runLoader = fLoader->GetRunLoader();
1563 // Display (current event_number + delta)
1564 // delta = 1 shown next event
1565 // delta = -1 show previous event
1567 //runLoader->CleanDetectors();
1568 //runLoader->CleanKinematics();
1569 Int_t currentEvent = runLoader->GetEventNumber();
1570 Int_t newEvent = currentEvent + delta;
1571 runLoader->GetEvent(newEvent);
1574 LoadDigits(fChamber, fCathode);
1579 //______________________________________________________________________________
1580 void AliMUONDisplay::UnZoom()
1583 if (fZooms <= 0) return;
1585 TPad *pad = (TPad*)gPad->GetPadSave();
1586 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1590 //_____________________________________________________________________________
1591 void AliMUONDisplay::ResetPoints()
1594 // Reset array of points
1602 //_____________________________________________________________________________
1603 void AliMUONDisplay::ResetPhits()
1606 // Reset array of points
1614 //_____________________________________________________________________________
1615 void AliMUONDisplay::ResetRpoints()
1618 // Reset array of points
1627 AliMUONDisplay & AliMUONDisplay::operator = (const AliMUONDisplay & rhs)
1629 // Protected assignement operator
1631 if (this == &rhs) return *this;
1633 AliFatal("Not implemented.");