]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSSurveyToAlign.cxx
Fixing coding violations. Method AliITSQADataMakerRec::AreEqual defined static
[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 "TGeoManager.h"
27 #include "TGeoPhysicalNode.h"
28 #include "TMatrixD.h"
29 #include "TMath.h"
30
31 #include "AliITSSurveyToAlign.h"
32 #include "AliSurveyPoint.h"
33 #include "AliAlignObjParams.h"
34 #include "AliGeomManager.h"
35
36 #include "AliLog.h"
37
38 #include "AliCDBManager.h"
39
40 #include "AliITSgeomTGeo.h"
41
42
43 //#include "dataSDDladder.h"
44 ClassImp(AliITSSurveyToAlign)
45
46 const Double_t AliITSSurveyToAlign::fgkLocR[6][3]={{ 3.24,0.21905,-2.4},
47                                                    { 3.58,0.21905, 0. },
48                                                    { 3.24,0.21905,+2.4},
49                                                    {-3.24,0.21905,+2.4},
50                                                    {-3.58,0.21905, 0. },
51                                                    {-3.24,0.21905,-2.4}};
52
53 const Double_t AliITSSurveyToAlign::fgkLocL[6][3]={{-3.24,0.21905, 2.4},
54                                                    {-3.58,0.21905, 0. },
55                                                    {-3.24,0.21905,-2.4},
56                                                    { 3.24,0.21905,-2.4},
57                                                    { 3.58,0.21905, 0. },
58                                                    { 3.24,0.21905, 2.4}};
59
60 const Double_t kRadToDeg = 180./TMath::Pi();
61
62 //________________________________________________________________________
63 AliITSSurveyToAlign::AliITSSurveyToAlign(Int_t run, Int_t repModSDD, Int_t repModVerSDD, Int_t repLadSDD, Int_t repLadVerSDD, Int_t repModSSD, Int_t repModVerSSD, Int_t repLaddSSD, Int_t repLaddVerSSD) :
64   AliSurveyToAlignObjs(),
65   fRun(run),
66   fSDDModuleRepNumber(repModSDD),
67   fSDDModuleRepVersion(repModVerSDD),
68   fSDDLadderRepNumber(repLadSDD),
69   fSDDLadderRepVersion(repLadVerSDD),
70   fSSDModuleRepNumber(repModSSD),
71   fSSDModuleRepVersion(repModVerSSD),
72   fSSDLadderRepNumber(repLaddSSD),
73   fSSDLadderRepVersion(repLaddVerSSD)
74  {
75    // Standard constructor
76   for(Int_t i=0; i<260; i++)
77     {
78       fuidSDDm[i]= 0;
79       fsymnameSDDm[i]=TString("");
80       fxSDDm[i]=0;
81       fySDDm[i]=0;
82       fzSDDm[i]=0;
83       fpsiSDDm[i]=0;
84       ftetSDDm[i]=0;
85       fphiSDDm[i]=0;
86       if(i>35) continue;
87       fuidSDDl[i]=0;
88       fsymnameSDDl[i]=TString("");
89       fxSDDl[i]=0;
90       fySDDl[i]=0;
91       fzSDDl[i]=0;
92       fpsiSDDl[i]=0;
93       ftetSDDl[i]=0;
94       fphiSDDl[i]=0;
95     }
96
97    //   ftypeSDDlad=0;
98   //
99   //  default constructor
100   //  Arguments are report numbers for survey data. 
101   //  The defaults point to reports from detector construction
102   // 
103 }
104
105 //_________________________________________________________________________
106 AliITSSurveyToAlign::AliITSSurveyToAlign(const AliITSSurveyToAlign &align) :
107   AliSurveyToAlignObjs(align),
108   fRun(align.fRun),
109   fSDDModuleRepNumber(align.fSDDModuleRepNumber),
110   fSDDModuleRepVersion(align.fSDDModuleRepVersion),
111   fSDDLadderRepNumber(align.fSDDLadderRepNumber),
112   fSDDLadderRepVersion(align.fSDDLadderRepVersion),
113   fSSDModuleRepNumber(align.fSSDModuleRepNumber),
114   fSSDModuleRepVersion(align.fSSDModuleRepVersion),
115   fSSDLadderRepNumber(align.fSSDLadderRepNumber),
116   fSSDLadderRepVersion(align.fSSDLadderRepVersion)
117 {
118   //
119   //  copy constructor 
120   //
121   for(Int_t i=0; i<260; i++)
122     {
123       fuidSDDm[i]= align.fuidSDDm[i];
124       fsymnameSDDm[i]=TString(align.fsymnameSDDm[i]);
125       fxSDDm[i]=align.fxSDDm[i];
126       fySDDm[i]=align.fySDDm[i];
127       fzSDDm[i]=align.fzSDDm[i];
128       fpsiSDDm[i]=align.fpsiSDDm[i];
129       ftetSDDm[i]=align.ftetSDDm[i];
130       fphiSDDm[i]=align.fphiSDDm[i];
131       if(i>35) continue;
132       fuidSDDl[i]=align.fuidSDDl[i];
133       fsymnameSDDl[i]=TString(align.fsymnameSDDl[i]);
134       fxSDDl[i]=align.fxSDDl[i];
135       fySDDl[i]=align.fySDDl[i];
136       fzSDDl[i]=align.fzSDDl[i];
137       fpsiSDDl[i]=align.fpsiSDDl[i];
138       ftetSDDl[i]=align.ftetSDDl[i];
139       fphiSDDl[i]=align.fphiSDDl[i];
140     }
141
142 }
143
144 //__________________________________________________________________________
145 AliITSSurveyToAlign & AliITSSurveyToAlign::operator =(const AliITSSurveyToAlign& align)  {
146   //
147   // assignment operator
148   //
149   this->~AliITSSurveyToAlign();
150   new(this) AliITSSurveyToAlign(align);
151   return *this;
152 }
153
154 //__________________________________________________________________________
155 AliITSSurveyToAlign::~AliITSSurveyToAlign() {
156   //
157   // destructor
158   //
159 }
160
161 //______________________________________________________________________
162 void AliITSSurveyToAlign::Run() { 
163   //
164   // Runs the full chain
165   // User should call StoreAlignObjToFile or StoreAlignObjToCDB afterwards to 
166   // store output (not included here to leave the choice between the two)
167   //
168
169   // Load ideal geometry from the OCDB
170   AliCDBManager *cdb = AliCDBManager::Instance();
171   cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
172   cdb->SetRun(fRun);
173   AliGeomManager::LoadGeometry();
174
175   if(!CreateAlignObjs()) AliError("Construction of alignment objects from survey failed!");
176
177
178 //______________________________________________________________________
179 Bool_t AliITSSurveyToAlign::CreateAlignObjs() { 
180   // Fill the array of alignment objects with alignment objects
181   // from survey for all three subdetectors
182   //
183
184   //for SPD
185   CreateAlignObjDummySPD();
186
187   ///////////////////////////
188   // for SDD
189   if(!LoadSurveyFromAlienFile("ITS", fSDDModuleRepNumber, fSDDModuleRepVersion)){
190       AliError("Loading of alignment objects from survey for SDD modules failed!");
191       return kFALSE;
192   }
193   CreateAlignObjSDDModules();
194
195   if(!LoadSurveyFromAlienFile("ITS", fSDDLadderRepNumber, fSDDLadderRepVersion)){
196       AliError("Loading of alignment objects from survey for SDD ladder failed!");
197       return kFALSE;
198   }
199   CreateAlignObjSDDLadders();
200   if(!ApplyAlignObjSDD()) return kFALSE;
201
202
203   // for SSD ladders
204   if(!LoadSurveyFromAlienFile("ITS", fSSDLadderRepNumber, fSSDLadderRepVersion)){
205       AliError("Loading of alignment objects from survey for SSD ladders failed!");
206       return kFALSE;
207   }
208   CreateAlignObjSSDLadders();
209
210   // for SSD modules
211   if(!ApplyAlignObjSSDLadders()) return kFALSE; // needed to build correctly the objects for SSD modules
212   if(!LoadSurveyFromAlienFile("ITS", fSSDModuleRepNumber, fSSDModuleRepVersion)){
213       AliError("Loading of alignment objects from survey for SSD modules failed!");
214       return kFALSE;
215   }
216   CreateAlignObjSSDModules();
217
218   return kTRUE;
219 }
220
221 //______________________________________________________________________
222 void AliITSSurveyToAlign::CreateAlignObjDummySPD(){
223   // 
224   // Create alignObjs for SPD
225   //    For the moment, uses 0,0,0,0,0,0
226   //
227   for(Int_t imod = 0; imod < 240; imod++) {
228     Int_t ilayer = (imod < 80) ? AliGeomManager::kSPD1 : AliGeomManager::kSPD2;
229     Int_t imodule = (imod < 80) ? imod : imod - 80;
230
231     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
232     const Char_t *symname = AliGeomManager::SymName(uid);
233
234     new((*fAlignObjArray)[imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
235   }//module loop
236
237 }
238 Bool_t AliITSSurveyToAlign::ApplyAlignObjSDD()
239 {
240   // Apply alignment for SDD
241   Int_t applied=0;
242
243   for(Int_t iLadd=0; iLadd<36; iLadd++)
244     {
245       new((*fAlignObjArray)[240+iLadd]) AliAlignObjParams(fsymnameSDDl[iLadd].Data(), fuidSDDl[iLadd], 
246                                                fxSDDl[iLadd]  , fySDDl[iLadd]  , fzSDDl[iLadd]  , 
247                                                fpsiSDDl[iLadd], ftetSDDl[iLadd], fphiSDDl[iLadd], kFALSE);
248       //      new((*fAlignObjArray)[240+iLadd]) AliAlignObjParams(fsymnameSDDl[iLadd].Data(), fuidSDDl[iLadd], 
249       //                               fxSDDl[iLadd]  , fySDDl[iLadd]  , fzSDDl[iLadd]  , 
250       //                               0, 0, 0, kFALSE);
251
252       AliAlignObjParams* ap = dynamic_cast<AliAlignObjParams*> (fAlignObjArray->UncheckedAt(240+iLadd));
253       //      printf("%s  %f  %f  %f\n",fxSDDl[iLadd], fsymnameSDDl[iLadd].Data(), fySDDl[iLadd]  , fzSDDl[iLadd]);
254       //printf("%d  %f\n", iLadd, fzSDDl[iLadd]);
255
256       if(fsymnameSDDl[iLadd].Contains("SDD") && fsymnameSDDl[iLadd].Contains("Ladder"))
257         {
258           //      printf("%d  %s  %d\n",240+iLadd, fsymnameSDDl[iLadd].Data(),fuidSDDl[iLadd] );
259
260           if(!ap->ApplyToGeometry()) return kFALSE;
261           applied++;
262         }
263       else
264         {
265           AliError("SDD Ladder array is not initialized correctly");
266           return kFALSE;
267         }
268     }
269
270   for(Int_t iMod=0; iMod<260; iMod++)
271     {
272       new((*fAlignObjArray)[240+36+iMod]) AliAlignObjParams(fsymnameSDDm[iMod].Data(), fuidSDDm[iMod], 
273                                              fxSDDm[iMod], fySDDm[iMod], fzSDDm[iMod], 
274                                              fpsiSDDm[iMod], ftetSDDm[iMod], fphiSDDm[iMod], kFALSE);
275
276       //           printf("%d %s  %f  %f  %f\n",240+36+iMod, fsymnameSDDm[iMod].Data(),fxSDDm[iMod], fySDDm[iMod], fzSDDm[iMod] );
277
278           if(!fsymnameSDDm[iMod].Contains("SDD") || !fsymnameSDDm[iMod].Contains("Sensor"))
279             {
280           AliError("SDD Module array is not initialized correctly\n");
281           return kFALSE;
282             }
283
284     }
285
286   AliInfo(Form(" %d alignment objects for SDD ladders applied to geometry.",applied));
287   return kTRUE;
288 }
289
290 //______________________________________________________________________
291 void AliITSSurveyToAlign::CreateAlignObjSDDModules(){
292   //
293   // Create alignment objects for SDD
294   // Called by Run()
295   //
296   Int_t uid = 0;
297   const char* symname = 0;
298   AliSurveyPoint* pt = 0;
299  
300   Int_t iModuleIndex=240;
301   Int_t iModule0=0;
302   Int_t iLadder0=0;
303   Int_t iLayer0=3;
304   Int_t nModules=0;
305
306   if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
307     AliWarning("SDD survey data are not available, using zero values");
308     CreateAlignObjDummySDDModules();
309     return;
310   }
311
312   for(Int_t imod = 1; imod < fSurveyPoints->GetEntries(); imod++) {
313     pt = (AliSurveyPoint*) fSurveyPoints->At(imod);
314     if(!pt) continue;
315
316     Int_t iLayer, iLadder, iModule, iPoint;
317     ReadPointNameSDD(pt->GetName(),iLayer, iLadder, iModule, iPoint);
318
319     if(iModule==iModule0)
320     {
321       fSDDmeP[iPoint][0]=pt->GetX();
322       fSDDmeP[iPoint][1]=pt->GetY();
323       fSDDmeP[iPoint][2]=pt->GetZ();
324       fSDDmeP[iPoint][3]=pt->GetPrecisionX();
325       fSDDmeP[iPoint][4]=pt->GetPrecisionY();
326       fSDDmeP[iPoint][5]=pt->GetPrecisionZ();
327       fSDDisMe[iPoint]=kTRUE;
328
329       if(iLayer==3) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-240);
330       if(iLayer==4) uid = AliGeomManager::LayerToVolUID(iLayer0,iModuleIndex-324);
331       symname = AliGeomManager::SymName(uid);
332       GetIdPosSDD(uid,iLayer0, iModule0, iPoint);
333       nModules++;
334     }
335     //    cout << "Points red module " << imod << endl;
336     if((iModule!=iModule0)||(imod==(fSurveyPoints->GetEntries()-1)))
337     {
338       ConvertToRSofModulesAndRotSDD(iLayer0, iModule0);
339
340       Double_t tet = 0.;
341       Double_t psi =0.;
342       Double_t phi = 0.;
343       Double_t x0  = 0.;
344       Double_t y0  =0.;
345       Double_t z0  = 0.;
346
347       if(nModules==2) CalcShiftSDD(x0,y0,z0);
348       if(nModules>2)   CalcShiftRotSDD(tet, psi, phi, x0, y0, z0);
349       tet*=kRadToDeg;
350       psi*=kRadToDeg;
351       phi*=kRadToDeg;
352
353 //    printf("%s  %d  %f  %f  %f  %f  %f  %f\n",symname, uid, x0/10., y0/10., z0/10., psi, tet, phi);
354 //      cout << "Allocate alignobjparams " << imod << endl;
355
356 //      new((*fAlignObjArray)[iModuleIndex]) AliAlignObjParams(symname, uid, x0/10., y0/10., z0/10., psi, tet, phi, kFALSE);
357 //       printf("INDEX:   Module: %d\n",iModuleIndex);
358
359       fsymnameSDDm[iModuleIndex-240]=TString(symname);
360       fuidSDDm[iModuleIndex-240]=uid;
361       fxSDDm[iModuleIndex-240]=x0/10.;
362       fySDDm[iModuleIndex-240]=y0/10.;
363       fzSDDm[iModuleIndex-240]=z0/10.;
364       fpsiSDDm[iModuleIndex-240]=psi;
365       ftetSDDm[iModuleIndex-240]=tet;
366       fphiSDDm[iModuleIndex-240]=phi;
367
368       //      new((*fAlignObjArray)[36+iModuleIndex]) AliAlignObjParams(fsymnameSDDm[iModuleIndex-240].Data(), fuidSDDm[iModuleIndex-240], 
369       //                                             fxSDDm[iModuleIndex-240], fySDDm[iModuleIndex-240], fzSDDm[iModuleIndex-240], 
370       //                                             fpsiSDDm[iModuleIndex-240], ftetSDDm[iModuleIndex-240], fphiSDDm[iModuleIndex-240], kFALSE);
371
372       iModule0=iModule;
373       iLayer0=iLayer;
374       iLadder0=iLadder;
375       nModules=0;
376       iModuleIndex = AliITSgeomTGeo::GetModuleIndex(iLayer,iLadder+1,iModule+1);
377       for(Int_t i=0; i<6;i++) fSDDisMe[i]=kFALSE;
378       if(imod!=(fSurveyPoints->GetEntries()-1)) imod--;
379     }
380   }//module loop
381 }
382
383 //______________________________________________________________________
384 void AliITSSurveyToAlign::CreateAlignObjDummySDDModules(){
385   // 
386   // Create empty alignment objects
387   // Used when fSurveySDD == 0
388   //
389   for(Int_t imod = 0; imod < 260; imod++) {
390
391     Int_t ilayer = (imod < 84) ? AliGeomManager::kSDD1 : AliGeomManager::kSDD2;
392     Int_t imodule = (imod < 84) ? imod : imod - 84;
393
394     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
395     const Char_t *symname = AliGeomManager::SymName(uid);
396
397       fsymnameSDDm[imod]=TString(symname);
398       fuidSDDm[imod]=uid;
399       fxSDDm[imod]=0.;
400       fySDDm[imod]=0.;
401       fzSDDm[imod]=0.;
402       fpsiSDDm[imod]=0.;
403       ftetSDDm[imod]=0.;
404       fphiSDDm[imod]=0.;
405
406     //    new((*fAlignObjArray)[imod+36+240]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
407   }//module loop
408 }
409
410 //______________________________________________________________________
411 void AliITSSurveyToAlign::CreateAlignObjSSDModules(){
412   //
413   // Create alignment objects for SSD modules
414   // Objects for SSD ladders must be applied to geometry first
415   //
416   Double_t sx, sz;
417   const Float_t kMu2Cm = 1e-4;
418   const Float_t kSensLength = 7.464;
419   const Int_t kSSDMODULES = 1698;
420
421   if (fSurveyPoints == 0 || fSurveyPoints->GetEntries() == 0) {
422     AliWarning("SSD module survey data not available; using dummy values");
423     CreateAlignObjDummySSDModules();
424     return;
425   }
426
427   // First do module-by-module
428
429   for(Int_t imod = 500; imod < kSSDMODULES + 500; imod++) {
430     Int_t iLayer, iLadder, iLaddMod;
431     AliITSgeomTGeo::GetModuleId(imod,iLayer,iLadder,iLaddMod);  // returns 1-based numbers
432  
433     TString pname="ITS/SSD";
434     pname += iLayer-1;
435     pname += "/Ladder";
436     pname += iLadder-1;
437     pname += "/Sensor";
438     pname += iLaddMod-1;
439     AliSurveyPoint *pt1 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point0");
440     AliSurveyPoint *pt2 = (AliSurveyPoint*) fSurveyPoints->FindObject(pname+"/Point1");
441     if(!pt1 || !pt2) {
442       AliWarning(Form("No Survey points for iladd %d imod %d",iLadder,iLaddMod));
443       continue;
444     }
445
446     sx = 0.5*(pt1->GetX() + pt2->GetX()) * kMu2Cm;
447     sz = 0.5*(pt1->GetZ() + pt2->GetZ()) * kMu2Cm;
448
449     // Minus sign to change local coordinate convention 
450     Float_t theta = -(pt2->GetZ() - pt1->GetZ())*kMu2Cm/kSensLength;
451
452     theta *= kRadToDeg;
453     Int_t iLayMod = imod - 500;
454     if (iLayer == 6)
455       iLayMod -= 748;
456     Int_t uid = AliGeomManager::LayerToVolUID(iLayer,iLayMod);
457
458     const Char_t *symname = AliGeomManager::SymName(uid);
459     if (pname.CompareTo(symname) != 0)
460       AliWarning(Form("Mapping mismatch survey point %s volume name %s",pname.Data(),symname));
461     /*
462     if (imod >= 676 && imod <= 697) {
463       cout << "ilayer " << iLayer << " imod " << imod 
464            << " uid " << uid << " name " << symname 
465            << " survey shift " << sx << " " << 0 << " " << sz << endl
466            << " theta " << theta << endl;
467     }
468     */
469     new((*fAlignObjArray)[imod+36]) AliAlignObjParams(symname, uid, sx, 0, sz, 0., theta, 0., kFALSE);
470   } //module loop
471 }
472
473 //______________________________________________________________________
474 Bool_t AliITSSurveyToAlign::ApplyAlignObjSSDLadders(){
475   //
476   //   Apply alignment objects for SSD ladders to geometry, needed to correctly
477   //   build alignment objects for SSD modules
478   // 
479   Int_t applied=0;
480
481   for(Int_t jj=0; jj<fAlignObjArray->GetEntriesFast(); jj++)
482   {
483       AliAlignObjParams* ap = dynamic_cast<AliAlignObjParams*> (fAlignObjArray->UncheckedAt(jj));
484       if(ap) 
485       {
486           TString sName(ap->GetSymName());
487           if(sName.Contains("SSD") && sName.Contains("Ladder"))
488           {
489               if(!ap->ApplyToGeometry()) return kFALSE;
490               applied++;
491           }
492       }
493   }
494   AliInfo(Form(" %d alignment objects for SSD ladders applied to geometry.",applied));
495
496   return kTRUE;
497 }
498 //______________________________________________________________________
499 /*
500 Bool_t AliITSSurveyToAlign::ApplyAlignObjSDDLadders(){
501   //
502   //   Apply alignment objects for SDD ladders to geometry, needed to correctly
503   //   build alignment objects for SDD modules
504   // 
505   Int_t applied=0;
506
507   for(Int_t jj=0; jj<fAlignObjArray->GetEntriesFast(); jj++)
508   {
509       AliAlignObjParams* ap = dynamic_cast<AliAlignObjParams*> (fAlignObjArray->UncheckedAt(jj));
510       if(ap) 
511       {
512           TString sName(ap->GetSymName());
513 //        printf("%s\n",sName.Data());
514           if(sName.Contains("SDD") && sName.Contains("Ladder"))
515           {
516               if(!ap->ApplyToGeometry()) return kFALSE;
517               applied++;
518           }
519       }
520   }
521   AliInfo(Form(" %d alignment objects for SDD ladders applied to geometry.",applied));
522
523   return kTRUE;
524 }
525 */
526 //______________________________________________________________________
527 void AliITSSurveyToAlign::CreateAlignObjDummySDDLadders()
528 {
529   // 
530   // Create empty alignment objects
531   TString sddName = "ITS/SDD";
532   Int_t iLadd = 0;
533
534   for (Int_t ilayer =  0; ilayer <  2; ilayer ++)
535     {
536     Int_t nLadder = 14;              // layer SDD1
537     if (ilayer == 1)  nLadder = 22;  // layer SDD2
538     for (Int_t iLadder = 0; iLadder < nLadder; iLadder++) {
539       TString ladName = sddName;
540       ladName += (ilayer+2);      
541       ladName += "/Ladder";
542       ladName += iLadder;
543       fsymnameSDDl[iLadd]=TString(ladName);
544       fuidSDDl[iLadd]=0;
545       fxSDDl[iLadd]=0;
546       fySDDl[iLadd]=10;
547       fzSDDl[iLadd]=0;
548       fpsiSDDl[iLadd]=0;
549       ftetSDDl[iLadd]=0;
550       fphiSDDl[iLadd]=0;
551       iLadd++;
552     }  // Ladder loop
553   }  // Layer loop
554
555 }
556
557
558
559 void AliITSSurveyToAlign::CreateAlignObjSDDLadders(){
560   //
561   //   Alignment objects from survey for SDD ladders
562   // 
563   const Float_t kLaddLenz[2] ={500, 650};  // Layer 2,3: distance between mouting points along z (mm)
564   const Float_t kLaddLenx    = 28 ;        // Layer 2,3: distance between mouting points along x (mm)
565
566   TString sddName = "ITS/SDD";
567
568  TObjArray *ladderPoints = fSurveyPoints;  
569   if (ladderPoints == 0 || ladderPoints->GetEntries() == 0) {
570     AliWarning("No SDD Ladder alignment points found. Skipping");
571     return;
572   }
573   if (ladderPoints->GetEntries()!= 72) {
574     AliWarning(Form("Unexpected number of survey points %d, should be 72",ladderPoints->GetEntries())); 
575   }
576  
577
578 /*
579 TAlien::Connect("alien://");
580 gSystem->Load("libXMLParser.so");
581 .x loadlibs.C 
582
583 AliCDBManager *cdb = AliCDBManager::Instance();
584 cdb->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
585 cdb->SetRun(0);
586 AliGeomManager::LoadGeometry();
587 AliITSSurveyToAlign *a = new AliITSSurveyToAlign(); 
588
589 a->CreateAlignObjSDDLadders()                     
590 a->ApplyAlignObjSDDLadders();
591
592 a->LoadSurveyFromAlienFile("ITS", 845069, 1);
593 a->CreateAlignObjSDD();
594
595 a->CreateAlignObjs();
596 */
597
598 //  Int_t type =0;
599   ////////////////////////////////////////////////////////////
600   // pRB2X[layer][ladder][dx,dy,dz]
601   Double_t pRB24[2][22][3]; //(mm)
602   Double_t pRB26[2][22][3]; //(mm)
603
604
605 pRB24[0][0][2] = 0.0228; pRB26[0][0][2] = 0.0284; pRB24[0][0][0] =-0.0377; pRB26[0][0][0] =-0.0008; pRB24[0][0][1] =-0.0336; pRB26[0][0][1] = 0.0540;
606 pRB24[0][1][2] =-0.0175; pRB26[0][1][2] =-0.0186; pRB24[0][1][0] =-0.0143; pRB26[0][1][0] = 0.0348; pRB24[0][1][1] =-0.0123; pRB26[0][1][1] = 0.0407;
607 pRB24[0][2][2] = 0.1313; pRB26[0][2][2] = 0.0302; pRB24[0][2][0] = 0.2408; pRB26[0][2][0] = 0.1113; pRB24[0][2][1] =-0.0150; pRB26[0][2][1] =-0.0194;
608 pRB24[0][3][2] = 0.0307; pRB26[0][3][2] = 0.0350; pRB24[0][3][0] = 0.1140; pRB26[0][3][0] = 0.0952; pRB24[0][3][1] =-0.0423; pRB26[0][3][1] =-0.0375;
609 pRB24[0][4][2] =-0.0159; pRB26[0][4][2] = 0.0185; pRB24[0][4][0] = 0.0713; pRB26[0][4][0] = 0.0803; pRB24[0][4][1] =-0.1311; pRB26[0][4][1] =-0.1201;
610 pRB24[0][5][2] =-0.0128; pRB26[0][5][2] =-0.0209; pRB24[0][5][0] = 0.0436; pRB26[0][5][0] = 0.0606; pRB24[0][5][1] =-0.1551; pRB26[0][5][1] =-0.0587;
611 pRB24[0][6][2] =-0.0257; pRB26[0][6][2] =-0.0044; pRB24[0][6][0] = 0.0380; pRB26[0][6][0] = 0.0264; pRB24[0][6][1] =-0.1913; pRB26[0][6][1] =-0.1851;
612 pRB24[0][7][2] = 0.0185; pRB26[0][7][2] = 0.1450; pRB24[0][7][0] =-0.0406; pRB26[0][7][0] =-0.0426; pRB24[0][7][1] = 0.1564; pRB26[0][7][1] = 0.2998;
613 pRB24[0][8][2] = 0.0048; pRB26[0][8][2] = 0.0077; pRB24[0][8][0] = 0.0290; pRB26[0][8][0] = 0.0361; pRB24[0][8][1] = 0.1321; pRB26[0][8][1] = 0.1679;
614 pRB24[0][9][2] = 0.0049; pRB26[0][9][2] = 0.0115; pRB24[0][9][0] = 0.0405; pRB26[0][9][0] = 0.1058; pRB24[0][9][1] = 0.0319; pRB26[0][9][1] = 0.2464;
615 pRB24[0][10][2]=-0.0017; pRB26[0][10][2]= 0.0100; pRB24[0][10][0]= 0.0724; pRB26[0][10][0]= 0.1035; pRB24[0][10][1]= 0.0561; pRB26[0][10][1]= 0.0713;
616 pRB24[0][11][2]=-0.0021; pRB26[0][11][2]= 0.0202; pRB24[0][11][0]= 0.0590; pRB26[0][11][0]= 0.0596; pRB24[0][11][1]=-0.0875; pRB26[0][11][1]=-0.0591;
617 pRB24[0][12][2]=-0.0200; pRB26[0][12][2]= 0.0242; pRB24[0][12][0]= 0.0664; pRB26[0][12][0]= 0.0780; pRB24[0][12][1]=-0.0197; pRB26[0][12][1]=-0.0227;
618 pRB24[0][13][2]=-0.0382; pRB26[0][13][2]=-0.0139; pRB24[0][13][0]=-0.0320; pRB26[0][13][0]=-0.0136; pRB24[0][13][1]=-0.0798; pRB26[0][13][1]=-0.0843;
619 pRB24[1][0][2] = 0.0191; pRB26[1][0][2] = 0.0295; pRB24[1][0][0] =-0.0776; pRB26[1][0][0] =-0.0585; pRB24[1][0][1] =-0.0148; pRB26[1][0][1] = 0.0123;
620 pRB24[1][1][2] = 0.0135; pRB26[1][1][2] = 0.0074; pRB24[1][1][0] =-0.0289; pRB26[1][1][0] =-0.0213; pRB24[1][1][1] = 0.0064; pRB26[1][1][1] = 0.0289;
621 pRB24[1][2][2] = 0.0096; pRB26[1][2][2] = 0.0011; pRB24[1][2][0] = 0.0123; pRB26[1][2][0] = 0.0404; pRB24[1][2][1] =-0.0246; pRB26[1][2][1] =-0.0064;
622 pRB24[1][3][2] = 0.0209; pRB26[1][3][2] =-0.0049; pRB24[1][3][0] = 0.0322; pRB26[1][3][0] = 0.0447; pRB24[1][3][1] = 0.0066; pRB26[1][3][1] = 0.0284;
623 pRB24[1][4][2] = 0.0149; pRB26[1][4][2] =-0.0033; pRB24[1][4][0] = 0.0530; pRB26[1][4][0] = 0.0822; pRB24[1][4][1] =-0.0455; pRB26[1][4][1] =-0.0365;
624 pRB24[1][5][2] = 0.0395; pRB26[1][5][2] = 0.0094; pRB24[1][5][0] = 0.0633; pRB26[1][5][0] = 0.0933; pRB24[1][5][1] =-0.1133; pRB26[1][5][1] =-0.1070;
625 pRB24[1][6][2] = 0.0288; pRB26[1][6][2] =-0.0002; pRB24[1][6][0] = 0.0692; pRB26[1][6][0] = 0.0916; pRB24[1][6][1] =-0.1670; pRB26[1][6][1] =-0.1670;
626 pRB24[1][7][2] = 0.0238; pRB26[1][7][2] =-0.0090; pRB24[1][7][0] = 0.0625; pRB26[1][7][0] = 0.0607; pRB24[1][7][1] =-0.1592; pRB26[1][7][1] =-0.1678;
627 pRB24[1][8][2] = 0.0196; pRB26[1][8][2] =-0.0738; pRB24[1][8][0] = 0.0639; pRB26[1][8][0] = 0.0686; pRB24[1][8][1] =-0.2050; pRB26[1][8][1] =-0.2056;
628 pRB24[1][9][2] = 0.0029; pRB26[1][9][2] =-0.0051; pRB24[1][9][0] = 0.0178; pRB26[1][9][0] = 0.0170; pRB24[1][9][1] =-0.1042; pRB26[1][9][1] =-0.1159;
629 pRB24[1][10][2]=-0.0108; pRB26[1][10][2]=-0.0023; pRB24[1][10][0]=-0.0005; pRB26[1][10][0]= 0.0119; pRB24[1][10][1]=-0.1353; pRB26[1][10][1]=-0.1461;
630 pRB24[1][11][2]=-0.0131; pRB26[1][11][2]=-0.0233; pRB24[1][11][0]=-0.0202; pRB26[1][11][0]=-0.0242; pRB24[1][11][1]=-0.1883; pRB26[1][11][1]=-0.2031;
631 pRB24[1][12][2]= 0.0028; pRB26[1][12][2]= 0.0036; pRB24[1][12][0]=-0.0066; pRB26[1][12][0]= 0.0011; pRB24[1][12][1]= 0.2024; pRB26[1][12][1]= 0.2382;
632 pRB24[1][13][2]= 0.0111; pRB26[1][13][2]= 0.0029; pRB24[1][13][0]= 0.0283; pRB26[1][13][0]= 0.0287; pRB24[1][13][1]= 0.2057; pRB26[1][13][1]= 0.2384;
633 pRB24[1][14][2]= 0.0140; pRB26[1][14][2]=-0.0657; pRB24[1][14][0]= 0.0682; pRB26[1][14][0]= 0.0825; pRB24[1][14][1]= 0.1650; pRB26[1][14][1]= 0.2545;
634 pRB24[1][15][2]= 0.0263; pRB26[1][15][2]=-0.0013; pRB24[1][15][0]= 0.0909; pRB26[1][15][0]= 0.0709; pRB24[1][15][1]= 0.1093; pRB26[1][15][1]= 0.1321;
635 pRB24[1][16][2]= 0.0025; pRB26[1][16][2]=-0.0045; pRB24[1][16][0]= 0.0672; pRB26[1][16][0]= 0.0955; pRB24[1][16][1]= 0.0745; pRB26[1][16][1]= 0.0901;
636 pRB24[1][17][2]= 0.0060; pRB26[1][17][2]= 0.0035; pRB24[1][17][0]= 0.0664; pRB26[1][17][0]= 0.0739; pRB24[1][17][1]= 0.0471; pRB26[1][17][1]= 0.0598;
637 pRB24[1][18][2]=-0.0124; pRB26[1][18][2]=-0.0168; pRB24[1][18][0]= 0.0710; pRB26[1][18][0]= 0.0866; pRB24[1][18][1]= 0.0123; pRB26[1][18][1]= 0.0237;
638 pRB24[1][19][2]=-0.0125; pRB26[1][19][2]=-0.0178; pRB24[1][19][0]= 0.0433; pRB26[1][19][0]= 0.0535; pRB24[1][19][1]= 0.0234; pRB26[1][19][1]= 0.0262;
639 pRB24[1][20][2]=-0.0021; pRB26[1][20][2]= 0.0100; pRB24[1][20][0]= 0.0213; pRB26[1][20][0]= 0.0394; pRB24[1][20][1]= 0.0734; pRB26[1][20][1]= 0.0677;
640 pRB24[1][21][2]=-0.0490; pRB26[1][21][2]=-0.0222; pRB24[1][21][0]=-0.0269; pRB26[1][21][0]=-0.0039; pRB24[1][21][1]= 0.0160; pRB26[1][21][1]= 0.0076;
641
642
643
644   Int_t iLadd = 0;
645
646   for (Int_t ilayer =  0; ilayer <  2; ilayer ++)
647     {
648     Int_t nLadder = 14;              // layer SDD1
649     if (ilayer == 1)  nLadder = 22;  // layer SDD2
650
651     for (Int_t iLadder = 0; iLadder < nLadder; iLadder++) {
652       TString ladName = sddName;
653       ladName += (ilayer+2);      
654       ladName += "/Ladder";
655       ladName += iLadder;
656       Double_t drLoc[3];
657       //      for(Int_t i=0; i<3; i++) drLoc[i]=(pRB26[ilayer][iLadder][i]+pRB24[ilayer][iLadder][i])/20.; // average 
658
659
660       /////////////////////////////////////////////////////////////////////////////
661       AliSurveyPoint *p24 =  (AliSurveyPoint*) ladderPoints->FindObject(ladName+"/RB24");
662       AliSurveyPoint *p26 =  (AliSurveyPoint*) ladderPoints->FindObject(ladName+"/RB26");
663       if (p24 == 0) {
664         AliWarning(Form("Cannot find RB24 side point for ladder %s",ladName.Data()));
665         continue;
666       }
667       if (p26 == 0) {
668         AliWarning(Form("Cannot find RB26 side point for ladder %s",ladName.Data()));
669         continue;
670       }
671
672       TString tmpStr;
673       tmpStr.Insert(0,p24->GetName(),3);
674       Int_t ladder = tmpStr.Atoi();
675       tmpStr="";
676       tmpStr.Insert(0,p26->GetName(),3);
677       if (tmpStr.Atoi() != ladder) 
678         AliError(Form("Survey data file error. Expect pairs of RB24, RB26 points. Got ladders %d %d",ladder,tmpStr.Atoi()));
679
680       for(Int_t i=0; i<3; i++) drLoc[i]=pRB24[ilayer][iLadder][i]/10.;
681
682       Double_t x24, y24, z24;
683       Double_t x26, y26, z26;
684
685       x24=p24->GetX();
686       y24=p24->GetY();
687       z24=p24->GetZ();
688       x26=p26->GetX();
689       y26=p26->GetY();
690       z26=p26->GetZ();
691
692       // for top ladders: RS(local) = RS(global) + Y_shift 
693       // rot around z-axis
694       Double_t phi = 0;  // Not measured
695       // rot around y-axis
696       Double_t theta = 0;
697       // rot around x-axis
698       Double_t psi = 0;
699
700
701       psi=TMath::ATan((y26-y24)/(kLaddLenz[ilayer]+z24-z26));
702       Double_t tgtet0 = kLaddLenx/kLaddLenz[ilayer];
703       Double_t tgtet1 = (x24-x26+kLaddLenx)/(kLaddLenz[ilayer]+z24-z26);
704       theta=TMath::ATan((tgtet1-tgtet0)/(1+tgtet1*tgtet0));
705
706       Double_t x0=x24-theta*kLaddLenz[ilayer]/2;
707       Double_t y0=y24+psi*kLaddLenz[ilayer]/2;
708       Double_t z0=z24+theta*kLaddLenx/2;
709
710       theta*= kRadToDeg;
711       psi*= kRadToDeg;
712
713       AliDebug(1,Form("ladname %f %f %f %f %f %f ",drLoc[0],drLoc[1],drLoc[2],psi,theta,phi));  
714       // local delta transformation by passing 3 shifts (in centimeters) and 3 angles (expressed in degrees)    
715       //      new((*fAlignObjArray)[500+1698+144+iLadd]) AliAlignObjParams(ladName,0,drLoc[0],drLoc[1],drLoc[2],psi,theta,phi,kFALSE);
716       fsymnameSDDl[iLadd]=TString(ladName);
717       fuidSDDl[iLadd]=0;
718       fxSDDl[iLadd]=x0/10.;
719       fySDDl[iLadd]=y0/10.;
720       fzSDDl[iLadd]=z0/10.;
721       fpsiSDDl[iLadd]=psi;
722       ftetSDDl[iLadd]=theta;
723       fphiSDDl[iLadd]=phi;
724       //      new((*fAlignObjArray)[240+iLadd]) AliAlignObjParams(fsymnameSDDl[iLadd].Data(), fuidSDDl[iLadd], 
725       //                                               fxSDDl[iLadd]  , fySDDl[iLadd]  , fzSDDl[iLadd]  , 
726       //                                               fpsiSDDl[iLadd], ftetSDDl[iLadd], fphiSDDl[iLadd], kFALSE);
727       //       printf("INDEX:   Ladder: %d\n",iLadd);
728       iLadd++;
729     }  // Ladder loop
730   }  // Layer loop
731 }
732 ////////////////////////////////////////////////////////////////////////////////////////
733
734
735
736 //______________________________________________________________________
737 void AliITSSurveyToAlign::CreateAlignObjSSDLadders(){
738   //
739   //   Alignment objects from survey for SSD ladders (Torino data)
740   // 
741   const Float_t kLaddLen5 = 90.27;  // Layer 5: distance between mouting points
742   const Float_t kLaddLen6 = 102.0;  // Layer 6: distance between mouting points
743   const Float_t zLag = 2.927;         // Distance between V mounting point and Zloc = 0
744                                     // = half ladder length - nom z-position of ladder from gGeoManager
745   const Float_t kMu2Cm = 1e-4;
746
747   TString ssdName = "ITS/SSD";
748
749   TObjArray *ladderPoints = fSurveyPoints;  
750   if (ladderPoints == 0 || ladderPoints->GetEntries() == 0) {
751     AliWarning("No SSD Ladder alignment points found. Skipping");
752     return;
753   }
754   if (ladderPoints->GetEntries()!= 2*(34+38)) {
755     AliWarning(Form("Unexpected number of survey points %d, should be 144",ladderPoints->GetEntries())); 
756   }
757   Int_t iLadd = 0;
758   for (Int_t ilayer =  4; ilayer <=  5; ilayer ++) {
759     Int_t nLadder = 34; // layer 5
760     if (ilayer == 5)
761       nLadder = 38;     // layer 6
762
763     for (Int_t iLadder = 0; iLadder < nLadder; iLadder++) {
764       TString ladName = ssdName;
765       ladName += ilayer;
766       ladName += "/Ladder";
767       ladName += iLadder;
768
769       AliSurveyPoint *vPoint =  (AliSurveyPoint*) ladderPoints->FindObject(ladName+"/V");
770       AliSurveyPoint *qPoint =  (AliSurveyPoint*) ladderPoints->FindObject(ladName+"/Q");
771       if (vPoint == 0) {
772         AliWarning(Form("Cannot find V side point for ladder %s",ladName.Data()));
773         continue;
774       }
775       if (qPoint == 0) {
776         AliWarning(Form("Cannot find Q side point for ladder %s",ladName.Data()));
777         continue;
778       }
779
780       TString tmpStr;
781       tmpStr.Insert(0,vPoint->GetName(),3);
782       Int_t ladder = tmpStr.Atoi();
783       tmpStr="";
784       tmpStr.Insert(0,qPoint->GetName(),3);
785       if (tmpStr.Atoi() != ladder) 
786         AliError(Form("Survey data file error. Expect pairs of V,Q points. Got ladders %d %d",ladder,tmpStr.Atoi()));
787
788       // Note: file gives meas-nom in local offline coordinates, 
789       // ie. local z = - global z and local x = - global x (for ladder 508, i.e. top ladder)
790       Double_t dxLoc = vPoint->GetX() * kMu2Cm;
791       Double_t dyLoc = vPoint->GetY() * kMu2Cm;
792       Double_t dzLoc = vPoint->GetZ() * kMu2Cm;
793
794       // rot around z-axis
795       Double_t phi = 0;  // Not measured
796       // rot around y-axis
797       Double_t theta = 0;
798       Double_t psi = 0;
799
800       // Note: local psi = -global psi, psi = atan(-(y(z1) - y(z0)) / (z1-z0))  
801       // local theta = global theta = atan(dx/dz) 
802       // V side is A side is large global z 
803       // Q side is C side is large local z
804
805
806       if (ladder >= 600) {
807         theta = TMath::ATan((qPoint->GetX() - vPoint->GetX())*kMu2Cm/kLaddLen6);
808         psi = TMath::ATan((vPoint->GetY() - qPoint->GetY())*kMu2Cm/kLaddLen6);
809       }
810       else {
811         theta = TMath::ATan((qPoint->GetX() - vPoint->GetX())*kMu2Cm/kLaddLen5);
812         psi = TMath::ATan((vPoint->GetY() - qPoint->GetY())*kMu2Cm/kLaddLen5);
813       } 
814
815       // Move along ladder to local Z = 0 point
816       dxLoc += zLag*theta;
817       dyLoc -= zLag*psi;
818
819       // Convert to degrees
820       theta *= kRadToDeg;
821       psi *= kRadToDeg;
822       AliDebug(1,Form("ladname %f %f %f %f %f %f ",dxLoc,dyLoc,dzLoc,psi,theta,phi));  
823       
824       new((*fAlignObjArray)[500+36+1698+iLadd]) AliAlignObjParams(ladName,0,dxLoc,dyLoc,dzLoc,psi,theta,phi,kFALSE);
825
826       iLadd++;
827     }  // Ladder loop
828   }  // Layer loop
829 }
830 ////////////////////////////////////////////////////////////////////////////////////////
831
832 //______________________________________________________________________
833 void AliITSSurveyToAlign::CreateAlignObjDummySSDModules(){
834   // 
835   // Create empty alignment objects
836   // Used when fSurveySSD == 0
837   //
838   for(Int_t imod = 0; imod < 1698; imod++) {
839     Int_t ilayer = (imod < 748) ? AliGeomManager::kSSD1 : AliGeomManager::kSSD2;
840     Int_t imodule = (imod < 748) ? imod : imod - 748;
841
842     Int_t uid = AliGeomManager::LayerToVolUID(ilayer,imodule);
843     const Char_t *symname = AliGeomManager::SymName(uid);
844
845     new((*fAlignObjArray)[500+36+imod]) AliAlignObjParams(symname, uid, 0., 0., 0., 0., 0., 0., kTRUE);
846   }//module loop
847 }
848
849
850 //______________________________________________________________________
851 void AliITSSurveyToAlign::GetIdPosSDD(Int_t uid, Int_t layer, Int_t module, Int_t iPoint)
852 {
853   // 
854   //    Utility function used by CreateAlignObjSDD
855   // 
856   TGeoHMatrix gMod = *AliGeomManager::GetMatrix(uid); //global matrix of sensor
857   TGeoPNEntry* pne = gGeoManager->GetAlignableEntryByUID(uid);
858   // TString ladderPath = AliGeomManager::SymName(uid);
859   TString ladderPath(pne->GetTitle());
860   if(ladderPath.EndsWith("/")) ladderPath.Remove(TString::kTrailing,'/');
861   ladderPath.Remove(ladderPath.Last('/'));
862   //  ladderPath.Remove(ladderPath.Last('/'));
863   gGeoManager->cd(ladderPath.Data());
864   TGeoHMatrix gLad = *gGeoManager->GetCurrentMatrix(); // global matrix of ladder
865
866
867   TGeoHMatrix rel = gMod; // to equal relative matrix ladder to sensor.
868   TGeoHMatrix invgLad = gLad.Inverse();
869   rel.MultiplyLeft(&invgLad);
870   TGeoRotation* rr = new TGeoRotation("rr",90,90,0,0,90,180);
871   TGeoCombiTrans* ct = 0;
872   if(layer==3) ct= new TGeoCombiTrans(25.,0.,0.,rr);
873   if(layer==4) ct= new TGeoCombiTrans(25.+7.5,0.,0.,rr);
874
875   rel.MultiplyLeft(ct);
876   
877   if((layer==3)&&(module<3)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
878   if((layer==3)&&(module>2)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
879   if((layer==4)&&(module<4)) rel.LocalToMaster(fgkLocR[iPoint],fSDDidP[iPoint]);
880   if((layer==4)&&(module>3)) rel.LocalToMaster(fgkLocL[iPoint],fSDDidP[iPoint]);
881
882   for(Int_t i=0; i<3; i++) fSDDidP[iPoint][i]*=10; 
883   fSDDidP[iPoint][2]-=0.5205;
884
885   //  if(ladderPath.Data(),"/ALIC_1/ITSV_1/ITSsddLayer3_1/ITSsddLadd_3");
886   //  if(ladderPath.Contains("ITSsddLayer3_1") && (ladderPath.Contains("ITSsddLadd_3")|| ladderPath.Contains("ITSsddLadd_10")))
887   //  {
888   ///ALIC_1/ITSV_1/ITSsddLayer4_1/ITSsddLadd_5
889   ///ALIC_1/ITSV_1/ITSsddLayer4_1/ITSsddLadd_16
890   //  gLad.Print();
891   //  printf("%s  : Module# %d  Point# %d\n",ladderPath.Data(), module, iPoint);
892   //  printf("ID   {%f, %f, %f}\n", fSDDidP[iPoint][0],fSDDidP[iPoint][1],fSDDidP[iPoint][2]);
893   //  printf("Me   {%f, %f, %f}\n", fSDDmeP[iPoint][0],fSDDmeP[iPoint][1],fSDDmeP[iPoint][2]);
894   //  }
895 }
896
897 //______________________________________________________________________
898 void AliITSSurveyToAlign::ReadPointNameSDD(const char str[], Int_t &iLayer, Int_t &iLader, Int_t &iModul, Int_t &iPoint) const
899 {
900   // 
901   //    Utility function used by CreateAlignObjSDD
902   // 
903   iLayer=-1;
904   iLader=-1;
905   iModul=-1;
906   iPoint=-1;
907
908   if(str[7]=='2') iLayer=3;
909   if(str[7]=='3') iLayer=4;
910
911   if(str[15]=='0') iLader=0;
912   if(str[15]=='1') iLader=1;
913   if(str[15]=='2') iLader=2;
914   if(str[15]=='3') iLader=3;
915   if(str[15]=='4') iLader=4;
916   if(str[15]=='5') iLader=5;
917   if(str[15]=='6') iLader=6;
918   if(str[15]=='7') iLader=7;
919   if(str[15]=='8') iLader=8;
920   if(str[15]=='9') iLader=9;
921
922   Int_t ord=0;
923   if(str[16]=='0') {iLader=10*iLader+0; ord=1;}
924   if(str[16]=='1') {iLader=10*iLader+1; ord=1;}
925   if(str[16]=='2') {iLader=10*iLader+2; ord=1;}
926   if(str[16]=='3') {iLader=10*iLader+3; ord=1;}
927   if(str[16]=='4') {iLader=10*iLader+4; ord=1;}
928   if(str[16]=='5') {iLader=10*iLader+5; ord=1;}
929   if(str[16]=='6') {iLader=10*iLader+6; ord=1;}
930   if(str[16]=='7') {iLader=10*iLader+7; ord=1;}
931   if(str[16]=='8') {iLader=10*iLader+8; ord=1;}
932   if(str[16]=='9') {iLader=10*iLader+9; ord=1;}
933
934   /*
935   //tmp solution
936   Int_t module=-1;
937   if(str[23+ord]=='0') module=0;
938   if(str[23+ord]=='1') module=1;
939   if(str[23+ord]=='2') module=2;
940   if(str[23+ord]=='3') module=3;
941   if(str[23+ord]=='4') module=4;
942   if(str[23+ord]=='5') module=5;
943   if(str[23+ord]=='6') module=6;
944   if(str[23+ord]=='7') module=7;
945   if(str[23+ord]=='8') module=8;
946   if(str[23+ord]=='9') module=9;
947
948   if(iLayer==3)
949     {
950       if(module==0) iModul= 5;
951       if(module==1) iModul= 4;
952       if(module==2) iModul= 3;
953       if(module==3) iModul= 2;
954       if(module==4) iModul= 1;
955       if(module==5) iModul= 0;
956     }
957
958
959   if(iLayer==4)
960     {
961       if(module==0) iModul= 7;
962       if(module==1) iModul= 6;
963       if(module==2) iModul= 5;
964       if(module==3) iModul= 4;
965       if(module==4) iModul= 3;
966       if(module==5) iModul= 2;
967       if(module==6) iModul= 1;
968       if(module==7) iModul= 0;
969     }
970
971   if(module<0) {printf("ERROR MOULE\n"); iModul=0;}
972   */
973
974   if(str[23+ord]=='0') iModul=0;
975   if(str[23+ord]=='1') iModul=1;
976   if(str[23+ord]=='2') iModul=2;
977   if(str[23+ord]=='3') iModul=3;
978   if(str[23+ord]=='4') iModul=4;
979   if(str[23+ord]=='5') iModul=5;
980   if(str[23+ord]=='6') iModul=6;
981   if(str[23+ord]=='7') iModul=7;
982   if(str[23+ord]=='8') iModul=8;
983   if(str[23+ord]=='9') iModul=9;
984
985
986   if((str[25+ord]=='R')&&(str[26+ord]=='D')) iPoint=0;
987   if((str[25+ord]=='R')&&(str[26+ord]=='C')) iPoint=1;
988   if((str[25+ord]=='R')&&(str[26+ord]=='U')) iPoint=2;
989   if((str[25+ord]=='L')&&(str[26+ord]=='U')) iPoint=3;
990   if((str[25+ord]=='L')&&(str[26+ord]=='C')) iPoint=4;
991   if((str[25+ord]=='L')&&(str[26+ord]=='D')) iPoint=5;
992
993
994   return;
995 }
996
997
998 //______________________________________________________________________
999 void AliITSSurveyToAlign::ConvertToRSofModulesAndRotSDD(Int_t Layer, Int_t Module)
1000 {
1001   // 
1002   //    Utility function used by CreateAlignObjSDD
1003   // 
1004
1005   Double_t ymId;
1006   Double_t zmId;
1007
1008   Double_t ymMe;
1009   Double_t zmMe;
1010   Double_t ymMeE;
1011   Double_t zmMeE;
1012
1013   Double_t x0=fSDDidP[1][0];
1014   Double_t z0=fSDDidP[1][2];//-0.5205;
1015   //  Double_t z0=fSDDidP[1][2]-0.5;
1016   for(Int_t i=0; i<6; i++)
1017     {
1018       //      fSDDidP[i][2]-=0.5205;
1019       //      fSDDidP[i][2]-=0.5;
1020
1021       if(!fSDDisMe[i]) continue; 
1022
1023       fSDDidP[i][0]-=x0;
1024       fSDDidP[i][2]-=z0;
1025       fSDDmeP[i][0]-=x0;
1026       fSDDmeP[i][2]-=z0;
1027                                 
1028       ymId=fSDDidP[i][1];
1029       zmId=fSDDidP[i][2];
1030                         
1031       fSDDidP[i][2]=fSDDidP[i][0];
1032       fSDDidP[i][0]=ymId;
1033       fSDDidP[i][1]=zmId;
1034                         
1035       ymMe=fSDDmeP[i][1];
1036       zmMe=fSDDmeP[i][2];
1037                         
1038       ymMeE=fSDDmeP[i][4];
1039       zmMeE=fSDDmeP[i][5];
1040                         
1041       fSDDmeP[i][2]=fSDDmeP[i][0];
1042       fSDDmeP[i][0]=ymMe;
1043       fSDDmeP[i][1]=zmMe;
1044       fSDDmeP[i][5]=fSDDmeP[i][3];
1045       fSDDmeP[i][3]=ymMeE;
1046       fSDDmeP[i][4]=zmMeE;
1047                         
1048
1049       if(((Layer==3)&&(Module>2))||((Layer==4)&&(Module>3)))
1050         {
1051           fSDDidP[i][0]*=(-1);
1052           fSDDidP[i][2]*=(-1);
1053           fSDDmeP[i][0]*=(-1);
1054           fSDDmeP[i][2]*=(-1);
1055         }
1056     }   
1057 }
1058
1059
1060 //______________________________________________________________________
1061 void AliITSSurveyToAlign::CalcShiftSDD(Double_t &x0,Double_t &y0,Double_t &z0) const
1062 {
1063     // Calculates the 3 shifts for the present SDD module
1064     // and sets the three reference arguments
1065     //
1066   Double_t xId, yId, zId;
1067   Double_t xMe, yMe, zMe, sX2, sY2, sZ2;
1068   Double_t aX=0., bX=0.;
1069   Double_t aY=0., bY=0.;
1070   Double_t aZ=0., bZ=0.;
1071   for(Int_t iP1=0; iP1<6; iP1++)
1072     {
1073       if(!fSDDisMe[iP1]) continue;
1074       xId=fSDDidP[iP1][0];
1075       yId=fSDDidP[iP1][1];
1076       zId=fSDDidP[iP1][2];
1077       xMe=fSDDmeP[iP1][0];
1078       yMe=fSDDmeP[iP1][1];
1079       zMe=fSDDmeP[iP1][2];
1080       sX2 =fSDDmeP[iP1][3]*fSDDmeP[iP1][3];
1081       sY2 =fSDDmeP[iP1][4]*fSDDmeP[iP1][4];
1082       sZ2 =fSDDmeP[iP1][5]*fSDDmeP[iP1][5];
1083       aX+=(1./sX2);
1084       bX+=((xMe-xId)/sX2); 
1085       aY+=(1./sY2);
1086       bY+=((yMe-yId)/sY2); 
1087       aZ+=(1./sZ2);
1088       bZ+=((zMe-zId)/sZ2); 
1089     }
1090   Double_t x1 = bX/aX;
1091   Double_t x2 = bY/aY;
1092   Double_t x3 = bZ/aZ;
1093   x0=x1;
1094   y0=x2;
1095   z0=x3;
1096   return;
1097 }
1098
1099
1100 //______________________________________________________________________
1101 void AliITSSurveyToAlign::CalcShiftRotSDD(Double_t &tet,Double_t &psi,Double_t &phi,Double_t &x0,Double_t &y0,Double_t &z0)
1102 {
1103     // Calculates the 3 shifts and 3 euler angles for the present SDD module
1104     // and sets the six reference arguments
1105     //
1106   TMatrixD pC(6,6);
1107
1108   Double_t a[6][6];
1109   for(Int_t ii=0; ii<6; ii++){
1110       for(Int_t jj=0; jj<6; jj++){
1111           a[ii][jj]=0.;
1112       }
1113   }
1114
1115   Double_t c[6];
1116   for(Int_t ii=0; ii<6; ii++)
1117       c[ii]=0.;
1118
1119   Double_t xId, yId, zId;
1120   Double_t xMe, yMe, zMe, sX2, sY2, sZ2;
1121
1122   //  printf("\n");
1123
1124   for(Int_t iP1=0; iP1<=6; iP1++)
1125     {
1126       if(!fSDDisMe[iP1]) continue;
1127
1128       //ideal x,y,z for fiducial mark iP1
1129       xId= fSDDidP[iP1][0];
1130       yId= fSDDidP[iP1][1];
1131       zId= fSDDidP[iP1][2];
1132
1133       //measured x,y,z for fiducial mark iP1
1134       xMe= fSDDmeP[iP1][0];
1135       yMe= fSDDmeP[iP1][1];
1136       zMe= fSDDmeP[iP1][2];
1137
1138
1139
1140       //squared precisions of measured x,y,z for fiducial mark iP1
1141       sX2 = fSDDmeP[iP1][3]* fSDDmeP[iP1][3];
1142       sY2 = fSDDmeP[iP1][4]* fSDDmeP[iP1][4];
1143       sZ2 = fSDDmeP[iP1][5]* fSDDmeP[iP1][5];
1144
1145       a[0][0]+=(zId*zId/sX2+xId*xId/sZ2);
1146       a[0][1]-=(zId*yId/sX2);
1147       a[0][2]-=(xId*yId/sZ2);
1148       a[0][3]-=(zId/sX2);
1149       a[0][4] =0.;
1150       a[0][5]+=(xId/sZ2);
1151       c[0]+=(xId*(zMe-zId)/sZ2-zId*(xMe-xId)/sX2); 
1152
1153       a[1][0]-=(yId*zId/sX2);
1154       a[1][1]+=(xId*xId/sY2+yId*yId/sX2);
1155       a[1][2]-=(xId*zId/sY2);
1156       a[1][3]+=(yId/sX2);
1157       a[1][4]-=(xId/sY2);
1158       a[1][5] =0.;
1159       c[1]+=(yId*(xMe-xId)/sX2-xId*(yMe-yId)/sY2); 
1160
1161       a[2][0]-=(yId*xId/sZ2);
1162       a[2][1]-=(xId*zId/sY2);
1163       a[2][2]+=(zId*zId/sY2+yId*yId/sZ2);
1164       a[2][3] =0.;
1165       a[2][4]+=(zId/sY2);
1166       a[2][5]-=(yId/sZ2);
1167       c[2]+=(zId*(yMe-yId)/sY2-yId*(zMe-zId)/sZ2); 
1168
1169       a[3][0]-=(zId/sX2);
1170       a[3][1]+=(yId/sX2);
1171       a[3][2] =0.;
1172       a[3][3]+=(1./sX2);
1173       a[3][4] =0.;
1174       a[3][5] =0.;
1175       c[3]+=((xMe-xId)/sX2); 
1176
1177       a[4][0] =0.;
1178       a[4][1]-=(xId/sY2);
1179       a[4][2]+=(zId/sY2);
1180       a[4][3] =0.;
1181       a[4][4]+=(1./sY2);
1182       a[4][5] =0.;
1183       c[4]+=((yMe-yId)/sY2); 
1184
1185       a[5][0]+=(xId/sZ2);
1186       a[5][1] =0.;
1187       a[5][2]-=(yId/sZ2);
1188       a[5][3] =0.;
1189       a[5][4] =0.;
1190       a[5][5]+=(1./sZ2);
1191       c[5]+=((zMe-zId)/sZ2); 
1192     }
1193
1194   ///////////////////////////////////////////////////////////////
1195
1196   pC.SetMatrixArray(&(a[0][0]));
1197   TMatrixD p1(pC);
1198   TMatrixD p2(pC);
1199   TMatrixD p3(pC);
1200   TMatrixD p4(pC);
1201   TMatrixD p5(pC);
1202   TMatrixD p6(pC);
1203
1204   for(Int_t raw=0; raw<6; raw++)
1205       p1[raw][0]=c[raw];
1206   for(Int_t raw=0; raw<6; raw++)
1207       p2[raw][1]=c[raw];
1208   for(Int_t raw=0; raw<6; raw++)
1209       p3[raw][2]=c[raw];
1210   for(Int_t raw=0; raw<6; raw++)
1211       p4[raw][3]=c[raw];
1212   for(Int_t raw=0; raw<6; raw++)
1213       p5[raw][4]=c[raw];
1214   for(Int_t raw=0; raw<6; raw++)
1215       p6[raw][5]=c[raw];
1216
1217   // cout << "calculating determinants" << endl;
1218   Double_t det0=pC.Determinant();
1219   Double_t x1 = p1.Determinant()/det0;
1220   Double_t x2 = p2.Determinant()/det0;
1221   Double_t x3 = p3.Determinant()/det0;
1222   Double_t x4 = p4.Determinant()/det0;
1223   Double_t x5 = p5.Determinant()/det0;
1224   Double_t x6 = p6.Determinant()/det0;
1225   //cout << "calculating determinants done" << endl;
1226   if (TMath::Abs(x1) < 1.e-10) {
1227     AliInfo("p1 singular ");
1228     p1.Print();
1229   }
1230   if (TMath::Abs(x2) < 1.e-10) {
1231     AliInfo("p2 singular ");
1232     p2.Print();
1233   }
1234   if (TMath::Abs(x3) < 1.e-10) {
1235     AliInfo("p3 singular ");
1236     p3.Print();
1237   }
1238   if (TMath::Abs(x4) < 1.e-10) {
1239     AliInfo("p4 singular ");
1240     p4.Print();
1241   }
1242   if (TMath::Abs(x5) < 1.e-10) {
1243     AliInfo("p5 singular ");
1244     p5.Print();
1245   }
1246   if (TMath::Abs(x6) < 1.e-10) {
1247     AliInfo("p6 singular ");
1248     p6.Print();
1249   }
1250
1251
1252   tet=x1;
1253   psi=x2;
1254   phi=x3;
1255   x0=x4;
1256   y0=x5;
1257   z0=x6;
1258   return;
1259 }