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