]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/DisplaySDDRawData.C
Implementation of air cooling on C side. New material - PVC - A. Pulvirenti, M. Sitta
[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(TString filename, Int_t firstEv=0, Int_t lastEv=5){
17
18
19   Bool_t writtenoutput=kFALSE;
20   AliITSDDLModuleMapSDD* ddlmap=new AliITSDDLModuleMapSDD();
21   ddlmap->SetJun09Map();
22
23   TH2F* hzphi3=new TH2F("hzphi3","Layer 3",1536,-0.5,1535.5,3584,-0.5,3584.5);
24   TH2F* hzphi4=new TH2F("hzphi4","Layer 4",2048,-0.5,2047.5,5632,-0.5,5631.5);
25
26   TLine** lA3=new TLine*[5];
27   for(Int_t ilin=0;ilin<5;ilin++){
28     lA3[ilin]=new TLine((ilin+1)*256,0,(ilin+1)*256,3584.5);
29     lA3[ilin]->SetLineColor(kGray);
30     lA3[ilin]->SetLineStyle(2);
31   }
32   TLine** lT3=new TLine*[13];
33   for(Int_t ilin=0;ilin<13;ilin++){
34     lT3[ilin]=new TLine(0,(ilin+1)*256,1535.5,(ilin+1)*256);
35     lT3[ilin]->SetLineColor(kGray);
36     lT3[ilin]->SetLineStyle(2);
37   }
38
39   TLine** lA4=new TLine*[7];
40   for(Int_t ilin=0;ilin<7;ilin++){
41     lA4[ilin]=new TLine((ilin+1)*256,0,(ilin+1)*256,5631.5);
42     lA4[ilin]->SetLineColor(kGray);
43     lA4[ilin]->SetLineStyle(2);
44   }
45   TLine** lT4=new TLine*[21];
46   for(Int_t ilin=0;ilin<21;ilin++){
47     lT4[ilin]=new TLine(0,(ilin+1)*256,2047.5,(ilin+1)*256);
48     lT4[ilin]->SetLineColor(kGray);
49     lT4[ilin]->SetLineStyle(2);
50   }
51
52   hzphi3->SetStats(0);
53   hzphi4->SetStats(0);
54
55   Int_t iev=firstEv;
56   AliRawReader *rd; 
57   if(filename.Contains(".root")){
58     rd=new AliRawReaderRoot(filename.Data(),iev);
59   }else{
60     rd=new AliRawReaderDate(filename.Data(),iev);
61   }
62
63   TStopwatch *evtime=new TStopwatch();
64   TCanvas* c0 = new TCanvas("cd0","c0",800,800);
65   gStyle->SetPalette(1);
66   do{
67     c0->Clear();                                
68     c0->Divide(1,2,0.001,0.001);
69
70     evtime->Start();
71     printf("Event # %d\n",iev);
72     rd->Reset();
73     hzphi3->Reset();
74     hzphi4->Reset();
75
76     UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rd);
77     UInt_t amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr);
78     AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rd,cdhAttr);
79     if(!writtenoutput){
80       printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq);
81       writtenoutput=kTRUE;
82
83     }
84
85     while(s->Next()){
86       
87       if(s->IsCompletedModule()==kFALSE && s->IsCompletedDDL()==kFALSE){
88         Int_t lay,lad,det;
89         Int_t modID=ddlmap->GetModuleNumber(rd->GetDDLID(),s->GetCarlosId());
90         AliITSgeomTGeo::GetModuleId(modID,lay,lad,det);
91         Int_t iz=s->GetCoord1()+256*(det-1);
92         Int_t iphi=s->GetCoord2()+256*(lad-1)+128*s->GetChannel();
93         if(lay==3){
94           hzphi3->SetBinContent(iz+1,iphi+1,s->GetSignal());
95         }else if(lay==4){
96           hzphi4->SetBinContent(iz+1,iphi+1,s->GetSignal());
97         }
98       }
99     }
100     evtime->Stop();
101     printf("**** Event=%d \n",iev);
102     evtime->Print("u");
103     evtime->Reset();
104     iev++;
105     
106     c0->cd(1);
107     hzphi3->Draw("colz");
108     for(Int_t ilin=0;ilin<5;ilin++) lA3[ilin]->Draw("same");
109     for(Int_t ilin=0;ilin<13;ilin++) lT3[ilin]->Draw("same");
110     hzphi3->GetXaxis()->SetTitle("Z (anode)");
111     hzphi3->GetYaxis()->SetTitle("PHI (time bin)");
112       
113     c0->cd(2);
114     hzphi4->Draw("colz");
115     for(Int_t ilin=0;ilin<7;ilin++) lA4[ilin]->Draw("same");
116     for(Int_t ilin=0;ilin<21;ilin++) lT4[ilin]->Draw("same");
117     hzphi4->GetXaxis()->SetTitle("Z (anode)");
118     hzphi4->GetYaxis()->SetTitle("PHI (time bin)");
119     c0->Update();
120   }while(rd->NextEvent()&&iev<=lastEv);
121
122 }
123
124 void DisplaySDDRawData(Int_t nrun, Int_t n2, Int_t chunk=10, Int_t year=2009, Char_t* dir="LHC09b", 
125                        Int_t firstEv=21*3, 
126                        Int_t lastEv=21*3+1){  
127   TGrid::Connect("alien:",0,0,"t");
128   TString filnam(Form("alien:///alice/data/%d/%s/%09d/raw/%02d%09d%03d.%02d.root",year,dir,nrun,year-2000,nrun,n2,chunk));
129   printf("Open file %s\n",filnam.Data());
130   DisplaySDDRawData(filnam,firstEv,lastEv);
131 }
132