]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONSurveyCh8L.C
added missing declaration in TriggerOnJet for AliHLTScalars
[u/mrichter/AliRoot.git] / MUON / MUONSurveyCh8L.C
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /// \ingroup macros
19 /// \file MUONSurveyCh8L.C
20 /// \brief Macro to process survey and photogrammetry data of chamber 8L
21 ///  
22 ///  Macro loads the survey data from .txt file using AliSurveyObj.
23 ///  Macro MUONSurveyUtil.C is then loaded.
24 ///
25 ///  The transformations of the slats are obatained in 2 steps:
26 ///  -  1. Fit a plane to the sticker targets -> psi, theta
27 ///  -  2. Using above psi in theta obtain xc, yc, zc and phi by solving 
28 ///        the equations from a local to global transformation of the
29 ///        fixed button targets
30 ///
31 ///  Various histograms are filled and printed for monitoring.
32 ///  MisAlignment object is then created.
33 /// 
34 /// \author Javier Castillo
35
36 #if !defined(__CINT__) || defined(__MAKECINT__)
37
38 #include "AliMUONGeometryTransformer.h"
39 #include "AliMUONGeometryMisAligner.h"
40
41 #include "AliSurveyObj.h"
42 #include "AliSurveyPoint.h"
43 #include "AliGeomManager.h"
44 #include "AliCDBManager.h"
45 #include "AliCDBMetaData.h"
46 #include "AliCDBId.h"
47
48 #include <TROOT.h>
49 #include <TGeoManager.h>
50 #include <TClonesArray.h>
51 #include <TObjArray.h>
52 #include <TObjString.h>
53 #include <TMath.h>
54 #include <TString.h>
55 #include <Riostream.h>
56 #include <TF2.h>
57 #include <TH2.h>
58 #include <TGraph2DErrors.h>
59 #include <TGraph.h>
60 #include <TVector3.h>
61 #include <TCanvas.h>
62 #include <TPad.h>
63 #include <TPostScript.h>
64 #include <TPaveLabel.h>
65 #include <TStyle.h>
66
67 #include <fstream>
68
69 #endif
70
71 Bool_t MatrixToAngles(const Double_t *rot, Double_t *angles);
72 Double_t eqPlane(Double_t *x, Double_t *par);
73 Double_t xpCenter(Double_t *x, Double_t *par);
74 Double_t xnCenter(Double_t *x, Double_t *par);
75 Double_t ypCenter(Double_t *x, Double_t *par);
76 Double_t ynCenter(Double_t *x, Double_t *par);
77 Double_t zpCenter(Double_t *x, Double_t *par);
78 Double_t znCenter(Double_t *x, Double_t *par);
79 Double_t phixpp(Double_t *x, Double_t *par);
80 Double_t phixpn(Double_t *x, Double_t *par);
81 Double_t phixnp(Double_t *x, Double_t *par);
82 Double_t phixnn(Double_t *x, Double_t *par);
83 Double_t phiypp(Double_t *x, Double_t *par);
84 Double_t phiypn(Double_t *x, Double_t *par);
85 Double_t phiynp(Double_t *x, Double_t *par);
86 Double_t phiynn(Double_t *x, Double_t *par);
87 AliMUONGeometryTransformer *ReAlign(const AliMUONGeometryTransformer * transformer, 
88                                     int rMod, TGeoCombiTrans deltaDetElemTransf[], Bool_t verbose);
89 void MUONSurveyCh8L() {
90   
91   char str[100];
92   char filename[100];
93   
94   int saveps = 1;
95   const int cWidth = (int)(700*(29./21.));
96   const int cHeight = 700;
97 //   const int   lineColor = 1;
98 //   const int   fillColor = 4;
99 //   const int   filetype  = 111; // portrait  
100   const int   filetype  = 112; // landscape  
101   
102   sprintf(filename,"surveyChamber8L.ps");
103   
104   Int_t nSlats = 13;
105
106   gROOT->LoadMacro("MUONSurveyUtil.C+");
107
108   AliSurveyObj *so = new AliSurveyObj();
109   
110   Int_t size = so->GetEntries();
111   printf("-> %d\n", size);
112   
113   so->FillFromLocalFile("Alice_MuonSystem_Chamber8LCavern_3561b.txt");
114   size = so->GetEntries();
115   printf("--> %d\n", size);
116
117   Printf("Title: \"%s\"", so->GetReportTitle().Data());
118   Printf("Date: \"%s\"", so->GetReportDate().Data());
119   Printf("Detector: \"%s\"", so->GetDetector().Data());
120   Printf("URL: \"%s\"", so->GetURL().Data());
121   Printf("Number: \"%d\"", so->GetReportNumber());
122   Printf("Version: \"%d\"", so->GetReportVersion());
123   Printf("Observations: \"%s\"", so->GetObservations().Data());
124   Printf("Coordinate System: \"%s\"", so->GetCoordSys().Data());
125   Printf("Measurement Units: \"%s\"", so->GetUnits().Data());
126   Printf("Nr Columns: \"%d\"", so->GetNrColumns());
127
128   TObjArray *colNames = so->GetColumnNames();
129   for (Int_t i = 0; i < colNames->GetEntries(); ++i)
130     Printf("  Column %d --> \"%s\"", i, ((TObjString *) colNames->At(i))->GetString().Data());
131
132   // Get Array of surveyed points
133   Printf("Points:");
134   TObjArray *points = so->GetData();
135   
136   for (Int_t i = 0; i < points->GetEntries(); ++i)
137     Printf("  Point %d --> \"%s\"  %s ", i, ((AliSurveyPoint *) points->At(i))->GetPointName().Data(), points->At(i)->GetName());
138
139
140   // Slats (#1 - #13) button targets local coordinates
141   Double_t lSBTLoc6[13][2][3] = {{{  -412.50,  0.0, -(11.75+ 8.20+20.00)},
142                                   {   412.50,  0.0, -(11.75+16.50+20.00)}},
143                                  {{  -612.50,  0.0,  (11.75+ 8.20+20.00)}, 
144                                   {   612.50,  0.0,  (11.75+16.50+20.00)}},
145                                  {{ - 812.50,  0.0, -(11.75+ 8.20+20.00)}, 
146                                   {   812.50,  0.0, -(11.75+16.50+20.00)}},
147                                  {{ -1012.50,  0.0,  (11.75+ 8.20+20.00)}, 
148                                   {  1012.50,  0.0,  (11.75+16.50+20.00)}},
149                                  {{ -1012.50,  0.0, -(11.75+ 8.20+20.00)}, 
150                                   {  1012.50,  0.0, -(11.75+16.50+20.00)}},
151                                  {{ -1212.50,  5.0,  (11.75+ 8.20+20.00)}, 
152                                   {  1212.50,  0.0,  (11.75+16.50+20.00)}},
153                                  {{ -1012.50,  0.0, -(11.75+ 8.20+20.00)}, 
154                                   {  1012.50,  0.0, -(11.75+16.50+20.00)}},
155                                  {{ -1212.50,  5.0,  (11.75+ 8.20+20.00)}, 
156                                   {  1212.50,  0.0,  (11.75+16.50+20.00)}},
157                                  {{ -1012.50,  0.0, -(11.75+ 8.20+20.00)}, 
158                                   {  1012.50,  0.0, -(11.75+16.50+20.00)}},
159                                  {{ -1012.50,  0.0,  (11.75+ 8.20+20.00)}, 
160                                   {  1012.50,  0.0,  (11.75+16.50+20.00)}},
161                                  {{  -812.50,  0.0, -(11.75+ 8.20+20.00)}, 
162                                   {   812.50,  0.0, -(11.75+16.50+20.00)}},
163                                  {{  -612.50,  0.0,  (11.75+ 8.20+20.00)}, 
164                                   {   612.50,  0.0,  (11.75+16.50+20.00)}},
165                                  {{  -412.50,  0.0, -(11.75+ 8.20+20.00)},
166                                   {   412.50,  0.0, -(11.75+16.50+20.00)}}};
167                                                     
168   
169   AliSurveyPoint *pointSST = 0;
170   AliSurveyPoint *pointCST = 0;
171   AliSurveyPoint *pointCPST = 0;
172   AliSurveyPoint **pointSBT = new AliSurveyPoint*[2];
173   
174   char sPointName[10] = "5000"; 
175   
176   // Print length of slats 
177   cout << "Slat lengths:" << endl;
178   TVector3 vTemp(0., 0., 0.);
179   TVector3 vSBT(0., 0., 0.);
180
181   for (Int_t iSlat=0; iSlat<nSlats; iSlat++){
182     // Get button targets survey points
183     vTemp.SetXYZ(0., 0., 0.);
184     if (iSlat+1<10) {
185       sprintf(sPointName,"60%d%d",iSlat+1,1);
186       pointSBT[0] = (AliSurveyPoint *)points->FindObject(sPointName);
187       if(!pointSBT[0]) {
188         cout << "Error! No button targets ... " << endl;
189         break;
190       }
191       vSBT.SetXYZ(pointSBT[0]->GetX(),pointSBT[0]->GetY(),pointSBT[0]->GetZ());
192       vTemp+=vSBT;
193       sprintf(sPointName,"60%d%d",iSlat+1,2);
194       pointSBT[1] = (AliSurveyPoint *)points->FindObject(sPointName);
195       if(!pointSBT[1]) {
196         cout << "Error! No button targets ... " << endl;
197         break;
198       }
199       vSBT.SetXYZ(pointSBT[1]->GetX(),pointSBT[1]->GetY(),pointSBT[1]->GetZ());
200       vTemp-=vSBT;
201     }
202     else {
203       sprintf(sPointName,"6%d%d",iSlat+1,1);
204       pointSBT[0] = (AliSurveyPoint *)points->FindObject(sPointName);
205       if(!pointSBT[0]) {
206         cout << "Error! No button targets ... " << endl;
207         break;
208       }
209       vSBT.SetXYZ(pointSBT[0]->GetX(),pointSBT[0]->GetY(),pointSBT[0]->GetZ());
210       vTemp+=vSBT;
211       sprintf(sPointName,"6%d%d",iSlat+1,2);
212       pointSBT[1] = (AliSurveyPoint *)points->FindObject(sPointName);
213       if(!pointSBT[1]) {
214         cout << "Error! No button targets ... " << endl;
215         break;
216       }
217       vSBT.SetXYZ(pointSBT[1]->GetX(),pointSBT[1]->GetY(),pointSBT[1]->GetZ());
218       vTemp-=vSBT;
219     }
220     cout << "Slat " << iSlat+1 << ": " << vTemp.Mag() << endl;
221   }
222  
223   // Histograms for monitoring
224   TH2F *hCPSTry = new TH2F("hCPSTry","hCPSTry",28,-600,2200,52,0,5200);
225   TH2F *hCPSTly = new TH2F("hCPSTly","hCPSTly",28,-600,2200,52,0,5200);
226   TH2F *hSSTry = new TH2F("hSSTry","hSSTry",70,-600,2200,130,0,5200);
227   TH2F *hSSTly = new TH2F("hSSTly","hSSTly",70,-600,2200,130,0,5200);
228   TH2F *hCSTy = new TH2F("hCSTy","hCSTy",70,-600,2200,130,0,5200);
229   TH2F *hSSTrpy = new TH2F("hSSTrpy","hSSTrpy",70,-600,2200,130,0,5200); 
230   TH2F *hSSTlpy = new TH2F("hSSTlpy","hSSTlpy",70,-600,2200,130,0,5200); 
231
232   // Chamber plane sticker targets
233   for (int iPoint=0; iPoint<9; iPoint++) {
234     sprintf(sPointName,"700%d",iPoint+1);
235     pointCPST = (AliSurveyPoint *)points->FindObject(sPointName);
236     if(!pointCPST) {
237       printf("Point %s is missing ...\n",sPointName);
238       break;       
239     }
240     hCPSTry->Fill(pointCPST->GetX(),pointCPST->GetZ(),pointCPST->GetY());
241   }
242   for (int iPoint=9; iPoint<18; iPoint++) {
243     sprintf(sPointName,"70%d",iPoint+1);
244     pointCPST = (AliSurveyPoint *)points->FindObject(sPointName);
245     if(!pointCPST) {
246       printf("Point %s is missing ...\n",sPointName);
247       break;       
248     }   
249     hCPSTly->Fill(pointCPST->GetX(),pointCPST->GetZ(),pointCPST->GetY());
250   }
251
252   // Chamber Side Targets
253   for (int iPoint=0; iPoint<25; iPoint++) {
254     if (iPoint+1<10) {   
255       sprintf(sPointName,"800%d",iPoint+1);
256     }
257     else {
258       sprintf(sPointName,"80%d",iPoint+1);
259     }
260     pointCST = (AliSurveyPoint *)points->FindObject(sPointName);
261     if(!pointCPST) {
262       printf("Point %s is missing ...\n",sPointName);
263       break;
264     }
265     hCSTy->Fill(pointCST->GetX(),pointCST->GetZ(),pointCST->GetY());
266   }
267
268   
269   // Graphs used for plane fitting
270   TGraph2DErrors **gSST5 = new TGraph2DErrors*[nSlats];
271   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
272     gSST5[iSlat] = new TGraph2DErrors();
273   }
274   
275   // Keep the id of slat sticker targets next to slat button targets
276   Int_t iSSBT[13][2] = {{0}};
277   
278   // Fill graph with sticker target positions
279   for (int iSlat=0; iSlat<nSlats; iSlat++){
280     for (int iPoint=0; iPoint<9; iPoint++) {
281       // Second sticker target is next to first button target
282       // Previous to last sticker target is next to second button target
283       iSSBT[iSlat][0] = 2;
284       iSSBT[iSlat][1] = 8;
285       if (iSlat+1<10) {
286         sprintf(sPointName,"50%d%d",iSlat+1,iPoint+1);
287       }
288       else {
289         sprintf(sPointName,"5%d%d",iSlat+1,iPoint+1);
290       }
291       pointSST = (AliSurveyPoint *)points->FindObject(sPointName);
292       if(!pointSST) {
293         printf("%s\n",sPointName);
294         cout << iSlat << " " << iPoint  << " " << pointSST << endl;
295         // Previous to last sticker target is next to second button target
296         iSSBT[iSlat][1] = iPoint+1-2;
297         break;
298       }
299      
300       gSST5[iSlat]->SetPoint(iPoint,-pointSST->GetX(),pointSST->GetZ(),pointSST->GetY());
301       gSST5[iSlat]->SetPointError(iPoint,pointSST->GetPrecisionX(),pointSST->GetPrecisionZ(),pointSST->GetPrecisionY());
302       
303       // Fill histograms of sticker targets. For monitoring purposes.
304       if((iSlat+1)%2==0){
305         hSSTly->Fill(pointSST->GetX(),pointSST->GetZ(),pointSST->GetY());
306       }
307       else {
308         hSSTry->Fill(pointSST->GetX(),pointSST->GetZ(),pointSST->GetY());
309       }
310     }
311   }
312
313   Float_t xMin = -2200.;
314   Float_t xMax =  2200.;
315   Float_t yMin = -5200.;
316   Float_t yMax =  5200.;
317   Float_t zMin = -200.;
318   Float_t zMax =  200.;
319
320   Double_t xMinSlat, xMaxSlat;
321   Double_t yMinSlat, yMaxSlat;  
322
323   // Slat plane function
324   char fsName[100] = "fSlat00"; 
325   TF2 **fSlat = new TF2*[nSlats];
326   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
327     sprintf(fsName,"fSlat%d",iSlat+1);
328     fSlat[iSlat] = new TF2(fsName,eqPlane,xMin,xMax,yMin,yMax,3);
329   }
330
331   // Xcenter functions
332   char fxcName[100] = "fXcnSlat00"; 
333   TF2 ***fXcSlat = new TF2**[nSlats];
334   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
335     fXcSlat[iSlat] = new TF2*[2];
336     sprintf(fxcName,"fXcnSlat%d",iSlat+1);
337     fXcSlat[iSlat][0] = new TF2(fxcName,xnCenter,xMin,xMax,yMin,yMax,7);
338     sprintf(fxcName,"fXcpSlat%d",iSlat+1);
339     fXcSlat[iSlat][1] = new TF2(fxcName,xpCenter,xMin,xMax,yMin,yMax,7);   
340   }
341
342   // Ycenter functions
343   char fycName[100] = "fYcnSlat00"; 
344   TF2 ***fYcSlat = new TF2**[nSlats];
345   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
346     fYcSlat[iSlat] = new TF2*[2];
347     sprintf(fycName,"fYcnSlat%d",iSlat+1);
348     fYcSlat[iSlat][0] = new TF2(fycName,ynCenter,yMin,yMax,yMin,yMax,8);
349     sprintf(fycName,"fYcpSlat%d",iSlat+1);
350     fYcSlat[iSlat][1] = new TF2(fycName,ypCenter,yMin,yMax,yMin,yMax,8);   
351   }
352
353   // Zcenter functions
354   char fzcName[100] = "fZcnSlat00"; 
355   TF2 ***fZcSlat = new TF2**[nSlats];
356   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
357     fZcSlat[iSlat] = new TF2*[2];
358     sprintf(fzcName,"fZcnSlat%d",iSlat+1);
359     fZcSlat[iSlat][0] = new TF2(fzcName,znCenter,zMin,zMax,zMin,zMax,8);
360     sprintf(fzcName,"fZcpSlat%d",iSlat+1);
361     fZcSlat[iSlat][1] = new TF2(fzcName,zpCenter,zMin,zMax,zMin,zMax,8);   
362   }
363
364   // Phi rotation using xglobal coords functions
365   char fphixName[100] = "fPhiXnnSlat00"; 
366   TF2 ****fPhiXSlat = new TF2***[nSlats];
367   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
368     fPhiXSlat[iSlat] = new TF2**[2];
369     for (Int_t iX =0; iX<2; iX++)
370       fPhiXSlat[iSlat][iX] = new TF2*[2];
371     sprintf(fphixName,"fPhiXnnSlat%d",iSlat+1);
372     fPhiXSlat[iSlat][0][0] = new TF2(fphixName,phixnn,xMin,xMax,xMin,xMax,7);
373     sprintf(fphixName,"fPhixnpSlat%d",iSlat+1);
374     fPhiXSlat[iSlat][0][1] = new TF2(fphixName,phixnp,xMin,xMax,xMin,xMax,7);   
375     sprintf(fphixName,"fPhiXpnSlat%d",iSlat+1);
376     fPhiXSlat[iSlat][1][0] = new TF2(fphixName,phixpn,xMin,xMax,xMin,xMax,7);
377     sprintf(fphixName,"fPhixppSlat%d",iSlat+1);
378     fPhiXSlat[iSlat][1][1] = new TF2(fphixName,phixpp,xMin,xMax,xMin,xMax,7);   
379   }
380
381   // Phi rotation using yglobal coords functions
382   char fphiyName[100] = "fPhiYnnSlat00"; 
383   TF2 ****fPhiYSlat = new TF2***[nSlats];
384   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
385     fPhiYSlat[iSlat] = new TF2**[2];
386     for (Int_t iY =0; iY<2; iY++)
387       fPhiYSlat[iSlat][iY] = new TF2*[2];
388     sprintf(fphiyName,"fPhiYnnSlat%d",iSlat+1);
389     fPhiYSlat[iSlat][0][0] = new TF2(fphiyName,phiynn,yMin,yMax,yMin,yMax,7);
390     sprintf(fphiyName,"fPhiYnpSlat%d",iSlat+1);
391     fPhiYSlat[iSlat][0][1] = new TF2(fphiyName,phiynp,yMin,yMax,yMin,yMax,7);   
392     sprintf(fphiyName,"fPhiYpnSlat%d",iSlat+1);
393     fPhiYSlat[iSlat][1][0] = new TF2(fphiyName,phiypn,yMin,yMax,yMin,yMax,7);
394     sprintf(fphiyName,"fPhiYppSlat%d",iSlat+1);
395     fPhiYSlat[iSlat][1][1] = new TF2(fphiyName,phiypp,yMin,yMax,yMin,yMax,7);   
396   }
397
398 //   Double_t *xce = new Double_t[nSlats]; 
399 //   Double_t *yce = new Double_t[nSlats]; 
400 //   Double_t *zce = new Double_t[nSlats]; 
401   Double_t *psi = new Double_t[nSlats]; 
402   Double_t *tht = new Double_t[nSlats]; 
403 //   Double_t *phi = new Double_t[nSlats]; 
404
405   Double_t **lCenSlat = new Double_t*[nSlats];
406   Double_t **lRotSlat = new Double_t*[nSlats];
407   Double_t **lDiffCenSlat0 = new Double_t*[nSlats];
408   Double_t **lDiffRotSlat0 = new Double_t*[nSlats];
409   Double_t **lDiffThCenSlat0 = new Double_t*[nSlats];
410   Double_t **lDiffThRotSlat0 = new Double_t*[nSlats];
411   Double_t **lDeltaDiffCenSlat0 = new Double_t*[nSlats];
412   Double_t **lDeltaDiffRotSlat0 = new Double_t*[nSlats];
413
414   for (int iSlat=0; iSlat<nSlats; iSlat++) {
415     lCenSlat[iSlat] = new Double_t[3];
416     lRotSlat[iSlat] = new Double_t[3];
417
418     lDiffCenSlat0[iSlat] = new Double_t[3];
419     lDiffRotSlat0[iSlat] = new Double_t[3];
420     lDiffThCenSlat0[iSlat] = new Double_t[3];
421     lDiffThRotSlat0[iSlat] = new Double_t[3];
422     lDeltaDiffCenSlat0[iSlat] = new Double_t[3];
423     lDeltaDiffRotSlat0[iSlat] = new Double_t[3];
424   }
425
426
427   TGeoTranslation transSlat[nSlats];
428   TGeoRotation rotSlat[nSlats];
429   TGeoCombiTrans trfSlat[nSlats];
430
431   TGeoTranslation dtransSlat[nSlats];
432   TGeoRotation drotSlat[nSlats];
433   TGeoCombiTrans dtrfSlat[nSlats];
434
435   TGeoTranslation transTemp;
436   TGeoRotation rotTemp;
437   TGeoCombiTrans trfTemp;
438
439
440   Double_t lCenTemp[3];
441   Double_t lRotTemp[3];
442
443   Double_t lDiffTemp[9];
444   Double_t lDiffMin[9];
445
446   double tempDiff = 0.;
447   double tempDiff1 = 0.;
448   double tempDiff2 = 0.;
449
450   AliSurveyPoint **pointSSBT = new AliSurveyPoint*[2];
451
452   //
453   // Get Slat transformation. 
454   // Psi and Theta are obtained by fitting a plane to the sticker targets.
455   // Then Xc, Yc, Zc and Phi are obtained by solving the equations to the ref.
456   // syst. transformation of the button targets
457   //
458   for (int iSlat=0; iSlat<nSlats; iSlat++) {
459     sprintf(fsName,"fSlat%d",iSlat+1);
460     cout << "Fitting Slat" << iSlat+1 << " ..." << endl;
461
462     // Fit a plane to the sticker targets
463     gSST5[iSlat]->Fit(fsName,"","same");
464
465     psi[iSlat] = TMath::ATan(fSlat[iSlat]->GetParameter(1));
466     tht[iSlat] = TMath::ATan(fSlat[iSlat]->GetParameter(0));
467     if (iSlat==5)
468       psi[iSlat] += TMath::Pi(); // Rotated slat
469
470     lRotSlat[iSlat][0] = psi[iSlat];
471     lRotSlat[iSlat][1] = tht[iSlat];
472
473     for(Int_t iS=0; iS<2; iS++){
474       fXcSlat[iSlat][iS]->SetParameters(lSBTLoc6[iSlat][0][0],lSBTLoc6[iSlat][0][1],lSBTLoc6[iSlat][0][2],lSBTLoc6[iSlat][1][0],lSBTLoc6[iSlat][1][1],lSBTLoc6[iSlat][1][2],tht[iSlat]);
475       fYcSlat[iSlat][iS]->SetParameters(lSBTLoc6[iSlat][0][0],lSBTLoc6[iSlat][0][1],lSBTLoc6[iSlat][0][2],lSBTLoc6[iSlat][1][0],lSBTLoc6[iSlat][1][1],lSBTLoc6[iSlat][1][2],psi[iSlat],tht[iSlat]);
476       fZcSlat[iSlat][iS]->SetParameters(lSBTLoc6[iSlat][0][0],lSBTLoc6[iSlat][0][1],lSBTLoc6[iSlat][0][2],lSBTLoc6[iSlat][1][0],lSBTLoc6[iSlat][1][1],lSBTLoc6[iSlat][1][2],psi[iSlat],tht[iSlat]);
477       for(Int_t jS=0; jS<2; jS++){
478         fPhiXSlat[iSlat][iS][jS]->SetParameters(lSBTLoc6[iSlat][0][0],lSBTLoc6[iSlat][0][1],lSBTLoc6[iSlat][0][2],lSBTLoc6[iSlat][1][0],lSBTLoc6[iSlat][1][1],lSBTLoc6[iSlat][1][2],tht[iSlat]);
479         fPhiYSlat[iSlat][iS][jS]->SetParameters(lSBTLoc6[iSlat][0][0],lSBTLoc6[iSlat][0][1],lSBTLoc6[iSlat][0][2],lSBTLoc6[iSlat][1][0],lSBTLoc6[iSlat][1][1],lSBTLoc6[iSlat][1][2],psi[iSlat],tht[iSlat]);
480       }
481     }
482
483     //
484     // Calculate Slat Center from button targets
485     //
486
487     // Get button targets survey points
488     for (Int_t iPoint=0; iPoint<2; iPoint++) {
489       if (iSlat+1<10) {
490         sprintf(sPointName,"60%d%d",iSlat+1,iPoint+1);
491         pointSBT[iPoint] = (AliSurveyPoint *)points->FindObject(sPointName);
492         if(!pointSBT[iPoint]) {
493           cout << "Error! No button targets ... " << endl;
494           break;       
495         }
496         sprintf(sPointName,"50%d%d",iSlat+1,iSSBT[iSlat][iPoint]);
497         pointSSBT[iPoint] = (AliSurveyPoint *)points->FindObject(sPointName);
498         if(!pointSSBT[iPoint]) {
499           cout << "Error! No sticker target ... " << sPointName << endl;
500           break;       
501         }
502       }
503       else {
504         sprintf(sPointName,"6%d%d",iSlat+1,iPoint+1);
505         pointSBT[iPoint] = (AliSurveyPoint *)points->FindObject(sPointName);
506         if(!pointSBT[iPoint]) {
507           cout << "Error! No button targets ... " << endl;
508           break;
509         }
510         sprintf(sPointName,"5%d%d",iSlat+1,iSSBT[iSlat][iPoint]);
511         pointSSBT[iPoint] = (AliSurveyPoint *)points->FindObject(sPointName);
512         if(!pointSSBT[iPoint]) {
513           cout << "Error! No sticker targets ... " << sPointName << endl;
514           break;
515         }
516       }
517     }
518
519     tempDiff += TMath::Power(-1,iSlat)*((pointSBT[1]->GetY() - pointSSBT[1]->GetY())-(pointSBT[0]->GetY() - pointSSBT[0]->GetY()));
520     tempDiff1 += TMath::Abs(pointSBT[0]->GetY() - pointSSBT[0]->GetY())-20;
521     tempDiff2 += TMath::Abs(pointSBT[1]->GetY() - pointSSBT[1]->GetY())-20;
522     cout << "BSdiff: " << TMath::Abs(pointSBT[0]->GetY() - pointSSBT[0]->GetY()) << " " << TMath::Abs(pointSBT[1]->GetY() - pointSSBT[1]->GetY()) << " " << tempDiff1/(iSlat+1) << " " << tempDiff2/(iSlat+1) << " " << tempDiff/(iSlat+1) << endl;
523
524
525     //    Double_t p0l[3] = {0};
526     Double_t p1l[3] = {0};
527     Double_t p2l[3] = {0};
528     //    Double_t p0g[3] = {0};
529     Double_t p1g[3] = {0};
530     Double_t p2g[3] = {0};
531
532     //    p0l[2] = lSBTLoc6[iSlat][0][2];
533     // Button targets local coordinates
534     for(Int_t iCor=0; iCor<3; iCor++){
535       p1l[iCor]= lSBTLoc6[iSlat][0][iCor];
536       p2l[iCor]= lSBTLoc6[iSlat][1][iCor];
537     }
538     for(Int_t i=0; i<9; i++){
539       lDiffMin[i]=1000000.;
540     }
541     // Trying 2x*2y*2z*2phi possibilities
542     for(Int_t iX=0; iX<2; iX++){
543       for(Int_t iY=0; iY<2; iY++){
544         for(Int_t iZ=0; iZ<2; iZ++){
545           lCenTemp[0] = fXcSlat[iSlat][iX]->Eval(-pointSBT[0]->GetX(),-pointSBT[1]->GetX());
546           lCenTemp[1] = fYcSlat[iSlat][iY]->Eval(pointSBT[0]->GetZ(),pointSBT[1]->GetZ());
547           lCenTemp[2] = fZcSlat[iSlat][iZ]->Eval(pointSBT[0]->GetY(),pointSBT[1]->GetY());
548           //     lCenTemp[2] = fZcSlat[iSlat][iZ]->Eval(pointSSBT[0]->GetY(),pointSSBT[1]->GetY());
549           lRotTemp[0] = psi[iSlat];
550           lRotTemp[1] = tht[iSlat];
551           for(Int_t iP=0; iP<2; iP++){
552             lRotTemp[2] = fPhiXSlat[iSlat][iX][iP]->Eval(-pointSBT[0]->GetX(),-pointSBT[1]->GetX());
553
554             trfTemp.SetTranslation(transTemp);
555             trfTemp.SetRotation(rotTemp);
556             trfTemp.Clear();
557             trfTemp.RotateZ(TMath::RadToDeg()*lRotTemp[2]);
558             trfTemp.RotateY(TMath::RadToDeg()*lRotTemp[1]);
559             trfTemp.RotateX(TMath::RadToDeg()*lRotTemp[0]);
560             trfTemp.SetTranslation(lCenTemp);
561            
562             //     trfTemp.LocalToMaster(p0l, p0g);
563             //     lCenTemp[2]= fSlat[iSlat]->Eval(p0g[0],p0g[1]) - trfTemp.GetRotationMatrix()[8]*p0l[2];
564             //     trfTemp.SetTranslation(lCenTemp);
565
566             trfTemp.LocalToMaster(p1l, p1g);
567             trfTemp.LocalToMaster(p2l, p2g);
568            
569             lDiffTemp[0] = (-pointSBT[0]->GetX()-p1g[0]);
570             lDiffTemp[1] = (pointSBT[0]->GetZ()-p1g[1]);
571             lDiffTemp[2] = (pointSBT[0]->GetY()-p1g[2]);
572             //     lDiffTemp[2] = (pointSSBT[0]->GetY()-p1g[2]);
573             lDiffTemp[3] = TMath::Sqrt(lDiffTemp[0]*lDiffTemp[0]+lDiffTemp[1]*lDiffTemp[1]+lDiffTemp[2]*lDiffTemp[2]);     
574             lDiffTemp[4] = (-pointSBT[1]->GetX()-p2g[0]);
575             lDiffTemp[5] = (pointSBT[1]->GetZ()-p2g[1]);
576             lDiffTemp[6] = (pointSBT[1]->GetY()-p2g[2]);
577             //     lDiffTemp[6] = (pointSSBT[1]->GetY()-p2g[2]);
578             lDiffTemp[7] = TMath::Sqrt(lDiffTemp[4]*lDiffTemp[4]+lDiffTemp[5]*lDiffTemp[5]+lDiffTemp[6]*lDiffTemp[6]);
579             lDiffTemp[8] = TMath::Sqrt(lDiffTemp[3]*lDiffTemp[3]+lDiffTemp[7]*lDiffTemp[7]);
580            
581             if(lDiffTemp[8]<lDiffMin[8]){
582               cout << "Diffs" ;
583               for(Int_t i=0; i<9; i++){
584                 cout << " " << lDiffTemp[i]; 
585               }
586               cout << endl;
587               cout << "Slat" << iSlat+1 << " : mycenX" << iX << "Y" << iY << "Z" << iZ << "(" << lCenTemp[0] << "," << lCenTemp[1] << "," << lCenTemp[2] << "); rotx" << iP << "(" << lRotTemp[0] << "," << lRotTemp[1] << "," << lRotTemp[2] << ")"  << endl;             
588               cout << p1g[0] << " " << p1g[1] << " " << p1g[2] << " " << p2g[0] << " " << p2g[1] << " " << p2g[2] << endl;
589               cout << "Transformation improved ..." << endl;
590               lCenSlat[iSlat][0] = lCenTemp[0]; lCenSlat[iSlat][1] = lCenTemp[1]; lCenSlat[iSlat][2] = lCenTemp[2]; 
591               lRotSlat[iSlat][2] = lRotTemp[2];
592               for(Int_t i=0; i<9; i++){
593                 lDiffMin[i]=lDiffTemp[i];
594               }
595               if((lDiffMin[3]*lDiffMin[3]<0.1*0.1+0.1*0.1+0.1*0.1)&&
596                  (lDiffMin[7]*lDiffMin[7]<0.1*0.1+0.1*0.1+0.1*0.1)){
597                 cout << "Correct Transformation found X " << iX  << " Y " << iY << " Z " << iZ << " xP " << iP << endl;
598                 lCenSlat[iSlat][0] = lCenTemp[0]; lCenSlat[iSlat][1] = lCenTemp[1]; lCenSlat[iSlat][2] = lCenTemp[2]; 
599                 lRotSlat[iSlat][2] = lRotTemp[2];
600               }
601             }
602           }
603           for(Int_t iP=0; iP<2; iP++){
604             lRotTemp[2] = fPhiYSlat[iSlat][iY][iP]->Eval(pointSBT[0]->GetZ(),pointSBT[1]->GetZ());
605             Double_t lPhi = TMath::ATan2((pointSBT[1]->GetZ()-pointSBT[0]->GetZ()),-(pointSBT[1]->GetX()-pointSBT[0]->GetX()));
606            
607             trfTemp.Clear();
608             trfTemp.RotateZ(TMath::RadToDeg()*lRotTemp[2]);
609             trfTemp.RotateY(TMath::RadToDeg()*lRotTemp[1]);
610             trfTemp.RotateX(TMath::RadToDeg()*lRotTemp[0]);
611             trfTemp.SetTranslation(lCenTemp);
612
613             //     trfTemp.LocalToMaster(p0l, p0g);
614             //     lCenTemp[2]= fSlat[iSlat]->Eval(p0g[0],p0g[1]) - trfTemp.GetRotationMatrix()[8]*p0l[2];
615             //     trfTemp.SetTranslation(lCenTemp);
616            
617             trfTemp.LocalToMaster(p1l, p1g);
618             trfTemp.LocalToMaster(p2l, p2g);
619
620             lDiffTemp[0] = (-pointSBT[0]->GetX()-p1g[0]);
621             lDiffTemp[1] = (pointSBT[0]->GetZ()-p1g[1]);        
622             lDiffTemp[2] = (pointSBT[0]->GetY()-p1g[2]);
623             //     lDiffTemp[2] = (pointSSBT[0]->GetY()-p1g[2]);
624             lDiffTemp[3] = TMath::Sqrt(lDiffTemp[0]*lDiffTemp[0]+lDiffTemp[1]*lDiffTemp[1]+lDiffTemp[2]*lDiffTemp[2]);     
625             lDiffTemp[4] = (-pointSBT[1]->GetX()-p2g[0]);
626             lDiffTemp[5] = (pointSBT[1]->GetZ()-p2g[1]);        
627             lDiffTemp[6] = (pointSBT[1]->GetY()-p2g[2]);
628             //     lDiffTemp[6] = (pointSSBT[1]->GetY()-p2g[2]);
629             lDiffTemp[7] = TMath::Sqrt(lDiffTemp[4]*lDiffTemp[4]+lDiffTemp[5]*lDiffTemp[5]+lDiffTemp[6]*lDiffTemp[6]);
630             lDiffTemp[8] = TMath::Sqrt(lDiffTemp[3]*lDiffTemp[3]+lDiffTemp[7]*lDiffTemp[7]);
631            
632             if(lDiffTemp[8]<lDiffMin[8]){
633               cout << "Diffs" ;
634               for(Int_t i=0; i<9; i++){
635                 cout << " " << lDiffTemp[i]; 
636               }
637               cout << endl;          
638               cout << "Slat" << iSlat+1 << " : mycenX" << iX << "Y" << iY << "Z" << iZ << "(" << lCenTemp[0] << "," << lCenTemp[1] << "," << lCenTemp[2] << "); roty" << iP << "(" << lRotTemp[0] << "," << lRotTemp[1] << "," << lRotTemp[2] << "(" << lPhi << "))"  << endl;
639               cout << p1g[0] << " " << p1g[1] << " " << p1g[2] << " " << p2g[0] << " " << p2g[1] << " " << p2g[2] << endl;           
640               cout << "Transformation improved ..." << endl;
641               lCenSlat[iSlat][0] = lCenTemp[0]; lCenSlat[iSlat][1] = lCenTemp[1]; lCenSlat[iSlat][2] = lCenTemp[2]; 
642               lRotSlat[iSlat][2] = lRotTemp[2];
643
644               for(Int_t i=0; i<9; i++){
645                 lDiffMin[i]=lDiffTemp[i];
646               }
647               if((lDiffMin[3]*lDiffMin[3]<0.1*0.1+0.1*0.1+0.1*0.1)&&
648                  (lDiffMin[7]*lDiffMin[7]<0.1*0.1+0.1*0.1+0.1*0.1)){
649                 cout << "Correct Transformation found X " << iX  << " Y " << iY << " Z " << iZ << " yP " << iP << endl;
650                 lCenSlat[iSlat][0] = lCenTemp[0]; lCenSlat[iSlat][1] = lCenTemp[1]; lCenSlat[iSlat][2] = lCenTemp[2]; 
651                 lRotSlat[iSlat][2] = lRotTemp[2];
652               }
653             }
654           } 
655         } 
656       } 
657     }
658
659     // Fill slat plane for fit monitor.
660     xMinSlat = hSSTrpy->GetXaxis()->FindBin(pointSBT[0]->GetX()); 
661     xMaxSlat = hSSTrpy->GetXaxis()->FindBin(pointSBT[1]->GetX()); 
662     yMinSlat = hSSTrpy->GetYaxis()->FindBin(pointSBT[0]->GetZ()-200.); 
663     yMaxSlat = hSSTrpy->GetYaxis()->FindBin(pointSBT[0]->GetZ()+200.); 
664    
665     for (int i=(int)xMinSlat; i<=(int)xMaxSlat; i++) {
666       for (int j=(int)yMinSlat; j<=(int)yMaxSlat; j++) {
667         Double_t zSlat = fSlat[iSlat]->Eval(-hSSTrpy->GetXaxis()->GetBinCenter(i),
668                                             hSSTrpy->GetYaxis()->GetBinCenter(j));
669         if((iSlat+1)%2==0){
670           hSSTlpy->SetBinContent(i,j,zSlat);
671         }
672         else {
673           hSSTrpy->SetBinContent(i,j,zSlat);
674         }
675       }
676     }   
677   }
678
679   //
680   // Compare transformations to expected ones 
681   //
682   Int_t iSlatToPos[13] = {0, 11, 2, 9, 4, 7, 6, 5, 8, 3, 10, 1, 12};
683  
684   // Theoretical differences with respect to Slat 513 which is here Slat07
685   lDiffThCenSlat0[0][0]  = (2-5)*400./2. +12.5 -375. -25.;
686   lDiffThCenSlat0[1][0]  = (3-5)*400./2. +12.5 -375. -25.;
687   lDiffThCenSlat0[2][0]  = (4-5)*400./2. +12.5 -375. -25.;
688   lDiffThCenSlat0[3][0]  = (5-5)*400./2. +12.5 -375. -25.;
689   lDiffThCenSlat0[4][0]  = (5-5)*400./2. +12.5 -375. -25.;
690   lDiffThCenSlat0[5][0]  = (6-5)*400./2. +12.5 -375. -25.;
691   lDiffThCenSlat0[6][0]  = (5-5)*400./2. +12.5 -375. -25. +375. +25. -12.5;
692   lDiffThCenSlat0[7][0]  = (6-5)*400./2. +12.5 -375. -25.;
693   lDiffThCenSlat0[8][0]  = (5-5)*400./2. +12.5 -375. -25.;
694   lDiffThCenSlat0[9][0]  = (5-5)*400./2. +12.5 -375. -25.;
695   lDiffThCenSlat0[10][0] = (4-5)*400./2. +12.5 -375. -25.;
696   lDiffThCenSlat0[11][0] = (3-5)*400./2. +12.5 -375. -25.;
697   lDiffThCenSlat0[12][0] = (2-5)*400./2. +12.5 -375. -25.;
698
699   lDiffThCenSlat0[12][1] = 382.0 +378.5 +375.5 +294.0 +370.0 +286.0; 
700   lDiffThCenSlat0[1][1]  = 382.0 +378.5 +375.5 +294.0 +370.0; 
701   lDiffThCenSlat0[10][1] = 382.0 +378.5 +375.5 +294.0; 
702   lDiffThCenSlat0[3][1]  = 382.0 +378.5 +375.5; 
703   lDiffThCenSlat0[8][1]  = 382.0 +378.5; 
704   lDiffThCenSlat0[5][1]  = 382.0; 
705   lDiffThCenSlat0[6][1]  = 0.0; 
706   lDiffThCenSlat0[7][1]  = -382.0; 
707   lDiffThCenSlat0[4][1]  = -382.0 -378.5; 
708   lDiffThCenSlat0[9][1]  = -382.0 -378.5 -375.5; 
709   lDiffThCenSlat0[2][1]  = -382.0 -378.5 -375.5 -294.0; 
710   lDiffThCenSlat0[11][1] = -382.0 -378.5 -375.5 -294.0 -370.0; 
711   lDiffThCenSlat0[0][1]  = -382.0 -378.5 -375.5 -294.0 -370.0 -286.0; 
712
713   lDiffThCenSlat0[12][2] =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
714   lDiffThCenSlat0[1][2]  = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
715   lDiffThCenSlat0[10][2] =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
716   lDiffThCenSlat0[3][2]  = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
717   lDiffThCenSlat0[8][2]  =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
718   lDiffThCenSlat0[5][2]  = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
719   lDiffThCenSlat0[6][2]  =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
720   lDiffThCenSlat0[7][2]  = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
721   lDiffThCenSlat0[4][2]  =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
722   lDiffThCenSlat0[9][2]  = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
723   lDiffThCenSlat0[2][2]  =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
724   lDiffThCenSlat0[11][2] = -(42.5)-(42.5); //-(42.5-1.175)-(42.5-1.175);
725   lDiffThCenSlat0[0][2]  =  (42.5)-(42.5); // (42.5-1.175)-(42.5-1.175);
726
727   TGraph *gDeltaDiffCenXSlat0 = new TGraph(nSlats);
728   TGraph *gDeltaDiffCenYSlat0 = new TGraph(nSlats);
729   TGraph *gDeltaDiffCenZSlat0 = new TGraph(nSlats);
730   TGraph *gDeltaDiffPsiSlat0 = new TGraph(nSlats);
731   TGraph *gDeltaDiffThtSlat0 = new TGraph(nSlats);
732   TGraph *gDeltaDiffPhiSlat0 = new TGraph(nSlats);
733
734   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
735     trfSlat[iSlat].SetTranslation(transSlat[iSlat]);
736     trfSlat[iSlat].SetRotation(rotSlat[iSlat]);
737     trfSlat[iSlat].Clear();
738     trfSlat[iSlat].RotateZ(TMath::RadToDeg()*lRotSlat[iSlat][2]);
739     trfSlat[iSlat].RotateY(TMath::RadToDeg()*lRotSlat[iSlat][1]);
740     trfSlat[iSlat].RotateX(TMath::RadToDeg()*lRotSlat[iSlat][0]);
741     trfSlat[iSlat].SetTranslation(lCenSlat[iSlat]); 
742   }
743
744   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
745     dtrfSlat[iSlat].SetTranslation(dtransSlat[iSlat]);
746     dtrfSlat[iSlat].SetRotation(drotSlat[iSlat]);
747     dtrfSlat[iSlat].Clear();
748     dtrfSlat[iSlat] = trfSlat[6].Inverse()*trfSlat[iSlat]; 
749     dtrfSlat[iSlat].Print();
750     lDiffCenSlat0[iSlat] = (Double_t*)dtrfSlat[iSlat].GetTranslation();
751     MatrixToAngles(dtrfSlat[iSlat].GetRotationMatrix(),lDiffRotSlat0[iSlat]);
752
753     lDeltaDiffCenSlat0[iSlat][0] = lDiffCenSlat0[iSlat][0]-lDiffThCenSlat0[iSlat][0];
754     lDeltaDiffCenSlat0[iSlat][1] =  -lDiffCenSlat0[iSlat][1]-lDiffThCenSlat0[iSlat][1];
755     lDeltaDiffCenSlat0[iSlat][2] =  -lDiffCenSlat0[iSlat][2]-lDiffThCenSlat0[iSlat][2];
756     lDeltaDiffRotSlat0[iSlat][0] = lDiffRotSlat0[iSlat][0];
757     lDeltaDiffRotSlat0[iSlat][1] = lDiffRotSlat0[iSlat][1];
758     lDeltaDiffRotSlat0[iSlat][2] = lDiffRotSlat0[iSlat][2];
759     gDeltaDiffCenXSlat0->SetPoint(iSlat,lDeltaDiffCenSlat0[iSlat][0],iSlatToPos[iSlat]+1);
760     gDeltaDiffCenYSlat0->SetPoint(iSlat,lDeltaDiffCenSlat0[iSlat][1],iSlatToPos[iSlat]+1);
761     gDeltaDiffCenZSlat0->SetPoint(iSlat,lDeltaDiffCenSlat0[iSlat][2],iSlatToPos[iSlat]+1);
762     gDeltaDiffPsiSlat0->SetPoint(iSlat,1e3*lDeltaDiffRotSlat0[iSlat][0],iSlatToPos[iSlat]+1);
763     gDeltaDiffThtSlat0->SetPoint(iSlat,1e3*lDeltaDiffRotSlat0[iSlat][1],iSlatToPos[iSlat]+1);
764     gDeltaDiffPhiSlat0->SetPoint(iSlat,1e3*lDeltaDiffRotSlat0[iSlat][2],iSlatToPos[iSlat]+1);
765   }
766
767   // Import TGeo geometry 
768   char* geoFilename = "geometry.root";
769   cout << "geometry imported" << endl;
770   if ( ! AliGeomManager::GetGeometry() ) {
771     AliGeomManager::LoadGeometry(geoFilename);
772     if (! AliGeomManager::GetGeometry() ) {
773       printf("MUONSurveyCh8L: getting geometry from file %s failed\n", geoFilename);
774       return;
775     }
776   }
777
778   AliMUONGeometryTransformer *transform = new AliMUONGeometryTransformer();
779 //   transform->ReadGeometryData("volpath.dat", gGeoManager);
780   transform->LoadGeometryData();
781   cout << "geometry data read" << endl;
782   AliMUONGeometryTransformer *newTransform = ReAlign(transform,11,dtrfSlat,true); 
783   newTransform->WriteTransformations("transform2ReAlign.dat");
784   cout << "newtransform read" << endl;
785   // Generate realigned data in local cdb
786   const TClonesArray* array = newTransform->GetMisAlignmentData();
787    
788   // CDB manager
789   AliCDBManager* cdbManager = AliCDBManager::Instance();
790   cdbManager->SetDefaultStorage("local://ReAlignCDB");
791   
792   AliCDBMetaData* cdbData = new AliCDBMetaData();
793   cdbData->SetResponsible("Dimuon Offline project");
794   cdbData->SetComment("MUON alignment objects with residual misalignment");
795   AliCDBId id("MUON/Align/Data", 0, 0); 
796   cdbManager->Put(const_cast<TClonesArray*>(array), id, cdbData);
797
798   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
799     cout << lDeltaDiffCenSlat0[iSlat][0] << " " << lDiffCenSlat0[iSlat][0] << " " << lDiffThCenSlat0[iSlat][0] << endl;
800   }
801   cout << endl;
802   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
803     //    lDiffCenSlat0[iSlat][1] = lCenSlat[iSlat][1]-lCenSlat[6][1];
804     //    lDeltaDiffCenSlat0[iSlat][1] = lDiffCenSlat0[iSlat][1]-lDiffThCenSlat0[iSlat][1];
805     cout << lDeltaDiffCenSlat0[iSlat][1] << " " << lDiffCenSlat0[iSlat][1] << " " << lDiffThCenSlat0[iSlat][1] << endl;
806   }
807   cout << endl;
808   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
809     cout << lDeltaDiffCenSlat0[iSlat][2] << " " << lDiffCenSlat0[iSlat][2] << " " << lDiffThCenSlat0[iSlat][2] << endl;
810   }
811   cout << endl;
812   cout << endl;
813   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
814     cout << lDeltaDiffRotSlat0[iSlat][0] << " " << lDiffRotSlat0[iSlat][0] << " " << lDiffThRotSlat0[iSlat][0] << endl;
815   }
816   cout << endl;
817   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
818     cout << lDeltaDiffRotSlat0[iSlat][1] << " " << lDiffRotSlat0[iSlat][1] << " " << lDiffThRotSlat0[iSlat][1] << endl;
819   }
820   cout << endl;
821   for(Int_t iSlat=0; iSlat<nSlats; iSlat++){
822     cout << lDeltaDiffRotSlat0[iSlat][2] << " " << lDiffRotSlat0[iSlat][2] << " " << lDiffThRotSlat0[iSlat][2] << endl;
823   }
824
825   TH1F *mySlatDeltaDiffCenX = new TH1F("mySlatDeltaDiffCenX","mySlatDeltaDiffCenX",100,-10,10);
826   mySlatDeltaDiffCenX->SetMaximum(15);
827   mySlatDeltaDiffCenX->SetMinimum(0);
828   TH1F *mySlatDeltaDiffCenY = new TH1F("mySlatDeltaDiffCenY","mySlatDeltaDiffCenY",100,-10,10);
829   mySlatDeltaDiffCenY->SetMaximum(15);
830   mySlatDeltaDiffCenY->SetMinimum(0);
831   TH1F *mySlatDeltaDiffCenZ = new TH1F("mySlatDeltaDiffCenZ","mySlatDeltaDiffCenZ",100,-20,20);
832   mySlatDeltaDiffCenZ->SetMaximum(15);
833   mySlatDeltaDiffCenZ->SetMinimum(0);
834
835   TH1F *mySlatDeltaDiffRotX = new TH1F("mySlatDeltaDiffRotX","mySlatDeltaDiffRotX",100,-10,10);
836   mySlatDeltaDiffRotX->SetMaximum(15);
837   mySlatDeltaDiffRotX->SetMinimum(0);
838   TH1F *mySlatDeltaDiffRotY = new TH1F("mySlatDeltaDiffRotY","mySlatDeltaDiffRotY",100,-10,10);
839   mySlatDeltaDiffRotY->SetMaximum(15);
840   mySlatDeltaDiffRotY->SetMinimum(0);
841   TH1F *mySlatDeltaDiffRotZ = new TH1F("mySlatDeltaDiffRotZ","mySlatDeltaDiffRotZ",100,-5,5);
842   mySlatDeltaDiffRotZ->SetMaximum(15);
843   mySlatDeltaDiffRotZ->SetMinimum(0);
844   //
845   // ******** Starting plots 
846   //
847   TCanvas *canvas;
848   TPad *pad;
849   TPaveLabel *theTitle;
850
851   TPostScript *ps = 0;
852
853   if( saveps ){
854     ps = new TPostScript(filename,filetype); 
855     ps->NewPage();
856   }
857  
858   // Inv Mass, Multiplicity
859   sprintf(str,"Chamber 8L");
860   TCanvas *cvn0 = new TCanvas("cvn0",str,cWidth,cHeight);
861   canvas = cvn0;
862   canvas->Range(0,0,21,29);
863  
864   TPaveLabel *theTitle0 = new TPaveLabel(3,27.0,18,28.5," Deformations of chamber 8L ","br");
865   theTitle = theTitle0;
866   theTitle->SetFillColor(18);
867   theTitle->SetTextFont(32);
868   theTitle->SetTextSize(0.4);
869   theTitle->SetTextColor(1);
870   theTitle->Draw();
871  
872   TPad *pad0 = new TPad("pad0","pad0",0.01,0.01,0.98,0.91,0);
873   pad = pad0;
874   pad->Draw();
875   pad->Divide(2,2);
876
877   pad->cd(1);
878   gStyle->SetPalette(1);
879   hCPSTry->SetMinimum(100);
880   hCPSTry->SetMaximum(120);
881   hCPSTry->Draw("lego2z");
882
883   pad->cd(2);
884   gStyle->SetPalette(1);
885   hSSTry->SetMinimum(60);
886   hSSTry->SetMaximum(80);
887   hSSTry->Draw("lego2z");
888
889   pad->cd(3);
890   gStyle->SetPalette(1);
891   hCSTy->SetMinimum(110);
892   hCSTy->SetMaximum(130);
893   hCSTy->Draw("lego2z");
894
895   pad->cd(4);
896   gStyle->SetPalette(1);
897   hSSTly->SetMinimum(165);
898   hSSTly->SetMaximum(185);
899   hSSTly->Draw("lego2z");
900
901   if(saveps){
902     ps->NewPage();
903   }
904
905   // Inv Mass, Multiplicity
906   sprintf(str,"Chamber 8L");
907   TCanvas *cvn1 = new TCanvas("cvn1",str,cWidth,cHeight);
908   canvas = cvn1;
909   canvas->Range(0,0,21,29);
910   
911   TPaveLabel *theTitle1 = new TPaveLabel(3,27.0,18,28.5," Deformations of chamber 8L ","br");
912   theTitle = theTitle1;
913   theTitle->SetFillColor(18);
914   theTitle->SetTextFont(32);
915   theTitle->SetTextSize(0.4);
916   theTitle->SetTextColor(1);
917   theTitle->Draw();
918  
919   TPad *pad1 = new TPad("pad1","pad1",0.01,0.01,0.98,0.91,0);
920   pad = pad1;
921   pad->Draw();
922   pad->Divide(2,2);
923
924   pad->cd(1);
925   gStyle->SetPalette(1);
926   hCPSTly->SetMinimum(120);
927   hCPSTly->SetMaximum(140);
928   hCPSTly->Draw("lego2z");
929
930   pad->cd(2);
931   gStyle->SetPalette(1);
932   hSSTrpy->SetMinimum(60);
933   hSSTrpy->SetMaximum(80);
934   hSSTrpy->Draw("surf2z");
935
936   pad->cd(3);
937   gStyle->SetPalette(1);
938   hCSTy->SetMinimum(110);
939   hCSTy->SetMaximum(130);
940   hCSTy->Draw("lego2z");
941
942   pad->cd(4);
943   gStyle->SetPalette(1);
944   hSSTlpy->SetMinimum(165);
945   hSSTlpy->SetMaximum(185);
946   hSSTlpy->Draw("surf2z");
947
948   // Inv Mass, Multiplicity
949   sprintf(str,"Chamber 8L");
950   TCanvas *cvn2 = new TCanvas("cvn2",str,cWidth,cHeight);
951   canvas = cvn2;
952   canvas->Range(0,0,21,29);
953   
954   TPaveLabel *theTitle2 = new TPaveLabel(3,27.0,18,28.5," Deformations of chamber 8L ","br");
955   theTitle = theTitle2;
956   theTitle->SetFillColor(18);
957   theTitle->SetTextFont(32);
958   theTitle->SetTextSize(0.4);
959   theTitle->SetTextColor(1);
960   theTitle->Draw();
961  
962   TPad *pad2 = new TPad("pad2","pad2",0.01,0.01,0.98,0.91,0);
963   pad = pad2;
964   pad->Draw();
965   pad->Divide(3,2);
966
967   pad->cd(1);
968   mySlatDeltaDiffCenX->Draw();
969   mySlatDeltaDiffCenX->SetXTitle("#Delta[(xc_{i}^{m}-xc_{0}^{m})-(xc_{i}^{th}-xc_{0}^{th})] (mm)");
970   mySlatDeltaDiffCenX->SetYTitle("Slat ascending vertical ordering");
971   gDeltaDiffCenXSlat0->SetMarkerStyle(20);
972   gDeltaDiffCenXSlat0->Draw("P");
973
974   pad->cd(2);
975   mySlatDeltaDiffCenY->Draw();
976   mySlatDeltaDiffCenY->SetXTitle("#Delta[(yc_{i}^{m}-yc_{0}^{m})-(yc_{i}^{th}-yc_{0}^{th})] (mm)");
977   mySlatDeltaDiffCenY->SetYTitle("Slat ascending vertical ordering");
978   gDeltaDiffCenYSlat0->SetMarkerStyle(20);
979   gDeltaDiffCenYSlat0->Draw("P");
980
981   pad->cd(3);
982   mySlatDeltaDiffCenZ->Draw();
983   mySlatDeltaDiffCenZ->SetXTitle("#Delta[(zc_{i}^{m}-zc_{0}^{m})-(zc_{i}^{th}-zc_{0}^{th})] (mm)");
984   mySlatDeltaDiffCenZ->SetYTitle("Slat ascending vertical ordering");
985   gDeltaDiffCenZSlat0->SetMarkerStyle(20);
986   gDeltaDiffCenZSlat0->Draw("P");
987
988   pad->cd(4);
989   mySlatDeltaDiffRotX->Draw();
990   mySlatDeltaDiffRotX->SetXTitle("#Delta[(#psi_{i}^{m}-#psi_{0}^{m})-(#psi_{i}^{th}-#psi_{0}^{th})] (mrad)");
991   mySlatDeltaDiffRotX->SetYTitle("Slat ascending vertical ordering");
992   gDeltaDiffPsiSlat0->SetMarkerStyle(20);
993   gDeltaDiffPsiSlat0->Draw("P");
994
995   pad->cd(5);
996   mySlatDeltaDiffRotY->Draw();
997   mySlatDeltaDiffRotY->SetXTitle("#Delta[(#theta_{i}^{m}-#theta_{0}^{m})-(#theta_{i}^{th}-#theta_{0}^{th})] (mrad)");
998   mySlatDeltaDiffRotY->SetYTitle("Slat ascending vertical ordering");
999   gDeltaDiffThtSlat0->SetMarkerStyle(20);
1000   gDeltaDiffThtSlat0->Draw("P");
1001
1002   pad->cd(6);
1003   mySlatDeltaDiffRotZ->Draw();
1004   mySlatDeltaDiffRotZ->SetXTitle("#Delta[(#phi_{i}^{m}-#phi_{0}^{m})-(#phi_{i}^{th}-#phi_{0}^{th})] (mrad)");
1005   mySlatDeltaDiffRotZ->SetYTitle("Slat ascending vertical ordering");
1006   gDeltaDiffPhiSlat0->SetMarkerStyle(20);
1007   gDeltaDiffPhiSlat0->Draw("P");
1008
1009   if( saveps ){
1010     ps->Close();
1011   }
1012 }