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