]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/macros/QA/check_radii.C
Starting a collection of QA/Comparison macros
[u/mrichter/AliRoot.git] / ITS / UPGRADE / macros / QA / check_radii.C
1
2 #if !defined(__CINT__) || defined(__MAKECINT__)
3   #include <TMath.h>
4   #include <TError.h>
5   #include <Riostream.h>
6   #include <TH1F.h>
7   #include <TH2F.h>
8   #include <TTree.h>
9   #include <TParticle.h>
10   #include <TCanvas.h>
11   #include <TLine.h>
12   #include <TText.h>
13   #include <TBenchmark.h>
14   #include <TStyle.h>
15   #include <TFile.h>
16   #include <TROOT.h>
17   #include <TNtuple.h>
18   #include <TEllipse.h>
19   #include <TGeoManager.h>
20
21   #include "AliStack.h"
22   #include "AliHeader.h"
23   #include "AliTrackReference.h"
24   #include "AliRunLoader.h"
25   #include "AliRun.h"
26   #include "AliESDEvent.h"
27   #include "AliESDtrack.h"
28
29   #include "UPGRADE/AliITSUClusterPix.h"
30   #include "UPGRADE/AliITSULoader.h"
31   #include "UPGRADE/AliITSUGeomTGeo.h"
32 #endif
33
34
35 Int_t check_radii(const Char_t *dir=".") {
36    TFile *f=TFile::Open("xyz.root","recreate");
37    TNtuple *nt=new TNtuple("nt","my ntuple","x:y:z");
38
39    if (gAlice) { 
40        delete AliRunLoader::Instance();
41        delete gAlice;//if everything was OK here it is already NULL
42        gAlice = 0x0;
43    }
44
45    Char_t fname[100];
46    sprintf(fname,"%s/galice.root",dir);
47
48    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
49    if (!rl) {
50       ::Error("GoodTracksITS","Can't start session !");
51       return 1;
52    }
53
54    rl->LoadgAlice();
55    rl->LoadHeader();
56    rl->LoadKinematics();
57
58    AliITSULoader* itsl = (AliITSULoader*)rl->GetLoader("ITSLoader");
59    if (itsl == 0x0) {
60        ::Error("GoodTracksITS","Can not find the ITSLoader");
61        delete rl;
62        return 4;
63    }
64    itsl->LoadRecPoints();
65   
66    TGeoManager::Import("geometry.root");
67    AliITSUGeomTGeo* gm = new AliITSUGeomTGeo(kTRUE,kTRUE);
68    AliITSUClusterPix::SetGeom(gm);
69
70    Int_t nev=rl->GetNumberOfEvents();
71    ::Info("GoodTracksITS","Number of events : %d\n",nev);  
72
73    //********  Loop over generated events 
74    for (Int_t e=0; e<nev; e++) {
75      Int_t k;
76
77      rl->GetEvent(e);
78
79      TTree *cTree=itsl->TreeR();
80      if (!cTree) {
81         ::Error("GoodTracksITS","Can't get the cluster tree !"); 
82         delete rl;
83         return 8;
84      }
85
86      const Int_t nLayers=7;
87      TBranch *branch[nLayers];
88      TClonesArray clusters[nLayers];
89      for (Int_t layer=0; layer<nLayers; layer++) {
90        TClonesArray *ptr = 
91        new(clusters+layer) TClonesArray("AliITSUClusterPix",1000);
92        Char_t bname[33];
93        sprintf(bname,"ITSRecPoints%d\0",layer);
94        branch[layer]=cTree->GetBranch(bname);
95        if (!branch[layer]) {
96           ::Error("GoodTracksITS","Can't get the clusters branch !"); 
97           delete rl;
98           return 9;
99        }
100        branch[layer]->SetAddress(&ptr);
101      }
102
103      Int_t entr=(Int_t)cTree->GetEntries();
104      for (k=0; k<entr; k++) {
105          cTree->GetEvent(k);
106          for (Int_t lay=0; lay<nLayers; lay++) {
107              Int_t ncl=clusters[lay].GetEntriesFast(); if (ncl==0) continue;
108              while (ncl--) {
109                 AliITSUClusterPix *pnt=
110                 (AliITSUClusterPix*)clusters[lay].UncheckedAt(ncl);
111                 Float_t g[3];
112                 pnt->GetGlobalXYZ(g);
113                 //pnt->GetLocalXYZ(g);
114                 //cout<<g[0]<<' '<<g[1]<<' '<<g[2]<<endl;
115                 nt->Fill(g);
116              }
117              clusters[lay].Clear();
118          }
119      }
120
121
122    } //*** end of the loop over generated events
123
124    nt->Draw("y:x");
125    Float_t rmin=19.44;
126    TEllipse *ellipse = new TEllipse(0,0,rmin,rmin,0,360,0);
127    ellipse->SetFillStyle(0);
128    ellipse->SetLineColor(2);
129    ellipse->Draw();
130    
131    Float_t rmax=19.77;
132    ellipse = new TEllipse(0,0,rmax,rmax,0,360,0);
133    ellipse->SetFillStyle(0);
134    ellipse->SetLineColor(4);
135    ellipse->Draw();
136
137    delete rl;
138    return 0;
139 }
140
141