]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDisplFast.cxx
Print removed
[u/mrichter/AliRoot.git] / RICH / AliRICHDisplFast.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15 #include "AliRICHDisplFast.h"
16 #include "AliRICH.h"
17 #include "AliRICHChamber.h"
18 #include "AliRICHParam.h"
19 #include <AliLoader.h>
20 #include <TCanvas.h>
21 #include <TPolyLine.h>
22 #include <TParticle.h>
23 #include <TStyle.h>
24 #include <TH2.h>
25 #include <TMath.h>
26 #include <TLatex.h>
27
28 ClassImp(AliRICHDisplFast)
29
30 //__________________________________________________________________________________________________
31 void AliRICHDisplFast::Exec(Option_t *)
32 {
33   AliRICH *pRich = (AliRICH*)gAlice->GetDetector("RICH");
34   Bool_t isHits    =!pRich->GetLoader()->LoadHits();
35   Bool_t isDigits  =!pRich->GetLoader()->LoadDigits();
36   Bool_t isClusters=!pRich->GetLoader()->LoadRecPoints();
37   
38   if(!isHits && !isDigits && !isClusters){Error("Exec","No hits digits and clusters. Nothing to display.");return;}
39   
40   TCanvas *Display = new TCanvas("Display","RICH Display",0,0,600,600);
41   gStyle->SetPalette(1);
42   
43   TH2F *pHitsH2=0,*pDigitsH2=0,*pClustersH2=0;
44   
45   if(isHits)     pHitsH2     = new TH2F("pHitsH2"  ,  "Event Display;x,cm;y,cm",165,0,AliRICHParam::PcSizeX(),
46                                                                       144,0,AliRICHParam::PcSizeY());
47   if(pHitsH2)    pHitsH2->SetStats(kFALSE);
48   
49   if(isDigits)   pDigitsH2   = new TH2F("pDigitsH2"  ,"Event Display",165,0,AliRICHParam::PcSizeX(),
50                                                                       144,0,AliRICHParam::PcSizeY());
51   if(isClusters) pClustersH2 = new TH2F("pClustersH2","Event Display",165,0,AliRICHParam::PcSizeX(),
52                                                                       144,0,AliRICHParam::PcSizeY());
53
54     
55
56   
57   for(Int_t iEventN=0;iEventN<gAlice->GetEventsPerRun();iEventN++){//events Loop
58     pRich->GetLoader()->GetRunLoader()->GetEvent(iEventN);
59
60   
61     Int_t nPrimaries = (Int_t)pRich->GetLoader()->TreeH()->GetEntries();
62     TObjArray * Hits = new TObjArray[nPrimaries];
63   
64     for(Int_t i=0;i<nPrimaries;i++) {
65       pRich->GetLoader()->TreeH()->GetEntry(i);
66       Int_t nHits = pRich->Hits()->GetEntries();
67       for(Int_t k=0;k<nHits;k++)         Hits[i].Add(pRich->Hits()->At(k));
68     
69     }
70 //display all the staff on chamber by chamber basis           
71     for(Int_t iChamber=1;iChamber<=7;iChamber++){//chambers loop       
72       if(isHits)     pHitsH2    ->Reset();     
73       if(isDigits)   pDigitsH2  ->Reset();     
74       if(isClusters) pClustersH2->Reset();
75 //deals with hits
76       for(Int_t i=0;i<nPrimaries;i++){//prims loop
77         pRich->GetLoader()->TreeH()->GetEntry(i);
78         Int_t nHits = pRich->Hits()->GetEntries();
79         for(Int_t j=0;j<nHits;j++){//hits loop
80           AliRICHhit *pHit = (AliRICHhit*)Hits[i].At(j);
81           if(pHit->C()==iChamber){
82             TVector3 hitGlobX3= pHit->OutX3();
83             TVector2 hitLocX2 = pRich->C(iChamber)->Mrs2Pc(hitGlobX3);
84             pHitsH2->Fill(hitLocX2.X(),hitLocX2.Y(),200);
85           }//if
86         }//hits loop         
87       }//prims loop
88       pHitsH2->SetTitle(Form("event %i chamber %2i",iEventN,iChamber));
89       pHitsH2->SetMarkerColor(kRed); pHitsH2->SetMarkerStyle(29); pHitsH2->SetMarkerSize(0.4);
90       pHitsH2->Draw();
91       DrawSectors();
92       TLatex l; l.SetNDC(); l.SetTextSize(0.02);
93       if(!isHits)     {l.SetTextColor(kRed)  ;l.DrawLatex(0.1,0.01,"No Hits"    );}
94       if(!isDigits)   {l.SetTextColor(kGreen);l.DrawLatex(0.4,0.01,"No DIGITS"  );}
95       if(!isClusters) {l.SetTextColor(kBlue) ;l.DrawLatex(0.8,0.01,"No CLUSTERS");}
96       Display->Update();
97       Display->Modified();
98       gPad->WaitPrimitive();
99 //deals with digits      
100       if(isDigits){
101         pRich->GetLoader()->TreeD()->GetEntry(0);
102         for(Int_t j=0;j<pRich->Digits(iChamber)->GetEntries();j++){//digits loop
103           AliRICHdigit *pDig = (AliRICHdigit*)pRich->Digits(iChamber)->At(j);
104           TVector2 x2=AliRICHParam::Pad2Loc(pDig->Pad());
105           pDigitsH2->Fill(x2.X(),x2.Y(),100);
106         }//digits loop
107         pDigitsH2->SetMarkerColor(kGreen); pDigitsH2->SetMarkerStyle(29); pDigitsH2->SetMarkerSize(0.4);
108         pDigitsH2->Draw("same");
109         Display->Update();
110         Display->Modified();       
111         gPad->WaitPrimitive();
112       }//if(isDigits)      
113 //deals with clusters      
114       if(isClusters){
115         pRich->GetLoader()->TreeR()->GetEntry(0);
116         for(Int_t j=0;j<pRich->Clusters(iChamber)->GetEntries();j++){//clusters loop
117           AliRICHcluster *pClus = (AliRICHcluster*)pRich->Clusters(iChamber)->At(j);
118           pClustersH2->Fill(pClus->X(),pClus->Y(),50);
119         }//clusters loop
120         pClustersH2->SetMarkerColor(kBlue); pClustersH2->SetMarkerStyle(29);  pClustersH2->SetMarkerSize(0.4);
121         pClustersH2->Draw("same");
122         Display->Update();
123         Display->Modified();
124         gPad->WaitPrimitive();
125       }//if(isClusters)
126     }//chambers loop
127     delete [] Hits;
128   }//events Loop
129   pRich->GetLoader()->UnloadHits();
130   if(isDigits)   pRich->GetLoader()->UnloadDigits();
131   if(isClusters) pRich->GetLoader()->UnloadRecPoints();
132 }//Exec()
133 //__________________________________________________________________________________________________
134 void AliRICHDisplFast::DrawSectors() 
135
136   AliRICHParam p;
137   Double_t xLeft[5]  = {0,0,p.SectorSizeX(),p.SectorSizeX(),0};
138   Double_t xRight[5] = {p.SectorSizeX()+p.DeadZone(),p.SectorSizeX()+p.DeadZone(),p.PcSizeX(),p.PcSizeX(),p.SectorSizeX()+p.DeadZone()};
139   
140   Double_t yDown[5]   = {0,p.SectorSizeY(),p.SectorSizeY(),0,0};
141   Double_t yCenter[5] = {  p.SectorSizeY()+p.DeadZone(),2*p.SectorSizeY()+p.DeadZone(),2*p.SectorSizeY()+p.DeadZone(),
142                            p.SectorSizeY()+p.DeadZone(),p.SectorSizeY()+p.DeadZone()};  
143   Double_t yUp[5]     = {2*p.SectorSizeY()+2*p.DeadZone(),p.PcSizeY(),p.PcSizeY(),2*p.SectorSizeY()+2*p.DeadZone(),2*p.SectorSizeY()+2*p.DeadZone()};
144   
145   TPolyLine *sec1 = new TPolyLine(5,xLeft ,yDown);    sec1->SetLineColor(21);  sec1->Draw();
146   TPolyLine *sec2 = new TPolyLine(5,xRight,yDown);    sec2->SetLineColor(21);  sec2->Draw();
147   TPolyLine *sec3 = new TPolyLine(5,xLeft ,yCenter);  sec3->SetLineColor(21);  sec3->Draw();
148   TPolyLine *sec4 = new TPolyLine(5,xRight,yCenter);  sec4->SetLineColor(21);  sec4->Draw();
149   TPolyLine *sec5 = new TPolyLine(5,xLeft, yUp);      sec5->SetLineColor(21);  sec5->Draw();
150   TPolyLine *sec6 = new TPolyLine(5,xRight,yUp);      sec6->SetLineColor(21);  sec6->Draw();
151 }//DrawSectors()