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