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