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