]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/DisplaySDDRawData.C
Changes to obey modified coding convention RC6 (R. Shaoyan)
[u/mrichter/AliRoot.git] / ITS / DisplaySDDRawData.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TH2F.h>
3 #include <TCanvas.h>
4 #include <TStopwatch.h>
5 #include <TStyle.h>
6 #include <TLine.h>
7 #include <TGrid.h>
8 #include "AliRawReaderDate.h"
9 #include "AliRawReaderRoot.h"
10 #include "AliITSRawStreamSDD.h"
11 #endif
12
13 // Macro for a z-phi event display of the SDD Raw Data
14 // Origin: F. Prino,   prino@to.infn.it
15
16 void DisplaySDDRawData(Char_t *datafil, Int_t firstEv=0, Int_t lastEv=5){
17
18
19   AliITSDDLModuleMapSDD* ddlmap=new AliITSDDLModuleMapSDD();
20   ddlmap->SetJun08Map();
21
22   TH2F* hzphi3=new TH2F("hzphi3","Layer 3",1536,-0.5,1535.5,3584,-0.5,3584.5);
23   TH2F* hzphi4=new TH2F("hzphi4","Layer 4",2048,-0.5,2047.5,5632,-0.5,5631.5);
24
25   TLine** lA3=new TLine*[5];
26   for(Int_t ilin=0;ilin<5;ilin++){
27     lA3[ilin]=new TLine((ilin+1)*256,0,(ilin+1)*256,3584.5);
28     lA3[ilin]->SetLineColor(kGray);
29     lA3[ilin]->SetLineStyle(2);
30   }
31   TLine** lT3=new TLine*[13];
32   for(Int_t ilin=0;ilin<13;ilin++){
33     lT3[ilin]=new TLine(0,(ilin+1)*256,1535.5,(ilin+1)*256);
34     lT3[ilin]->SetLineColor(kGray);
35     lT3[ilin]->SetLineStyle(2);
36   }
37
38   TLine** lA4=new TLine*[7];
39   for(Int_t ilin=0;ilin<7;ilin++){
40     lA4[ilin]=new TLine((ilin+1)*256,0,(ilin+1)*256,5631.5);
41     lA4[ilin]->SetLineColor(kGray);
42     lA4[ilin]->SetLineStyle(2);
43   }
44   TLine** lT4=new TLine*[21];
45   for(Int_t ilin=0;ilin<21;ilin++){
46     lT4[ilin]=new TLine(0,(ilin+1)*256,2047.5,(ilin+1)*256);
47     lT4[ilin]->SetLineColor(kGray);
48     lT4[ilin]->SetLineStyle(2);
49   }
50
51   hzphi3->SetStats(0);
52   hzphi4->SetStats(0);
53
54   Int_t iev=firstEv;
55   AliRawReader *rd; 
56   if(strstr(datafil,".root")!=0){
57     rd=new AliRawReaderRoot(datafil,iev);
58   }else{
59     rd=new AliRawReaderDate(datafil,iev);
60   }
61   TStopwatch *evtime=new TStopwatch();
62   TCanvas* c0 = new TCanvas("cd0","c0",800,800);
63   gStyle->SetPalette(1);
64   Int_t idev;
65   do{
66     c0->Clear();                                
67     c0->Divide(1,2,0.001,0.001);
68
69     evtime->Start();
70     printf("Event # %d\n",iev);
71     rd->Reset();
72     hzphi3->Reset();
73     hzphi4->Reset();
74     AliITSRawStreamSDD s(rd);
75     Int_t iCountNext=0;    
76     while(s.Next()){
77       iCountNext++;
78       
79       if(s.IsCompletedModule()==kFALSE && s.IsCompletedDDL()==kFALSE){
80         Int_t lay,lad,det;
81         Int_t modID=ddlmap->GetModuleNumber(rd->GetDDLID(),s.GetCarlosId());
82         AliITSgeomTGeo::GetModuleId(modID,lay,lad,det);
83         Int_t iz=s.GetCoord1()+256*(det-1);
84         Int_t iphi=s.GetCoord2()+256*(lad-1)+128*s.GetChannel();
85         if(lay==3){
86           hzphi3->SetBinContent(iz+1,iphi+1,s.GetSignal());
87         }else if(lay==4){
88           hzphi4->SetBinContent(iz+1,iphi+1,s.GetSignal());
89         }
90       }
91     }
92     idev=s.GetEventId();
93     evtime->Stop();
94     printf("**** Event=%d  ID=%d\n",iev,idev);
95     evtime->Print("u");
96     evtime->Reset();
97     iev++;
98     
99     c0->cd(1);
100     hzphi3->Draw("colz");
101     for(Int_t ilin=0;ilin<5;ilin++) lA3[ilin]->Draw("same");
102     for(Int_t ilin=0;ilin<13;ilin++) lT3[ilin]->Draw("same");
103     hzphi3->GetXaxis()->SetTitle("Z (anode)");
104     hzphi3->GetYaxis()->SetTitle("PHI (time bin)");
105       
106     c0->cd(2);
107     hzphi4->Draw("colz");
108     for(Int_t ilin=0;ilin<7;ilin++) lA4[ilin]->Draw("same");
109     for(Int_t ilin=0;ilin<21;ilin++) lT4[ilin]->Draw("same");
110     hzphi4->GetXaxis()->SetTitle("Z (anode)");
111     hzphi4->GetYaxis()->SetTitle("PHI (time bin)");
112     c0->Update();
113   }while(rd->NextEvent()&&iev<=lastEv);
114
115 }
116
117 void DisplaySDDRawData(Int_t nrun, Int_t n2, Int_t nchunk=10, Int_t firstEv=1, Int_t lastEv=20){  
118   TGrid::Connect("alien:",0,0,"t");
119   Char_t filnam[200];
120   sprintf(filnam,"alien:///alice/data/2008/LHC08c/%09d/raw/08%09d%03d.%02d.root",nrun,nrun,n2,nchunk);
121   printf("Open file %s\n",filnam);
122   DisplaySDDRawData(filnam,firstEv,lastEv);
123 }
124