]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSSurveyToAlign.cxx
Warning removal (Marian)
[u/mrichter/AliRoot.git] / ITS / AliITSSurveyToAlign.cxx
1 /**************************************************************************
2  * Copyright(c) 2008-2010, 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 //////////////////////////////////////////////////////////////////////////
19 //   Class to convert survey tables in alignment objects
20 //   for SSD and SDD
21 //   origin: Marco Van Leeuwen (m.vanleeuwen1@uu.nl)
22 //           Panos.Christakoglou (Panos.Christakoglou@cern.ch)
23 //           Martin Poghosyan (Martin.Poghosyan@to.infn.it)
24 //////////////////////////////////////////////////////////////////////////
25
26 #include "Riostream.h"
27 #include "TFile.h"
28 #include "TSystem.h"
29 #include "TClonesArray.h"
30 #include "TGeoManager.h"
31 #include "TGeoMatrix.h"
32 #include "TGeoPhysicalNode.h"
33 #include "TMatrixD.h"
34 #include "TMath.h"
35
36 #include "AliITSSurveyToAlign.h"
37 #include "AliSurveyObj.h"
38 #include "AliSurveyPoint.h"
39 #include "AliAlignObjParams.h"
40 #include "AliGeomManager.h"
41
42 #include "AliLog.h"
43
44 #include "AliCDBManager.h"
45 #include "AliCDBEntry.h"
46 #include "AliCDBStorage.h"
47
48 #include "AliITSgeomTGeo.h"
49
50 ClassImp(AliITSSurveyToAlign)
51
52 const Double_t AliITSSurveyToAlign::fgkLocR[6][3]={{ 3.24,0.21905,-2.4},
53                                                    { 3.58,0.21905, 0. },
54                                                    { 3.24,0.21905,+2.4},
55                                                    {-3.24,0.21905,+2.4},
56                                                    {-3.58,0.21905, 0. },
57                                                    {-3.24,0.21905,-2.4}};
58
59 const Double_t AliITSSurveyToAlign::fgkLocL[6][3]={{-3.24,0.21905, 2.4},
60                                                    {-3.58,0.21905, 0. },
61                                                    {-3.24,0.21905,-2.4},
62                                                    { 3.24,0.21905,-2.4},
63                                                    { 3.58,0.21905, 0. },
64                                                    { 3.24,0.21905, 2.4}};
65
66
67 //________________________________________________________________________
68 AliITSSurveyToAlign::AliITSSurveyToAlign(Int_t run, Int_t repSDD, Int_t repVerSDD, Int_t repModSSD, Int_t repModVerSSD, Int_t repLaddSSD, Int_t repLaddVerSSD) :
69   //TObject(),
70   AliSurveyToAlignObjs(),
71   fRun(run),
72   fSDDrepNumber(repSDD),
73   fSDDrepVersion(repVerSDD),
74   fSSDModuleRepNumber(repModSSD),
75   fSSDModuleRepVersion(repModVerSSD),
76   fSSDLadderRepNumber(repLaddSSD),
77   fSSDLadderRepVersion(repLaddVerSSD)
78  {
79   //
80   //  default constructor
81   //  Arguments are report numbers for survey data. 
82   //  The defaults point to reports from detector construction
83   // 
84 }
85
86 //_________________________________________________________________________
87 AliITSSurveyToAlign::AliITSSurveyToAlign(const AliITSSurveyToAlign &align) :
88   //  TObject(),
89   AliSurveyToAlignObjs(align),
90   fRun(align.fRun),
91   fSDDrepNumber(align.fSDDrepNumber),
92   fSDDrepVersion(align.fSDDrepVersion),
93   fSSDModuleRepNumber(align.fSSDModuleRepNumber),
94   fSSDModuleRepVersion(align.fSSDModuleRepVersion),
95   fSSDLadderRepNumber(align.fSSDLadderRepNumber),
96   fSSDLadderRepVersion(align.fSSDLadderRepVersion)
97 {
98   //
99   //  copy constructor 
100   //
101 }
102
103 //__________________________________________________________________________
104 AliITSSurveyToAlign & AliITSSurveyToAlign::operator =(const AliITSSurveyToAlign & align) {
105   //
106   // assignment operator - dummy
107   //
108
109   return (*this);
110 }
111
112 //__________________________________________________________________________
113 AliITSSurveyToAlign::~AliITSSurveyToAlign() {
114   //
115   // destructor
116   //
117 }
118
119 //______________________________________________________________________
120 void AliITSSurveyToAlign::Run() { 
121   //
122   // Runs the full chain
123   // User should call StoreAlignObj() or StoreAlignObjCDB() afterwards to 
124   // store output
125   //
126
127   // Load ideal geometry from the OCDB
128   AliCDBManager *cdb = AliCDBManager::Instance();
129   cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
130   cdb->SetRun(fRun);
131   AliGeomManager::LoadGeometry();
132   
133   CreateAlignObjSPD();
134   LoadSurveyFromAlienFile("ITS", fSDDrepNumber, fSDDrepVersion);
135   CreateAlignObjSDD();
136   LoadSurveyFromAlienFile("ITS", fSSDModuleRepNumber, fSSDModuleRepVersion);
137   CreateAlignObjSSDModules();
138   LoadSurveyFromAlienFile("ITS", fSSDLadderRepNumber, fSSDLadderRepVersion);
139   CreateAlignObjSSDLadders();
140 }
141
142 void AliITSSurveyToAlign::CreateAlignObjSPD(){
143   // 
144   // Create alignObjs for SPD
145   //    For the moment, uses 0,0,0,0,0,0
146   //
147   for(Int_t imod = 0; imod < 240; imod++) {
148     Int_t ilayer = (imod < 80) ? AliGeomManager::kSPD1 : AliGeomManager::kSPD2;
149     Int_t imodule = (imod < 80) ? imod : imod - 80;
150
151     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
152     const Char_t *symname = AliGeomManager::SymName(uid);
153
154     new((*fAlignObjArray)[imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
155   }//module loop
156
157 }
158
159 void AliITSSurveyToAlign::CreateAlignObjSDD(){
160   //
161   // Create alignment objects for SDD
162   // Called by Run()
163   //
164   Int_t uid = 0;
165   const char* symname = 0;
166   AliSurveyPoint* pt = 0;
167  
168   Int_t iModuleIndex=240;
169   Int_t iModule0=0;
170   Int_t iLadder0=0;
171   Int_t iLayer0=3;
172   Int_t nModules=0;
173
174   if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
175     AliWarning("SDD survey data are not available, using zero values");
176     CreateAlignObjDummySDD();
177     return;
178   }
179
180   for(Int_t imod = 1; imod < fSurveyPoints->GetEntries(); imod++) {
181     pt = (AliSurveyPoint*) fSurveyPoints->At(imod);
182     if(!pt) continue;
183
184     Int_t iLayer, iLadder, iModule, iPoint;
185     ReadPointNameSDD(pt->GetName(),iLayer, iLadder, iModule, iPoint);
186
187     if(iModule==iModule0)
188     {
189       fSDDmeP[iPoint][0]=pt->GetX();
190       fSDDmeP[iPoint][1]=pt->GetY();
191       fSDDmeP[iPoint][2]=pt->GetZ();
192       fSDDmeP[iPoint][3]=pt->GetPrecisionX();
193       fSDDmeP[iPoint][4]=pt->GetPrecisionY();
194       fSDDmeP[iPoint][5]=pt->GetPrecisionZ();
195       fSDDisMe[iPoint]=kTRUE;
196
197       if(iLayer==3) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-240);
198       if(iLayer==4) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-324);
199       symname = AliGeomManager::SymName(uid);
200       GetIdPosSDD(uid,iLayer0, iModule0, iPoint);
201       nModules++;
202     }
203     //    cout << "Points red module " << imod << endl;
204     if((iModule!=iModule0)||(imod==(fSurveyPoints->GetEntries()-1)))
205     {
206       ConvertToRSofModulesAndRotSDD(iLayer0, iModule0);
207
208       Double_t tet = 0.;
209       Double_t psi =0.;
210       Double_t phi = 0.;
211       Double_t x0  = 0.;
212       Double_t y0  =0.;
213       Double_t z0  = 0.;
214
215       if(nModules==2) CalcShiftSDD(x0,y0,z0);
216       if(nModules>2)   CalcShiftRotSDD(tet, psi, phi, x0, y0, z0);
217       tet*=(180/TMath::Pi());
218       psi*=(180/TMath::Pi());
219       phi*=(180/TMath::Pi());
220
221 //    printf("%s  %d  %f  %f  %f  %f  %f  %f\n",symname, uid, x0/10., y0/10., z0/10., psi, tet, phi);
222 //      cout << "Allocate alignobjparams " << imod << endl;
223       new((*fAlignObjArray)[iModuleIndex]) AliAlignObjParams(symname, uid, x0/10., y0/10., z0/10., psi, tet, phi, kFALSE);
224
225       iModule0=iModule;
226       iLayer0=iLayer;
227       iLadder0=iLadder;
228       nModules=0;
229       iModuleIndex = AliITSgeomTGeo::GetModuleIndex(iLayer,iLadder+1,iModule+1);
230       for(Int_t i=0; i<6;i++) fSDDisMe[i]=kFALSE;
231       if(imod!=(fSurveyPoints->GetEntries()-1)) imod--;
232     }
233   }//module loop
234 }
235
236 void AliITSSurveyToAlign::CreateAlignObjDummySDD(){
237   // 
238   // Create empty alignment objects
239   // Used when fSurveySDD == 0
240   //
241   for(Int_t imod = 0; imod < 260; imod++) {
242
243     Int_t ilayer = (imod < 84) ? AliGeomManager::kSDD1 : AliGeomManager::kSDD2;
244     Int_t imodule = (imod < 84) ? imod : imod - 84;
245
246     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
247     const Char_t *symname = AliGeomManager::SymName(uid);
248
249     new((*fAlignObjArray)[imod+240]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
250   }//module loop
251 }
252
253 void AliITSSurveyToAlign::CreateAlignObjSSDModules(){
254   //
255   // Create alignment objects for SDD
256   // Called by Run()
257   // Two sets of objects are made: 1 per module + 1 per ladder
258   //
259   Double_t sx, sz;
260   const Float_t kMu2Cm = 1e-4;
261   const Float_t kSensLength = 7.464;
262   const Int_t kSSDMODULES = 1698;
263
264   if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
265     AliWarning("SSD module survey data not available; using dummy values");
266     CreateAlignObjDummySSDModules();
267     return;
268   }
269
270   // First do module-by-module
271
272   for(Int_t imod = 500; imod < kSSDMODULES + 500; imod++) {
273     Int_t iLayer, iLadder, iLaddMod;
274     AliITSgeomTGeo::GetModuleId(imod,iLayer,iLadder,iLaddMod);  // returns 1-based numbers
275  
276     TString pname="ITS/SSD";
277     pname += iLayer-1;
278     pname += "/Ladder";
279     pname += iLadder-1;
280     pname += "/Sensor";
281     pname += iLaddMod-1;
282     AliSurveyPoint *pt1 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point0");
283     AliSurveyPoint *pt2 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point1");
284     if(!pt1 || !pt2) {
285       AliWarning(Form("No Survey points for iladd %d imod %d",iLadder,iLaddMod));
286       continue;
287     }
288
289     sx = 0.5*(pt1->GetX() + pt2->GetX()) * kMu2Cm;
290     sz = 0.5*(pt1->GetZ() + pt2->GetZ()) * kMu2Cm;
291
292     // Minus sign to change local coordinate convention 
293     Float_t theta = -(pt2->GetZ() - pt1->GetZ())*kMu2Cm/kSensLength;
294
295     theta *= 180./TMath::Pi();
296     Int_t iLayMod = imod - 500;
297     if (iLayer == 6)
298       iLayMod -= 748;
299     Int_t uid = AliGeomManager::LayerToVolUID(iLayer,iLayMod);
300
301     const Char_t *symname = AliGeomManager::SymName(uid);
302     if (pname.CompareTo(symname) != 0)
303       AliWarning(Form("Mapping mismatch survey point %s volume name %s",pname.Data(),symname));
304     /*
305     if (imod >= 676 && imod <= 697) {
306       cout << "ilayer " << iLayer << " imod " << imod 
307            << " uid " << uid << " name " << symname 
308            << " survey shift " << sx << " " << 0 << " " << sz << endl
309            << " theta " << theta << endl;
310     }
311     */
312     new((*fAlignObjArray)[imod]) AliAlignObjParams(symname, uid, sx, 0, sz, 0., theta, 0., kFALSE);
313   } //module loop
314 }
315
316 void AliITSSurveyToAlign::CreateAlignObjSSDLadders(){
317   //
318   //   Now do ladders  (Torino data)
319   // 
320   const Float_t kLaddLen5 = 90.27;  // Layer 5: distance between mouting points
321   const Float_t kLaddLen6 = 102.0;  // Layer 6: distance between mouting points
322   const Float_t Z0 = 2.927;         // Distance between V mounting point and Zloc = 0
323                                     // = half ladder length - nom z-position of ladder from gGeoManager
324   const Float_t kMu2Cm = 1e-4;
325
326   TString ssdName = "ITS/SSD";
327
328   TObjArray *ladd_points = fSurveyPoints;  
329   if (ladd_points == 0 || ladd_points->GetEntries() == 0) {
330     AliWarning("No SSD Ladder alignment points found. Skipping");
331     return;
332   }
333   if (ladd_points->GetEntries()!= 2*(34+38)) {
334     AliWarning(Form("Unexpected number of survey points %d, should be 144",ladd_points->GetEntries())); 
335   }
336   Int_t iLadd = 0;
337   for (Int_t ilayer =  4; ilayer <=  5; ilayer ++) {
338     Int_t nLadder = 34; // layer 5
339     if (ilayer == 5)
340       nLadder = 38;     // layer 6
341
342     for (Int_t iLadder = 0; iLadder < nLadder; iLadder++) {
343       TString ladName = ssdName;
344       ladName += ilayer;
345       ladName += "/Ladder";
346       ladName += iLadder;
347
348       AliSurveyPoint *v_point =  (AliSurveyPoint*) ladd_points->FindObject(ladName+"/V");
349       AliSurveyPoint *q_point =  (AliSurveyPoint*) ladd_points->FindObject(ladName+"/Q");
350       if (v_point == 0) {
351         AliWarning(Form("Cannot find V side point for ladder %s",ladName.Data()));
352         continue;
353       }
354       if (q_point == 0) {
355         AliWarning(Form("Cannot find Q side point for ladder %s",ladName.Data()));
356         continue;
357       }
358
359       TString tmp_str;
360       tmp_str.Insert(0,v_point->GetName(),3);
361       Int_t ladder = tmp_str.Atoi();
362       tmp_str="";
363       tmp_str.Insert(0,q_point->GetName(),3);
364       if (tmp_str.Atoi() != ladder) 
365         AliError(Form("Survey data file error. Expect pairs of V,Q points. Got ladders %d %d",ladder,tmp_str.Atoi()));
366
367       // Note: file gives meas-nom in local offline coordinates, 
368       // ie. local z = - global z and local x = - global x (for ladder 508, i.e. top ladder)
369       Double_t dx_loc = v_point->GetX() * kMu2Cm;
370       Double_t dy_loc = v_point->GetY() * kMu2Cm;
371       Double_t dz_loc = v_point->GetZ() * kMu2Cm;
372
373       // rot around z-axis
374       Double_t phi = 0;  // Not measured
375       // rot around y-axis
376       Double_t theta = 0;
377       Double_t psi = 0;
378
379       // Note: local psi = -global psi, psi = atan(-(y(z1) - y(z0)) / (z1-z0))  
380       // local theta = global theta = atan(dx/dz) 
381       // V side is A side is large global z 
382       // Q side is C side is large local z
383
384       if (ladder >= 600) {
385         theta = TMath::ATan((q_point->GetX() - v_point->GetX())*kMu2Cm/kLaddLen6);
386         psi = TMath::ATan((v_point->GetY() - q_point->GetY())*kMu2Cm/kLaddLen6);
387       }
388       else {
389         theta = TMath::ATan((q_point->GetX() - v_point->GetX())*kMu2Cm/kLaddLen5);
390         psi = TMath::ATan((v_point->GetY() - q_point->GetY())*kMu2Cm/kLaddLen5);
391       } 
392
393       // Move along ladder to local Z = 0 point
394       dx_loc += Z0*theta;
395       dy_loc -= Z0*psi;
396
397       // Convert to degrees
398       theta *= 180./TMath::Pi();
399       psi *= 180./TMath::Pi();
400       AliDebug(1,Form("ladname %f %f %f %f %f %f ",dx_loc,dy_loc,dz_loc,psi,theta,phi));  
401       
402       new((*fAlignObjArray)[500+1698+iLadd]) AliAlignObjParams(ladName,0,dx_loc,dy_loc,dz_loc,psi,theta,phi,kFALSE);
403
404       iLadd++;
405     }  // Ladder loop
406   }  // Layer loop
407 }
408
409
410 void AliITSSurveyToAlign::CreateAlignObjDummySSDModules(){
411   // 
412   // Create empty alignment objects
413   // Used when fSurveySSD == 0
414   //
415   for(Int_t imod = 0; imod < 1698; imod++) {
416     Int_t ilayer = (imod < 748) ? AliGeomManager::kSSD1 : AliGeomManager::kSSD2;
417     Int_t imodule = (imod < 748) ? imod : imod - 748;
418
419     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
420     const Char_t *symname = AliGeomManager::SymName(uid);
421
422     new((*fAlignObjArray)[500+imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
423   }//module loop
424 }
425
426
427 void AliITSSurveyToAlign::GetIdPosSDD(Int_t uid, Int_t layer, Int_t module, Int_t iPoint)
428 {
429   // 
430   //    Utility function used by CreateAlignObjSDD
431   // 
432   TGeoHMatrix gMod = *AliGeomManager::GetMatrix(uid); //global matrix of sensor
433   TGeoPNEntry* pne = gGeoManager->GetAlignableEntryByUID(uid);
434   // TString ladderPath = AliGeomManager::SymName(uid);
435   TString ladderPath(pne->GetTitle());
436   if(ladderPath.EndsWith("/")) ladderPath.Remove(TString::kTrailing,'/');
437   ladderPath.Remove(ladderPath.Last('/'));
438   ladderPath.Remove(ladderPath.Last('/'));
439   gGeoManager->cd(ladderPath.Data());
440   TGeoHMatrix gLad = *gGeoManager->GetCurrentMatrix(); // global matrix of ladder
441   TGeoHMatrix rel = gMod; // to equal relative matrix ladder to sensor.
442   TGeoHMatrix invgLad = gLad.Inverse();
443   rel.MultiplyLeft(&invgLad);
444   TGeoRotation* rr = new TGeoRotation("rr",90,90,0,0,90,180);
445   TGeoCombiTrans* ct = 0;
446   if(layer==3) ct= new TGeoCombiTrans(25.,0.,0.,rr);
447   if(layer==4) ct= new TGeoCombiTrans(25.+7.5,0.,0.,rr);
448
449   rel.MultiplyLeft(ct);
450   
451   if((layer==3)&&(module<3)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
452   if((layer==3)&&(module>2)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
453   if((layer==4)&&(module<4)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
454   if((layer==4)&&(module>3)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
455
456   for(Int_t i=0; i<3; i++) fSDDidP[iPoint][i]*=10;
457
458 }
459
460 void AliITSSurveyToAlign::ReadPointNameSDD(const char str[], Int_t &iLayer, Int_t &iLader, Int_t &iModul, Int_t &iPoint)
461 {
462   // 
463   //    Utility function used by CreateAlignObjSDD
464   // 
465   iLayer=-1;
466   iLader=-1;
467   iModul=-1;
468   iPoint=-1;
469
470   if(str[7]=='2') iLayer=3;
471   if(str[7]=='3') iLayer=4;
472
473   if(str[15]=='0') iLader=0;
474   if(str[15]=='1') iLader=1;
475   if(str[15]=='2') iLader=2;
476   if(str[15]=='3') iLader=3;
477   if(str[15]=='4') iLader=4;
478   if(str[15]=='5') iLader=5;
479   if(str[15]=='6') iLader=6;
480   if(str[15]=='7') iLader=7;
481   if(str[15]=='8') iLader=8;
482   if(str[15]=='9') iLader=9;
483
484   Int_t ord=0;
485   if(str[16]=='0') {iLader=10*iLader+0; ord=1;}
486   if(str[16]=='1') {iLader=10*iLader+1; ord=1;}
487   if(str[16]=='2') {iLader=10*iLader+2; ord=1;}
488   if(str[16]=='3') {iLader=10*iLader+3; ord=1;}
489   if(str[16]=='4') {iLader=10*iLader+4; ord=1;}
490   if(str[16]=='5') {iLader=10*iLader+5; ord=1;}
491   if(str[16]=='6') {iLader=10*iLader+6; ord=1;}
492   if(str[16]=='7') {iLader=10*iLader+7; ord=1;}
493   if(str[16]=='8') {iLader=10*iLader+8; ord=1;}
494   if(str[16]=='9') {iLader=10*iLader+9; ord=1;}
495
496   if(str[23+ord]=='0') iModul=0;
497   if(str[23+ord]=='1') iModul=1;
498   if(str[23+ord]=='2') iModul=2;
499   if(str[23+ord]=='3') iModul=3;
500   if(str[23+ord]=='4') iModul=4;
501   if(str[23+ord]=='5') iModul=5;
502   if(str[23+ord]=='6') iModul=6;
503   if(str[23+ord]=='7') iModul=7;
504   if(str[23+ord]=='8') iModul=8;
505   if(str[23+ord]=='9') iModul=9;
506
507   if((str[25+ord]=='R')&&(str[26+ord]=='D')) iPoint=0;
508   if((str[25+ord]=='R')&&(str[26+ord]=='C')) iPoint=1;
509   if((str[25+ord]=='R')&&(str[26+ord]=='U')) iPoint=2;
510   if((str[25+ord]=='L')&&(str[26+ord]=='U')) iPoint=3;
511   if((str[25+ord]=='L')&&(str[26+ord]=='C')) iPoint=4;
512   if((str[25+ord]=='L')&&(str[26+ord]=='D')) iPoint=5;
513   return;
514 }
515
516
517 void AliITSSurveyToAlign::ConvertToRSofModulesAndRotSDD(Int_t Layer, Int_t Module)
518 {
519   // 
520   //    Utility function used by CreateAlignObjSDD
521   // 
522
523   Double_t YmId;
524   Double_t ZmId;
525
526   Double_t YmMe;
527   Double_t ZmMe;
528   Double_t YmMeE;
529   Double_t ZmMeE;
530
531   Double_t x0=fSDDidP[1][0];
532   Double_t z0=fSDDidP[1][2]-0.52;
533   for(Int_t i=0; i<6; i++)
534     {
535       fSDDidP[i][2]-=0.52;
536
537       if(!fSDDisMe[i]) continue; 
538
539       fSDDidP[i][0]-=x0;
540       fSDDidP[i][2]-=z0;
541       fSDDmeP[i][0]-=x0;
542       fSDDmeP[i][2]-=z0;
543                                 
544       YmId=fSDDidP[i][1];
545       ZmId=fSDDidP[i][2];
546                         
547       fSDDidP[i][2]=fSDDidP[i][0];
548       fSDDidP[i][0]=YmId;
549       fSDDidP[i][1]=ZmId;
550                         
551       YmMe=fSDDmeP[i][1];
552       ZmMe=fSDDmeP[i][2];
553                         
554       YmMeE=fSDDmeP[i][4];
555       ZmMeE=fSDDmeP[i][5];
556                         
557       fSDDmeP[i][2]=fSDDmeP[i][0];
558       fSDDmeP[i][0]=YmMe;
559       fSDDmeP[i][1]=ZmMe;
560       fSDDmeP[i][5]=fSDDmeP[i][3];
561       fSDDmeP[i][3]=YmMeE;
562       fSDDmeP[i][4]=ZmMeE;
563                         
564
565       if(((Layer==3)&&(Module>2))||((Layer==4)&&(Module>3)))
566         {
567           fSDDidP[i][0]*=(-1);
568           fSDDidP[i][2]*=(-1);
569           fSDDmeP[i][0]*=(-1);
570           fSDDmeP[i][2]*=(-1);
571         }
572     }   
573 }
574
575
576 void AliITSSurveyToAlign::CalcShiftSDD(Double_t &x0,Double_t &y0,Double_t &z0)
577 {
578   Double_t Xid, Yid, Zid;
579   Double_t Xme, Yme, Zme, sX2, sY2, sZ2;
580   Double_t aX=0., bX=0.;
581   Double_t aY=0., bY=0.;
582   Double_t aZ=0., bZ=0.;
583   for(Int_t iP1=0; iP1<6; iP1++)
584     {
585       if(!fSDDisMe[iP1]) continue;
586       Xid=fSDDidP[iP1][0];
587       Yid=fSDDidP[iP1][1];
588       Zid=fSDDidP[iP1][2];
589       Xme=fSDDmeP[iP1][0];
590       Yme=fSDDmeP[iP1][1];
591       Zme=fSDDmeP[iP1][2];
592       sX2 =fSDDmeP[iP1][3]*fSDDmeP[iP1][3];
593       sY2 =fSDDmeP[iP1][4]*fSDDmeP[iP1][4];
594       sZ2 =fSDDmeP[iP1][5]*fSDDmeP[iP1][5];
595       aX+=(1./sX2);
596       bX+=((Xme-Xid)/sX2); 
597       aY+=(1./sY2);
598       bY+=((Yme-Yid)/sY2); 
599       aZ+=(1./sZ2);
600       bZ+=((Zme-Zid)/sZ2); 
601     }
602   Double_t x1 = bX/aX;
603   Double_t x2 = bY/aY;
604   Double_t x3 = bZ/aZ;
605   x0=x1;
606   y0=x2;
607   z0=x3;
608   return;
609 }
610
611
612 void AliITSSurveyToAlign::CalcShiftRotSDD(Double_t &tet,Double_t &psi,Double_t &phi,Double_t &x0,Double_t &y0,Double_t &z0)
613 {
614   TMatrixD p1(6,6);
615   TMatrixD p2(6,6);
616   TMatrixD p3(6,6);
617   TMatrixD p4(6,6);
618   TMatrixD p5(6,6);
619   TMatrixD p6(6,6);
620   TMatrixD pC(6,6);
621
622   Double_t a11 =0.;
623   Double_t a12 =0.;
624   Double_t a13 =0.;
625   Double_t a14 =0.;
626   Double_t a15 =0.;
627   Double_t a16 =0.;
628
629   Double_t a21 =0.;
630   Double_t a22 =0.;
631   Double_t a23 =0.;
632   Double_t a24 =0.;
633   Double_t a25 =0.;
634   Double_t a26 =0.;
635
636   Double_t a31 =0.;
637   Double_t a32 =0.;
638   Double_t a33 =0.;
639   Double_t a34 =0.;
640   Double_t a35 =0.;
641   Double_t a36 =0.;
642
643   Double_t a41 =0.;
644   Double_t a42 =0.;
645   Double_t a43 =0.;
646   Double_t a44 =0.;
647   Double_t a45 =0.;
648   Double_t a46 =0.;
649
650   Double_t a51 =0.;
651   Double_t a52 =0.;
652   Double_t a53 =0.;
653   Double_t a54 =0.;
654   Double_t a55 =0.;
655   Double_t a56 =0.;
656
657   Double_t a61 =0.;
658   Double_t a62 =0.;
659   Double_t a63 =0.;
660   Double_t a64 =0.;
661   Double_t a65 =0.;
662   Double_t a66 =0.;
663
664   Double_t c1 =0.;
665   Double_t c2 =0.;
666   Double_t c3 =0.;
667   Double_t c4 =0.;
668   Double_t c5 =0.;
669   Double_t c6 =0.;
670
671   Double_t Xid, Yid, Zid;
672   Double_t Xme, Yme, Zme, sX2, sY2, sZ2;
673
674   for(Int_t iP1=0; iP1<=6; iP1++)
675     {
676       if(!fSDDisMe[iP1]) continue;
677
678       Xid= fSDDidP[iP1][0];
679       Yid= fSDDidP[iP1][1];
680       Zid= fSDDidP[iP1][2];
681
682       Xme= fSDDmeP[iP1][0];
683       Yme= fSDDmeP[iP1][1];
684       Zme= fSDDmeP[iP1][2];
685
686       sX2 = fSDDmeP[iP1][3]* fSDDmeP[iP1][3];
687       sY2 = fSDDmeP[iP1][4]* fSDDmeP[iP1][4];
688       sZ2 = fSDDmeP[iP1][5]* fSDDmeP[iP1][5];
689
690       a11+=(Zid*Zid/sX2+Xid*Xid/sZ2);
691       a12-=(Zid*Yid/sX2);
692       a13-=(Xid*Yid/sZ2);
693       a14-=(Zid/sX2);
694       a15 =0.;
695       a16+=(Xid/sZ2);
696       c1+=(Xid*(Zme-Zid)/sZ2-Zid*(Xme-Xid)/sX2); 
697
698       a21-=(Yid*Zid/sX2);
699       a22+=(Xid*Xid/sY2+Yid*Yid/sX2);
700       a23-=(Xid*Zid/sY2);
701       a24+=(Yid/sX2);
702       a25-=(Xid/sY2);
703       a26 =0.;
704       c2+=(Yid*(Xme-Xid)/sX2-Xid*(Yme-Yid)/sY2); 
705
706       a31-=(Yid*Xid/sZ2);
707       a32-=(Xid*Zid/sY2);
708       a33+=(Zid*Zid/sY2+Yid*Yid/sZ2);
709       a34 =0.;
710       a35+=(Zid/sY2);
711       a36-=(Yid/sZ2);
712       c3+=(Zid*(Yme-Yid)/sY2-Yid*(Zme-Zid)/sZ2); 
713
714       a41-=(Zid/sX2);
715       a42+=(Yid/sX2);
716       a43 =0.;
717       a44+=(1./sX2);
718       a45 =0.;
719       a46 =0.;
720       c4+=((Xme-Xid)/sX2); 
721
722       a51 =0.;
723       a52-=(Xid/sY2);
724       a53+=(Zid/sY2);
725       a54 =0.;
726       a55+=(1./sY2);
727       a56 =0.;
728       c5+=((Yme-Yid)/sY2); 
729
730       a61+=(Xid/sZ2);
731       a62 =0.;
732       a63-=(Yid/sZ2);
733       a64 =0.;
734       a65 =0.;
735       a66+=(1./sZ2);
736       c6+=((Zme-Zid)/sZ2); 
737     }
738
739   ///////////////////////////////////////////////////////////////
740
741   Double_t dataC[]={a11, a12, a13, a14, a15, a16, a21, a22, a23, a24, a25, a26,
742                     a31, a32, a33, a34, a35, a36, a41, a42, a43, a44, a45, a46,
743                     a51, a52, a53, a54, a55, a56, a61, a62, a63, a64, a65, a66};
744   Double_t data1[]={c1, a12, a13, a14, a15, a16, c2, a22, a23, a24, a25, a26,
745                     c3, a32, a33, a34, a35, a36, c4, a42, a43, a44, a45, a46,
746                     c5, a52, a53, a54, a55, a56, c6, a62, a63, a64, a65, a66};
747   Double_t data2[]={a11, c1, a13, a14, a15, a16, a21, c2, a23, a24, a25, a26,
748                     a31, c3, a33, a34, a35, a36, a41, c4, a43, a44, a45, a46,
749                     a51, c5, a53, a54, a55, a56, a61, c6, a63, a64, a65, a66};
750   Double_t data3[]={a11, a12, c1, a14, a15, a16, a21, a22, c2, a24, a25, a26,
751                     a31, a32, c3, a34, a35, a36, a41, a42, c4, a44, a45, a46,
752                     a51, a52, c5, a54, a55, a56, a61, a62, c6, a64, a65, a66};
753   Double_t data4[]={a11, a12, a13, c1, a15, a16, a21, a22, a23, c2, a25, a26,
754                     a31, a32, a33, c3, a35, a36, a41, a42, a43, c4, a45, a46,
755                     a51, a52, a53, c5, a55, a56, a61, a62, a63, c6, a65, a66};
756   Double_t data5[]={a11, a12, a13, a14, c1, a16, a21, a22, a23, a24, c2, a26,
757                     a31, a32, a33, a34, c3, a36, a41, a42, a43, a44, c4, a46,
758                     a51, a52, a53, a54, c5, a56, a61, a62, a63, a64, c6, a66};
759   Double_t data6[]={a11, a12, a13, a14, a15, c1, a21, a22, a23, a24, a25, c2,
760                     a31, a32, a33, a34, a35, c3, a41, a42, a43, a44, a45, c4,
761                     a51, a52, a53, a54, a55, c5, a61, a62, a63, a64, a65, c6};
762
763   p1.SetMatrixArray(data1);
764   p2.SetMatrixArray(data2);
765   p3.SetMatrixArray(data3);
766   p4.SetMatrixArray(data4);
767   p5.SetMatrixArray(data5);
768   p6.SetMatrixArray(data6);
769   pC.SetMatrixArray(dataC);
770
771   // cout << "calculating determinants" << endl;
772   Double_t det0=pC.Determinant();
773   Double_t x1 = p1.Determinant()/det0;
774   Double_t x2 = p2.Determinant()/det0;
775   Double_t x3 = p3.Determinant()/det0;
776   Double_t x4 = p4.Determinant()/det0;
777   Double_t x5 = p5.Determinant()/det0;
778   Double_t x6 = p6.Determinant()/det0;
779   //cout << "calculating determinants done" << endl;
780   if (x1 == 0) {
781     AliInfo("p1 singular ");
782     p1.Print();
783   }
784   if (x2 == 0) {
785     AliInfo("p2 singular ");
786     p2.Print();
787   }
788   if (x3 == 0) {
789     AliInfo("p3 singular ");
790     p3.Print();
791   }
792   if (x4 == 0) {
793     AliInfo("p4 singular ");
794     p4.Print();
795   }
796   if (x5 == 0) {
797     AliInfo("p5 singular ");
798     p5.Print();
799   }
800   if (x6 == 0) {
801     AliInfo("p6 singular ");
802     p6.Print();
803   }
804
805
806   tet=x1;
807   psi=x2;
808   phi=x3;
809   x0=x4;
810   y0=x5;
811   z0=x6;
812   return;
813 }
814
815
816