]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFAlignment.cxx
Include Methods to derive TOF AlignObjs from Survey Data
[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.12  2007/02/28 18:09:23  arcelli
19 Add protection against failed retrieval of the CDB cal object
20
21 Revision 1.11  2006/09/19 14:31:26  cvetan
22 Bugfixes and clean-up of alignment object classes. Introduction of so called symbolic names used to identify the alignable volumes (Raffaele and Cvetan)
23
24 Revision 1.10  2006/08/22 13:26:05  arcelli
25 removal of effective c++ warnings (C.Zampolli)
26
27 Revision 1.9  2006/08/10 14:46:54  decaro
28 TOF raw data format: updated version
29
30 Revision 1.8  2006/05/04 19:41:42  hristov
31 Possibility for partial TOF geometry (S.Arcelli)
32
33 Revision 1.7  2006/04/27 13:13:29  hristov
34 Moving the destructor to the implementation file
35
36 Revision 1.6  2006/04/20 22:30:49  hristov
37 Coding conventions (Annalisa)
38
39 Revision 1.5  2006/04/16 22:29:05  hristov
40 Coding conventions (Annalisa)
41
42 Revision 1.4  2006/04/05 08:35:38  hristov
43 Coding conventions (S.Arcelli, C.Zampolli)
44
45 Revision 1.3  2006/03/31 13:49:07  arcelli
46 Removing some junk printout
47
48 Revision 1.2  2006/03/31 11:26:30  arcelli
49  changing CDB Ids according to standard convention
50
51 Revision 1.1  2006/03/28 14:54:48  arcelli
52 class for TOF alignment
53
54 author: Silvia Arcelli, arcelli@bo.infn.it
55 */  
56
57 /////////////////////////////////////////////////////////
58 //                                                     //
59 //            Class for alignment procedure            //
60 //                                                     //
61 //                                                     //
62 //                                                     //
63 /////////////////////////////////////////////////////////
64
65 #include <Rtypes.h>
66
67 #include "TMath.h"
68 #include "TFile.h"
69 #include "TRandom.h"
70
71 #include "AliLog.h"
72 #include "AliAlignObj.h"
73 #include "AliAlignObjAngles.h"
74 #include "AliAlignObjMatrix.h"
75 #include "AliCDBManager.h"
76 #include "AliCDBMetaData.h"
77 #include "AliCDBId.h"
78 #include "AliCDBEntry.h"
79 #include "AliTOFAlignment.h"
80
81
82 ClassImp(AliTOFAlignment)
83 const Double_t AliTOFAlignment::fgkXsizeTOF  = 124.5; // x size of the TOF ext. volume, cm
84 const Double_t AliTOFAlignment::fgkYsizeTOF  = 29.0;  // y size of the TOF ext. volume, cm
85 const Double_t AliTOFAlignment::fgkZsizeTOF  = 913.8; // z size of the TOF ext. volume, cm
86 const Double_t AliTOFAlignment::fgkRorigTOF  = 384.5; // Mean Radius of the TOF ext. volume, cm
87 const Double_t AliTOFAlignment::fgkXFM = 38.0; //x pos, cm 
88 const Double_t AliTOFAlignment::fgkYFM = 11.2; //y pos, cm
89 const Double_t AliTOFAlignment::fgkZFM = 457.3;//z pos, cm
90 const Double_t AliTOFAlignment::fgkZsizeTOFSens=741.2; //z size of the TOF sensitive volume, cm
91
92 //_____________________________________________________________________________
93 AliTOFAlignment::AliTOFAlignment():
94   TTask("AliTOFAlignment",""),
95   fNTOFAlignObj(0),
96   fTOFmgr(0x0),
97   fTOFAlignObjArray(0x0)
98  { 
99    //AliTOFalignment main Ctor
100    for(Int_t ism=0;ism<18;ism++){
101      for(Int_t iFM=0;iFM<4;iFM++){
102        for(Int_t iFMc=0;iFMc<3;iFMc++){
103          fTOFSurveyFM[ism][iFM][iFMc]=-1.;
104        }
105      }
106    }
107 }
108 //_____________________________________________________________________________
109 AliTOFAlignment::AliTOFAlignment(const AliTOFAlignment &t):
110   TTask("AliTOFAlignment",""),
111   fNTOFAlignObj(0),
112   fTOFmgr(0x0),
113   fTOFAlignObjArray(0x0)
114
115   //AliTOFAlignment copy Ctor
116
117   fNTOFAlignObj=t.fNTOFAlignObj;
118   fTOFAlignObjArray=t.fTOFAlignObjArray;
119   //AliTOFalignment main Ctor
120   for(Int_t iSM=0;iSM<18;iSM++){
121      for(Int_t iFM=0;iFM<4;iFM++){
122        for(Int_t iFMc=0;iFMc<3;iFMc++){
123          fTOFSurveyFM[iSM][iFM][iFMc]=-1.;
124        }
125      }
126   }  
127 }
128 //_____________________________________________________________________________
129 AliTOFAlignment& AliTOFAlignment::operator=(const AliTOFAlignment &t){ 
130   //AliTOFAlignment assignment operator
131
132   this->fNTOFAlignObj=t.fNTOFAlignObj;
133   this->fTOFmgr=t.fTOFmgr;
134   this->fTOFAlignObjArray=t.fTOFAlignObjArray;
135   return *this;
136
137 }
138 //_____________________________________________________________________________
139 AliTOFAlignment::~AliTOFAlignment() {delete fTOFAlignObjArray;}
140
141 //_____________________________________________________________________________
142 void AliTOFAlignment::Smear( Float_t *tr, Float_t *rot)
143 {
144   //Introduce Random Offset/Tilts
145   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
146   Float_t dx, dy, dz;  // shifts
147   Float_t dpsi, dtheta, dphi; // angular displacements
148   TRandom *rnd   = new TRandom(1567);
149  
150   Int_t nSMTOF = 18;
151   AliAlignObj::ELayerID iLayer = AliAlignObj::kInvalidLayer;
152   UShort_t iIndex=0; //dummy volume index
153   //  AliAlignObj::ELayerID iLayer = AliAlignObj::kTOF;
154   //  Int_t iIndex=1; //dummy volume index
155   UShort_t dvoluid = AliAlignObj::LayerToVolUID(iLayer,iIndex); //dummy volume identity 
156   Int_t i;
157   for (i = 0; i<nSMTOF ; i++) {
158     Char_t  path[100];
159     sprintf(path,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
160
161     dx = (rnd->Gaus(0.,1.))*tr[0];
162     dy = (rnd->Gaus(0.,1.))*tr[1];
163     dz = (rnd->Gaus(0.,1.))*tr[2];
164     dpsi   = rot[0];
165     dtheta = rot[1];
166     dphi   = rot[2];
167     AliAlignObjAngles *o =new AliAlignObjAngles(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
168     fTOFAlignObjArray->Add(o);
169   }
170
171   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
172   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
173   delete rnd;
174 }
175
176 //_____________________________________________________________________________
177 void AliTOFAlignment::Align( Float_t *tr, Float_t *rot)
178 {
179   //Introduce Offset/Tilts
180
181   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
182   Float_t dx, dy, dz;  // shifts
183   Float_t dpsi, dtheta, dphi; // angular displacements
184
185
186   Int_t nSMTOF = 18;
187   AliAlignObj::ELayerID iLayer = AliAlignObj::kInvalidLayer;
188   UShort_t iIndex=0; //dummy volume index
189   UShort_t dvoluid = AliAlignObj::LayerToVolUID(iLayer,iIndex); //dummy volume identity 
190   Int_t i;
191   for (i = 0; i<nSMTOF ; i++) {
192
193     Char_t  path[100];
194     sprintf(path,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
195     dx = tr[0];
196     dy = tr[1];
197     dz = tr[2];
198     dpsi   = rot[0];
199     dtheta = rot[1];
200     dphi   = rot[2];
201     
202     AliAlignObjAngles *o =new AliAlignObjAngles(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
203     fTOFAlignObjArray->Add(o);
204   }
205   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
206   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
207 }
208 //_____________________________________________________________________________
209 void AliTOFAlignment::WriteParOnCDB(Char_t *sel, Int_t minrun, Int_t maxrun)
210 {
211   //Write Align Par on CDB
212   AliCDBManager *man = AliCDBManager::Instance();
213   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
214   Char_t *sel1 = "AlignPar" ;
215   Char_t  out[100];
216   sprintf(out,"%s/%s",sel,sel1); 
217   AliCDBId idTOFAlign(out,minrun,maxrun);
218   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
219   mdTOFAlign->SetResponsible("TOF");
220   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
221   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
222 }
223 //_____________________________________________________________________________
224 void AliTOFAlignment::ReadParFromCDB(Char_t *sel, Int_t nrun)
225 {
226   //Read Align Par from CDB
227   AliCDBManager *man = AliCDBManager::Instance();
228   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
229   Char_t *sel1 = "AlignPar" ;
230   Char_t  out[100];
231   sprintf(out,"%s/%s",sel,sel1); 
232   AliCDBEntry *entry = man->Get(out,nrun);
233   if (!entry) { 
234     AliError(Form("Failed to get entry: %s",out));
235     return; 
236   }
237   fTOFAlignObjArray=(TObjArray*)entry->GetObject();
238   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
239   AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
240
241 }
242 //_____________________________________________________________________________
243 void AliTOFAlignment::WriteSimParOnCDB(Char_t *sel, Int_t minrun, Int_t maxrun)
244 {
245   //Write Sim Align Par on CDB
246   AliCDBManager *man = AliCDBManager::Instance();
247   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
248   Char_t *sel1 = "AlignSimPar" ;
249   Char_t  out[100];
250   sprintf(out,"%s/%s",sel,sel1); 
251   AliCDBId idTOFAlign(out,minrun,maxrun);
252   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
253   mdTOFAlign->SetResponsible("TOF");
254   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
255   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
256 }
257 //_____________________________________________________________________________
258 void AliTOFAlignment::ReadSimParFromCDB(Char_t *sel, Int_t nrun){
259   //Read Sim Align Par from CDB
260   AliCDBManager *man = AliCDBManager::Instance();
261   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
262   Char_t *sel1 = "AlignSimPar" ;
263   Char_t  out[100];
264   sprintf(out,"%s/%s",sel,sel1); 
265   AliCDBEntry *entry = man->Get(out,nrun);
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::WriteOnCDBforDC()
273 {
274   //Write Align Par on CDB for DC06
275   AliCDBManager *man = AliCDBManager::Instance();
276   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
277   AliCDBId idTOFAlign("TOF/Align/Data",0,0);
278   AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
279   mdTOFAlign->SetComment("Alignment objects for ideal geometry, i.e. applying them to TGeo has to leave geometry unchanged");
280   mdTOFAlign->SetResponsible("TOF");
281   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
282   man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
283 }
284 //_____________________________________________________________________________
285 void AliTOFAlignment::ReadFromCDBforDC()
286 {
287   //Read Sim Align Par from CDB for DC06
288   AliCDBManager *man = AliCDBManager::Instance();
289   if(!man->IsDefaultStorageSet())man->SetDefaultStorage("local://$ALICE_ROOT");
290   AliCDBEntry *entry = man->Get("TOF/Align/Data",0);
291   fTOFAlignObjArray=(TObjArray*)entry->GetObject();
292   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
293   AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
294
295 }
296 //_____________________________________________________________________________
297 void AliTOFAlignment::BuildGeomForSurvey()
298 {
299
300   //Simulates a Misalignment and generates a list of FM coordinates in the 
301   //global RS, to be passed to the survey-to-alignment algo. 
302   //Highly inspired to Raffaele's example... 
303
304   fTOFmgr = new TGeoManager("Geom","survey to alignment for TOF");
305   TGeoMedium *medium = 0;
306   TGeoVolume *top = fTOFmgr->MakeBox("TOP",medium,1000,1000,1000);
307   fTOFmgr->SetTopVolume(top);
308   // make shape components:  
309   // This is the big box containing the TOF master sensitive volume+services  
310   TGeoBBox *sbox0  = new TGeoBBox(fgkXsizeTOF*0.5,fgkYsizeTOF*0.5,fgkZsizeTOF*0.5);
311   TGeoVolume* box0[18];
312   // This is the big box containing the TOF master sensitive volume  
313   TGeoBBox *sbox1  = new TGeoBBox(fgkXsizeTOF*0.5,fgkYsizeTOF*0.5,fgkZsizeTOFSens*0.5);
314   TGeoVolume* box1 = new TGeoVolume("B1",sbox1);
315   box1->SetLineColor(3);//green
316
317   // Now four fiducial marks on SM, expressed in local coordinates
318   // They are positioned at x=+/- 38 cm, y=11.2, z=+/- 456.94 cm
319
320   TGeoBBox *fmbox  = new TGeoBBox(1,1,1);
321   TGeoVolume* fm = new TGeoVolume("FM",fmbox);
322   fm->SetLineColor(2);//color
323
324   TGeoTranslation* Atr = new TGeoTranslation("Atr",-fgkXFM, fgkYFM ,fgkZFM);
325   TGeoTranslation* Btr = new TGeoTranslation("Btr", fgkXFM, fgkYFM, fgkZFM);
326   TGeoTranslation* Ctr = new TGeoTranslation("Ctr", fgkXFM, fgkYFM,-fgkZFM);
327   TGeoTranslation* Dtr = new TGeoTranslation("Dtr",-fgkXFM, fgkYFM,-fgkZFM);
328
329   // position all this stuff in the global ALICE frame
330
331   char name[16];
332   Double_t smX = 0.;
333   Double_t smY = 0.;
334   Double_t smZ = 0.;
335   Float_t  smR = fgkRorigTOF;
336
337   for (Int_t iSM = 0; iSM < 18; iSM++) {
338     Int_t mod = iSM + 13;
339     if (mod > 17) mod -= 18;
340     sprintf(name, "BTOF%d",mod);
341     box0[iSM] = new TGeoVolume(name,sbox0);
342     Float_t phi  = iSM * 20.;
343     Float_t phirot = 180 + phi;    
344     smX =  TMath::Sin(phi*TMath::Pi()/180.)*smR;
345     smY = -TMath::Cos(phi*TMath::Pi()/180.)*smR;
346     smZ = 0.;
347     TGeoRotation* smRot = new TGeoRotation("smRot",phirot,0,0.);    
348     TGeoCombiTrans trans = *(new TGeoCombiTrans(smX,smY,smZ, smRot));
349     TGeoMatrix* id = new TGeoHMatrix();
350     TGeoHMatrix  transMat = *id * trans;
351     TGeoHMatrix  *smTrans = new TGeoHMatrix(transMat);
352     box0[iSM]->SetVisDaughters();
353     box0[iSM]->SetLineColor(1); //black
354     top->AddNode(box0[iSM],1,smTrans); //place the extended SM volume
355     box0[iSM]->AddNode(box1,1); //place the inner SM volume
356     box0[iSM]->AddNode(fm,1,Atr);
357     box0[iSM]->AddNode(fm,2,Btr);
358     box0[iSM]->AddNode(fm,3,Ctr);
359     box0[iSM]->AddNode(fm,4,Dtr);
360   }  
361
362   fTOFmgr->CloseGeometry();
363   fTOFmgr->GetTopVolume()->Draw();
364   fTOFmgr->SetVisOption(0);
365   fTOFmgr->SetVisLevel(6);
366
367   // Now Store the ideal Matrices for later use....
368
369   for (Int_t iSM = 0; iSM < 18; iSM++) {
370
371     sprintf(name, "TOP_1/BTOF%d_1", iSM);
372     printf("\n\n************  SuperModule N.r: ************** %s \n",name);
373     TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
374     fTOFMatrixId[iSM] = pn3->GetMatrix(); //save "ideal" global matrix 
375     printf("\n\n************  Ideal global matrix, 1 **************\n");
376     fTOFMatrixId[iSM]->Print();
377
378   }
379 }
380 //_____________________________________________________________________________
381 void AliTOFAlignment::TestAlignFromSurvey( Float_t *mis)
382 {
383   // Now Apply the Displacements and store the misaligned FM positions...
384
385   Double_t A[3]={-fgkXFM,fgkYFM, fgkZFM};
386   Double_t B[3]={ fgkXFM,fgkYFM, fgkZFM};
387   Double_t C[3]={ fgkXFM,fgkYFM,-fgkZFM};
388   Double_t D[3]={-fgkXFM,fgkYFM,-fgkZFM};
389
390   for(Int_t iSM=0;iSM<18;iSM++){
391   // ************* get ideal global matrix *******************
392     char name[16];
393     sprintf(name, "TOP_1/BTOF%d_1", iSM);
394     fTOFmgr->cd(name);
395     printf("\n\n************  SuperModule N.r: ************** %s \n",name);
396
397   // ************* get ideal local matrix *******************
398     TGeoHMatrix g3 = *fTOFmgr->GetCurrentMatrix(); 
399     TGeoNode* n3 = fTOFmgr->GetCurrentNode();
400     TGeoMatrix* l3 = n3->GetMatrix(); 
401
402     Double_t gA[3], gB[3], gC[3], gD[3]; // ideal FM point coord., global RS
403     g3.LocalToMaster(A,gA);
404     g3.LocalToMaster(B,gB);
405     g3.LocalToMaster(C,gC);
406     g3.LocalToMaster(D,gD);
407
408     // We apply a delta transformation to the surveyed vol to represent
409     // its real position, given below by ng3 nl3, which differs from its
410     // ideal position saved above in g3 and l3
411
412
413     Double_t dx     = mis[0]; // shift along x
414     Double_t dy     = mis[1]; // shift along y
415     Double_t dz     = mis[2]; // shift along z
416     Double_t dphi   = mis[3]; // rot around z
417     Double_t dtheta = mis[4]; // rot around x'
418     Double_t dpsi   = mis[5]; // rot around z'
419
420     TGeoRotation* rrot = new TGeoRotation("rot",dphi,dtheta,dpsi);
421     TGeoCombiTrans localdelta = *(new TGeoCombiTrans(dx,dy,dz, rrot));
422   // new local matrix, representing real position
423     TGeoHMatrix nlocal = *l3 * localdelta;
424     TGeoHMatrix* nl3 = new TGeoHMatrix(nlocal);
425
426     TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
427     TGeoHMatrix* ng2 = pn3->GetMatrix(); //"real" global matrix, what survey sees 
428     printf("\n\n************  Ideal global matrix, 2 **************\n");
429     ng2->Print();
430
431     pn3->Align(nl3);    //Align....    
432     
433     TGeoHMatrix* ng3 = pn3->GetMatrix(); //"real" global matrix, what survey sees 
434     printf("\n\n************  Misaligned global matrix **************\n");
435     ng3->Print();
436     Double_t ngA[3], ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
437     ng3->LocalToMaster(A,ngA);
438     ng3->LocalToMaster(B,ngB);
439     ng3->LocalToMaster(C,ngC);
440     ng3->LocalToMaster(D,ngD);    
441
442     for(Int_t iFM=0;iFM<3;iFM++){
443       fTOFSurveyFM[iSM][0][iFM]=ngA[iFM];
444       fTOFSurveyFM[iSM][1][iFM]=ngB[iFM];
445       fTOFSurveyFM[iSM][2][iFM]=ngC[iFM];
446       fTOFSurveyFM[iSM][3][iFM]=ngD[iFM];
447     }
448   }
449 }
450
451 //_____________________________________________________________________________
452 void AliTOFAlignment::AlignFromSurvey()
453 {
454   //From survey Data, derive the needed transformations to get the 
455   //Alignment Objects. 
456   //Again, highly "inspired" to Raffaele's example... 
457
458   fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
459   Int_t index=0; //let all SM modules have index=0
460   AliAlignObj::ELayerID layer = AliAlignObj::kInvalidLayer;
461   UShort_t dvoluid = AliAlignObj::LayerToVolUID(layer,index); //dummy vol id 
462   
463   for(Int_t iSM=0;iSM<18;iSM++){
464     Double_t ngA[3], ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
465  
466    // Get the input from the Survey Matrix
467     for(Int_t iFM=0;iFM<3;iFM++){
468       ngA[iFM]=   fTOFSurveyFM[iSM][0][iFM];
469       ngB[iFM]=   fTOFSurveyFM[iSM][1][iFM];
470       ngD[iFM]=   fTOFSurveyFM[iSM][2][iFM];
471       ngC[iFM]=   fTOFSurveyFM[iSM][3][iFM];
472     }
473
474     // From the new fiducial marks coordinates derive back the
475     // new global position of the surveyed volume
476     //*** What follows is the actual survey-to-alignment procedure
477     
478     Double_t ab[3], bc[3], n[3];
479     Double_t plane[4], s=1.;
480     
481     // first vector on the plane of the fiducial marks
482     for(Int_t i=0;i<3;i++){
483       ab[i] = (ngB[i] - ngA[i]);
484     }
485     
486     // second vector on the plane of the fiducial marks
487     for(Int_t i=0;i<3;i++){
488       bc[i] = (ngC[i] - ngB[i]);
489     }
490     
491     // vector normal to the plane of the fiducial marks obtained
492     // as cross product of the two vectors on the plane d0^d1
493     n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
494     n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
495     n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
496     
497     Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
498     if(sizen>1.e-8){
499       s = Double_t(1.)/sizen ; //normalization factor
500     }else{
501       AliInfo("Problem in normalizing the vector");
502     }
503     
504     // plane expressed in the hessian normal form, see:
505     // http://mathworld.wolfram.com/HessianNormalForm.html
506     // the first three are the coordinates of the orthonormal vector
507     // the fourth coordinate is equal to the distance from the origin
508   
509     for(Int_t i=0;i<3;i++){
510       plane[i] = n[i] * s;
511     }
512     plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
513     
514     // The center of the square with fiducial marks as corners
515     // as the middle point of one diagonal - md
516     // Used below to get the center - orig - of the surveyed box
517
518     Double_t orig[3], md[3];
519     for(Int_t i=0;i<3;i++){
520       md[i] = (ngA[i] + ngC[i]) * 0.5;
521     }
522     
523     // The center of the box, gives the global translation
524
525     for(Int_t i=0;i<3;i++){
526       orig[i] = md[i] - plane[i]*fgkYFM;
527     }
528     
529     // get local directions needed to write the global rotation matrix
530     // for the surveyed volume by normalising vectors ab and bc
531
532     Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
533     if(sx>1.e-8){
534       for(Int_t i=0;i<3;i++){
535         ab[i] /= sx;
536       }
537     }
538     Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
539     if(sy>1.e-8){
540       for(Int_t i=0;i<3;i++){
541         bc[i] /= sy;
542       }
543     }
544         
545     Double_t rot[9] = {ab[0],plane[0],bc[0],ab[1],plane[1],-bc[1],ab[2],plane[2],-bc[2]}; // the rotation matrix
546     TGeoHMatrix ng;
547     ng.SetTranslation(orig);
548     ng.SetRotation(rot);
549     ng.Print();
550     // Calculate the delta transformation wrt Ideal
551     TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
552     gdelta.MultiplyLeft(&ng);
553
554   // Now Get the alignment Objects....
555      TString symname(Form("TOF/sm%02d",iSM));
556 //  // if the volume is in the look-up table use something like this instead:
557      AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
558      fTOFAlignObjArray->Add(o);
559   }
560
561   // saving TOF AligObjs from survey on a file, for the moment.. 
562   fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
563   AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
564   TFile f("TOFAlignFromSurvey.root","RECREATE");
565   f.cd();
566   f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
567   f.Close();
568 }