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