]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/readDigits.C
f35342bc431ffb7895aa37f938b3a070f5c8d564
[u/mrichter/AliRoot.git] / ITS / UPGRADE / readDigits.C
1 void readDigits(){
2
3   gSystem->Load("libITSUpgradeBase");
4   gSystem->Load("libITSUpgradeSim");
5   gROOT->SetStyle("Plain");
6
7
8   gAlice=NULL;
9   AliRunLoader* runLoader = AliRunLoader::Open("galice.root");
10   runLoader->LoadgAlice();
11
12   gAlice = runLoader->GetAliRun();
13
14   runLoader->LoadHeader();
15   runLoader->LoadKinematics();
16   runLoader->LoadSDigits();
17   runLoader->LoadDigits();
18
19   AliGeomManager::LoadGeometry("geometry.root");
20   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE);
21   TObjArray segmArr;
22   AliITSUSegmentationPix::LoadSegmentations(&segmArr, AliITSUGeomTGeo::GetITSsegmentationFileName());
23   //
24   Int_t nLayers = gm->GetNLayers();
25   Int_t nModules = gm->GetNModules();
26   Int_t nbins=100;
27   Int_t xmin=0;
28   Int_t xmax=50000;//00*1e-09;
29
30
31   TH1D **hNelSDig = new TH1D*[nLayers], **hNelDig = new TH1D*[nLayers];
32   //
33   for(Int_t i=0; i< nLayers; i++ ) {
34     hNelSDig[i] = new TH1D(Form("hNelSDig%i",i),Form("electron distribution in SDigits [ Layer %i] ",i),nbins,xmin,xmax);
35     hNelSDig[i]->SetXTitle("N electrons");
36     hNelDig[i] = new TH1D(Form("hNelDig%i",i),Form("electron distribution in Digits [ Layer %i] ",i),nbins,xmin,xmax);
37     hNelDig[i]->SetXTitle("N electrons");
38   }
39   
40   AliLoader *dl = runLoader->GetDetectorLoader("ITS");
41
42
43   //SDIGITS INIT
44   TTree * sDigTree = 0x0;
45   TClonesArray *sDigArr= new TClonesArray("AliITSUSDigit");
46
47   //DIGITS INIT
48   TTree * digTree = 0x0;
49   TClonesArray *digArr= new TClonesArray("AliITSUDigitPix");
50
51   printf("N Events : %i \n",(Int_t)runLoader->GetNumberOfEvents());
52
53   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
54     printf("\n Event %i \n",iEvent);
55     runLoader->GetEvent(iEvent);
56     AliStack *stack = runLoader->Stack();
57     sDigTree=dl->TreeS();
58     digTree=dl->TreeD();
59     //
60     sDigTree->SetBranchAddress("ITS",&sDigArr);
61     digTree->SetBranchAddress("ITSDigitsPix",&digArr);
62
63     for (int imod=0;imod<nModules;imod++) {
64       sDigTree->GetEntry(imod);
65       digTree->GetEntry(imod);      
66       int detType = gm->GetModuleDetTypeID(imod);
67       AliITSUSegmentationPix* segm = (AliITSUSegmentationPix*)segmArr.At(detType);
68       int lay,lad,det;
69       int nsdig = sDigArr->GetEntries();
70       int ndig  = digArr->GetEntries();
71       if (nsdig<1) continue;
72       gm->GetModuleId(imod, lay,lad,det);
73       printf("\nModule %3d: (det %2d in ladder %2d of Layer %d) |NSDigits: %4d NDigits: %4d\n",imod,det,lad,lay,nsdig,ndig);
74       //
75       for (int isdig=0;isdig<nsdig;isdig++) {
76         AliITSUSDigit *pSdig = (AliITSUSDigit*)sDigArr->At(isdig);
77         int sdinfo = pSdig->GetUniqueID();
78         printf("#%3d Sdigit col:%4d/row:%4d generated by track %5d (%s)\t",isdig, sdinfo/segm->Npx(),sdinfo%segm->Npx(),
79                pSdig->GetTrack(0),stack->Particle(pSdig->GetTrack(0))->GetName());
80         pSdig->Print();
81         hNelSDig[lay]->Fill(pSdig->GetSignal()); 
82       }
83       //
84       for (int idig=0;idig<ndig;idig++) {
85         AliITSUDigitPix *pDig = (AliITSUDigitPix*)digArr->At(idig);
86         printf("#%3d digit, col:%4d/row:%4d signal: %.2e,  generated by tracks ",idig,pDig->GetCoord1(),pDig->GetCoord2(),pDig->GetSignalPix()); 
87         for (int itr=0;itr<pDig->GetNTracks();itr++) if (pDig->GetTrack(itr)>=0) printf(" %5d",pDig->GetTrack(itr)); printf("\n");
88         //
89         hNelDig[lay]->Fill(pDig->GetSignalPix()); 
90       }
91       //
92     }
93
94   }//event loop
95
96   TCanvas *cSd = new TCanvas("cSd","Summable Digits",1000,800);
97   cSd->Divide(nLayers&0x1 ?  (nLayers/2+1) : nLayers/2,2); 
98   TCanvas *cD = new TCanvas("cD","Digits",1000,800);
99   cD->Divide(nLayers&0x1 ?  (nLayers/2+1) : nLayers/2,2); 
100
101   for(Int_t ip =1; ip<=nLayers; ip++){
102     cSd->cd(ip);
103     hNelSDig[ip-1]->Draw();
104     cD->cd(ip);
105     hNelDig[ip-1]->Draw();
106   }
107
108
109 }