]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/readDigits.C
update for the NUA
[u/mrichter/AliRoot.git] / ITS / UPGRADE / readDigits.C
1 void readDigits(int nev=-1,int evStart=0)
2 {
3
4   gSystem->Load("libITSUpgradeBase");
5   gSystem->Load("libITSUpgradeSim");
6   gROOT->SetStyle("Plain");
7   const Int_t kMaxROCycleAccept=126;
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,kTRUE);
21   //
22   Int_t nLayers = gm->GetNLayers();
23   Int_t nChips = gm->GetNChips();
24   Int_t nbins=100;
25   Int_t xmin=0;
26   Int_t xmax=50000;//00*1e-09;
27
28
29   TH1D **hNelSDig = new TH1D*[nLayers], **hNelDig = new TH1D*[nLayers];
30   //
31   for(Int_t i=0; i< nLayers; i++ ) {
32     hNelSDig[i] = new TH1D(Form("hNelSDig%i",i),Form("electron distribution in SDigits [ Layer %i] ",i),nbins,xmin,xmax);
33     hNelSDig[i]->SetXTitle("N electrons");
34     hNelDig[i] = new TH1D(Form("hNelDig%i",i),Form("electron distribution in Digits [ Layer %i] ",i),nbins,xmin,xmax);
35     hNelDig[i]->SetXTitle("N electrons");
36   }
37   
38   AliLoader *dl = runLoader->GetDetectorLoader("ITS");
39
40
41   //SDIGITS INIT
42   TTree * sDigTree = 0x0;
43   TClonesArray *sDigArr= new TClonesArray("AliITSUSDigit");
44
45   //DIGITS INIT
46   TTree * digTree = 0x0;
47   TClonesArray *digArr= new TClonesArray("AliITSUDigitPix");
48
49   int nevTot = (Int_t)runLoader->GetNumberOfEvents();
50   printf("N Events : %i \n",nevTot);
51   evStart = evStart<nevTot ? evStart : nevTot-1;
52   if (evStart<0) evStart = 0;
53   //
54   int lastEv = nev<0 ? nevTot : evStart+nev;
55   if (lastEv > nevTot) lastEv = nevTot;
56   //
57   printf("N Events : %i \n",(Int_t)nevTot);
58
59   for (Int_t iEvent = evStart; iEvent < lastEv; iEvent++) {
60     printf("\n Event %i \n",iEvent);
61     runLoader->GetEvent(iEvent);
62     AliStack *stack = runLoader->Stack();
63     sDigTree=dl->TreeS();
64     digTree=dl->TreeD();
65     //
66     if (sDigTree) sDigTree->SetBranchAddress("ITS",&sDigArr);
67     digTree->SetBranchAddress("ITSDigitsPix",&digArr);
68
69     for (int imod=0;imod<nChips;imod++) {
70       if (sDigTree) sDigTree->GetEntry(imod);
71       digTree->GetEntry(imod);      
72       int detType = gm->GetChipChipTypeID(imod);
73       AliITSUSegmentationPix* segm = (AliITSUSegmentationPix*)gm->GetSegmentationByID(detType);
74       int lay,sta,ssta,mod,chip;
75       int nsdig = sDigArr->GetEntries();
76       int ndig  = digArr->GetEntries();
77       if (ndig<1) continue;
78       gm->GetChipId(imod, lay,sta,ssta,mod,chip);
79       printf("\nChip %3d: (chip %2d in module:%d/substave:%1d/stave:%2d/Layer:%d) |NSDigits: %4d NDigits: %4d\n",imod,chip,mod,ssta,sta,lay,nsdig,ndig);
80       //
81       for (int isdig=0;isdig<nsdig;isdig++) {
82         AliITSUSDigit *pSdig = (AliITSUSDigit*)sDigArr->At(isdig);
83         int sdinfo = pSdig->GetUniqueID();
84         UInt_t row,col;
85         Int_t cycle;
86         AliITSUSensMap::GetCell(sdinfo,segm->Npz(),segm->Npx(),kMaxROCycleAccept,col,row,cycle);
87         printf("#%3d Sdigit col:%4d/row:%4d/cycle:%d generated by track %5d (%s)\t",isdig, col,row,cycle,
88                pSdig->GetTrack(0),stack->Particle(pSdig->GetTrack(0))->GetName());
89         pSdig->Print();
90         hNelSDig[lay]->Fill(pSdig->GetSignal()); 
91       }
92       //
93       for (int idig=0;idig<ndig;idig++) {
94         AliITSUDigitPix *pDig = (AliITSUDigitPix*)digArr->At(idig);
95         printf("#%3d digit, col:%4d/row:%4d ROCycle:%d signal: %.2e,  generated by tracks ",idig,pDig->GetCoord1(),pDig->GetCoord2(),
96                pDig->GetROCycle(),pDig->GetSignalPix()); 
97         for (int itr=0;itr<pDig->GetNTracks();itr++) if (pDig->GetTrack(itr)>=0) printf(" %5d",pDig->GetTrack(itr)); printf("\n");
98         //
99         hNelDig[lay]->Fill(pDig->GetSignalPix()); 
100       }
101       //
102     }
103
104   }//event loop
105
106   TCanvas *cSd = new TCanvas("cSd","Summable Digits",1000,800);
107   cSd->Divide(nLayers&0x1 ?  (nLayers/2+1) : nLayers/2,2); 
108   TCanvas *cD = new TCanvas("cD","Digits",1000,800);
109   cD->Divide(nLayers&0x1 ?  (nLayers/2+1) : nLayers/2,2); 
110
111   for(Int_t ip =1; ip<=nLayers; ip++){
112     cSd->cd(ip);
113     hNelSDig[ip-1]->Draw();
114     cD->cd(ip);
115     hNelDig[ip-1]->Draw();
116   }
117
118
119 }