]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ShowITSRecPoints.C
example on how to run alignment with MillePede - ITSAlignMille.C - and a configuratio...
[u/mrichter/AliRoot.git] / ITS / ShowITSRecPoints.C
CommitLineData
29d63bb9 1#if !defined(__CINT__) || defined(__MAKECINT__)
2#include <TCanvas.h>
3#include <TClassTable.h>
4#include <TGraph.h>
5#include <TGraph2D.h>
6#include <TGeoManager.h>
7#include <TH1.h>
8#include <TInterpreter.h>
9#include "AliCDBManager.h"
10#include "AliGeomManager.h"
11#include "AliHeader.h"
12#include "AliITS.h"
13#include "AliITSDetTypeRec.h"
14#include "AliITSgeom.h"
15#include "AliITSRecPoint.h"
16#include "AliRun.h"
17#endif
18
19Int_t ShowITSRecPoints(){
20 ///////////////////////////////////////////////////////////////////////
21 // Macro to check clusters in the 6 ITS layers //
22 // Provides: //
23 // 6 canvases with 9 plots each (1 canvas for each layer) //
24 // 3 canvases with cluster XY coordinates for the first 3 events //
25 ///////////////////////////////////////////////////////////////////////
26
27 if (gClassTable->GetID("AliRun") < 0) {
28 gInterpreter->ExecuteMacro("loadlibs.C");
29 }
29d63bb9 30 // Set OCDB if needed
31 AliCDBManager* man = AliCDBManager::Instance();
32 if (!man->IsDefaultStorageSet()) {
33 printf("Setting a local default storage and run number 0\n");
34 man->SetDefaultStorage("local://$ALICE_ROOT");
35 man->SetRun(0);
36 }
37 else {
38 printf("Using deafult storage \n");
39 }
40
41 // retrives geometry
42 if(!gGeoManager){
43 AliGeomManager::LoadGeometry("geometry.root");
44 }
45
46 AliRunLoader* rl = AliRunLoader::Open("galice.root");
47 if (rl == 0x0){
48 cerr<<"Can not open session RL=NULL"<< endl;
49 return -1;
50 }
682085e2 51 Int_t retval = rl->LoadHeader();
29d63bb9 52 if (retval){
53 cerr<<"LoadHeader returned error"<<endl;
54 return -1;
55 }
56
57
58 AliITSLoader* ITSloader = (AliITSLoader*) rl->GetLoader("ITSLoader");
59 if(!ITSloader){
60 cerr<<"ITS loader not found"<<endl;
61 return -1;
62 }
63 ITSloader->LoadRecPoints("read");
682085e2 64 AliITSgeom *geom = ITSloader->GetITSgeom();
29d63bb9 65
66 Float_t cluglo[3]={0.,0.,0.};
29d63bb9 67 AliITSDetTypeRec* detTypeRec = new AliITSDetTypeRec();
682085e2 68 detTypeRec->SetITSgeom(geom);
29d63bb9 69 detTypeRec->SetDefaults();
70
71 Int_t modmin=geom->GetStartDet(0);
72 Int_t modmax=geom->GetLastDet(2);
73 Int_t totmod=modmax-modmin;
74 Float_t xlim[6]={4.5,7.5,16.,26.,40.,45.};
75 Float_t zlim[6]={15.,15.,22.,30.,45.,55.};
76
77 TH1F* hlayer=new TH1F("hlayer","",6,0.5,6.5);
78 TH1F** hmod=new TH1F*[6];
79 TH1F** hxl=new TH1F*[6];
80 TH1F** hzl=new TH1F*[6];
81 TH1F** hxg=new TH1F*[6];
82 TH1F** hyg=new TH1F*[6];
83 TH1F** hzg=new TH1F*[6];
84 TH1F** hr=new TH1F*[6];
85 TH1F** hphi=new TH1F*[6];
86 TH1F** hq=new TH1F*[6];
87
88 Char_t name[10];
89 for(Int_t iLay=0;iLay<6;iLay++){
90 sprintf(name,"hmod%d",iLay+1);
91 hmod[iLay]=new TH1F(name,"",totmod,modmin-0.5,modmax+0.5);
92 hmod[iLay]->GetXaxis()->SetTitle("Module");
93 hmod[iLay]->GetXaxis()->CenterTitle();
94 sprintf(name,"hxloc%d",iLay+1);
95 hxl[iLay]=new TH1F(name,"",100,-4.,4.);
96 hxl[iLay]->GetXaxis()->SetTitle("Xloc (cm)");
97 hxl[iLay]->GetXaxis()->CenterTitle();
98 sprintf(name,"hzloc%d",iLay+1);
99 hzl[iLay]=new TH1F(name,"",100,-4.,4.);
100 hzl[iLay]->GetXaxis()->SetTitle("Zloc (cm)");
101 hzl[iLay]->GetXaxis()->CenterTitle();
102 sprintf(name,"hxgl%d",iLay+1);
103 hxg[iLay]=new TH1F(name,"",100,-xlim[iLay],xlim[iLay]);
104 hxg[iLay]->GetXaxis()->SetTitle("Xglob (cm)");
105 hxg[iLay]->GetXaxis()->CenterTitle();
106 sprintf(name,"hygl%d",iLay+1);
107 hyg[iLay]=new TH1F(name,"",100,-xlim[iLay],xlim[iLay]);
108 hyg[iLay]->GetXaxis()->SetTitle("Yglob (cm)");
109 hyg[iLay]->GetXaxis()->CenterTitle();
110 sprintf(name,"hzgl%d",iLay+1);
111 hzg[iLay]=new TH1F(name,"",100,-zlim[iLay],zlim[iLay]);
112 hzg[iLay]->GetXaxis()->SetTitle("Zglob (cm)");
113 hzg[iLay]->GetXaxis()->CenterTitle();
114 sprintf(name,"hr%d",iLay+1);
115 hr[iLay]=new TH1F(name,"",100,0.,50.);
116 hr[iLay]->GetXaxis()->SetTitle("r (cm)");
117 hr[iLay]->GetXaxis()->CenterTitle();
118 sprintf(name,"hphi%d",iLay+1);
119 hphi[iLay]=new TH1F(name,"",100,-TMath::Pi(),TMath::Pi());
120 hphi[iLay]->GetXaxis()->SetTitle("#varphi (rad)");
121 hphi[iLay]->GetXaxis()->CenterTitle();
122 sprintf(name,"hq%d",iLay+1);
123 hq[iLay]=new TH1F(name,"",100,0.,300.);
124 hq[iLay]->GetXaxis()->SetTitle("Charge (keV)");
125 hq[iLay]->GetXaxis()->CenterTitle();
126 }
127
128 TGraph **gpts=new TGraph*[3];
129 TGraph2D **gpts3d=new TGraph2D*[3];
130 for(Int_t i=0;i<3;i++){
131 gpts[i]=new TGraph(0);
132 gpts3d[i]=new TGraph2D(0);
133 }
134 Int_t totev=rl->GetNumberOfEvents();
135 printf("Total Number of events = %d\n",totev);
136
137 for(Int_t iev=0;iev<totev;iev++){
138 rl->GetEvent(iev);
139 TTree *TR = ITSloader->TreeR();
140 TClonesArray *ITSrec = detTypeRec->RecPoints();
141 TBranch *branch = 0;
142 if(TR && ITSrec){
143 branch = ITSloader->TreeR()->GetBranch("ITSRecPoints");
144 if(branch)branch->SetAddress(&ITSrec);
145 }
146 Int_t nparticles = rl->GetHeader()->GetNtrack();
147 cout<<"Event #"<<iev<<" #Particles="<<nparticles<<endl;
148
149
150 Int_t ipt=0;
151 for(Int_t subd=0;subd<3;subd++){
152
153 Int_t first = geom->GetStartDet(subd);
154 Int_t last = geom->GetLastDet(subd);
155
156 for (Int_t mod=first; mod<=last; mod++){
157 detTypeRec->ResetRecPoints();
158 branch->GetEvent(mod);
159 Int_t nrecp = ITSrec->GetEntries();
160 if(nrecp>0){
161 for(Int_t irec=0;irec<nrecp;irec++) {
162 AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(irec);
163 Int_t lay=recp->GetLayer();
164 hlayer->Fill(lay);
165 recp->GetGlobalXYZ(cluglo);
166 Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]);
167 Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);
168 if(iev<3){
169 gpts[iev]->SetPoint(ipt,cluglo[0],cluglo[1]);
170 gpts3d[iev]->SetPoint(ipt,cluglo[0],cluglo[1],cluglo[2]);
171 ipt++;
172 }
173 hmod[lay]->Fill(mod);
174 hzl[lay]->Fill(recp->GetDetLocalZ());
175 hxl[lay]->Fill(recp->GetDetLocalX());
176 hzg[lay]->Fill(cluglo[2]);
177 hyg[lay]->Fill(cluglo[1]);
178 hxg[lay]->Fill(cluglo[0]);
179 hr[lay]->Fill(rad);
180 hphi[lay]->Fill(phi);
181 hq[lay]->Fill(recp->GetQ());
182 }
183 }
184 }
185 }
186 }
187 TCanvas **c=new TCanvas*[6];
188 Char_t ctit[30];
189 for(Int_t iLay=0;iLay<6;iLay++){
190 sprintf(name,"can%d",iLay+1);
191 sprintf(ctit,"Layer %d",iLay+1);
192 c[iLay]=new TCanvas(name,ctit,1200,900);
193 c[iLay]->Divide(3,3,0.001,0.001);
194 c[iLay]->cd(1);
195 hmod[iLay]->Draw();
196 c[iLay]->cd(2);
197 hxl[iLay]->Draw();
198 c[iLay]->cd(3);
199 hzl[iLay]->Draw();
200 c[iLay]->cd(4);
201 hxg[iLay]->Draw();
202 c[iLay]->cd(5);
203 hyg[iLay]->Draw();
204 c[iLay]->cd(6);
205 hzg[iLay]->Draw();
206 c[iLay]->cd(7);
207 hr[iLay]->Draw();
208 c[iLay]->cd(8);
209 hphi[iLay]->Draw();
210 c[iLay]->cd(9);
211 hq[iLay]->Draw();
212 }
213
214 TCanvas *cev0;
215 cev0=new TCanvas("cev0","Event 0",600,600);
216 gpts[0]->SetMarkerStyle(7);
217 gpts[0]->SetTitle(0);
218 gpts[0]->Draw("AP");
219
220 TCanvas *cev1;
221 cev1=new TCanvas("cev1","Event 1",600,600);
222 gpts[1]->SetMarkerStyle(7);
223 gpts[1]->SetTitle(0);
224 gpts[1]->Draw("AP");
225
226 TCanvas *cev2;
227 cev2=new TCanvas("cev2","Event 2",600,600);
228 gpts[2]->SetMarkerStyle(7);
229 gpts[2]->SetTitle(0);
230 gpts[2]->Draw("AP");
231
232
233 return 0;
234}
235
236