]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/DrawITSAlignObj.C
Add -Wconversion to Mac warnings
[u/mrichter/AliRoot.git] / ITS / DrawITSAlignObj.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TCanvas.h>
3 #include <TClonesArray.h>
4 #include <TFile.h>
5 #include <TGeoManager.h>
6 #include <TH1F.h>
7 #include <TString.h>
8 #include <Riostream.h>
9 #include "AliAlignObj.h"
10 #endif
11 void DrawITSAlignObj(Bool_t local=kFALSE) {  //
12   // Draw distribution of alignment constants
13   // Set TOCDB and STORAGE environment variables to run the macro
14   // on an AliCDBEntry instead of on a file
15   //
16
17   const char* filenameOut="ITSfullModuleMisalignment.root";
18   TClonesArray* ar = 0;
19   
20   AliCDBManager* cdb = AliCDBManager::Instance();
21
22   // Activate CDB storage and load geometry from CDB
23   if( TString(gSystem->Getenv("TOCDB")) == TString("kTRUE") ){
24     TString Storage = gSystem->Getenv("STORAGE");
25     if(!Storage.BeginsWith("local://") && !Storage.BeginsWith("alien://")) {
26       Error(macroname,"STORAGE variable set to %s is not valid. Exiting\n",Storage.Data());
27       return;
28     }
29     cdb->SetDefaultStorage(Storage.Data());
30     cdb->SetRun(0);
31     AliCDBEntry *entry = cdb->Get("GRP/Geometry/Data");
32     if(!entry) Fatal("Could not get the specified CDB entry!");
33     entry->SetOwner(0);
34     TGeoManager* geom = (TGeoManager*) entry->GetObject();
35     AliGeomManager::SetGeometry(geom);
36     AliCDBEntry* eItsAlign = cdb->Get("ITS/Align/Data");
37     ar = (TClonesArray*)eItsAlign->GetObject();
38     if(!ar) {
39       Fatal("Could not get the alignment-objects array from the CDB entry!");
40       return;
41     }
42   }else{
43     cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
44     cdb->SetRun(0);
45     AliGeomManager::LoadGeometry("geometry.root"); //load geom from default CDB storage
46     const char* filename="ITSfullMisalignment.root";
47     TFile* f = TFile::Open(filename);
48     if(!f) {cerr<<"cannot open input file\n"; return;}
49     ar = (TClonesArray*)f->Get("ITSAlignObjs");
50     if(!ar) {
51       Fatal("Could not get the alignment-objects array from the file %s!", filename);
52       return;
53     }
54   }    
55                   
56   TCanvas  *cSPDsector = new TCanvas("cSPDsector","SPD sectors alignobjs",0,0,800,600);
57   cSPDsector->Divide(3,2);                  
58   TH1F* hxSPDsector = new TH1F("hxSPDsector","x in SPD sectors",100,-0.1,0.1);
59   hxSPDsector->SetXTitle("[cm]");
60   TH1F* hySPDsector = new TH1F("hySPDsector","y in SPD sectors",100,-0.1,0.1);
61   hySPDsector->SetXTitle("[cm]");
62   TH1F* hzSPDsector = new TH1F("hzSPDsector","z in SPD sectors",100,-0.1,0.1);
63   hzSPDsector->SetXTitle("[cm]");
64   TH1F* hrxSPDsector = new TH1F("hrxSPDsector","rx in SPD sectors",100,-0.5,0.5);
65   hrxSPDsector->SetXTitle("[deg]");
66   TH1F* hrySPDsector = new TH1F("hrySPDsector","ry in SPD sectors",100,-0.5,0.5);
67   hrySPDsector->SetXTitle("[deg]");
68   TH1F* hrzSPDsector = new TH1F("hrzSPDsector","rz in SPD sectors",100,-0.5,0.5);
69   hrzSPDsector->SetXTitle("[deg]");
70
71   TCanvas  *cSPDhalfstave = new TCanvas("cSPDhalfstave","SPD halfstaves alignobjs",0,0,800,600);
72   cSPDhalfstave->Divide(3,2);               
73   TH1F* hxSPDhalfstave = new TH1F("hxSPDhalfstave","x in SPD halfstaves",100,-0.01,0.01);
74   hxSPDhalfstave->SetXTitle("[cm]");
75   TH1F* hySPDhalfstave = new TH1F("hySPDhalfstave","y in SPD halfstaves",100,-0.01,0.01);
76   hySPDhalfstave->SetXTitle("[cm]");
77   TH1F* hzSPDhalfstave = new TH1F("hzSPDhalfstave","z in SPD halfstaves",100,-0.01,0.01);
78   hzSPDhalfstave->SetXTitle("[cm]");
79   TH1F* hrxSPDhalfstave = new TH1F("hrxSPDhalfstave","rx in SPD halfstaves",100,-0.1,0.1);
80   hrxSPDhalfstave->SetXTitle("[deg]");
81   TH1F* hrySPDhalfstave = new TH1F("hrySPDhalfstave","ry in SPD halfstaves",100,-0.1,0.1);
82   hrySPDhalfstave->SetXTitle("[deg]");
83   TH1F* hrzSPDhalfstave = new TH1F("hrzSPDhalfstave","rz in SPD halfstaves",100,-0.5,0.5);
84   hrzSPDhalfstave->SetXTitle("[deg]");
85
86   TCanvas  *cSPDladder = new TCanvas("cSPDladder","SPD ladders alignobjs",0,0,800,600);
87   cSPDladder->Divide(3,2);                  
88   TH1F* hxSPDladder = new TH1F("hxSPDladder","x in SPD ladders",100,-0.01,0.01);
89   hxSPDladder->SetXTitle("[cm]");
90   TH1F* hySPDladder = new TH1F("hySPDladder","y in SPD ladders",100,-0.01,0.01);
91   hySPDladder->SetXTitle("[cm]");
92   TH1F* hzSPDladder = new TH1F("hzSPDladder","z in SPD ladders",100,-0.01,0.01);
93   hzSPDladder->SetXTitle("[cm]");
94   TH1F* hrxSPDladder = new TH1F("hrxSPDladder","rx in SPD ladders",100,-0.01,0.01);
95   hrxSPDladder->SetXTitle("[deg]");
96   TH1F* hrySPDladder = new TH1F("hrySPDladder","ry in SPD ladders",100,-0.01,0.01);
97   hrySPDladder->SetXTitle("[deg]");
98   TH1F* hrzSPDladder = new TH1F("hrzSPDladder","rz in SPD ladders",100,-0.01,0.01);
99   hrzSPDladder->SetXTitle("[deg]");
100
101   TCanvas  *cSDDladder = new TCanvas("cSDDladder","SDD ladders alignobjs",0,0,800,600);
102   cSDDladder->Divide(3,2);                  
103   TH1F* hxSDDladder = new TH1F("hxSDDladder","x in SDD ladders",100,-0.01,0.01);
104   hxSDDladder->SetXTitle("[cm]");
105   TH1F* hySDDladder = new TH1F("hySDDladder","y in SDD ladders",100,-0.01,0.01);
106   hySDDladder->SetXTitle("[cm]");
107   TH1F* hzSDDladder = new TH1F("hzSDDladder","z in SDD ladders",100,-0.01,0.01);
108   hzSDDladder->SetXTitle("[cm]");
109   TH1F* hrxSDDladder = new TH1F("hrxSDDladder","rx in SDD ladders",100,-0.01,0.01);
110   hrxSDDladder->SetXTitle("[deg]");
111   TH1F* hrySDDladder = new TH1F("hrySDDladder","ry in SDD ladders",100,-0.01,0.01);
112   hrySDDladder->SetXTitle("[deg]");
113   TH1F* hrzSDDladder = new TH1F("hrzSDDladder","rz in SDD ladders",100,-0.01,0.01);
114   hrzSDDladder->SetXTitle("[deg]");
115
116   TCanvas  *cSDDsensor = new TCanvas("cSDDsensor","SDD sensors alignobjs",0,0,800,600);
117   cSDDsensor->Divide(3,2);                  
118   TH1F* hxSDDsensor = new TH1F("hxSDDsensor","x in SDD sensors",100,-0.01,0.01);
119   hxSDDsensor->SetXTitle("[cm]");
120   TH1F* hySDDsensor = new TH1F("hySDDsensor","y in SDD sensors",100,-0.01,0.01);
121   hySDDsensor->SetXTitle("[cm]");
122   TH1F* hzSDDsensor = new TH1F("hzSDDsensor","z in SDD sensors",100,-0.01,0.01);
123   hzSDDsensor->SetXTitle("[cm]");
124   TH1F* hrxSDDsensor = new TH1F("hrxSDDsensor","rx in SDD sensors",100,-0.01,0.01);
125   hrxSDDsensor->SetXTitle("[deg]");
126   TH1F* hrySDDsensor = new TH1F("hrySDDsensor","ry in SDD sensors",100,-0.01,0.01);
127   hrySDDsensor->SetXTitle("[deg]");
128   TH1F* hrzSDDsensor = new TH1F("hrzSDDsensor","rz in SDD sensors",100,-0.01,0.01);
129   hrzSDDsensor->SetXTitle("[deg]");
130
131
132   TCanvas  *cSSDladder = new TCanvas("cSSDladder","SSD ladders alignobjs",0,0,800,600);
133   cSSDladder->Divide(3,2);                  
134   TH1F* hxSSDladder = new TH1F("hxSSDladder","x in SSD ladders",100,-0.01,0.01);
135   hxSSDladder->SetXTitle("[cm]");
136   TH1F* hySSDladder = new TH1F("hySSDladder","y in SSD ladders",100,-0.01,0.01);
137   hySSDladder->SetXTitle("[cm]");
138   TH1F* hzSSDladder = new TH1F("hzSSDladder","z in SSD ladders",100,-0.01,0.01);
139   hzSSDladder->SetXTitle("[cm]");
140   TH1F* hrxSSDladder = new TH1F("hrxSSDladder","rx in SSD ladders",100,-0.01,0.01);
141   hrxSSDladder->SetXTitle("[deg]");
142   TH1F* hrySSDladder = new TH1F("hrySSDladder","ry in SSD ladders",100,-0.01,0.01);
143   hrySSDladder->SetXTitle("[deg]");
144   TH1F* hrzSSDladder = new TH1F("hrzSSDladder","rz in SSD ladders",100,-0.01,0.01);
145   hrzSSDladder->SetXTitle("[deg]");
146
147   TCanvas  *cSSDsensor = new TCanvas("cSSDsensor","SSD sensors alignobjs",0,0,800,600);
148   cSSDsensor->Divide(3,2);                  
149   TH1F* hxSSDsensor = new TH1F("hxSSDsensor","x in SSD sensors",100,-0.01,0.01);
150   hxSSDsensor->SetXTitle("[cm]");
151   TH1F* hySSDsensor = new TH1F("hySSDsensor","y in SSD sensors",100,-0.01,0.01);
152   hySSDsensor->SetXTitle("[cm]");
153   TH1F* hzSSDsensor = new TH1F("hzSSDsensor","z in SSD sensors",100,-0.01,0.01);
154   hzSSDsensor->SetXTitle("[cm]");
155   TH1F* hrxSSDsensor = new TH1F("hrxSSDsensor","rx in SSD sensors",100,-0.01,0.01);
156   hrxSSDsensor->SetXTitle("[deg]");
157   TH1F* hrySSDsensor = new TH1F("hrySSDsensor","ry in SSD sensors",100,-0.01,0.01);
158   hrySSDsensor->SetXTitle("[deg]");
159   TH1F* hrzSSDsensor = new TH1F("hrzSSDsensor","rz in SSD sensors",100,-0.01,0.01);
160   hrzSSDsensor->SetXTitle("[deg]");
161
162
163   Int_t entries = ar->GetEntriesFast();
164   Printf("number of alignment objects: %d",entries);
165   //Bool_t overlaps;
166   if(!AliGeomManager::GetGeometry()) return;
167   AliGeomManager::ApplyAlignObjsToGeom(*ar); //,overlaps);
168
169   AliAlignObj* a;
170   Option_t *opt = NULL;
171   Double_t tr[3];
172   Double_t rot[3];
173   TGeoHMatrix delta;
174
175
176   TClonesArray *arrayOut = new TClonesArray("AliAlignObjParams",4000);
177   TClonesArray &alobjOut = *arrayOut;
178   Int_t j=0;
179
180   for(Int_t i=0; i<entries; i++){
181     a=(AliAlignObj*)ar->UncheckedAt(i);
182     TString symName = a->GetSymName();
183     UShort_t volUID = a->GetVolUID();
184     //printf("VolId %d    %s\n",volUID,symName.Data());
185     
186     if(local) {
187       a->GetLocalPars(tr,rot);
188     } else {
189       a->GetPars(tr,rot);
190     }
191
192     AliGeomManager::GetDeltaForBranch(*a,delta);
193     //delta.Print();
194
195     // create AliAlignObjParam with total delta
196     if(i==0 || volUID!=0) {
197       new(alobjOut[j]) AliAlignObjParams(symName.Data(),volUID,delta,kTRUE);
198       j++;
199     }
200
201
202     // plots
203     //
204     if(!symName.Contains("SPD") && !symName.Contains("SDD") && !symName.Contains("SSD")) {
205       a->Print(opt);
206     }
207     if(symName.Contains("SPD") && symName.Contains("Sector") && !symName.Contains("HalfStave")) {
208       hxSPDsector->Fill(tr[0]);
209       hySPDsector->Fill(tr[1]);
210       hzSPDsector->Fill(tr[2]);
211       hrxSPDsector->Fill(rot[0]);
212       hrySPDsector->Fill(rot[1]);
213       hrzSPDsector->Fill(rot[2]);
214     }
215     if(symName.Contains("SPD") && symName.Contains("Sector") && symName.Contains("HalfStave") && !symName.Contains("Ladder")) {
216       hxSPDhalfstave->Fill(tr[0]);
217       hySPDhalfstave->Fill(tr[1]);
218       hzSPDhalfstave->Fill(tr[2]);
219       hrxSPDhalfstave->Fill(rot[0]);
220       hrySPDhalfstave->Fill(rot[1]);
221       hrzSPDhalfstave->Fill(rot[2]);
222     }
223     if(symName.Contains("SPD") && symName.Contains("Sector") && symName.Contains("HalfStave") && symName.Contains("Ladder")) {
224       hxSPDladder->Fill(tr[0]);
225       hySPDladder->Fill(tr[1]);
226       hzSPDladder->Fill(tr[2]);
227       hrxSPDladder->Fill(rot[0]);
228       hrySPDladder->Fill(rot[1]);
229       hrzSPDladder->Fill(rot[2]);
230     }
231     if(symName.Contains("SDD") && symName.Contains("Ladder") && !symName.Contains("Sensor")) {
232       hxSDDladder->Fill(tr[0]);
233       hySDDladder->Fill(tr[1]);
234       hzSDDladder->Fill(tr[2]);
235       hrxSDDladder->Fill(rot[0]);
236       hrySDDladder->Fill(rot[1]);
237       hrzSDDladder->Fill(rot[2]);
238     }
239     if(symName.Contains("SDD") && symName.Contains("Ladder") && symName.Contains("Sensor")) {
240       hxSDDsensor->Fill(tr[0]);
241       hySDDsensor->Fill(tr[1]);
242       hzSDDsensor->Fill(tr[2]);
243       hrxSDDsensor->Fill(rot[0]);
244       hrySDDsensor->Fill(rot[1]);
245       hrzSDDsensor->Fill(rot[2]);
246     }
247     if(symName.Contains("SSD") && symName.Contains("Ladder") && !symName.Contains("Sensor")) {
248       hxSSDladder->Fill(tr[0]);
249       hySSDladder->Fill(tr[1]);
250       hzSSDladder->Fill(tr[2]);
251       hrxSSDladder->Fill(rot[0]);
252       hrySSDladder->Fill(rot[1]);
253       hrzSSDladder->Fill(rot[2]);
254     }
255     if(symName.Contains("SSD") && symName.Contains("Ladder") && symName.Contains("Sensor")) {
256       hxSSDsensor->Fill(tr[0]);
257       hySSDsensor->Fill(tr[1]);
258       hzSSDsensor->Fill(tr[2]);
259       hrxSSDsensor->Fill(rot[0]);
260       hrySSDsensor->Fill(rot[1]);
261       hrzSSDsensor->Fill(rot[2]);
262     }
263
264
265   }
266
267
268   cSPDsector->cd(1);
269   hxSPDsector->Draw();
270   cSPDsector->cd(2);
271   hySPDsector->Draw();
272   cSPDsector->cd(3);
273   hzSPDsector->Draw();
274   cSPDsector->cd(4);
275   hrxSPDsector->Draw();
276   cSPDsector->cd(5);
277   hrySPDsector->Draw();
278   cSPDsector->cd(6);
279   hrzSPDsector->Draw();
280
281   cSPDhalfstave->cd(1);
282   hxSPDhalfstave->Draw();
283   cSPDhalfstave->cd(2);
284   hySPDhalfstave->Draw();
285   cSPDhalfstave->cd(3);
286   hzSPDhalfstave->Draw();
287   cSPDhalfstave->cd(4);
288   hrxSPDhalfstave->Draw();
289   cSPDhalfstave->cd(5);
290   hrySPDhalfstave->Draw();
291   cSPDhalfstave->cd(6);
292   hrzSPDhalfstave->Draw();
293
294   cSPDladder->cd(1);
295   hxSPDladder->Draw();
296   cSPDladder->cd(2);
297   hySPDladder->Draw();
298   cSPDladder->cd(3);
299   hzSPDladder->Draw();
300   cSPDladder->cd(4);
301   hrxSPDladder->Draw();
302   cSPDladder->cd(5);
303   hrySPDladder->Draw();
304   cSPDladder->cd(6);
305   hrzSPDladder->Draw();
306
307
308   cSDDladder->cd(1);
309   hxSDDladder->Draw();
310   cSDDladder->cd(2);
311   hySDDladder->Draw();
312   cSDDladder->cd(3);
313   hzSDDladder->Draw();
314   cSDDladder->cd(4);
315   hrxSDDladder->Draw();
316   cSDDladder->cd(5);
317   hrySDDladder->Draw();
318   cSDDladder->cd(6);
319   hrzSDDladder->Draw();
320  
321   cSDDsensor->cd(1);
322   hxSDDsensor->Draw();
323   cSDDsensor->cd(2);
324   hySDDsensor->Draw();
325   cSDDsensor->cd(3);
326   hzSDDsensor->Draw();
327   cSDDsensor->cd(4);
328   hrxSDDsensor->Draw();
329   cSDDsensor->cd(5);
330   hrySDDsensor->Draw();
331   cSDDsensor->cd(6);
332   hrzSDDsensor->Draw();
333
334   cSSDladder->cd(1);
335   hxSSDladder->Draw();
336   cSSDladder->cd(2);
337   hySSDladder->Draw();
338   cSSDladder->cd(3);
339   hzSSDladder->Draw();
340   cSSDladder->cd(4);
341   hrxSSDladder->Draw();
342   cSSDladder->cd(5);
343   hrySSDladder->Draw();
344   cSSDladder->cd(6);
345   hrzSSDladder->Draw();
346  
347   cSSDsensor->cd(1);
348   hxSSDsensor->Draw();
349   cSSDsensor->cd(2);
350   hySSDsensor->Draw();
351   cSSDsensor->cd(3);
352   hzSSDsensor->Draw();
353   cSSDsensor->cd(4);
354   hrxSSDsensor->Draw();
355   cSSDsensor->cd(5);
356   hrySSDsensor->Draw();
357   cSSDsensor->cd(6);
358   hrzSSDsensor->Draw();
359
360
361   TFile fout(filenameOut,"RECREATE");
362   fout.cd();
363   fout.WriteObject(arrayOut,"ITSAlignObjs","kSingleKey");
364   fout.Close();
365   
366   return;
367 }
368