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