]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/readHit.C
7b43122973e37558f796fded9a816e8f6b40409b
[u/mrichter/AliRoot.git] / ITS / UPGRADE / readHit.C
1
2 void readHit(){
3   gROOT->SetStyle("Plain");
4
5   gSystem->Load("libITSUpgradeBase");
6   gSystem->Load("libITSUpgradeSim");
7
8   TH2F *xyGlob = new TH2F("xyGlob"," X - Y Global coordinates ",100,-50,50,100,-50,50);
9   xyGlob->SetXTitle("cm");
10   xyGlob->SetMarkerStyle(7);
11   TH1F *zGlob  = new TH1F("zGlob", " Z Global coordinates ",200, -100,100 );
12   zGlob->SetXTitle("cm");
13
14   Int_t nbins=100;
15   Double_t xmin=0;
16   Double_t xmax=1e-04;
17
18   gAlice=NULL;
19   AliRunLoader* runLoader = AliRunLoader::Open("galice.root");
20   runLoader->LoadgAlice();
21
22   gAlice = runLoader->GetAliRun();
23
24   runLoader->LoadHeader();
25   runLoader->LoadKinematics();
26   runLoader->LoadSDigits();
27   runLoader->LoadHits();
28
29   AliGeomManager::LoadGeometry("geometry.root");
30   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE);
31
32   AliLoader *dl = runLoader->GetDetectorLoader("ITS");
33
34
35
36   Int_t nLayers = gm->GetNLayers();
37
38   TH1D **hDeLoss = new TH1D*[nLayers];
39   for(Int_t i=0; i< nLayers; i++ ) {
40     hDeLoss[i] = new TH1D(Form("hDeLossl%i",i),Form("E loss distribution [ Layer %i] ",i),nbins,xmin,xmax);
41     hDeLoss[i]->SetXTitle("GeV");
42   }
43
44
45   //HIT INIT
46
47   TTree *hitTree = 0x0;
48   TClonesArray *hitList=new TClonesArray("AliITSUHit");
49
50   for (Int_t iEvent = 0; iEvent < runLoader->GetNumberOfEvents(); iEvent++) {
51  
52     runLoader->GetEvent(iEvent);
53     
54     hitTree=dl->TreeH();
55     hitTree->SetBranchAddress("ITS",&hitList);
56     for(Int_t iEnt=0;iEnt<hitTree->GetEntries();iEnt++){//entries loop degli hits
57       hitTree->GetEntry(iEnt);
58       for(Int_t iHit=0; iHit<hitList->GetEntries();iHit++){
59               
60         AliITSUHit *pHit = (AliITSUHit*)hitList->At(iHit);
61         int id = pHit->GetModule();
62         int lr = gm->GetLayer(id);
63         int ld = gm->GetLadder(id);
64         //
65         //      if(pHit->GetParticle()->IsPrimary()){
66         Double_t xg,yg,zg=0.;
67         pHit->GetPositionG(xg,yg,zg);
68         xyGlob->Fill(xg,yg);
69         zGlob->Fill(zg);
70         printf("Module %d | Lr:%d Ladder: %d, X:%+.3e Y:%+.3e Z:%+.3e TrackID: %d\n",id,lr,ld,xg,yg,zg,pHit->GetTrack());
71         hDeLoss[lr]->Fill(pHit->GetIonization());
72         //      } // is primary
73       }//loop hit 
74     }//entryloopHitList
75         
76   }//event loop
77
78   TCanvas *xyCanv =  new TCanvas("xyCanv","Hit X-Y positions",500,500);
79   xyCanv->cd();
80   xyGlob->Draw();
81   TCanvas *zCanv =  new TCanvas("zCanv","Hit Z positions",500,500);
82   zCanv->cd();
83   zGlob->Draw();
84
85   TCanvas *c = new TCanvas("c","E loss  distribution per layer",1000,800);
86   c->Divide(nLayers&0x1 ?  (nLayers/2+1) : nLayers/2,2);
87   for(Int_t ip =1; ip<=nLayers; ip++){
88     c->cd(ip);
89     hDeLoss[ip-1]->Draw();
90   }
91
92 }