]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/readClusters.C
Updated hydro macro
[u/mrichter/AliRoot.git] / ITS / UPGRADE / readClusters.C
1 void readClusters(int nev=-1,int evStart=0)
2 {
3
4   gSystem->Load("libITSUpgradeBase");
5   gSystem->Load("libITSUpgradeRec");
6   gROOT->SetStyle("Plain");
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->LoadRecPoints();
17
18   AliLoader *dl = runLoader->GetDetectorLoader("ITS");
19
20   AliGeomManager::LoadGeometry("geometry.root");
21   AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE);
22   Int_t nLayers = gm->GetNLayers();
23   AliITSUClusterPix::SetGeom(gm);
24
25   TH2F *xyGlob = new TH2F("xyGlob"," X - Y Global coordinates ",500,-50,50,500,-50,50);
26   xyGlob->SetXTitle("cm"); 
27   xyGlob->SetMarkerStyle(7); 
28   TH1F *zGlob  = new TH1F("zGlob", " Z Global coordinates ",200, -50,50 );
29   zGlob->SetXTitle("cm"); 
30   //
31   TH1F* rGlob = new TH1F("rGlob","R global", 5000, 0,50.);
32
33   TTree * cluTree = 0x0;
34   TObjArray layerClus;
35   
36   int nevTot = (Int_t)runLoader->GetNumberOfEvents();
37   printf("N Events : %i \n",nevTot);
38   evStart = evStart<nevTot ? evStart : nevTot-1;
39   if (evStart<0) evStart = 0;
40   //
41   int lastEv = nev<0 ? nevTot : evStart+nev;
42   if (lastEv > nevTot) lastEv = nevTot;
43   //
44   for (Int_t iEvent = evStart; iEvent < lastEv; iEvent++) {
45     printf("\n Event %i \n",iEvent);
46     runLoader->GetEvent(iEvent);
47     //   AliStack *stack = runLoader->Stack();
48     cluTree=dl->TreeR();
49     int nlr=0;
50     while(1) {
51       TBranch* br = cluTree->GetBranch(Form("ITSRecPoints%d",nlr));
52       if (!br) break;
53       TClonesArray* clr = 0;
54       br->SetAddress(&clr);
55       layerClus.AddLast(clr);
56       nlr++;
57     }
58
59     printf(" tree entries: %d\n",cluTree->GetEntries());
60     cluTree->GetEntry(0);      
61     //
62     for (int ilr=0;ilr<nlr;ilr++) {
63       TClonesArray* clr = (TClonesArray*)layerClus.At(ilr);
64       int nClu = clr->GetEntries();
65       printf("Layer %d : %d clusters\n",ilr,nClu);
66       //
67       for (int icl=0;icl<nClu;icl++) {
68         AliITSUClusterPix *cl = (AliITSUClusterPix*)clr->At(icl);
69         cl->Print("glo");
70         Double_t loc[3]={cl->GetX(),cl->GetY(),cl->GetZ()}; 
71         Double_t glob[3]; 
72         gm->LocalToGlobal(cl->GetVolumeId(),loc,glob);
73         //printf("%d: mod %d: loc(%.4lf,%.4lf,%.4lf); glob(%.4lf,%.4lf,%.4lf); \n",icl,cl->GetVolumeId(),
74         //       loc[0],loc[1],loc[2],glob[0],glob[1],glob[2]);
75         xyGlob->Fill(glob[0],glob[1]);
76         zGlob->Fill(glob[2]);
77         rGlob->Fill(TMath::Sqrt(glob[0]*glob[0]+glob[1]*glob[1]));
78       }
79     }
80     layerClus.Clear();
81   }//event loop
82
83   Int_t size = 400;
84
85   TCanvas *xyCanv =  new TCanvas("xvCanvClus","RecPoint X-Y Positions",10,10,size,size);
86   xyCanv->cd();
87   xyGlob->Draw("P");
88
89   TCanvas *zCanv =  new TCanvas("zCanvClus","RecPoint Z Positions",size+20,10,size,size);
90   zCanv->cd();
91   zGlob->Draw();
92
93   TCanvas *rCanv =  new TCanvas("rCanvClus","RecPoint R Positions",size+20,10,size,size);
94   rCanv->cd();
95   rGlob->Draw();
96
97 }