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