]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFAlignment.cxx
Reducing logs for FDR
[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(Char_t *sel, Int_t minrun, Int_t maxrun)
239 {
240   //Write Align Par on CDB
241   AliCDBManager *man = AliCDBManager::Instance();
242   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(Char_t *sel, Int_t nrun)
253 {
254   //Read Align Par from CDB
255   AliCDBManager *man = AliCDBManager::Instance();
256   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(Char_t *sel, Int_t minrun, Int_t maxrun)
271 {
272   //Write Sim Align Par on CDB
273   AliCDBManager *man = AliCDBManager::Instance();
274   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(Char_t *sel, Int_t nrun){
285   //Read Sim Align Par from CDB
286   AliCDBManager *man = AliCDBManager::Instance();
287   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 i=0;i<nf; i++)
571     {
572       AliSurveyObj *so = new AliSurveyObj();
573       const Char_t *nome=namefiles[i];
574       so->FillFromLocalFile(nome);
575       TObjArray *points = so->GetData();
576       Int_t nSurveyPoint=points->GetEntries();
577       for(Int_t i=0;i<nSurveyPoint;i++){
578         const char* pointName= ((AliSurveyPoint *) points->At(i))->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(i))->GetX();
582         data[nsm*4+nfm][2][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(i))->GetY();
583         data[nsm*4+nfm][4][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(i))->GetZ();
584         data[nsm*4+nfm][1][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(i))->GetPrecisionX();
585         data[nsm*4+nfm][3][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(i))->GetPrecisionY();
586         data[nsm*4+nfm][5][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(i))->GetPrecisionZ();
587         totdata[nsm*4+nfm]=totdata[nsm*4+nfm]+1;
588       } 
589       delete so;
590     }
591
592   
593
594   for(Int_t i=0; i<72 ;i++){
595     Float_t numx=0, numy=0,numz=0, comodox=0, comodoy=0, comodoz=0,denx=0, deny=0, denz=0;
596     if(totdata[i]!=0){    
597       for(Int_t j=0; j<totdata[i]; j++){
598         comodox=1/(data[i][1][j]/10*data[i][1][j]/10);//precision in mm, position in cm
599         numx=numx+data[i][0][j]*comodox;
600         denx=denx+comodox;
601         comodoy=1/(data[i][3][j]/10*data[i][3][j]/10);
602         numy=numy+data[i][2][j]*comodoy;
603         deny=deny+comodoy;
604         comodoz=1/(data[i][5][j]/10*data[i][5][j]/10);
605         numz=numz+data[i][4][j]*comodoz;
606         denz=denz+comodoz;
607         }
608       fCombFMData[i][1]=TMath::Sqrt(1/denx); //error for x position
609       fCombFMData[i][3]=TMath::Sqrt(1/deny); //error for y position
610       fCombFMData[i][5]=TMath::Sqrt(1/denz); //error for z position
611       fCombFMData[i][0]=numx/denx;           //combined survey data for x position of FM
612       fCombFMData[i][2]=numy/deny;           //combined survey data for y position of FM
613       fCombFMData[i][4]=numz/denz;           //combined survey data for z position of FM
614       } else continue;
615     }
616
617   for(Int_t i=0;i<72;i++)
618     if (fCombFMData[i][0]!=0){
619       fNFMforSM[(i-i%4)/4][i%4]=1;
620       fNFMforSM[(i-i%4)/4][4]=fNFMforSM[(i-i%4)/4][4]+1;
621     }
622 }
623
624 //_____________________________________________________________________________
625 void AliTOFAlignment::ReadSurveyDataAndAlign(){
626   //
627   // read the survey data and, if we know the positions of at least 3 FM 
628   //for a SM, call the right Alignement procedure  
629
630   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
631
632   Float_t deltaFM0=0, deltaFM1=0, deltaFM2=0, deltaFM3=0;
633
634   for(Int_t i=0; i<18; i++){
635     switch(fNFMforSM[i][4]){
636     case 0:
637       printf("we don't know the position of any FM of SM %i\n",i);
638       break;
639     case 1:
640       printf("we know the position of only one FM for SM %i\n",i);
641      
642       break;
643     case 2:
644       printf("we know the position of only 2 FM for SM %i\n",i);
645       
646       break;
647     case 3:
648       if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1){
649         printf("we know the position of FM A B C for SM %i\n",i);
650         AliTOFAlignment::AlignFromSurveyABC(i);};
651
652         
653       if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][3]==1){
654         printf("we know the position of FM A B D for SM %i\n",i);
655         AliTOFAlignment::AlignFromSurveyABD(i);};
656
657         
658       if (fNFMforSM[i][0]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
659         printf("we know the position of FM A C D for SM %i\n",i);
660         AliTOFAlignment::AlignFromSurveyACD(i);};
661
662         
663       if (fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
664         printf("we know the position of FM B C D for SM %i\n",i);
665         AliTOFAlignment::AlignFromSurveyBCD(i);};
666
667         
668       break;
669     case 4:
670       printf("we know the position of all the 4 FM for SM %i\n",i);
671       //check the precision of the measurement
672
673       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]);
674       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]);
675       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]);
676       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]);
677
678       //to AlignFromSurvey we use the 3 FM whose positions are known with greatest precision
679       if(deltaFM0>=deltaFM1 && deltaFM0>=deltaFM2 && deltaFM0>=deltaFM3){
680         printf("to Align we use FM B,C,D");
681         AliTOFAlignment::AlignFromSurveyBCD(i);} else
682           if(deltaFM1>=deltaFM0 && deltaFM1>=deltaFM2 && deltaFM1>=deltaFM3){
683            printf("to Align we use FM A,C,D");
684            AliTOFAlignment::AlignFromSurveyACD(i);} else
685              if(deltaFM2>=deltaFM0 && deltaFM2>=deltaFM1 && deltaFM2>=deltaFM3){
686                printf("to Align we use FM A,B,D");
687                AliTOFAlignment::AlignFromSurveyABD(i);} else{
688                  printf("to Align we use FM A,B,C");
689                  AliTOFAlignment::AlignFromSurveyABC(i);}
690      
691       break;
692     }
693   
694   }
695
696     // saving TOF AligObjs from survey on a file, for the moment.. 
697   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
698   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
699   TFile f("TOFAlignFromSurvey.root","RECREATE");
700   f.cd();
701   f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
702   f.Close();
703   
704
705 }
706
707 //_____________________________________________________________________________
708 void AliTOFAlignment::AlignFromSurveyABC(Int_t iSM)
709 {
710
711   //From Survey data, derive the needed transformations to get the 
712   //Alignment Objects. 
713   //Again, highly "inspired" to Raffaele's example... 
714   //we use FM A,B,C
715     
716     Double_t ngA[3], ngB[3], ngC[3]; // real FM point coord., global RS
717     // Get the 'realistic' input from the Survey Matrix
718       for(Int_t coord=0;coord<3;coord++){
719       ngA[coord]=   fCombFMData[iSM*4][coord*2];
720       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
721       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
722       }
723
724     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
725
726     // From the real fiducial marks coordinates derive back the
727     // new global position of the surveyed volume
728     //*** What follows is the actual survey-to-alignment procedure
729     
730     Double_t ab[3], bc[3], n[3];
731     Double_t plane[4], s=1.;
732     
733     // first vector on the plane of the fiducial marks
734     for(Int_t i=0;i<3;i++){
735       ab[i] = (ngB[i] - ngA[i]);
736     }
737     
738     // second vector on the plane of the fiducial marks
739     for(Int_t i=0;i<3;i++){
740       bc[i] = (ngC[i] - ngB[i]);
741     }
742     
743     // vector normal to the plane of the fiducial marks obtained
744     // as cross product of the two vectors on the plane d0^d1
745     n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
746     n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
747     n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
748     
749     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
750     if(sizen>1.e-8){
751       s = Double_t(1.)/sizen ; //normalization factor
752     }else{
753       AliInfo("Problem in normalizing the vector");
754     }
755     
756     // plane expressed in the hessian normal form, see:
757     // http://mathworld.wolfram.com/HessianNormalForm.html
758     // the first three are the coordinates of the orthonormal vector
759     // the fourth coordinate is equal to the distance from the origin
760   
761     for(Int_t i=0;i<3;i++){
762       plane[i] = n[i] * s;
763     }
764     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
765     
766     // The center of the square with fiducial marks as corners
767     // as the middle point of one diagonal - md
768     // Used below to get the center - orig - of the surveyed box
769
770     Double_t orig[3], md[3];
771     for(Int_t i=0;i<3;i++){
772       md[i] = (ngA[i] + ngC[i]) * 0.5;
773     }
774     
775     // The center of the box, gives the global translation
776     for(Int_t i=0;i<3;i++){
777       orig[i] = md[i] - plane[i]*fgkZFM;
778     }
779     
780     // get local directions needed to write the global rotation matrix
781     // for the surveyed volume by normalising vectors ab and bc
782     Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
783
784
785     if(sx>1.e-8){
786       for(Int_t i=0;i<3;i++){
787         ab[i] /= sx;
788       }
789     }
790     Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
791     if(sy>1.e-8){
792       for(Int_t i=0;i<3;i++){
793         bc[i] /= sy;
794       }
795     }
796     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
797     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey
798     TGeoHMatrix ng;              
799     ng.SetTranslation(orig);
800     ng.SetRotation(rot);
801     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
802     ng.Print();    
803
804     // Calculate the delta transformation wrt Ideal geometry
805     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
806     
807     printf("\n\n**** The ideal matrix ***\n"); 
808     fTOFMatrixId[iSM]->Print();   
809     
810     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
811     printf("\n\n**** The inverse of the ideal matrix ***\n");
812     gdelta.Print();
813  
814     gdelta.MultiplyLeft(&ng);
815     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
816     gdelta.Print(); //this is the global delta trasformation
817     
818     // Now Write the Alignment Objects....
819     Int_t index=0; //let all SM modules have index=0
820     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
821     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
822     TString symname(Form("TOF/sm%02d",iSM));
823     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
824     fTOFAlignObjArray->Add(o);
825
826   }
827
828
829 //_____________________________________________________________________________
830 void AliTOFAlignment::AlignFromSurveyABD(Int_t iSM)
831 {
832   
833   //From Survey data, derive the needed transformations to get the 
834   //Alignment Objects. 
835   //Again, highly "inspired" to Raffaele's example... 
836   //we use FM A,B,D
837     
838   Double_t ngA[3], ngB[3], ngD[3];// real FM point coord., global RS
839     
840    // Get the 'realistic' input from the Survey Matrix
841       for(Int_t coord=0;coord<3;coord++){
842       ngA[coord]=   fCombFMData[iSM*4][coord*2];
843       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
844       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
845       }
846
847     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
848
849     // From the new fiducial marks coordinates derive back the
850     // new global position of the surveyed volume
851     //*** What follows is the actual survey-to-alignment procedure
852     
853     Double_t ab[3], ad[3], n[3];
854     Double_t plane[4], s=1.;
855     
856     // first vector on the plane of the fiducial marks
857     for(Int_t i=0;i<3;i++){
858       ab[i] = (ngB[i] - ngA[i]);
859     }
860     
861     // second vector on the plane of the fiducial marks
862     for(Int_t i=0;i<3;i++){
863       ad[i] = (ngD[i] - ngA[i]);
864     }
865     
866     // vector normal to the plane of the fiducial marks obtained
867     // as cross product of the two vectors on the plane d0^d1
868     n[0] = (ab[1] * ad[2] - ab[2] * ad[1]);
869     n[1] = (ab[2] * ad[0] - ab[0] * ad[2]);
870     n[2] = (ab[0] * ad[1] - ab[1] * ad[0]);
871     
872     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
873     if(sizen>1.e-8){
874       s = Double_t(1.)/sizen ; //normalization factor
875     }else{
876       AliInfo("Problem in normalizing the vector");
877     }
878     
879     // plane expressed in the hessian normal form, see:
880     // http://mathworld.wolfram.com/HessianNormalForm.html
881     // the first three are the coordinates of the orthonormal vector
882     // the fourth coordinate is equal to the distance from the origin
883   
884     for(Int_t i=0;i<3;i++){
885       plane[i] = n[i] * s;
886     }
887     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
888     
889     // The center of the square with fiducial marks as corners
890     // as the middle point of one diagonal - md
891     // Used below to get the center - orig - of the surveyed box
892
893     Double_t orig[3], md[3];
894     for(Int_t i=0;i<3;i++){
895       md[i] = (ngB[i] + ngD[i]) * 0.5;
896     }
897     
898     // The center of the box, gives the global translation
899     for(Int_t i=0;i<3;i++){
900       orig[i] = md[i] - plane[i]*fgkZFM;
901     }
902     
903     // get local directions needed to write the global rotation matrix
904     // for the surveyed volume by normalising vectors ab and bc
905     Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
906     if(sx>1.e-8){
907       for(Int_t i=0;i<3;i++){
908         ab[i] /= sx;
909       }
910     }
911     Double_t sy = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
912     if(sy>1.e-8){
913       for(Int_t i=0;i<3;i++){
914         ad[i] /= sy;
915       }
916     }
917     Double_t rot[9] = {ab[0],ad[0],plane[0],ab[1],ad[1],plane[1],ab[2],ad[2],plane[2]};
918     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
919     TGeoHMatrix ng;              
920     ng.SetTranslation(orig);
921     ng.SetRotation(rot);
922     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
923     ng.Print();    
924
925     // Calculate the delta transformation wrt Ideal geometry
926     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
927     
928     printf("\n\n**** The ideal matrix ***\n"); 
929     fTOFMatrixId[iSM]->Print();   
930     
931     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
932     printf("\n\n**** The inverse of the ideal matrix ***\n");
933     gdelta.Print();
934  
935     gdelta.MultiplyLeft(&ng);
936     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
937     gdelta.Print();  //global delta trasformation
938     
939     // Now Write the Alignment Objects....
940     Int_t index=0; //let all SM modules have index=0
941     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
942     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
943     TString symname(Form("TOF/sm%02d",iSM));
944     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
945     fTOFAlignObjArray->Add(o);
946
947   }
948 //_____________________________________________________________________________
949 void AliTOFAlignment::AlignFromSurveyACD(Int_t iSM)
950 {
951   //From Survey data, derive the needed transformations to get the 
952   //Alignment Objects. 
953   //Again, highly "inspired" to Raffaele's example... 
954   //we use FM A,C,D
955   
956     
957     Double_t ngA[3], ngC[3], ngD[3];// real FM point coord., global RS
958     
959    // Get the 'realistic' input from the Survey Matrix
960       for(Int_t coord=0;coord<3;coord++){
961       ngA[coord]=   fCombFMData[iSM*4][coord*2];
962       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
963       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
964       }
965
966     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
967
968     // From the new fiducial marks coordinates derive back the
969     // new global position of the surveyed volume
970     //*** What follows is the actual survey-to-alignment procedure
971     
972     Double_t cd[3], ad[3], n[3];
973     Double_t plane[4], s=1.;
974     
975     // first vector on the plane of the fiducial marks
976     for(Int_t i=0;i<3;i++){
977       cd[i] = (ngC[i] - ngD[i]);
978     }
979     
980     // second vector on the plane of the fiducial marks
981     for(Int_t i=0;i<3;i++){
982       ad[i] = (ngD[i] - ngA[i]);
983     }
984     
985     // vector normal to the plane of the fiducial marks obtained
986     // as cross product of the two vectors on the plane d0^d1
987     n[0] = (ad[1] * cd[2] - ad[2] * cd[1]);
988     n[1] = (ad[2] * cd[0] - ad[0] * cd[2]);
989     n[2] = (ad[0] * cd[1] - ad[1] * cd[0]);
990     
991     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
992     if(sizen>1.e-8){
993       s = Double_t(1.)/sizen ; //normalization factor
994     }else{
995       AliInfo("Problem in normalizing the vector");
996     }
997     
998     // plane expressed in the hessian normal form, see:
999     // http://mathworld.wolfram.com/HessianNormalForm.html
1000     // the first three are the coordinates of the orthonormal vector
1001     // the fourth coordinate is equal to the distance from the origin
1002   
1003     for(Int_t i=0;i<3;i++){
1004       plane[i] = n[i] * s;
1005     }
1006     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
1007     
1008     // The center of the square with fiducial marks as corners
1009     // as the middle point of one diagonal - md
1010     // Used below to get the center - orig - of the surveyed box
1011
1012     Double_t orig[3], md[3];
1013     for(Int_t i=0;i<3;i++){
1014       md[i] = (ngA[i] + ngC[i]) * 0.5;
1015     }
1016     
1017     // The center of the box, gives the global translation
1018     for(Int_t i=0;i<3;i++){
1019       orig[i] = md[i] + plane[i]*fgkZFM;
1020     }
1021     
1022     // get local directions needed to write the global rotation matrix
1023     // for the surveyed volume by normalising vectors ab and bc
1024     Double_t sx = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
1025     if(sx>1.e-8){
1026       for(Int_t i=0;i<3;i++){
1027         ad[i] /= sx;
1028       }
1029     }
1030     Double_t sy = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1031     if(sy>1.e-8){
1032       for(Int_t i=0;i<3;i++){
1033         cd[i] /= sy;
1034       }
1035     }
1036     Double_t rot[9] = {cd[0],ad[0],-plane[0],cd[1],ad[1],-plane[1],cd[2],ad[2],-plane[2]};
1037     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1038     TGeoHMatrix ng;              
1039     ng.SetTranslation(orig);
1040     ng.SetRotation(rot);
1041     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1042     ng.Print();    
1043
1044     // Calculate the delta transformation wrt Ideal geometry
1045     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1046     
1047     printf("\n\n**** The ideal matrix ***\n"); 
1048     fTOFMatrixId[iSM]->Print();
1049     
1050     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1051     printf("\n\n**** The inverse of the ideal matrix ***\n");
1052     gdelta.Print();
1053  
1054     gdelta.MultiplyLeft(&ng);
1055     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1056     gdelta.Print(); //global delta trasformation
1057     
1058     // Now Write the Alignment Objects....
1059     Int_t index=0; //let all SM modules have index=0
1060     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1061     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
1062     TString symname(Form("TOF/sm%02d",iSM));
1063     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1064     fTOFAlignObjArray->Add(o);
1065   }
1066
1067 //___________________________________________________________________________
1068 void AliTOFAlignment::AlignFromSurveyBCD(Int_t iSM)
1069 {
1070   //From Survey data, derive the needed transformations to get the 
1071   //Alignment Objects. 
1072   //Again, highly "inspired" to Raffaele's example... 
1073   //we use FM B,C,D
1074
1075     Double_t ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
1076     
1077     
1078    // Get the 'realistic' input from the Survey Matrix
1079       for(Int_t coord=0;coord<3;coord++){
1080       ngB[coord]=   fCombFMData[iSM*4+1][coord*2];
1081       ngC[coord]=   fCombFMData[iSM*4+2][coord*2];
1082       ngD[coord]=   fCombFMData[iSM*4+3][coord*2];
1083       }
1084
1085     printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
1086
1087     // From the new fiducial marks coordinates derive back the
1088     // new global position of the surveyed volume
1089     //*** What follows is the actual survey-to-alignment procedure
1090     
1091     Double_t cd[3], bc[3], n[3];
1092     Double_t plane[4], s=1.;
1093     
1094     // first vector on the plane of the fiducial marks
1095     for(Int_t i=0;i<3;i++){
1096       cd[i] = (ngC[i] - ngD[i]);
1097     }
1098     
1099     // second vector on the plane of the fiducial marks
1100     for(Int_t i=0;i<3;i++){
1101       bc[i] = (ngC[i] - ngB[i]);
1102     }
1103     
1104     // vector normal to the plane of the fiducial marks obtained
1105     // as cross product of the two vectors on the plane d0^d1
1106     n[0] = (bc[1] * cd[2] - bc[2] * cd[1]);
1107     n[1] = (bc[2] * cd[0] - bc[0] * cd[2]);
1108     n[2] = (bc[0] * cd[1] - bc[1] * cd[0]);
1109     
1110     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1111     if(sizen>1.e-8){
1112       s = Double_t(1.)/sizen ; //normalization factor
1113     }else{
1114       AliInfo("Problem in normalizing the vector");
1115     }
1116     
1117     // plane expressed in the hessian normal form, see:
1118     // http://mathworld.wolfram.com/HessianNormalForm.html
1119     // the first three are the coordinates of the orthonormal vector
1120     // the fourth coordinate is equal to the distance from the origin
1121   
1122     for(Int_t i=0;i<3;i++){
1123       plane[i] = n[i] * s;
1124     }
1125     plane[3] = ( plane[0] * ngB[0] + plane[1] * ngB[1] + plane[2] * ngB[2] );
1126     
1127     // The center of the square with fiducial marks as corners
1128     // as the middle point of one diagonal - md
1129     // Used below to get the center - orig - of the surveyed box
1130
1131     Double_t orig[3], md[3];
1132     for(Int_t i=0;i<3;i++){
1133       md[i] = (ngB[i] + ngD[i]) * 0.5;
1134     }
1135     
1136     // The center of the box, gives the global translation
1137     for(Int_t i=0;i<3;i++){
1138       orig[i] = md[i] + plane[i]*fgkZFM;
1139     }
1140     
1141     // get local directions needed to write the global rotation matrix
1142     // for the surveyed volume by normalising vectors ab and bc
1143     Double_t sx = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1144     if(sx>1.e-8){
1145       for(Int_t i=0;i<3;i++){
1146         cd[i] /= sx;
1147       }
1148     }
1149     Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
1150     if(sy>1.e-8){
1151       for(Int_t i=0;i<3;i++){
1152         bc[i] /= sy;
1153       }
1154     }
1155     Double_t rot[9] = {cd[0],bc[0],-plane[0],cd[1],bc[1],-plane[1],cd[2],bc[2],-plane[2]};
1156     // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1157     TGeoHMatrix ng;              
1158     ng.SetTranslation(orig);
1159     ng.SetRotation(rot);
1160     printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1161     ng.Print();    
1162
1163     // Calculate the delta transformation wrt Ideal geometry
1164     // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1165     
1166     printf("\n\n**** The ideal matrix ***\n"); 
1167     fTOFMatrixId[iSM]->Print();   
1168     
1169     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1170     printf("\n\n**** The inverse of the ideal matrix ***\n");
1171     gdelta.Print();
1172  
1173     gdelta.MultiplyLeft(&ng);
1174     printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1175     gdelta.Print(); //global delta trasformation
1176     
1177     // Now Write the Alignment Objects....
1178     Int_t index=0; //let all SM modules have index=0
1179     AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1180     UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id 
1181     TString symname(Form("TOF/sm%02d",iSM));
1182     AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1183     fTOFAlignObjArray->Add(o);
1184   }
1185
1186