]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/readClusters.C
Install macros
[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   AliITSURecoDet *its = new AliITSURecoDet(gm, "ITSinterface");
25   its->CreateClusterArrays();
26   //
27   TH2F *xyGlob = new TH2F("xyGlob"," X - Y Global coordinates ",500,-50,50,500,-50,50);
28   xyGlob->SetXTitle("cm"); 
29   xyGlob->SetMarkerStyle(7); 
30   TH1F *zGlob  = new TH1F("zGlob", " Z Global coordinates ",200, -50,50 );
31   zGlob->SetXTitle("cm"); 
32   //
33   TH1F* rGlob = new TH1F("rGlob","R global", 5000, 0,50.);
34
35   TTree * cluTree = 0x0;
36   
37   int nevTot = (Int_t)runLoader->GetNumberOfEvents();
38   printf("N Events : %i \n",nevTot);
39   evStart = evStart<nevTot ? evStart : nevTot-1;
40   if (evStart<0) evStart = 0;
41   //
42   int lastEv = nev<0 ? nevTot : evStart+nev;
43   if (lastEv > nevTot) lastEv = nevTot;
44   //
45   for (Int_t iEvent = evStart; iEvent < lastEv; iEvent++) {
46     printf("\n Event %i \n",iEvent);
47     runLoader->GetEvent(iEvent);
48     //   AliStack *stack = runLoader->Stack();
49     cluTree=dl->TreeR();
50     int nlr=0;
51     while(1) {
52       TBranch* br = cluTree->GetBranch(Form("ITSRecPoints%d",nlr));
53       if (!br) break;
54       br->SetAddress(its->GetLayerActive(nlr)->GetClustersAddress());
55       nlr++;
56     }    
57     printf(" tree entries: %d\n",cluTree->GetEntries());
58     cluTree->GetEntry(0);      
59     AliITSUClusterPix::SetSortMode( AliITSUClusterPix::SortModeIdTrkYZ());
60     for (int ilr=0;ilr<nlr;ilr++) {its->GetLayerActive(ilr)->GetClusters()->Sort();}
61     its->ProcessClusters();
62     //
63     for (int ilr=0;ilr<nlr;ilr++) {
64       AliITSURecoLayer* lr = its->GetLayerActive(ilr);
65       TClonesArray* clr = lr->GetClusters();
66       
67       int nClu = clr->GetEntries();
68       printf("Layer %d : %d clusters\n",ilr,nClu);
69       //
70       for (int icl=0;icl<nClu;icl++) {
71         AliITSUClusterPix *cl = (AliITSUClusterPix*)clr->At(icl);
72         printf("#%4d | ",icl); cl->Print("glo p");
73         Float_t loc[3];
74         cl->GetLocalXYZ(loc);
75         Float_t glob[3]; 
76         cl->GetGlobalXYZ(glob);
77         //printf("%d: mod %d: loc(%.4lf,%.4lf,%.4lf); glob(%.4lf,%.4lf,%.4lf); \n",icl,cl->GetVolumeId(),
78         //       loc[0],loc[1],loc[2],glob[0],glob[1],glob[2]);
79         xyGlob->Fill(glob[0],glob[1]);
80         zGlob->Fill(glob[2]);
81         rGlob->Fill(TMath::Sqrt(glob[0]*glob[0]+glob[1]*glob[1]));
82       }
83     }
84   }//event loop
85
86   Int_t size = 400;
87
88   TCanvas *xyCanv =  new TCanvas("xvCanvClus","RecPoint X-Y Positions",10,10,size,size);
89   xyCanv->cd();
90   xyGlob->Draw("P");
91
92   TCanvas *zCanv =  new TCanvas("zCanvClus","RecPoint Z Positions",size+20,10,size,size);
93   zCanv->cd();
94   zGlob->Draw();
95
96   TCanvas *rCanv =  new TCanvas("rCanvClus","RecPoint R Positions",size+20,10,size,size);
97   rCanv->cd();
98   rGlob->Draw();
99
100 }