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>
44 #include "AliMUONDisplay.h"
47 #include "AliMUONPoints.h"
48 #include "AliMUONGlobalTrigger.h"
49 #include "AliHeader.h"
51 #include "AliMUONHit.h"
52 #include "AliMUONDigit.h"
53 #include "AliMUONRawCluster.h"
54 #include "AliMUONTrack.h"
55 #include "AliMUONTrackParam.h"
57 #include "AliSegmentation.h"
58 #include "AliMUONResponse.h"
59 #include "AliMUONChamber.h"
60 #include "AliMUONConstants.h"
62 // to manage the same zoom on both cathodes
66 ClassImp(AliMUONDisplay)
69 //_____________________________________________________________________________
70 AliMUONDisplay::AliMUONDisplay()
78 fNextCathode = kFALSE;
82 //_____________________________________________________________________________
83 AliMUONDisplay::AliMUONDisplay(Int_t size, AliLoader * loader)
86 // Create an event display object.
87 // A canvas named "edisplay" is created with a vertical size in pixels
89 // A QUICK Overview of the Event Display functions
90 // ===============================================
92 // The event display can ve invoked by executing the macro "display.C"
93 // A canvas like in the picture below will appear.
95 // On the left side of the canvas, the following buttons appear:
96 // *Next* to move to the next event
97 // *Previous* to move to the previous event
99 // *Pick* Select this option to be able to point on a track with the
100 // mouse. Once on the track, use the right button to select
101 // an action. For example, select SetMarkerAttributes to
102 // change the marker type/color/size for the track.
103 // *Zoom* Select this option (default) if you want to zoom.
104 // To zoom, simply select the selected area with the left button.
105 // *UnZoom* To revert to the previous picture size.
107 // slider R On the left side, the vertical slider can be used to
108 // set the default picture size.
110 // When you are in Zoom mode, you can click on the black part of the canvas
111 // to select special options with the right mouse button.
114 // When you are in pick mode, you can "Inspect" the object pointed by the mouse.
115 // When you are on a track, select the menu item "InspectParticle"
116 // to display the current particle attributes.
118 // You can activate the Root browser by selecting the Inspect menu
119 // in the canvas tool bar menu. Then select "Start Browser"
120 // This will open a new canvas with the browser. At this point, you may want
121 // to display some histograms (from the Trees). Go to the "File" menu
122 // of the browser and click on "New canvas".
123 // In the browser, click on item "ROOT files" in the left pane.
124 // Click on galice.root.
126 // Click on TPC for example
127 // Click on any variable (eg TPC.fX) to histogram the variable.
129 // If you are lost, you can click on HELP in any Root canvas or browser.
132 <img src="gif/AliMUONDisplay.gif">
139 gAlice->SetDisplay(this);
141 // Initialize display default parameters
143 // Set front view by default
150 fDrawClusters = kTRUE;
152 fDrawTracks = kFALSE;
161 // Create display canvas
163 if (ysize < 100) ysize = 750;
164 Int_t xsize = Int_t(size*830./ysize);
165 fCanvas = new TCanvas("Canvas", "MUON Display",14,47,xsize,ysize);
166 fCanvas->ToggleEventStatus();
168 // Create main display pad
169 fPad = new TPad("viewpad", "MUON display",0.15,0,0.9,1);
172 fPad->SetFillColor(30);
173 fPad->SetBorderSize(2);
178 fColPad = new TPad("colpad", "Colors pad",0.9,0,1,1);
180 fColPad->SetFillColor(17);
181 fColPad->SetBorderSize(2);
186 // Create user interface control pad
190 // Create Range and mode pad
193 fTrigPad = new TPad("trigger", "range and mode pad",0,0,dxtr,dytr);
196 fTrigPad->SetFillColor(22);
197 fTrigPad->SetBorderSize(2);
198 fRangeSlider = new TSlider("range","range",0.7,0.42,0.9,0.98);
199 fRangeSlider->SetObject(this);
200 char pickmode[] = "gAlice->Display()->SetPickMode()";
202 fPickButton = new TButton("Pick",pickmode,0.05,0.32,0.65,0.32+db);
203 fPickButton->SetFillColor(38);
205 char zoommode[] = "gAlice->Display()->SetZoomMode()";
206 fZoomButton = new TButton("Zoom",zoommode,0.05,0.21,0.65,0.21+db);
207 fZoomButton->SetFillColor(38);
209 fArcButton = new TArc(.8,fZoomButton->GetYlowNDC()+0.5*db,0.33*db);
210 fArcButton->SetFillColor(kGreen);
212 char butUnzoom[] = "gAlice->Display()->UnZoom()";
213 TButton *button = new TButton("UnZoom",butUnzoom,0.05,0.05,0.95,0.15);
214 button->SetFillColor(38);
216 AppendPad(); // append display object as last object to force selection
219 fTrigPad->SetEditable(kFALSE);
220 fButtons->SetEditable(kFALSE);
222 fNextCathode = kFALSE;
224 // initialize container
226 fMUONData = new AliMUONData(fLoader,"MUON","MUON");
231 AliMUONDisplay::AliMUONDisplay(const AliMUONDisplay & display)
232 : AliDisplay(display)
234 // Protected copy constructor
236 Fatal("AliMUONDisplay", "Not implemented.");
241 //_____________________________________________________________________________
242 AliMUONDisplay::~AliMUONDisplay()
244 // Delete space point structure
245 if (fPoints) fPoints->Delete();
249 if (fPhits) fPhits->Delete();
253 if (fRpoints) fRpoints->Delete();
258 //_____________________________________________________________________________
259 void AliMUONDisplay::Clear(Option_t *)
261 // Delete graphics temporary objects
264 //_____________________________________________________________________________
265 void AliMUONDisplay::DisplayButtons()
267 // Create the user interface buttons
270 fButtons = new TPad("buttons", "newpad",0,0.45,0.15,1);
272 fButtons->SetFillColor(38);
273 fButtons->SetBorderSize(2);
276 // Int_t butcolor = 33;
277 Float_t dbutton = 0.08;
284 char but1[] = "gAlice->Display()->ShowNextEvent(1)";
285 button = new TButton("Event +", but1, x0, y - dbutton, x1, y);
286 button->SetFillColor(38);
290 char but2[] = "gAlice->Display()->ShowNextEvent(-1)";
291 button = new TButton("Event -", but2, x0, y - dbutton, x1, y);
292 button->SetFillColor(38);
296 char but3[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(1)";
297 button = new TButton("Chamber +", but3, x0, y - dbutton, x1, y);
298 button->SetFillColor(38);
302 char but4[] = "((AliMUONDisplay*)(gAlice->Display()))->NextChamber(-1)";
303 button = new TButton("Chamber -", but4, x0, y - dbutton, x1, y);
304 button->SetFillColor(38);
308 char but5[] = "((AliMUONDisplay*)(gAlice->Display()))->SetChamberAndCathode(1,1)";
309 button = new TButton("Chamber 1", but5, x0, y - dbutton, x1, y);
310 button->SetFillColor(38);
314 char but6[] = "((AliMUONDisplay*)(gAlice->Display()))->NextCathode()";
315 button = new TButton("Cathode <>", but6, x0, y - dbutton, x1, y);
316 button->SetFillColor(38);
320 char but7[] = "((AliMUONDisplay*)(gAlice->Display()))->Trigger()";
321 button = new TButton("Trigger", but7, x0, y - dbutton, x1, y);
322 button->SetFillColor(38);
326 char but8[] = "((AliMUONDisplay*)(gAlice->Display()))->DrawReco()";
327 button = new TButton("Tracking", but8, x0, y - dbutton, x1, y);
328 button->SetFillColor(38);
332 TDiamond *diamond = new TDiamond(0.05,0.015,0.95,0.22);
333 diamond->SetFillColor(50);
334 diamond->SetTextAlign(22);
335 diamond->SetTextColor(5);
336 diamond->SetTextSize(0.11);
338 diamond->AddText(".. ");
339 diamond->AddText("ROOT");
340 diamond->AddText("MUON");
341 diamond->AddText("... ");
342 diamond->AddText(" ");
345 //_____________________________________________________________________________
346 void AliMUONDisplay::CreateColors() const
348 // Create the colors palette used to display clusters
363 new TColor(color,r,g,b);
373 new TColor(color,r,g,b);
383 new TColor(color,r,g,b);
393 new TColor(color,r,g,b);
403 new TColor(color,r,g,b);
410 //_____________________________________________________________________________
411 void AliMUONDisplay::DisplayColorScale()
413 // Display pulse height color scale
416 Float_t xlow, ylow, xup, yup, hs;
417 Float_t x1, y1, x2, y2;
421 TText *text = new TText(0,0,"");
422 text->SetTextFont(61);
423 text->SetTextSize(0.2);
424 text->SetTextAlign(22);
426 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
427 AliMUONChamber *iChamber = &(pMUON->Chamber(fChamber-1));
428 AliMUONResponse * response=iChamber->ResponseModel();
430 if (response) adcmax = (Int_t) response->MaxAdc();
435 //*-* draw colortable boxes
436 hs = (y2-y1)/Float_t(22);
440 ylow = y1 + hs*(Float_t(i));
441 yup = y1 + hs*(Float_t(i+1));
443 Double_t logscale=Double_t(i+1)*(TMath::Log(adcmax)/22);
444 Int_t scale=(Int_t)TMath::Exp(logscale);
445 sprintf(label,"%d",scale);
446 box = new TBox(xlow, ylow, xup, yup);
448 box->SetFillColor(color);
449 text->DrawText(xlow+0.7, 0.5*(ylow+yup),label);
453 //______________________________________________________________________________
454 Int_t AliMUONDisplay::DistancetoPrimitive(Int_t px, Int_t)
456 // Compute distance from point px,py to objects in event
458 gPad->SetCursor(kCross);
460 if (gPad == fTrigPad) return 9999;
462 const Int_t kBig = 9999;
464 Float_t xmin = gPad->GetX1();
465 Float_t xmax = gPad->GetX2();
466 Float_t dx = 0.02*(xmax - xmin);
467 Float_t x = gPad->AbsPixeltoX(px);
468 if (x < xmin+dx || x > xmax-dx) return dist;
470 if (fZoomMode) return 0;
474 //_____________________________________________________________________________
475 void AliMUONDisplay::Draw(Option_t *)
477 // Display current event
485 //_____________________________________________________________________________
486 void AliMUONDisplay::DrawChamber()
489 fDrawTracks = kFALSE;
491 DrawView(fTheta, fPhi, fPsi);
492 // Display the event number and title
497 //_____________________________________________________________________________
498 void AliMUONDisplay::DrawReco(Option_t *)
500 // Display current event
503 // print kinematics of generated particles
505 // Draw global view of muon system
507 DrawGlobalView(135, -50, -140);
509 // Display the event number and title
514 //_____________________________________________________________________________
515 void AliMUONDisplay::PrintKinematics()
517 AliRunLoader * runLoader;
518 TParticle *particle = new TParticle();
520 Float_t vertex[3], momentum[3];
523 runLoader = fLoader->GetRunLoader();
527 printf("****** Event # %d ******\n",runLoader->GetEventNumber());
528 runLoader->TreeK()->GetBranch("Particles")->SetAddress(&particle);
529 nPart = (Int_t)runLoader->TreeK()->GetEntries();
530 for(Int_t iPart = 0; iPart < nPart; iPart++) {
531 runLoader->TreeK()->GetEvent(iPart);
532 vertex[0] = particle->Vx();
533 vertex[1] = particle->Vy();
534 vertex[2] = particle->Vz();
535 momentum[0] = particle->Px();
536 momentum[1] = particle->Py();
537 momentum[2] = particle->Pz();
539 printf("===================================================\n");
540 printf(" Generated particle # %d \n",iPart);
541 printf(" name: %s \n",particle->GetName());
542 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
543 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
548 void AliMUONDisplay::DrawSegmentation()
550 // Draw graphical representation of segmenatation
551 // Attention: still experimental code
554 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
555 AliMUONChamber* iChamber;
556 AliSegmentation* seg;
557 iChamber = &(pMUON->Chamber(fChamber));
558 seg=iChamber->SegmentationModel(icat);
559 Float_t zpos=iChamber->Z();
560 Float_t r=iChamber->ROuter();
562 TMarker3DBox *marker;
564 for (Int_t j=0; j<seg->Npy(); j++) {
566 y0=j*seg->Dpy()-seg->Dpy()/2.;
567 for (seg->FirstPad(0.,y0,0,300,0.);
571 if (seg->ISector()==0) continue;
573 seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
574 Float_t dpx=seg->Dpx(seg->ISector())/2;
575 Float_t dpy=seg->Dpy(seg->ISector())/2;
576 marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
577 marker->SetLineColor(seg->ISector()+1);
578 marker->SetFillStyle(1001);
579 marker->SetFillColor(0);
584 for (Int_t j=0; j<250; j++) {
585 Float_t x0=j*seg->Dpx();
586 Float_t y0=TMath::Sqrt(r*r-x0*x0);
588 for (seg->FirstPad(x0,0,0,0,y0);
592 if (seg->ISector()==0) continue;
595 seg->GetPadC(seg->Ix(), seg->Iy(), x, y, z);
596 Float_t dpx=seg->Dpx(seg->ISector())/2;
597 Float_t dpy=seg->Dpy(seg->ISector())/2;
598 marker=new TMarker3DBox(x,y,zpos,dpx,dpy,0,0,0);
599 marker->SetLineColor(seg->ISector()+1);
600 marker->SetFillStyle(1001);
601 marker->SetFillColor(0);
608 //_____________________________________________________________________________
609 void AliMUONDisplay::DrawClusters()
611 // Draw clusters for MUON chambers
613 Int_t ndigits, digit;
621 ndigits = points->GetEntriesFast();
622 for (digit=0;digit<ndigits;digit++){
623 pm = (AliMUONPoints*)points->UncheckedAt(digit);
627 for (Int_t im=0;im<3;im++) {
628 TMarker3DBox *marker=pm->GetMarker(im);
633 fClustersCuts +=pm->GetN();
637 //_____________________________________________________________________________
638 void AliMUONDisplay::DrawHits()
640 // Draw hits for MUON chambers
644 Int_t ntracks, track;
651 ntracks = points->GetEntriesFast();
652 for (track=0;track<ntracks;track++) {
653 pm = (AliMUONPoints*)points->UncheckedAt(track);
656 fHitsCuts += pm->GetN();
661 //_____________________________________________________________________________
662 void AliMUONDisplay::DrawCoG()
664 // Draw hits for MUON chambers
665 if (!fDrawCoG) return;
666 if (fChamber > 10) return;
667 LoadCoG(fChamber,fCathode);
675 ncog = points->GetEntriesFast();
676 for (icog=0;icog<ncog;icog++) {
677 pm = (AliMUONPoints*)points->UncheckedAt(icog);
682 //_____________________________________________________________________________
683 void AliMUONDisplay::DrawTracks()
686 if (!fDrawTracks) return;
689 Int_t nTrack, iTrack;
695 nTrack = points->GetEntriesFast();
696 for ( iTrack = 0; iTrack < nTrack; iTrack++) {
697 pm = (TPolyLine3D*)points->UncheckedAt(iTrack);
702 //_____________________________________________________________________________
704 void AliMUONDisplay::DrawTitle(Option_t *option)
706 // Draw the event title
708 Float_t xmin = gPad->GetX1();
709 Float_t xmax = gPad->GetX2();
710 Float_t ymin = gPad->GetY1();
711 Float_t ymax = gPad->GetY2();
712 Float_t dx = xmax-xmin;
713 Float_t dy = ymax-ymin;
715 AliRunLoader * runLoader;
717 runLoader = fLoader->GetRunLoader();
722 if (strlen(option) == 0) {
723 TPaveText *title = new TPaveText(xmin +0.01*dx, ymax-0.09*dy, xmin +0.5*dx, ymax-0.01*dy);
724 // title->SetTextSize(0.023932);
725 title->SetTextSize(0.02);
726 title->SetBit(kCanDelete);
727 title->SetFillColor(42);
730 sprintf(ptitle, "Alice event:%d Run:%d Chamber:%d Cathode:%d",
731 runLoader->GetEventNumber(),
732 gAlice->GetHeader()->GetRun(),
735 title->AddText(ptitle);
736 Int_t nparticles = gAlice->GetMCApp()->Particles()->GetEntriesFast();
737 sprintf(ptitle,"Nparticles = %d Nhits = %d Npads fired = %d",
738 nparticles, fHitsCuts,fClustersCuts);
739 title->AddText(ptitle);
741 TPaveLabel *label = new TPaveLabel(xmin +0.01*dx, ymax-0.07*dy, xmin +0.2*dx, ymax-0.01*dy,option);
742 label->SetBit(kCanDelete);
743 label->SetFillColor(42);
748 //_____________________________________________________________________________
749 void AliMUONDisplay::DrawView(Float_t theta, Float_t phi, Float_t psi)
751 // Draw a view of MUON clusters
752 printf("\n Draw View\n");
754 gPad->SetCursor(kWatch);
755 // gPad->SetFillColor(39);
756 gPad->SetFillColor(1);
758 // gPad->SetFillColor(39);
759 gPad->SetFillColor(1);
763 TView *view = new TView(1);
765 Float_t range = fRrange*fRangeSlider->GetMaximum();
766 view->SetRange(-range,-range,-range,range, range, range);
767 // zoom back to full scale only if DrawView not called from NextCathode
776 // Display MUON Chamber Geometry
778 sprintf(nodeName,"MUON%d",100+fChamber);
780 TNode *node1=gAlice->GetGeometry()->GetNode(nodeName);
781 if (node1) node1->Draw("same");
782 //add clusters to the pad
786 // DrawSegmentation();
787 // add itself to the list (must be last)
789 view->SetView(phi, theta, psi, iret);
792 //_____________________________________________________________________________
793 void AliMUONDisplay::DrawGlobalView(Float_t theta, Float_t phi, Float_t psi)
795 // Draw a view of muons chambers with tracks
797 gPad->SetCursor(kWatch);
798 // gPad->SetFillColor(39);
799 gPad->SetFillColor(1);
801 // gPad->SetFillColor(39);
802 gPad->SetFillColor(1);
806 TView *view = new TView(1);
808 Float_t range = fRrange*fRangeSlider->GetMaximum()*3.;
809 view->SetRange(-range,-range,-range,range,range,range);
811 // Display all MUON Chambers segmentation
814 sprintf(nodeName,"alice");
816 node1=gAlice->GetGeometry()->GetNode(nodeName);
817 if (node1) node1->Draw("same");
820 // Draw clusters for all chambers
821 Int_t chamberSave = fChamber;
822 for (fChamber = 1; fChamber <= 10; fChamber++){
825 fChamber = chamberSave;
826 // Draw reconstructed tracks
833 Float_t x0 = (-1+shift)/zoom;
834 Float_t y0 = (-1+shift)/zoom;
835 Float_t x1 = (1+shift)/zoom;
836 Float_t y1 = (1+shift)/zoom;
837 gPad->Range(x0,y0,x1,y1);
838 view->SetView(phi, theta, psi, iret);
842 //______________________________________________________________________________
843 void AliMUONDisplay::ExecuteEvent(Int_t event, Int_t px, Int_t py)
845 // Execute action corresponding to the mouse event
847 static Float_t x0, y0, x1, y1;
849 static Int_t pxold, pyold;
850 static Int_t px0, py0;
851 static Int_t linedrawn;
854 if (px == 0 && py == 0) { //when called by sliders
855 if (event == kButton1Up) {
860 if (!fZoomMode && gPad->GetView()) {
861 gPad->GetView()->ExecuteRotateView(event, px, py);
865 // something to zoom ?
866 gPad->SetCursor(kCross);
871 gVirtualX->SetLineColor(-1);
872 gPad->TAttLine::Modify(); //Change line attributes only if necessary
873 x0 = gPad->AbsPixeltoX(px);
874 y0 = gPad->AbsPixeltoY(py);
876 pxold = px; pyold = py;
881 if (linedrawn) gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
885 gVirtualX->DrawBox(px0, py0, pxold, pyold, TVirtualX::kHollow);
889 gPad->GetCanvas()->FeedbackMode(kFALSE);
890 if (px == px0) return;
891 if (py == py0) return;
892 x1 = gPad->AbsPixeltoX(px);
893 y1 = gPad->AbsPixeltoY(py);
895 if (x1 < x0) {temp = x0; x0 = x1; x1 = temp;}
896 if (y1 < y0) {temp = y0; y0 = y1; y1 = temp;}
897 gPad->Range(x0,y0,x1,y1);
898 if (fZooms < AliMUONConstants::MaxZoom()-1) {
900 fZoomX0[fZooms] = x0;
901 fZoomY0[fZooms] = y0;
902 fZoomX1[fZooms] = x1;
903 fZoomY1[fZooms] = y1;
905 gPad->Modified(kTRUE);
910 //___________________________________________
911 void AliMUONDisplay::LoadDigits(Int_t chamber, Int_t cathode)
913 // Read digits info and store x,y,z info in arrays fPoints
914 // Loop on all detectors
916 if (chamber > 14) return;
922 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
923 AliMUONChamber* iChamber;
924 AliSegmentation* segmentation;
925 AliMUONResponse* response;
927 GetMUONData()->SetTreeAddress("D");
929 TClonesArray *muonDigits = GetMUONData()->Digits(chamber-1);
930 if (muonDigits == 0) return;
932 gAlice->ResetDigits();
935 if (GetLoader()->TreeD()) {
936 nent = (Int_t) GetLoader()->TreeD()->GetEntries();
937 // gAlice->TreeD()->GetEvent(nent-2+cathode-1);
938 GetMUONData()->GetCathode(cathode-1);
941 Int_t ndigits = muonDigits->GetEntriesFast();
942 if (ndigits == 0) return;
943 if (fPoints == 0) fPoints = new TObjArray(ndigits);
945 iChamber = &(pMUON->Chamber(chamber-1));
947 segmentation = iChamber->SegmentationModel(cathode);
948 response = iChamber->ResponseModel();
949 Float_t zpos = iChamber->Z();
952 AliMUONPoints *points = 0;
953 TMarker3DBox *marker = 0;
955 //loop over all digits and store their position
958 Float_t adcmax = 1024;
959 if (response&&chamber<11) adcmax = response->MaxAdc();
961 for (Int_t digit = 0; digit < ndigits; digit++) {
962 mdig = (AliMUONDigit*)muonDigits->UncheckedAt(digit);
963 if (mdig->Cathode() != cathode-1) continue;
966 // First get all needed parameters
968 Int_t charge = mdig->Signal();
969 Int_t index = Int_t(TMath::Log(charge)/(TMath::Log(adcmax)/22));
970 Int_t color = 261+index;
971 Int_t colorTrigger = 2;
972 if (color > 282) color = 282;
974 if (chamber > 10) { // trigger chamber
977 for (Int_t icharge = 0; icharge < 10; icharge++) {
978 sumCharge = sumCharge+mdig->TrackCharge(icharge);
980 Int_t testCharge = sumCharge-(Int_t(sumCharge/10))*10;
981 if(sumCharge <= 10 || testCharge > 0) {
982 colorTrigger = color;
988 // get the center of the pad - add on x and y half of pad size
989 Float_t xpad, ypad, zpad;
990 segmentation->GetPadC(mdig->PadX(), mdig->PadY(), xpad, ypad, zpad);
992 Int_t isec = segmentation->Sector(mdig->PadX(), mdig->PadY());
993 Float_t dpx = segmentation->Dpx(isec)/2;
994 Float_t dpy = segmentation->Dpy(isec)/2;
996 // segmentation->Dump();
999 // Then set the objects
1001 points = new AliMUONPoints(npoints);
1002 fPoints->AddAt(points,digit);
1004 points->SetMarkerColor(colorTrigger);
1006 points->SetMarkerColor(color);
1008 points->SetMarkerStyle(21);
1009 points->SetMarkerSize(0.5);
1010 points->SetParticle(-1);
1011 points->SetHitIndex(-1);
1012 points->SetTrackIndex(-1);
1013 points->SetDigitIndex(digit);
1014 points->SetPoint(0,xpad,ypad,zpos);
1016 Int_t lineColor = (zpad-zpos > 0) ? 2:3;
1017 marker=new TMarker3DBox(xpad,ypad,zpos,dpx,dpy,0,0,0);
1020 marker->SetLineColor(lineColor);
1021 marker->SetFillStyle(1001);
1022 marker->SetFillColor(color);
1023 marker->SetRefObject((TObject*)points);
1024 points->Set3DMarker(0, marker);
1027 //___________________________________________
1028 void AliMUONDisplay::LoadCoG(Int_t chamber, Int_t /*cathode*/)
1030 // Read raw clusters info and store x,y,z info in arrays fRpoints
1031 // Loop on all detectors
1033 if (chamber > 10) return;
1037 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1038 AliMUONChamber* iChamber;
1040 GetMUONData()->SetTreeAddress("RC");
1041 TClonesArray *muonRawClusters = GetMUONData()->RawClusters(chamber-1);
1043 if (muonRawClusters == 0) return;
1046 if (GetMUONData()->TreeR()) {
1047 nent=(Int_t) GetMUONData()->TreeR()->GetEntries();
1048 GetMUONData()->TreeR()->GetEvent(0);
1051 Int_t nrawcl = muonRawClusters->GetEntriesFast();
1052 if (nrawcl == 0) return;
1053 if (fRpoints == 0) fRpoints = new TObjArray(nrawcl);
1055 iChamber = &(pMUON->Chamber(chamber-1));
1056 Float_t zpos=iChamber->Z();
1057 AliMUONRawCluster *mRaw;
1058 AliMUONPoints *points = 0;
1060 //loop over all raw clusters and store their position
1061 points = new AliMUONPoints(nrawcl);
1062 for (Int_t iraw=0;iraw<nrawcl;iraw++) {
1063 mRaw = (AliMUONRawCluster*)muonRawClusters->UncheckedAt(iraw);
1064 points->SetMarkerColor(51);
1065 points->SetMarkerStyle(2);
1066 points->SetMarkerSize(1.);
1067 points->SetParticle(-1);
1068 points->SetHitIndex(-1);
1069 points->SetTrackIndex(-1);
1070 points->SetDigitIndex(-1);
1071 points->SetPoint(iraw,mRaw->GetX(0),mRaw->GetY(0),zpos);
1072 fRpoints->AddAt(points,iraw);
1073 // printf("%f and %f and %f\n",mRaw->GetX(0),mRaw->GetY(0),mRaw->GetZ(0));
1077 //___________________________________________
1078 void AliMUONDisplay::LoadTracks()
1081 AliMUONTrack* recTrack = 0;
1082 AliMUONTrackParam* trackParam = 0;
1083 TClonesArray * trackParamAtHit = 0;
1087 GetMUONData()->SetTreeAddress("RT");
1088 TClonesArray* recTracksArray = GetMUONData()->RecTracks();
1089 if (recTracksArray == NULL) return;
1090 GetMUONData()->GetRecTracks();
1092 Int_t nRecTracks = 0;
1094 nRecTracks = (Int_t) recTracksArray->GetEntriesFast();
1097 if (fRpoints == 0) fRpoints = new TObjArray(nRecTracks);
1099 for (Int_t iRecTracks = 0; iRecTracks < nRecTracks; iRecTracks++) {
1100 // reading info from tracks
1101 recTrack = (AliMUONTrack*) recTracksArray->At(iRecTracks);
1103 Int_t nTrackHits = recTrack->GetNTrackHits();
1105 if (nTrackHits == 0) continue;
1108 TPolyLine3D *points = new TPolyLine3D(nTrackHits+1);
1109 points->SetLineColor(6);
1110 points->SetLineWidth(1);
1111 fRpoints->AddAt(points,iRecTracks);
1117 trackParam = recTrack->GetTrackParamAtVertex();
1118 xRec = trackParam->GetNonBendingCoor();
1119 yRec = trackParam->GetBendingCoor();
1120 zRec = trackParam->GetZ();
1121 points->SetPoint(iPoint,xRec,yRec,zRec);
1124 for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
1125 trackParamAtHit = recTrack->GetTrackParamAtHit();
1126 trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit);
1127 xRec = trackParam->GetNonBendingCoor();
1128 yRec = trackParam->GetBendingCoor();
1129 zRec = trackParam->GetZ();
1130 points->SetPoint(iPoint,xRec,yRec,zRec);
1132 } // end loop rec. hits
1133 PrintTrack(iRecTracks,recTrack);
1134 } // end loop tracks
1139 void AliMUONDisplay::PrintTrack(Int_t iRecTracks, AliMUONTrack *recTrack)
1141 AliMUONTrackParam *trackParam;
1142 Float_t vertex[3], momentum[3];
1143 Float_t pYZ, bendingSlope, nonBendingSlope, chi2dof;
1146 trackParam = recTrack->GetTrackParamAtVertex();
1147 vertex[0] = trackParam->GetNonBendingCoor();
1148 vertex[1] = trackParam->GetBendingCoor();
1149 vertex[2] = trackParam->GetZ();
1150 pYZ = 1./TMath::Abs(trackParam->GetInverseBendingMomentum());
1151 bendingSlope = trackParam->GetBendingSlope();
1152 nonBendingSlope = trackParam->GetNonBendingSlope();
1153 momentum[2] = -pYZ / TMath::Sqrt(1.0 + bendingSlope*bendingSlope);
1154 momentum[0] = momentum[2] * nonBendingSlope;
1155 momentum[1] = momentum[2] * bendingSlope;
1156 charge = Int_t(TMath::Sign(1.,trackParam->GetInverseBendingMomentum()));
1157 chi2dof = recTrack->GetFitFMin()/(2.0 * recTrack->GetNTrackHits() - 5.);
1159 printf("===================================================\n");
1160 printf(" Reconstructed track # %d \n",iRecTracks);
1161 printf(" charge: %d \n",charge);
1162 printf(" vertex x,y,z (cm): %f %f %f \n",vertex[0],vertex[1],vertex[2]);
1163 printf(" momentum Px,Py,Pz (GeV/c): %f %f %f \n",momentum[0],momentum[1],momentum[2]);
1164 printf(" track chi2/dof: %f \n",chi2dof);
1168 //___________________________________________
1169 void AliMUONDisplay::LoadHits(Int_t chamber)
1171 // Read hits info and store x,y,z info in arrays fPhits
1172 // Loop on all detectors
1174 if (chamber > 14) return;
1181 AliMUON *pMUON = (AliMUON*)gAlice->GetModule("MUON");
1182 AliMUONChamber* iChamber;
1184 iChamber = &(pMUON->Chamber(chamber-1));
1185 Float_t zpos=iChamber->Z();
1188 if (GetMUONData()->TreeH()) {
1189 GetMUONData()->SetTreeAddress("H");
1190 Int_t ntracks = (Int_t)GetMUONData()->TreeH()->GetEntries(); //skowron
1192 for (track = 0; track < ntracks; track++) {
1193 GetMUONData()->ResetHits();
1194 GetMUONData()->GetTrack(track);//skowron
1195 TClonesArray *muonHits = GetMUONData()->Hits();
1196 if (muonHits == 0) return;
1197 nthits += muonHits->GetEntriesFast();
1199 if (fPhits == 0) fPhits = new TObjArray(nthits);
1201 for (track=0; track<ntracks;track++) {
1202 GetMUONData()->ResetHits();
1203 GetMUONData()->GetTrack(track);//skowron
1204 TClonesArray *muonHits = GetMUONData()->Hits();
1205 if (muonHits == 0) return;
1206 Int_t nhits = muonHits->GetEntriesFast();
1207 if (nhits == 0) continue;
1209 AliMUONPoints *points = 0;
1211 for (Int_t hit=0;hit<nhits;hit++) {
1212 mHit = (AliMUONHit*)muonHits->UncheckedAt(hit);
1213 Int_t nch = mHit->Chamber(); // chamber number
1214 if (nch != chamber) continue;
1216 // Retrieve info and set the objects
1218 points = new AliMUONPoints(npoints);
1219 fPhits->AddAt(points,nhold+hit);
1220 points->SetMarkerColor(kRed);
1221 points->SetMarkerStyle(5);
1222 points->SetMarkerSize(1.);
1223 points->SetParticle(mHit->Track());
1224 points->SetHitIndex(hit);
1225 points->SetTrackIndex(track);
1226 points->SetDigitIndex(-1);
1227 points->SetPoint(0,mHit->X(),mHit->Y(),zpos);
1228 // printf("%f and %f and %f\n",mHit->X(),mHit->Y(),mHit->Z());
1235 //_____________________________________________________________________________
1236 void AliMUONDisplay::Paint(Option_t *)
1238 // Paint miscellaneous items
1241 //_____________________________________________________________________________
1242 void AliMUONDisplay::SetPickMode()
1244 // Set parameters for pick mode.
1248 fArcButton->SetY1(fPickButton->GetYlowNDC()+0.5*fPickButton->GetHNDC());
1249 fTrigPad->Modified();
1252 //_____________________________________________________________________________
1253 void AliMUONDisplay::SetZoomMode()
1255 // Set parameters for zoom mode
1258 fArcButton->SetY1(fZoomButton->GetYlowNDC()+0.5*fZoomButton->GetHNDC());
1259 fTrigPad->Modified();
1262 //_____________________________________________________________________________
1263 void AliMUONDisplay::NextChamber(Int_t delta)
1265 // to go from chamber to next chamber if delta = 1
1266 // or previous chamber otherwise
1268 if (fChamber < AliMUONConstants::NCh()) fChamber++;
1270 if (fChamber > 1) fChamber--;
1274 LoadDigits(fChamber, fCathode);
1278 //_____________________________________________________________________________
1279 void AliMUONDisplay::NextCathode()
1281 // to switch to other cathode plane
1284 if (fCathode == 1) {
1285 LoadDigits(fChamber, 2);
1287 LoadDigits(fChamber, 1);
1289 fNextCathode = kTRUE; // to keep the same zoom
1291 fNextCathode = kFALSE;
1292 TPad *pad = (TPad*)gPad->GetPadSave();
1293 pad->Range(fZoomX0[fZooms], fZoomY0[fZooms],
1294 fZoomX1[fZooms], fZoomY1[fZooms]);
1300 //_____________________________________________________________________________
1301 void AliMUONDisplay::Trigger()
1304 AliMUONGlobalTrigger* globalTrig;
1306 GetMUONData()->SetTreeAddress("GLT");
1307 GetMUONData()->GetTriggerD();
1309 globalTrig = (AliMUONGlobalTrigger*)GetMUONData()->GlobalTrigger()->UncheckedAt(0);
1310 if (globalTrig == 0) return;
1312 printf("===================================================\n");
1313 printf(" Global Trigger output Low pt High pt All\n");
1315 printf(" number of Single Plus :\t");
1316 printf("%i\t",globalTrig->SinglePlusLpt());
1317 printf("%i\t",globalTrig->SinglePlusHpt());
1318 printf("%i\t",globalTrig->SinglePlusApt());
1321 printf(" number of Single Minus :\t");
1322 printf("%i\t",globalTrig->SingleMinusLpt());
1323 printf("%i\t",globalTrig->SingleMinusHpt());
1324 printf("%i\t",globalTrig->SingleMinusApt());
1327 printf(" number of Single Undefined :\t");
1328 printf("%i\t",globalTrig->SingleUndefLpt());
1329 printf("%i\t",globalTrig->SingleUndefHpt());
1330 printf("%i\t",globalTrig->SingleUndefApt());
1333 printf(" number of UnlikeSign pair :\t");
1334 printf("%i\t",globalTrig->PairUnlikeLpt());
1335 printf("%i\t",globalTrig->PairUnlikeHpt());
1336 printf("%i\t",globalTrig->PairUnlikeApt());
1339 printf(" number of LikeSign pair :\t");
1340 printf("%i\t",globalTrig->PairLikeLpt());
1341 printf("%i\t",globalTrig->PairLikeHpt());
1342 printf("%i\t",globalTrig->PairLikeApt());
1344 printf("===================================================\n");
1348 // // returns Trigger Decision for current event
1349 // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(GetLoader(),1);
1351 // // AliMUONTriggerDecision* decision= new AliMUONTriggerDecision(1);
1352 // AliMUONData* muonData = decision->GetMUONData();
1353 // muonData->SetTreeAddress("D");
1354 // decision->Trigger();
1356 //_____________________________________________________________________________
1357 void AliMUONDisplay::SetChamberAndCathode(Int_t chamber, Int_t cathode)
1359 // Set chamber and cathode number
1365 LoadDigits(chamber,cathode);
1369 void AliMUONDisplay::SetEvent(Int_t newevent)
1372 gAlice->GetEvent(newevent);
1374 if (!gAlice->TreeD()) return;
1377 LoadDigits(fChamber,fCathode);
1381 //_____________________________________________________________________________
1382 void AliMUONDisplay::SetRange(Float_t rrange, Float_t zrange)
1384 // Set view range along R and Z
1393 //_____________________________________________________________________________
1394 void AliMUONDisplay::SetView(Float_t theta, Float_t phi, Float_t psi)
1396 // change viewing angles for current event
1404 TView *view = gPad->GetView();
1405 if (view) view->SetView(fPhi, fTheta, fPsi, iret);
1410 //_____________________________________________________________________________
1411 void AliMUONDisplay::ShowNextEvent(Int_t delta)
1413 AliRunLoader * runLoader;
1415 runLoader = fLoader->GetRunLoader();
1419 // Display (current event_number + delta)
1420 // delta = 1 shown next event
1421 // delta = -1 show previous event
1423 //runLoader->CleanDetectors();
1424 //runLoader->CleanKinematics();
1425 Int_t currentEvent = runLoader->GetEventNumber();
1426 Int_t newEvent = currentEvent + delta;
1427 runLoader->GetEvent(newEvent);
1430 LoadDigits(fChamber, fCathode);
1435 //______________________________________________________________________________
1436 void AliMUONDisplay::UnZoom()
1439 if (fZooms <= 0) return;
1441 TPad *pad = (TPad*)gPad->GetPadSave();
1442 pad->Range(fZoomX0[fZooms],fZoomY0[fZooms], fZoomX1[fZooms],fZoomY1[fZooms]);
1446 //_____________________________________________________________________________
1447 void AliMUONDisplay::ResetPoints()
1450 // Reset array of points
1458 //_____________________________________________________________________________
1459 void AliMUONDisplay::ResetPhits()
1462 // Reset array of points
1470 //_____________________________________________________________________________
1471 void AliMUONDisplay::ResetRpoints()
1474 // Reset array of points
1483 AliMUONDisplay & AliMUONDisplay::operator = (const AliMUONDisplay & rhs)
1485 // Protected assignement operator
1487 if (this == &rhs) return *this;
1489 Fatal("operator=", "Not implemented.");