]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFAlignment.cxx
Correct character constantness
[u/mrichter/AliRoot.git] / TOF / AliTOFAlignment.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14 ***************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.19  2007/10/02 09:46:08  arcelli
19 add methods to retrieve real survey data, and make some analysis (by B. Guerzoni)
20
21 Revision 1.17  2007/06/06 16:26:46  arcelli
22 remove fall-back call to local CDB storage
23
24 Revision 1.16  2007/05/15 16:25:44  cvetan
25 Moving the alignment-related static methods from AliAlignObj to the new geometry steering class AliGeomManager (macro from Raffaele)
26
27 Revision 1.15  2007/05/03 09:25:10  decaro
28 Coding convention: RN13 violation -> suppression
29
30 Revision 1.14  2007/04/18 14:49:54  arcelli
31 Some code cleanup, added more debug info
32
33 Revision 1.13  2007/04/17 16:38:36  arcelli
34 Include Methods to derive TOF AlignObjs from Survey Data
35
36 Revision 1.12  2007/02/28 18:09:23  arcelli
37 Add protection against failed retrieval of the CDB cal object
38
39 Revision 1.11  2006/09/19 14:31:26  cvetan
40 Bugfixes and clean-up of alignment object classes. Introduction of so called symbolic names used to identify the alignable volumes (Raffaele and Cvetan)
41
42 Revision 1.10  2006/08/22 13:26:05  arcelli
43 removal of effective c++ warnings (C.Zampolli)
44
45 Revision 1.9  2006/08/10 14:46:54  decaro
46 TOF raw data format: updated version
47
48 Revision 1.8  2006/05/04 19:41:42  hristov
49 Possibility for partial TOF geometry (S.Arcelli)
50
51 Revision 1.7  2006/04/27 13:13:29  hristov
52 Moving the destructor to the implementation file
53
54 Revision 1.6  2006/04/20 22:30:49  hristov
55 Coding conventions (Annalisa)
56
57 Revision 1.5  2006/04/16 22:29:05  hristov
58 Coding conventions (Annalisa)
59
60 Revision 1.4  2006/04/05 08:35:38  hristov
61 Coding conventions (S.Arcelli, C.Zampolli)
62
63 Revision 1.3  2006/03/31 13:49:07  arcelli
64 Removing some junk printout
65
66 Revision 1.2  2006/03/31 11:26:30  arcelli
67  changing CDB Ids according to standard convention
68
69 Revision 1.1  2006/03/28 14:54:48  arcelli
70 class for TOF alignment
71
72 author: Silvia Arcelli, arcelli@bo.infn.it
73 */  
74
75 /////////////////////////////////////////////////////////
76 //                                                     //
77 //            Class for alignment procedure            //
78 //                                                     //
79 //                                                     //
80 //                                                     //
81 /////////////////////////////////////////////////////////
82
83 #include <Rtypes.h>
84
85 #include "TGeoMatrix.h"
86 #include "TMath.h"
87 #include "TFile.h"
88 #include "TRandom.h"
89 #include "TGeoManager.h"
90 #include "TGeoVolume.h"
91 #include "TGeoBBox.h"
92 #include "TGeoTrd1.h"
93 #include "TGeoPhysicalNode.h"
94 #include "TGeoNode.h"
95 #include "TObjString.h"
96
97 #include "AliLog.h"
98 //#include "AliAlignObj.h"
99 #include "AliAlignObjParams.h"
100 #include "AliAlignObjMatrix.h"
101 #include "AliCDBManager.h"
102 #include "AliCDBMetaData.h"
103 #include "AliCDBId.h"
104 #include "AliCDBEntry.h"
105 #include "AliTOFAlignment.h"
106 #include "AliSurveyObj.h"
107 #include "AliSurveyPoint.h"
108
109 ClassImp(AliTOFAlignment)
110
111 const Double_t AliTOFAlignment::fgkRorigTOF  = 384.5; // Mean Radius of the TOF ext. volume, cm
112 const Double_t AliTOFAlignment::fgkX1BTOF = 124.5;    //x1 size of BTOF
113 const Double_t AliTOFAlignment::fgkX2BTOF = 134.7262; //x2 size of BTOF
114 const Double_t AliTOFAlignment::fgkYBTOF = 747.2;     //y size of BTOF
115 const Double_t AliTOFAlignment::fgkZBTOF = 29.0;      //z size of BTOF
116 const Double_t AliTOFAlignment::fgkXFM = 38.0;     //x pos of FM in BTOF, cm 
117 const Double_t AliTOFAlignment::fgkYFM = 457.3;    //y pos of FM in BTOF, cm
118 const Double_t AliTOFAlignment::fgkZFM = 11.2;     //z pos of FM in BTOF, cm
119
120 //_____________________________________________________________________________
121 AliTOFAlignment::AliTOFAlignment():
122   TTask("AliTOFAlignment",""),
123   fNTOFAlignObj(0),
124   fTOFmgr(0x0),
125   fTOFAlignObjArray(0x0)
126  { 
127    //AliTOFalignment main Ctor
128    for(Int_t i=0; i<18;i++)
129      for(Int_t j=0; j<5; j++)
130        fNFMforSM[i][j]=0;
131    for(Int_t i=0; i<72; i++)
132     for (Int_t j=0; j<6; j++)
133       fCombFMData[i][j]=0;
134 }
135 //_____________________________________________________________________________
136 AliTOFAlignment::AliTOFAlignment(const AliTOFAlignment &t):
137   TTask("AliTOFAlignment",""),
138   fNTOFAlignObj(0),
139   fTOFmgr(0x0),
140   fTOFAlignObjArray(0x0)
141
142   //AliTOFAlignment copy Ctor
143
144   fNTOFAlignObj=t.fNTOFAlignObj;
145   fTOFAlignObjArray=t.fTOFAlignObjArray;
146   //AliTOFalignment main Ctor
147   for(Int_t i=0; i<18;i++)
148      for(Int_t j=0; j<5; j++)
149        fNFMforSM[i][j]=t.fNFMforSM[i][j];
150   for(Int_t i=0; i<72; i++)
151     for (Int_t j=0; j<6; j++)
152       fCombFMData[i][j]=t.fCombFMData[i][j]; 
153 }
154 //_____________________________________________________________________________
155 AliTOFAlignment& AliTOFAlignment::operator=(const AliTOFAlignment &t){ 
156   //AliTOFAlignment assignment operator
157
158   this->fNTOFAlignObj=t.fNTOFAlignObj;
159   this->fTOFmgr=t.fTOFmgr;
160   this->fTOFAlignObjArray=t.fTOFAlignObjArray;
161   return *this;
162
163 }
164 //_____________________________________________________________________________
165 AliTOFAlignment::~AliTOFAlignment() {
166   delete fTOFAlignObjArray;
167   delete fTOFmgr;
168 }
169
170 //_____________________________________________________________________________
171 void AliTOFAlignment::Smear( Float_t *tr, Float_t *rot)
172 {
173   //Introduce Random Offset/Tilts
174   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
175   Float_t dx, dy, dz;  // shifts
176   Float_t dpsi, dtheta, dphi; // angular displacements
177   TRandom *rnd   = new TRandom(1567);
178  
179   Int_t nSMTOF = 18;
180   AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
181   UShort_t iIndex=0; //dummy volume index
182   //  AliGeomManager::ELayerID iLayer = AliGeomManager::kTOF;
183   //  Int_t iIndex=1; //dummy volume index
184   UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity 
185   Int_t i;
186   for (i = 0; i<nSMTOF ; i++) {
187     Char_t  path[100];
188     sprintf(path,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
189
190     dx = (rnd->Gaus(0.,1.))*tr[0];
191     dy = (rnd->Gaus(0.,1.))*tr[1];
192     dz = (rnd->Gaus(0.,1.))*tr[2];
193     dpsi   = rot[0];
194     dtheta = rot[1];
195     dphi   = rot[2];
196     AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
197     fTOFAlignObjArray->Add(o);
198   }
199
200   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
201   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
202   delete rnd;
203 }
204
205 //_____________________________________________________________________________
206 void AliTOFAlignment::Align( Float_t *tr, Float_t *rot)
207 {
208   //Introduce Offset/Tilts
209
210   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
211   Float_t dx, dy, dz;  // shifts
212   Float_t dpsi, dtheta, dphi; // angular displacements
213
214
215   Int_t nSMTOF = 18;
216   AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
217   UShort_t iIndex=0; //dummy volume index
218   UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity 
219   Int_t i;
220   for (i = 0; i<nSMTOF ; i++) {
221
222     Char_t  path[100];
223     sprintf(path,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
224     dx = tr[0];
225     dy = tr[1];
226     dz = tr[2];
227     dpsi   = rot[0];
228     dtheta = rot[1];
229     dphi   = rot[2];
230     
231     AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
232     fTOFAlignObjArray->Add(o);
233   }
234   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
235   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
236 }
237 //_____________________________________________________________________________
238 void AliTOFAlignment::WriteParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
239 {
240   //Write Align Par on CDB
241   AliCDBManager *man = AliCDBManager::Instance();
242   const Char_t *sel1 = "AlignPar" ;
243   Char_t  out[100];
244   sprintf(out,"%s/%s",sel,sel1); 
245   AliCDBId idTOFAlign(out,minrun,maxrun);
246   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
247   mdTOFAlign->SetResponsible("TOF");
248   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
249   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
250 }
251 //_____________________________________________________________________________
252 void AliTOFAlignment::ReadParFromCDB(const Char_t *sel, Int_t nrun)
253 {
254   //Read Align Par from CDB
255   AliCDBManager *man = AliCDBManager::Instance();
256   const Char_t *sel1 = "AlignPar" ;
257   Char_t  out[100];
258   sprintf(out,"%s/%s",sel,sel1); 
259   AliCDBEntry *entry = man->Get(out,nrun);
260   if (!entry) { 
261     AliError(Form("Failed to get entry: %s",out));
262     return; 
263   }
264   fTOFAlignObjArray=(TObjArray*)entry->GetObject();
265   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
266   AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
267
268 }
269 //_____________________________________________________________________________
270 void AliTOFAlignment::WriteSimParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
271 {
272   //Write Sim Align Par on CDB
273   AliCDBManager *man = AliCDBManager::Instance();
274   const Char_t *sel1 = "AlignSimPar" ;
275   Char_t  out[100];
276   sprintf(out,"%s/%s",sel,sel1); 
277   AliCDBId idTOFAlign(out,minrun,maxrun);
278   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
279   mdTOFAlign->SetResponsible("TOF");
280   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
281   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
282 }
283 //_____________________________________________________________________________
284 void AliTOFAlignment::ReadSimParFromCDB(const Char_t *sel, Int_t nrun){
285   //Read Sim Align Par from CDB
286   AliCDBManager *man = AliCDBManager::Instance();
287   const Char_t *sel1 = "AlignSimPar" ;
288   Char_t  out[100];
289   sprintf(out,"%s/%s",sel,sel1); 
290   AliCDBEntry *entry = man->Get(out,nrun);
291   fTOFAlignObjArray=(TObjArray*)entry->GetObject();
292   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
293   AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
294
295 }
296 //_____________________________________________________________________________
297 void AliTOFAlignment::WriteOnCDBforDC()
298 {
299   //Write Align Par on CDB for DC06
300   AliCDBManager *man = AliCDBManager::Instance();
301   AliCDBId idTOFAlign("TOF/Align/Data",0,0);
302   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
303   mdTOFAlign->SetComment("Alignment objects for ideal geometry, i.e. applying them to TGeo has to leave geometry unchanged");
304   mdTOFAlign->SetResponsible("TOF");
305   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
306   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
307 }
308 //_____________________________________________________________________________
309 void AliTOFAlignment::ReadFromCDBforDC()
310 {
311   //Read Sim Align Par from CDB for DC06
312   AliCDBManager *man = AliCDBManager::Instance();
313   AliCDBEntry *entry = man->Get("TOF/Align/Data",0);
314   fTOFAlignObjArray=(TObjArray*)entry->GetObject();
315   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
316   AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
317
318 }
319
320 //_____________________________________________________________________________
321 void AliTOFAlignment::BuildGeomForSurvey()
322 {
323
324   //Generates the ideal TOF structure with four Fiducial Marks in each 
325   //supermodule (two on each z side) in their expected position. 
326   //Make BTOF
327
328   fTOFmgr = new TGeoManager("Geom","survey to alignment for TOF");
329   TGeoMedium *medium = 0;
330   TGeoVolume *top = fTOFmgr->MakeBox("TOP",medium,1000,1000,1000);
331   fTOFmgr->SetTopVolume(top);
332   // make shape components:  
333   // This is the BTOF containing the FTOA  
334   TGeoTrd1 *strd1  = new TGeoTrd1(fgkX1BTOF*0.5,fgkX2BTOF*0.5, fgkYBTOF*0.5,fgkZBTOF*0.5);
335   TGeoVolume* trd1[18];
336
337   // Now four fiducial marks on SM, expressed in local coordinates
338   // They are positioned at x=+/- 38 cm, y=+/- 457.3 cm, z=11.2 cm
339   
340   TGeoBBox *fmbox  = new TGeoBBox(1,1,1);
341   TGeoVolume* fm = new TGeoVolume("FM",fmbox);
342   fm->SetLineColor(2);
343   
344
345   TGeoTranslation* mAtr = new TGeoTranslation("mAtr",-fgkXFM, -fgkYFM ,fgkZFM);
346   TGeoTranslation* mBtr = new TGeoTranslation("mBtr",fgkXFM, -fgkYFM ,fgkZFM );
347   TGeoTranslation* mCtr = new TGeoTranslation("mCtr",fgkXFM, fgkYFM ,fgkZFM );
348   TGeoTranslation* mDtr = new TGeoTranslation("mDtr",-fgkXFM, fgkYFM ,fgkZFM );
349
350   // position all this stuff in the global ALICE frame
351
352   char name[16];
353   Double_t smX = 0.;
354   Double_t smY = 0.;
355   Double_t smZ = 0.;
356   Float_t  smR = fgkRorigTOF;
357   for (Int_t iSM = 0; iSM < 18; iSM++) {
358     Int_t mod = iSM + 13;
359     if (mod > 17) mod -= 18;
360     sprintf(name, "BTOF%d",mod);
361     trd1[iSM] = new TGeoVolume(name,strd1);
362     Float_t phi  = iSM * 20.;
363     Float_t phi2 = 270 + phi;
364     if (phi2 >= 360.) phi2 -= 360.;
365     smX =  TMath::Sin(phi*TMath::Pi()/180.)*smR;
366     smY = -TMath::Cos(phi*TMath::Pi()/180.)*smR;
367     smZ = 0.;  
368     TGeoRotation* bTOFRot = new TGeoRotation("bTOFRot",phi,90,0.);
369     TGeoCombiTrans trans = *(new TGeoCombiTrans(smX,smY,smZ, bTOFRot));
370     TGeoMatrix* id = new TGeoHMatrix();
371     TGeoHMatrix  transMat = *id * trans;
372     TGeoHMatrix  *smTrans = new TGeoHMatrix(transMat);
373     
374     trd1[iSM]->AddNode(fm,1,mAtr);        //place FM in BTOF
375     trd1[iSM]->AddNode(fm,2,mBtr);
376     trd1[iSM]->AddNode(fm,3,mCtr);
377     trd1[iSM]->AddNode(fm,4,mDtr);
378     top->AddNode(trd1[iSM],1,smTrans);    //place BTOF_iSM in ALICE
379     trd1[iSM]->SetVisDaughters();
380     trd1[iSM]->SetLineColor(iSM);         //black
381     
382   }  
383
384   fTOFmgr->CloseGeometry();
385   fTOFmgr->GetTopVolume()->Draw();
386   fTOFmgr->SetVisOption(0);
387   fTOFmgr->SetVisLevel(6);
388
389   // Now Store the "Ideal"  Global Matrices (local to global) for later use
390   
391   for (Int_t iSM = 0; iSM < 18; iSM++) {
392
393     sprintf(name, "TOP_1/BTOF%d_1", iSM);
394     printf("\n\n*****************  TOF SuperModule:  %s ****************** \n",name);
395     TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
396     fTOFMatrixId[iSM] = pn3->GetMatrix(); //save "ideal" global matrix
397     printf("\n\n***************  The Ideal Matrix in GRS *****************\n");
398     fTOFMatrixId[iSM]->Print();
399
400   }
401 }
402
403 //_____________________________________________________________________________
404 void AliTOFAlignment::InsertMisAlignment(Float_t *mis)
405 {
406   // Now Apply the Displacements and store the misaligned FM positions...
407   //
408   //
409
410   Double_t lA[3]={-fgkXFM, -fgkYFM ,fgkZFM};
411   Double_t lB[3]={fgkXFM, -fgkYFM ,fgkZFM};
412   Double_t lC[3]={fgkXFM, fgkYFM ,fgkZFM};
413   Double_t lD[3]={-fgkXFM, fgkYFM ,fgkZFM};
414
415   for(Int_t iSM=0;iSM<18;iSM++){
416      char name[16];
417      sprintf(name, "TOP_1/BTOF%d_1", iSM);
418      fTOFmgr->cd(name);
419      printf("\n\n******Misaligning TOF SuperModule ************** %s \n",name);
420
421     // ************* get ideal global matrix *******************
422     TGeoHMatrix g3 = *fTOFmgr->GetCurrentMatrix(); 
423     AliInfo(Form("This is the ideal global trasformation of SM %i",iSM));
424     g3.Print(); // g3 is the local(BTOF) to global (ALICE) matrix and is the same of fTOFMatrixId
425     TGeoNode* n3 = fTOFmgr->GetCurrentNode(); 
426     TGeoMatrix* l3 = n3->GetMatrix(); 
427     
428     Double_t gA[3], gB[3], gC[3], gD[3]; // ideal global FM point coord.
429     g3.LocalToMaster(lA,gA);
430     g3.LocalToMaster(lB,gB);
431     g3.LocalToMaster(lC,gC);
432     g3.LocalToMaster(lD,gD);
433    
434     //  We apply a delta transformation to the surveyed vol to represent
435     //  its real position, given below by ng3 nl3, which differs from its
436     //  ideal position saved above in g3 and l3
437
438     //we have to express the displacements as regards the old local RS (non misaligned BTOF)
439     Double_t dx     = mis[0]; // shift along x 
440     Double_t dy     = mis[1]; // shift along y 
441     Double_t dz     = mis[2]; // shift along z 
442     Double_t dphi   = mis[3]; // rot around z 
443     Double_t dtheta = mis[4]; // rot around x' 
444     Double_t dpsi   = mis[5]; // rot around z''
445
446     TGeoRotation* rrot = new TGeoRotation("rot",dphi,dtheta,dpsi);
447     TGeoCombiTrans localdelta = *(new TGeoCombiTrans(dx,dy,dz, rrot));
448     AliInfo(Form("This is the local delta trasformation for SM %i \n",iSM));
449     localdelta.Print();
450     TGeoHMatrix nlocal = *l3 * localdelta;
451     TGeoHMatrix* nl3 = new TGeoHMatrix(nlocal); // new matrix, representing real position (from new local mis RS to the global one)
452    
453     TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
454
455     pn3->Align(nl3);   
456     
457     TGeoHMatrix* ng3 = pn3->GetMatrix(); //"real" global matrix, what survey sees 
458     printf("\n\n*************  The Misaligned Matrix in GRS **************\n");
459     ng3->Print();
460     Double_t ngA[3], ngB[3], ngC[3], ngD[3];// real FM point coord., global RS 
461     ng3->LocalToMaster(lA,ngA);
462     ng3->LocalToMaster(lB,ngB);
463     ng3->LocalToMaster(lC,ngC);
464     ng3->LocalToMaster(lD,ngD);    
465
466     for(Int_t coord=0;coord<3;coord++){
467       fCombFMData[iSM*4][2*coord]=ngA[coord];
468       fCombFMData[iSM*4][2*coord+1]=1;
469       fCombFMData[iSM*4+1][2*coord]=ngB[coord];
470       fCombFMData[iSM*4+1][2*coord+1]=1;
471       fCombFMData[iSM*4+2][2*coord]=ngC[coord];
472       fCombFMData[iSM*4+2][2*coord+1]=1;
473       fCombFMData[iSM*4+3][2*coord]=ngD[coord];
474       fCombFMData[iSM*4+3][2*coord+1]=1;
475       }
476     }
477
478 }
479
480 //____________________________________________________________________________
481 void AliTOFAlignment::WriteCombData(const Char_t *nomefile, Int_t option)
482 {
483   // 1 for simulated data; 0 for data from survey file
484   // write combined data on a file
485   //
486
487   FILE *data;
488   /* Open file in text mode: */
489   if( (data = fopen( nomefile, "w+t" )) != NULL ){
490     if (option==1){
491       fprintf( data, "simulated data\n" );} else {
492         fprintf( data, "survey data\n" );}
493     if (option==1){
494       fprintf( data, "data from InsertMisAlignmentBTOF method\n");}
495     else {fprintf( data, "real survey data from text file (coordinate in global RS)\n");}
496     fprintf( data, "Point Name,XPH,YPH,ZPH,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
497     fprintf( data, "> Data:\n");
498     for(Int_t i=0;i<72;i++){
499       if (fCombFMData[i][0]!=0){
500         fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1]*10,fCombFMData[i][3]*10,fCombFMData[i][5]*10); 
501       }
502     }
503     fclose( data );
504    }
505   else{
506     printf(  "Problem opening the file\n" );
507   }
508   
509   return;  
510 }
511
512 //____________________________________________________________________________
513 void AliTOFAlignment::WriteSimSurveyData(const Char_t *nomefile)
514 {
515   // write sim data in standard format
516   //
517   //
518
519   FILE *data;
520   /* Open file in text mode: */
521   if( (data = fopen( nomefile, "w+t" )) != NULL )
522    {
523       fprintf( data, "> Title:\n" );
524       fprintf( data, "simulated data\n" );
525       fprintf( data, "> Date:\n" );
526       fprintf( data, "24.09.2007\n" );
527       fprintf( data, "> Subdetector:\n" );
528       fprintf( data, "TOF\n" );
529       fprintf( data, "> Report URL:\n" );
530       fprintf( data, "https://edms.cern.ch/document/835615\n" );
531       fprintf( data, "> Version:\n" );
532       fprintf( data, "1\n");
533       fprintf( data, "> General Observations:\n"); 
534       fprintf( data, "data from InsertMisAlignmentBTOF method\n");
535       fprintf( data, "> Coordinate System:\n");
536       fprintf( data, "\\ALICEPH\n");
537       fprintf( data, "> Units:\n");
538       fprintf( data, "cm\n");
539       fprintf( data, "> Nr Columns:\n");
540       fprintf( data, "9\n");
541       fprintf( data, "> Column Names:\n");
542       fprintf( data, "Point Name,XPH,YPH,ZPH,Point Type,Target Used,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
543       fprintf( data, "> Data:\n");
544       for(Int_t i=0;i<72;i++)
545         if (fCombFMData[i][0]!=0)
546           fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1],fCombFMData[i][3],fCombFMData[i][5]); 
547       
548        fclose( data );
549    }
550    else
551      printf(  "Problem opening the file\n" );
552 }
553
554 //____________________________________________________________________________
555 void AliTOFAlignment::MakeDefData(const Int_t nf,TString namefiles[])
556 {
557   //this method combines survey data from different files (namefiles[]) 
558   //
559   // 
560  
561   Float_t data[72][6][100];
562   for (Int_t i=0;i<72;i++)
563     for (Int_t j=0; j<6; j++)
564       for(Int_t k=0; k<100; k++)
565         data[i][j][k]=0;
566   Int_t nfm=0;
567   Int_t nsm=0;
568   Long64_t totdata[72]={0};
569
570   for (Int_t ii=0;ii<nf; ii++)
571     {
572       AliSurveyObj *so = new AliSurveyObj();
573       const Char_t *nome=namefiles[ii];
574       so->FillFromLocalFile(nome);
575       TObjArray *points = so->GetData();
576       Int_t nSurveyPoint=points->GetEntries();
577       for(Int_t jj=0;jj<nSurveyPoint;jj++){
578         const char* pointName= ((AliSurveyPoint *) points->At(jj))->GetPointName().Data();
579         nfm=atoi(&pointName[6]);
580         nsm=atoi(&pointName[2]);
581         data[nsm*4+nfm][0][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetX();
582         data[nsm*4+nfm][2][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetY();
583         data[nsm*4+nfm][4][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetZ();
584         data[nsm*4+nfm][1][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionX();
585         data[nsm*4+nfm][3][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionY();
586         data[nsm*4+nfm][5][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionZ();
587         totdata[nsm*4+nfm]=totdata[nsm*4+nfm]+1;
588       } 
589       delete so;
590     }
591
592   
593   for(Int_t i=0; i<72 ;i++){
594     Float_t numx=0, numy=0,numz=0, comodox=0, comodoy=0, comodoz=0,denx=0, deny=0, denz=0;
595     if(totdata[i]!=0){    
596       for(Int_t j=0; j<totdata[i]; j++){
597         comodox=1/(data[i][1][j]/10*data[i][1][j]/10);//precision in mm, position in cm
598         numx=numx+data[i][0][j]*comodox;
599         denx=denx+comodox;
600         comodoy=1/(data[i][3][j]/10*data[i][3][j]/10);
601         numy=numy+data[i][2][j]*comodoy;
602         deny=deny+comodoy;
603         comodoz=1/(data[i][5][j]/10*data[i][5][j]/10);
604         numz=numz+data[i][4][j]*comodoz;
605         denz=denz+comodoz;
606         }
607       fCombFMData[i][1]=TMath::Sqrt(1/denx); //error for x position
608       fCombFMData[i][3]=TMath::Sqrt(1/deny); //error for y position
609       fCombFMData[i][5]=TMath::Sqrt(1/denz); //error for z position
610       fCombFMData[i][0]=numx/denx;           //combined survey data for x position of FM
611       fCombFMData[i][2]=numy/deny;           //combined survey data for y position of FM
612       fCombFMData[i][4]=numz/denz;           //combined survey data for z position of FM
613       } else continue;
614     }
615
616   for(Int_t i=0;i<72;i++)
617     if (fCombFMData[i][0]!=0){
618       fNFMforSM[(i-i%4)/4][i%4]=1;
619       fNFMforSM[(i-i%4)/4][4]=fNFMforSM[(i-i%4)/4][4]+1;
620     }
621 }
622
623 //_____________________________________________________________________________
624 void AliTOFAlignment::ReadSurveyDataAndAlign(){
625   //
626   // read the survey data and, if we know the positions of at least 3 FM 
627   //for a SM, call the right Alignement procedure  
628
629   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
630
631   Float_t deltaFM0=0, deltaFM1=0, deltaFM2=0, deltaFM3=0;
632
633   for(Int_t i=0; i<18; i++){
634     switch(fNFMforSM[i][4]){
635     case 0:
636       printf("we don't know the position of any FM of SM %i\n",i);
637       break;
638     case 1:
639       printf("we know the position of only one FM for SM %i\n",i);
640      
641       break;
642     case 2:
643       printf("we know the position of only 2 FM for SM %i\n",i);
644       
645       break;
646     case 3:
647       if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1){
648         printf("we know the position of FM A B C for SM %i\n",i);
649         AliTOFAlignment::AlignFromSurveyABC(i);};
650
651         
652       if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][3]==1){
653         printf("we know the position of FM A B D for SM %i\n",i);
654         AliTOFAlignment::AlignFromSurveyABD(i);};
655
656         
657       if (fNFMforSM[i][0]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
658         printf("we know the position of FM A C D for SM %i\n",i);
659         AliTOFAlignment::AlignFromSurveyACD(i);};
660
661         
662       if (fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
663         printf("we know the position of FM B C D for SM %i\n",i);
664         AliTOFAlignment::AlignFromSurveyBCD(i);};
665
666         
667       break;
668     case 4:
669       printf("we know the position of all the 4 FM for SM %i\n",i);
670       //check the precision of the measurement
671
672       deltaFM0=fCombFMData[i*4][1]/TMath::Abs(fCombFMData[i*4][0])+fCombFMData[i*4][3]/TMath::Abs(fCombFMData[i*4][2])+fCombFMData[i*4][5]/TMath::Abs(fCombFMData[i*4][4]);
673       deltaFM1=fCombFMData[i*4+1][1]/TMath::Abs(fCombFMData[i*4+1][0])+fCombFMData[i*4+1][3]/TMath::Abs(fCombFMData[i*4+1][2])+fCombFMData[i*4+1][5]/TMath::Abs(fCombFMData[i*4+1][4]);
674       deltaFM2=fCombFMData[i*4+2][1]/TMath::Abs(fCombFMData[i*4+2][0])+fCombFMData[i*4+2][3]/TMath::Abs(fCombFMData[i*4+2][2])+fCombFMData[i*4+2][5]/TMath::Abs(fCombFMData[i*4+2][4]);
675       deltaFM3=fCombFMData[i*4+3][1]/TMath::Abs(fCombFMData[i*4+3][0])+fCombFMData[i*4+3][3]/TMath::Abs(fCombFMData[i*4+3][2])+fCombFMData[i*4+3][5]/TMath::Abs(fCombFMData[i*4+3][4]);
676
677       //to AlignFromSurvey we use the 3 FM whose positions are known with greatest precision
678       if(deltaFM0>=deltaFM1 && deltaFM0>=deltaFM2 && deltaFM0>=deltaFM3){
679         printf("to Align we use FM B,C,D");
680         AliTOFAlignment::AlignFromSurveyBCD(i);} else
681           if(deltaFM1>=deltaFM0 && deltaFM1>=deltaFM2 && deltaFM1>=deltaFM3){
682            printf("to Align we use FM A,C,D");
683            AliTOFAlignment::AlignFromSurveyACD(i);} else
684              if(deltaFM2>=deltaFM0 && deltaFM2>=deltaFM1 && deltaFM2>=deltaFM3){
685                printf("to Align we use FM A,B,D");
686                AliTOFAlignment::AlignFromSurveyABD(i);} else{
687                  printf("to Align we use FM A,B,C");
688                  AliTOFAlignment::AlignFromSurveyABC(i);}
689      
690       break;
691     }
692   
693   }
694
695     // saving TOF AligObjs from survey on a file, for the moment.. 
696   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
697   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
698   TFile f("TOFAlignFromSurvey.root","RECREATE");
699   f.cd();
700   f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
701   f.Close();
702   
703
704 }
705
706 //_____________________________________________________________________________
707 void AliTOFAlignment::AlignFromSurveyABC(Int_t iSM)
708 {
709
710   //From Survey data, derive the needed transformations to get the 
711   //Alignment Objects. 
712   //Again, highly "inspired" to Raffaele's example... 
713   //we use FM A,B,C
714     
715     Double_t ngA[3], ngB[3], ngC[3]; // real FM point coord., global RS
716     // Get the 'realistic' input from the Survey Matrix
717       for(Int_t coord=0;coord<3;coord++){
718       ngA[coord]=   fCombFMData[iSM*4][coord*2];
719       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
720       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
721       }
722
723     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
724
725     // From the real fiducial marks coordinates derive back the
726     // new global position of the surveyed volume
727     //*** What follows is the actual survey-to-alignment procedure
728     
729     Double_t ab[3], bc[3], n[3];
730     Double_t plane[4], s=1.;
731     
732     // first vector on the plane of the fiducial marks
733     for(Int_t i=0;i<3;i++){
734       ab[i] = (ngB[i] - ngA[i]);
735     }
736     
737     // second vector on the plane of the fiducial marks
738     for(Int_t i=0;i<3;i++){
739       bc[i] = (ngC[i] - ngB[i]);
740     }
741     
742     // vector normal to the plane of the fiducial marks obtained
743     // as cross product of the two vectors on the plane d0^d1
744     n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
745     n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
746     n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
747     
748     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
749     if(sizen>1.e-8){
750       s = Double_t(1.)/sizen ; //normalization factor
751     }else{
752       AliInfo("Problem in normalizing the vector");
753     }
754     
755     // plane expressed in the hessian normal form, see:
756     // http://mathworld.wolfram.com/HessianNormalForm.html
757     // the first three are the coordinates of the orthonormal vector
758     // the fourth coordinate is equal to the distance from the origin
759   
760     for(Int_t i=0;i<3;i++){
761       plane[i] = n[i] * s;
762     }
763     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
764     
765     // The center of the square with fiducial marks as corners
766     // as the middle point of one diagonal - md
767     // Used below to get the center - orig - of the surveyed box
768
769     Double_t orig[3], md[3];
770     for(Int_t i=0;i<3;i++){
771       md[i] = (ngA[i] + ngC[i]) * 0.5;
772     }
773     
774     // The center of the box, gives the global translation
775     for(Int_t i=0;i<3;i++){
776       orig[i] = md[i] - plane[i]*fgkZFM;
777     }
778     
779     // get local directions needed to write the global rotation matrix
780     // for the surveyed volume by normalising vectors ab and bc
781     Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
782
783
784     if(sx>1.e-8){
785       for(Int_t i=0;i<3;i++){
786         ab[i] /= sx;
787       }
788     }
789     Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
790     if(sy>1.e-8){
791       for(Int_t i=0;i<3;i++){
792         bc[i] /= sy;
793       }
794     }
795     Double_t rot[9] = {ab[0],bc[0],plane[0],ab[1],bc[1],plane[1],ab[2],bc[2],plane[2]}; // the rotation matrix
796     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey
797     TGeoHMatrix ng;              
798     ng.SetTranslation(orig);
799     ng.SetRotation(rot);
800     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
801     ng.Print();    
802
803     // Calculate the delta transformation wrt Ideal geometry
804     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
805     
806     printf("\n\n**** The ideal matrix ***\n"); 
807     fTOFMatrixId[iSM]->Print();   
808     
809     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
810     printf("\n\n**** The inverse of the ideal matrix ***\n");
811     gdelta.Print();
812  
813     gdelta.MultiplyLeft(&ng);
814     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
815     gdelta.Print(); //this is the global delta trasformation
816     
817     // Now Write the Alignment Objects....
818     Int_t index=0; //let all SM modules have index=0
819     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
820     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
821     TString symname(Form("TOF/sm%02d",iSM));
822     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
823     fTOFAlignObjArray->Add(o);
824
825   }
826
827
828 //_____________________________________________________________________________
829 void AliTOFAlignment::AlignFromSurveyABD(Int_t iSM)
830 {
831   
832   //From Survey data, derive the needed transformations to get the 
833   //Alignment Objects. 
834   //Again, highly "inspired" to Raffaele's example... 
835   //we use FM A,B,D
836     
837   Double_t ngA[3], ngB[3], ngD[3];// real FM point coord., global RS
838     
839    // Get the 'realistic' input from the Survey Matrix
840       for(Int_t coord=0;coord<3;coord++){
841       ngA[coord]=   fCombFMData[iSM*4][coord*2];
842       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
843       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
844       }
845
846     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
847
848     // From the new fiducial marks coordinates derive back the
849     // new global position of the surveyed volume
850     //*** What follows is the actual survey-to-alignment procedure
851     
852     Double_t ab[3], ad[3], n[3];
853     Double_t plane[4], s=1.;
854     
855     // first vector on the plane of the fiducial marks
856     for(Int_t i=0;i<3;i++){
857       ab[i] = (ngB[i] - ngA[i]);
858     }
859     
860     // second vector on the plane of the fiducial marks
861     for(Int_t i=0;i<3;i++){
862       ad[i] = (ngD[i] - ngA[i]);
863     }
864     
865     // vector normal to the plane of the fiducial marks obtained
866     // as cross product of the two vectors on the plane d0^d1
867     n[0] = (ab[1] * ad[2] - ab[2] * ad[1]);
868     n[1] = (ab[2] * ad[0] - ab[0] * ad[2]);
869     n[2] = (ab[0] * ad[1] - ab[1] * ad[0]);
870     
871     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
872     if(sizen>1.e-8){
873       s = Double_t(1.)/sizen ; //normalization factor
874     }else{
875       AliInfo("Problem in normalizing the vector");
876     }
877     
878     // plane expressed in the hessian normal form, see:
879     // http://mathworld.wolfram.com/HessianNormalForm.html
880     // the first three are the coordinates of the orthonormal vector
881     // the fourth coordinate is equal to the distance from the origin
882   
883     for(Int_t i=0;i<3;i++){
884       plane[i] = n[i] * s;
885     }
886     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
887     
888     // The center of the square with fiducial marks as corners
889     // as the middle point of one diagonal - md
890     // Used below to get the center - orig - of the surveyed box
891
892     Double_t orig[3], md[3];
893     for(Int_t i=0;i<3;i++){
894       md[i] = (ngB[i] + ngD[i]) * 0.5;
895     }
896     
897     // The center of the box, gives the global translation
898     for(Int_t i=0;i<3;i++){
899       orig[i] = md[i] - plane[i]*fgkZFM;
900     }
901     
902     // get local directions needed to write the global rotation matrix
903     // for the surveyed volume by normalising vectors ab and bc
904     Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
905     if(sx>1.e-8){
906       for(Int_t i=0;i<3;i++){
907         ab[i] /= sx;
908       }
909     }
910     Double_t sy = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
911     if(sy>1.e-8){
912       for(Int_t i=0;i<3;i++){
913         ad[i] /= sy;
914       }
915     }
916     Double_t rot[9] = {ab[0],ad[0],plane[0],ab[1],ad[1],plane[1],ab[2],ad[2],plane[2]};
917     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
918     TGeoHMatrix ng;              
919     ng.SetTranslation(orig);
920     ng.SetRotation(rot);
921     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
922     ng.Print();    
923
924     // Calculate the delta transformation wrt Ideal geometry
925     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
926     
927     printf("\n\n**** The ideal matrix ***\n"); 
928     fTOFMatrixId[iSM]->Print();   
929     
930     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
931     printf("\n\n**** The inverse of the ideal matrix ***\n");
932     gdelta.Print();
933  
934     gdelta.MultiplyLeft(&ng);
935     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
936     gdelta.Print();  //global delta trasformation
937     
938     // Now Write the Alignment Objects....
939     Int_t index=0; //let all SM modules have index=0
940     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
941     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
942     TString symname(Form("TOF/sm%02d",iSM));
943     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
944     fTOFAlignObjArray->Add(o);
945
946   }
947 //_____________________________________________________________________________
948 void AliTOFAlignment::AlignFromSurveyACD(Int_t iSM)
949 {
950   //From Survey data, derive the needed transformations to get the 
951   //Alignment Objects. 
952   //Again, highly "inspired" to Raffaele's example... 
953   //we use FM A,C,D
954   
955     
956     Double_t ngA[3], ngC[3], ngD[3];// real FM point coord., global RS
957     
958    // Get the 'realistic' input from the Survey Matrix
959       for(Int_t coord=0;coord<3;coord++){
960       ngA[coord]=   fCombFMData[iSM*4][coord*2];
961       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
962       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
963       }
964
965     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
966
967     // From the new fiducial marks coordinates derive back the
968     // new global position of the surveyed volume
969     //*** What follows is the actual survey-to-alignment procedure
970     
971     Double_t cd[3], ad[3], n[3];
972     Double_t plane[4], s=1.;
973     
974     // first vector on the plane of the fiducial marks
975     for(Int_t i=0;i<3;i++){
976       cd[i] = (ngC[i] - ngD[i]);
977     }
978     
979     // second vector on the plane of the fiducial marks
980     for(Int_t i=0;i<3;i++){
981       ad[i] = (ngD[i] - ngA[i]);
982     }
983     
984     // vector normal to the plane of the fiducial marks obtained
985     // as cross product of the two vectors on the plane d0^d1
986     n[0] = (ad[1] * cd[2] - ad[2] * cd[1]);
987     n[1] = (ad[2] * cd[0] - ad[0] * cd[2]);
988     n[2] = (ad[0] * cd[1] - ad[1] * cd[0]);
989     
990     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
991     if(sizen>1.e-8){
992       s = Double_t(1.)/sizen ; //normalization factor
993     }else{
994       AliInfo("Problem in normalizing the vector");
995     }
996     
997     // plane expressed in the hessian normal form, see:
998     // http://mathworld.wolfram.com/HessianNormalForm.html
999     // the first three are the coordinates of the orthonormal vector
1000     // the fourth coordinate is equal to the distance from the origin
1001   
1002     for(Int_t i=0;i<3;i++){
1003       plane[i] = n[i] * s;
1004     }
1005     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
1006     
1007     // The center of the square with fiducial marks as corners
1008     // as the middle point of one diagonal - md
1009     // Used below to get the center - orig - of the surveyed box
1010
1011     Double_t orig[3], md[3];
1012     for(Int_t i=0;i<3;i++){
1013       md[i] = (ngA[i] + ngC[i]) * 0.5;
1014     }
1015     
1016     // The center of the box, gives the global translation
1017     for(Int_t i=0;i<3;i++){
1018       orig[i] = md[i] + plane[i]*fgkZFM;
1019     }
1020     
1021     // get local directions needed to write the global rotation matrix
1022     // for the surveyed volume by normalising vectors ab and bc
1023     Double_t sx = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
1024     if(sx>1.e-8){
1025       for(Int_t i=0;i<3;i++){
1026         ad[i] /= sx;
1027       }
1028     }
1029     Double_t sy = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1030     if(sy>1.e-8){
1031       for(Int_t i=0;i<3;i++){
1032         cd[i] /= sy;
1033       }
1034     }
1035     Double_t rot[9] = {cd[0],ad[0],-plane[0],cd[1],ad[1],-plane[1],cd[2],ad[2],-plane[2]};
1036     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1037     TGeoHMatrix ng;              
1038     ng.SetTranslation(orig);
1039     ng.SetRotation(rot);
1040     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1041     ng.Print();    
1042
1043     // Calculate the delta transformation wrt Ideal geometry
1044     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1045     
1046     printf("\n\n**** The ideal matrix ***\n"); 
1047     fTOFMatrixId[iSM]->Print();
1048     
1049     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1050     printf("\n\n**** The inverse of the ideal matrix ***\n");
1051     gdelta.Print();
1052  
1053     gdelta.MultiplyLeft(&ng);
1054     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1055     gdelta.Print(); //global delta trasformation
1056     
1057     // Now Write the Alignment Objects....
1058     Int_t index=0; //let all SM modules have index=0
1059     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1060     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
1061     TString symname(Form("TOF/sm%02d",iSM));
1062     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1063     fTOFAlignObjArray->Add(o);
1064   }
1065
1066 //___________________________________________________________________________
1067 void AliTOFAlignment::AlignFromSurveyBCD(Int_t iSM)
1068 {
1069   //From Survey data, derive the needed transformations to get the 
1070   //Alignment Objects. 
1071   //Again, highly "inspired" to Raffaele's example... 
1072   //we use FM B,C,D
1073
1074     Double_t ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
1075     
1076     
1077    // Get the 'realistic' input from the Survey Matrix
1078       for(Int_t coord=0;coord<3;coord++){
1079       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
1080       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
1081       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
1082       }
1083
1084     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
1085
1086     // From the new fiducial marks coordinates derive back the
1087     // new global position of the surveyed volume
1088     //*** What follows is the actual survey-to-alignment procedure
1089     
1090     Double_t cd[3], bc[3], n[3];
1091     Double_t plane[4], s=1.;
1092     
1093     // first vector on the plane of the fiducial marks
1094     for(Int_t i=0;i<3;i++){
1095       cd[i] = (ngC[i] - ngD[i]);
1096     }
1097     
1098     // second vector on the plane of the fiducial marks
1099     for(Int_t i=0;i<3;i++){
1100       bc[i] = (ngC[i] - ngB[i]);
1101     }
1102     
1103     // vector normal to the plane of the fiducial marks obtained
1104     // as cross product of the two vectors on the plane d0^d1
1105     n[0] = (bc[1] * cd[2] - bc[2] * cd[1]);
1106     n[1] = (bc[2] * cd[0] - bc[0] * cd[2]);
1107     n[2] = (bc[0] * cd[1] - bc[1] * cd[0]);
1108     
1109     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1110     if(sizen>1.e-8){
1111       s = Double_t(1.)/sizen ; //normalization factor
1112     }else{
1113       AliInfo("Problem in normalizing the vector");
1114     }
1115     
1116     // plane expressed in the hessian normal form, see:
1117     // http://mathworld.wolfram.com/HessianNormalForm.html
1118     // the first three are the coordinates of the orthonormal vector
1119     // the fourth coordinate is equal to the distance from the origin
1120   
1121     for(Int_t i=0;i<3;i++){
1122       plane[i] = n[i] * s;
1123     }
1124     plane[3] = ( plane[0] * ngB[0] + plane[1] * ngB[1] + plane[2] * ngB[2] );
1125     
1126     // The center of the square with fiducial marks as corners
1127     // as the middle point of one diagonal - md
1128     // Used below to get the center - orig - of the surveyed box
1129
1130     Double_t orig[3], md[3];
1131     for(Int_t i=0;i<3;i++){
1132       md[i] = (ngB[i] + ngD[i]) * 0.5;
1133     }
1134     
1135     // The center of the box, gives the global translation
1136     for(Int_t i=0;i<3;i++){
1137       orig[i] = md[i] + plane[i]*fgkZFM;
1138     }
1139     
1140     // get local directions needed to write the global rotation matrix
1141     // for the surveyed volume by normalising vectors ab and bc
1142     Double_t sx = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1143     if(sx>1.e-8){
1144       for(Int_t i=0;i<3;i++){
1145         cd[i] /= sx;
1146       }
1147     }
1148     Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
1149     if(sy>1.e-8){
1150       for(Int_t i=0;i<3;i++){
1151         bc[i] /= sy;
1152       }
1153     }
1154     Double_t rot[9] = {cd[0],bc[0],-plane[0],cd[1],bc[1],-plane[1],cd[2],bc[2],-plane[2]};
1155     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1156     TGeoHMatrix ng;              
1157     ng.SetTranslation(orig);
1158     ng.SetRotation(rot);
1159     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1160     ng.Print();    
1161
1162     // Calculate the delta transformation wrt Ideal geometry
1163     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1164     
1165     printf("\n\n**** The ideal matrix ***\n"); 
1166     fTOFMatrixId[iSM]->Print();   
1167     
1168     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1169     printf("\n\n**** The inverse of the ideal matrix ***\n");
1170     gdelta.Print();
1171  
1172     gdelta.MultiplyLeft(&ng);
1173     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1174     gdelta.Print(); //global delta trasformation
1175     
1176     // Now Write the Alignment Objects....
1177     Int_t index=0; //let all SM modules have index=0
1178     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1179     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
1180     TString symname(Form("TOF/sm%02d",iSM));
1181     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1182     fTOFAlignObjArray->Add(o);
1183   }
1184
1185