]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/macrosSDD/CheckDataSizeSDD.C
New macro to plot SDD occupancy per module/ddl
[u/mrichter/AliRoot.git] / ITS / macrosSDD / CheckDataSizeSDD.C
CommitLineData
3f55c67a 1#if !defined(__CINT__) || defined(__MAKECINT__)
2#include <TH2F.h>
3#include <TGrid.h>
4#include <TCanvas.h>
5#include <TStopwatch.h>
6#include <TStyle.h>
7#include "AliRawReaderDate.h"
8#include "AliRawReaderRoot.h"
9#include "AliITSRawStreamSDD.h"
10#include "AliITSRawStreamSDDCompressed.h"
11#endif
12
13
14void CheckDataSizeSDD(TString datafil="12000188359075.10.root",
15 Int_t firstEv=0,
16 Int_t lastEv=123456){
17
18 gStyle->SetOptStat(0);
19 gStyle->SetTitleFont(42,"XY");
20 gStyle->SetLabelFont(42,"XYZ");
21
22 Double_t maxOcc=2500.;
23 TH2F* hCellsOnMod=new TH2F("hCellsOn","",260,239.5,499.5,500,0.,maxOcc);
24 TH2F* hCellsOnDDL=new TH2F("hCellsOnDDL","",24,-0.5,23.5,500,0.,maxOcc*12);
25 TH1F* hMaxOccMod=new TH1F("hMaxOccMod","",260,239.5,499.5);
26 TH1F* hMaxOccDDL=new TH1F("hMaxOccDDL","",24,-0.5,23.5);
27
28 AliITSDDLModuleMapSDD* dmap=new AliITSDDLModuleMapSDD();
29 dmap->SetJun09Map();
30
31
32 Int_t iev=firstEv;
33 AliRawReader *rd;
34 if(datafil.Contains(".root")){
35 rd=new AliRawReaderRoot(datafil.Data(),iev);
36 }else{
37 rd=new AliRawReaderDate(datafil.Data(),iev);
38 }
39
40 Bool_t writtenoutput=kFALSE;
41 Int_t countMod[260],countDDL[24];
42 Int_t maxOccMod[260],maxOccDDL[24];
43 for(Int_t im=0; im<260; im++) maxOccMod[im]=0;
44 for(Int_t id=0; id<24; id++) maxOccDDL[id]=0;
45
46 do{
47
48 printf("Event # %d\n",iev);
49 for(Int_t im=0; im<260; im++) countMod[im]=0;
50 for(Int_t id=0; id<24; id++) countDDL[id]=0;
51 rd->Reset();
52 UChar_t cdhAttr=AliITSRawStreamSDD::ReadBlockAttributes(rd);
53 UInt_t amSamplFreq=AliITSRawStreamSDD::ReadAMSamplFreqFromCDH(cdhAttr);
54 AliITSRawStream* s=AliITSRawStreamSDD::CreateRawStreamSDD(rd,cdhAttr);
55 if(!writtenoutput){
56 printf("Use %s raw stream, sampling frequency %d MHz\n",s->ClassName(),amSamplFreq);
57 writtenoutput=kTRUE;
58 }
59
60 while(s->Next()){
61 if(s->IsCompletedModule()==kFALSE && s->IsCompletedDDL()==kFALSE){
62 Int_t counts=s->GetSignal();
63 if(counts>0){
64 Int_t iDDL=rd->GetDDLID();
65 ++countDDL[iDDL];
66 Int_t iMod=s->GetCarlosId();
67 Int_t modInd=dmap->GetModuleNumber(iDDL,iMod)-240;
68 if(modInd>=0 && modInd<260) ++countMod[modInd];
69 }
70 }
71 }
72 for(Int_t im=0; im<260; im++){
73 hCellsOnMod->Fill(im+240,countMod[im]);
74 if(countMod[im]>maxOccMod[im]) maxOccMod[im]=countMod[im];
75 }
76 for(Int_t id=0; id<24; id++){\
77 hCellsOnDDL->Fill(id,countDDL[id]);
78 if(countDDL[id]>maxOccDDL[id]) maxOccDDL[id]=countDDL[id];
79 }
80 iev++;
81
82 }while(rd->NextEvent()&&iev<=lastEv);
83
84 for(Int_t im=0; im<260; im++) hMaxOccMod->SetBinContent(im+1,maxOccMod[im]);
85 for(Int_t id=0; id<24; id++) hMaxOccDDL->SetBinContent(id+1,maxOccDDL[id]);
86
87 TCanvas* c1= new TCanvas("c1","Module occupancy",600,750);
88 c1->Divide(1,2);
89 c1->cd(1);
90 gPad->SetLogz();
91 hCellsOnMod->GetXaxis()->SetTitle("Module Id");
92 hCellsOnMod->GetYaxis()->SetTitle("Number of cells on per event");
93 hCellsOnMod->GetYaxis()->SetTitleOffset(1.3);
94 hCellsOnMod->Draw("colz");
95 c1->cd(2);
96 hMaxOccMod->GetXaxis()->SetTitle("Module Id");
97 hMaxOccMod->GetYaxis()->SetTitle("Maximum Occupancy");
98 hMaxOccMod->GetYaxis()->SetTitleOffset(1.3);
99 hMaxOccMod->Draw();
100
101 TCanvas* c2= new TCanvas("c2","DDL occupancy",600,750);
102 c2->Divide(1,2);
103 c2->cd(1);
104 gPad->SetLogz();
105 hCellsOnDDL->GetXaxis()->SetTitle("DDL Number");
106 hCellsOnDDL->GetYaxis()->SetTitle("Number of cells on per event");
107 hCellsOnDDL->GetYaxis()->SetTitleOffset(1.3);
108 hCellsOnDDL->Draw("colz");
109 c2->cd(2);
110 hMaxOccDDL->GetXaxis()->SetTitle("DDL Number");
111 hMaxOccDDL->GetYaxis()->SetTitle("Maximum Occupancy");
112 hMaxOccDDL->GetYaxis()->SetTitleOffset(1.3);
113 hMaxOccDDL->Draw();
114
115}