1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 ***************************************************************************/
18 Revision 1.19 2007/10/02 09:46:08 arcelli
19 add methods to retrieve real survey data, and make some analysis (by B. Guerzoni)
21 Revision 1.17 2007/06/06 16:26:46 arcelli
22 remove fall-back call to local CDB storage
24 Revision 1.16 2007/05/15 16:25:44 cvetan
25 Moving the alignment-related static methods from AliAlignObj to the new geometry steering class AliGeomManager (macro from Raffaele)
27 Revision 1.15 2007/05/03 09:25:10 decaro
28 Coding convention: RN13 violation -> suppression
30 Revision 1.14 2007/04/18 14:49:54 arcelli
31 Some code cleanup, added more debug info
33 Revision 1.13 2007/04/17 16:38:36 arcelli
34 Include Methods to derive TOF AlignObjs from Survey Data
36 Revision 1.12 2007/02/28 18:09:23 arcelli
37 Add protection against failed retrieval of the CDB cal object
39 Revision 1.11 2006/09/19 14:31:26 cvetan
40 Bugfixes and clean-up of alignment object classes. Introduction of so called symbolic names used to identify the alignable volumes (Raffaele and Cvetan)
42 Revision 1.10 2006/08/22 13:26:05 arcelli
43 removal of effective c++ warnings (C.Zampolli)
45 Revision 1.9 2006/08/10 14:46:54 decaro
46 TOF raw data format: updated version
48 Revision 1.8 2006/05/04 19:41:42 hristov
49 Possibility for partial TOF geometry (S.Arcelli)
51 Revision 1.7 2006/04/27 13:13:29 hristov
52 Moving the destructor to the implementation file
54 Revision 1.6 2006/04/20 22:30:49 hristov
55 Coding conventions (Annalisa)
57 Revision 1.5 2006/04/16 22:29:05 hristov
58 Coding conventions (Annalisa)
60 Revision 1.4 2006/04/05 08:35:38 hristov
61 Coding conventions (S.Arcelli, C.Zampolli)
63 Revision 1.3 2006/03/31 13:49:07 arcelli
64 Removing some junk printout
66 Revision 1.2 2006/03/31 11:26:30 arcelli
67 changing CDB Ids according to standard convention
69 Revision 1.1 2006/03/28 14:54:48 arcelli
70 class for TOF alignment
72 author: Silvia Arcelli, arcelli@bo.infn.it
75 /////////////////////////////////////////////////////////
77 // Class for alignment procedure //
81 /////////////////////////////////////////////////////////
85 #include "TGeoMatrix.h"
89 #include "TGeoManager.h"
90 #include "TGeoVolume.h"
93 #include "TGeoPhysicalNode.h"
95 #include "TObjString.h"
98 //#include "AliAlignObj.h"
99 #include "AliAlignObjParams.h"
100 #include "AliAlignObjMatrix.h"
101 #include "AliCDBManager.h"
102 #include "AliCDBMetaData.h"
103 #include "AliCDBId.h"
104 #include "AliCDBEntry.h"
105 #include "AliTOFAlignment.h"
106 #include "AliSurveyObj.h"
107 #include "AliSurveyPoint.h"
109 ClassImp(AliTOFAlignment)
111 const Double_t AliTOFAlignment::fgkRorigTOF = 384.5; // Mean Radius of the TOF ext. volume, cm
112 const Double_t AliTOFAlignment::fgkX1BTOF = 124.5; //x1 size of BTOF
113 const Double_t AliTOFAlignment::fgkX2BTOF = 134.7262; //x2 size of BTOF
114 const Double_t AliTOFAlignment::fgkYBTOF = 747.2; //y size of BTOF
115 const Double_t AliTOFAlignment::fgkZBTOF = 29.0; //z size of BTOF
116 const Double_t AliTOFAlignment::fgkXFM = 38.0; //x pos of FM in BTOF, cm
117 const Double_t AliTOFAlignment::fgkYFM = 457.3; //y pos of FM in BTOF, cm
118 const Double_t AliTOFAlignment::fgkZFM = 11.2; //z pos of FM in BTOF, cm
120 //_____________________________________________________________________________
121 AliTOFAlignment::AliTOFAlignment():
122 TTask("AliTOFAlignment",""),
125 fTOFAlignObjArray(0x0)
127 //AliTOFalignment main Ctor
128 for(Int_t i=0; i<18;i++)
129 for(Int_t j=0; j<5; j++)
131 for(Int_t i=0; i<72; i++)
132 for (Int_t j=0; j<6; j++)
135 for(Int_t i=0; i<18;i++)
139 //_____________________________________________________________________________
140 AliTOFAlignment::AliTOFAlignment(const AliTOFAlignment &t):
142 fNTOFAlignObj(t.fNTOFAlignObj),
144 fTOFAlignObjArray(t.fTOFAlignObjArray)
146 //AliTOFAlignment copy Ctor
148 //AliTOFalignment main Ctor
149 for(Int_t i=0; i<18;i++)
150 for(Int_t j=0; j<5; j++)
151 fNFMforSM[i][j]=t.fNFMforSM[i][j];
152 for(Int_t i=0; i<72; i++)
153 for (Int_t j=0; j<6; j++)
154 fCombFMData[i][j]=t.fCombFMData[i][j];
156 for(Int_t i=0; i<18;i++)
157 fTOFMatrixId[i]=t.fTOFMatrixId[i];
160 //_____________________________________________________________________________
161 AliTOFAlignment& AliTOFAlignment::operator=(const AliTOFAlignment &t){
162 //AliTOFAlignment assignment operator
168 fNTOFAlignObj=t.fNTOFAlignObj;
170 fTOFAlignObjArray=t.fTOFAlignObjArray;
171 for(Int_t i=0; i<18;i++)
172 fTOFMatrixId[i]=t.fTOFMatrixId[i];
177 //_____________________________________________________________________________
178 AliTOFAlignment::~AliTOFAlignment() {
179 delete fTOFAlignObjArray;
183 //_____________________________________________________________________________
184 void AliTOFAlignment::Smear(Float_t * const tr, Float_t * const rot)
186 //Introduce Random Offset/Tilts
187 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
188 Float_t dx, dy, dz; // shifts
189 Float_t dpsi, dtheta, dphi; // angular displacements
190 TRandom *rnd = new TRandom(1567);
193 AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
194 UShort_t iIndex=0; //dummy volume index
195 // AliGeomManager::ELayerID iLayer = AliGeomManager::kTOF;
196 // Int_t iIndex=1; //dummy volume index
197 UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity
200 const Int_t kSize=100;
202 for (i = 0; i<nSMTOF ; i++) {
203 snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
205 dx = (rnd->Gaus(0.,1.))*tr[0];
206 dy = (rnd->Gaus(0.,1.))*tr[1];
207 dz = (rnd->Gaus(0.,1.))*tr[2];
211 AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
212 fTOFAlignObjArray->Add(o);
215 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
216 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
220 //_____________________________________________________________________________
221 void AliTOFAlignment::Align(Float_t * const tr, Float_t * const rot)
223 //Introduce Offset/Tilts
225 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
226 Float_t dx, dy, dz; // shifts
227 Float_t dpsi, dtheta, dphi; // angular displacements
231 AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
232 UShort_t iIndex=0; //dummy volume index
233 UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity
235 const Int_t kSize=100;
238 for (i = 0; i<nSMTOF ; i++) {
240 snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
248 AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
249 fTOFAlignObjArray->Add(o);
251 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
252 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
254 //_____________________________________________________________________________
255 void AliTOFAlignment::WriteParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
257 //Write Align Par on CDB
258 AliCDBManager *man = AliCDBManager::Instance();
259 const Char_t *sel1 = "AlignPar" ;
260 const Int_t kSize=100;
262 snprintf(out,kSize,"%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);
269 //_____________________________________________________________________________
270 void AliTOFAlignment::ReadParFromCDB(const Char_t *sel, Int_t nrun)
272 //Read Align Par from CDB
273 AliCDBManager *man = AliCDBManager::Instance();
274 const Char_t *sel1 = "AlignPar" ;
275 const Int_t kSize=100;
278 snprintf(out,kSize,"%s/%s",sel,sel1);
279 AliCDBEntry *entry = man->Get(out,nrun);
281 AliError(Form("Failed to get entry: %s",out));
284 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
285 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
286 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
289 //_____________________________________________________________________________
290 void AliTOFAlignment::WriteSimParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
292 //Write Sim Align Par on CDB
293 AliCDBManager *man = AliCDBManager::Instance();
294 const Char_t *sel1 = "AlignSimPar" ;
295 const Int_t kSize=100;
297 snprintf(out,kSize,"%s/%s",sel,sel1);
298 AliCDBId idTOFAlign(out,minrun,maxrun);
299 AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
300 mdTOFAlign->SetResponsible("TOF");
301 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
302 man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
304 //_____________________________________________________________________________
305 void AliTOFAlignment::ReadSimParFromCDB(const Char_t *sel, Int_t nrun){
306 //Read Sim Align Par from CDB
307 AliCDBManager *man = AliCDBManager::Instance();
308 const Char_t *sel1 = "AlignSimPar" ;
309 const Int_t kSize=100;
311 snprintf(out,kSize,"%s/%s",sel,sel1);
312 AliCDBEntry *entry = man->Get(out,nrun);
314 AliError(Form("Failed to get entry: %s",out));
317 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
318 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
319 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
322 //_____________________________________________________________________________
323 void AliTOFAlignment::WriteOnCDBforDC()
325 //Write Align Par on CDB for DC06
326 AliCDBManager *man = AliCDBManager::Instance();
327 AliCDBId idTOFAlign("TOF/Align/Data",0,0);
328 AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
329 mdTOFAlign->SetComment("Alignment objects for ideal geometry, i.e. applying them to TGeo has to leave geometry unchanged");
330 mdTOFAlign->SetResponsible("TOF");
331 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
332 man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
334 //_____________________________________________________________________________
335 void AliTOFAlignment::ReadFromCDBforDC()
337 //Read Sim Align Par from CDB for DC06
338 AliCDBManager *man = AliCDBManager::Instance();
339 AliCDBEntry *entry = man->Get("TOF/Align/Data",0);
340 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
341 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
342 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
346 //_____________________________________________________________________________
347 void AliTOFAlignment::BuildGeomForSurvey()
350 //Generates the ideal TOF structure with four Fiducial Marks in each
351 //supermodule (two on each z side) in their expected position.
354 fTOFmgr = new TGeoManager("Geom","survey to alignment for TOF");
355 TGeoMedium *medium = 0;
356 TGeoVolume *top = fTOFmgr->MakeBox("TOP",medium,1000,1000,1000);
357 fTOFmgr->SetTopVolume(top);
358 // make shape components:
359 // This is the BTOF containing the FTOA
360 TGeoTrd1 *strd1 = new TGeoTrd1(fgkX1BTOF*0.5,fgkX2BTOF*0.5, fgkYBTOF*0.5,fgkZBTOF*0.5);
361 TGeoVolume* trd1[18];
363 // Now four fiducial marks on SM, expressed in local coordinates
364 // They are positioned at x=+/- 38 cm, y=+/- 457.3 cm, z=11.2 cm
366 TGeoBBox *fmbox = new TGeoBBox(1,1,1);
367 TGeoVolume* fm = new TGeoVolume("FM",fmbox);
371 TGeoTranslation* mAtr = new TGeoTranslation("mAtr",-fgkXFM, -fgkYFM ,fgkZFM);
372 TGeoTranslation* mBtr = new TGeoTranslation("mBtr",fgkXFM, -fgkYFM ,fgkZFM );
373 TGeoTranslation* mCtr = new TGeoTranslation("mCtr",fgkXFM, fgkYFM ,fgkZFM );
374 TGeoTranslation* mDtr = new TGeoTranslation("mDtr",-fgkXFM, fgkYFM ,fgkZFM );
376 // position all this stuff in the global ALICE frame
378 const Int_t kSize=100;
383 Float_t smR = fgkRorigTOF;
384 for (Int_t iSM = 0; iSM < 18; iSM++) {
385 Int_t mod = iSM + 13;
386 if (mod > 17) mod -= 18;
387 snprintf(name,kSize, "BTOF%d",mod);
388 trd1[iSM] = new TGeoVolume(name,strd1);
389 Float_t phi = iSM * 20.;
390 Float_t phi2 = 270 + phi;
391 if (phi2 >= 360.) phi2 -= 360.;
392 smX = TMath::Sin(phi*TMath::Pi()/180.)*smR;
393 smY = -TMath::Cos(phi*TMath::Pi()/180.)*smR;
395 TGeoRotation* bTOFRot = new TGeoRotation("bTOFRot",phi,90,0.);
396 TGeoCombiTrans trans = *(new TGeoCombiTrans(smX,smY,smZ, bTOFRot));
397 TGeoMatrix* id = new TGeoHMatrix();
398 TGeoHMatrix transMat = *id * trans;
399 TGeoHMatrix *smTrans = new TGeoHMatrix(transMat);
401 trd1[iSM]->AddNode(fm,1,mAtr); //place FM in BTOF
402 trd1[iSM]->AddNode(fm,2,mBtr);
403 trd1[iSM]->AddNode(fm,3,mCtr);
404 trd1[iSM]->AddNode(fm,4,mDtr);
405 top->AddNode(trd1[iSM],1,smTrans); //place BTOF_iSM in ALICE
406 trd1[iSM]->SetVisDaughters();
407 trd1[iSM]->SetLineColor(iSM); //black
411 fTOFmgr->CloseGeometry();
412 fTOFmgr->GetTopVolume()->Draw();
413 fTOFmgr->SetVisOption(0);
414 fTOFmgr->SetVisLevel(6);
416 // Now Store the "Ideal" Global Matrices (local to global) for later use
418 for (Int_t iSM = 0; iSM < 18; iSM++) {
420 snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
421 printf("\n\n***************** TOF SuperModule: %s ****************** \n",name);
422 TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
423 fTOFMatrixId[iSM] = pn3->GetMatrix(); //save "ideal" global matrix
424 printf("\n\n*************** The Ideal Matrix in GRS *****************\n");
425 fTOFMatrixId[iSM]->Print();
430 //_____________________________________________________________________________
431 void AliTOFAlignment::InsertMisAlignment(Float_t * const mis)
433 // Now Apply the Displacements and store the misaligned FM positions...
437 Double_t lA[3]={-fgkXFM, -fgkYFM ,fgkZFM};
438 Double_t lB[3]={fgkXFM, -fgkYFM ,fgkZFM};
439 Double_t lC[3]={fgkXFM, fgkYFM ,fgkZFM};
440 Double_t lD[3]={-fgkXFM, fgkYFM ,fgkZFM};
442 const Int_t kSize=16;
445 for(Int_t iSM=0;iSM<18;iSM++){
446 snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
448 printf("\n\n******Misaligning TOF SuperModule ************** %s \n",name);
450 // ************* get ideal global matrix *******************
451 TGeoHMatrix g3 = *fTOFmgr->GetCurrentMatrix();
452 AliInfo(Form("This is the ideal global trasformation of SM %i",iSM));
453 g3.Print(); // g3 is the local(BTOF) to global (ALICE) matrix and is the same of fTOFMatrixId
454 TGeoNode* n3 = fTOFmgr->GetCurrentNode();
455 TGeoMatrix* l3 = n3->GetMatrix();
457 Double_t gA[3], gB[3], gC[3], gD[3]; // ideal global FM point coord.
458 g3.LocalToMaster(lA,gA);
459 g3.LocalToMaster(lB,gB);
460 g3.LocalToMaster(lC,gC);
461 g3.LocalToMaster(lD,gD);
463 // We apply a delta transformation to the surveyed vol to represent
464 // its real position, given below by ng3 nl3, which differs from its
465 // ideal position saved above in g3 and l3
467 //we have to express the displacements as regards the old local RS (non misaligned BTOF)
468 Double_t dx = mis[0]; // shift along x
469 Double_t dy = mis[1]; // shift along y
470 Double_t dz = mis[2]; // shift along z
471 Double_t dphi = mis[3]; // rot around z
472 Double_t dtheta = mis[4]; // rot around x'
473 Double_t dpsi = mis[5]; // rot around z''
475 TGeoRotation* rrot = new TGeoRotation("rot",dphi,dtheta,dpsi);
476 TGeoCombiTrans localdelta = *(new TGeoCombiTrans(dx,dy,dz, rrot));
477 AliInfo(Form("This is the local delta trasformation for SM %i \n",iSM));
479 TGeoHMatrix nlocal = *l3 * localdelta;
480 TGeoHMatrix* nl3 = new TGeoHMatrix(nlocal); // new matrix, representing real position (from new local mis RS to the global one)
482 TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
486 TGeoHMatrix* ng3 = pn3->GetMatrix(); //"real" global matrix, what survey sees
487 printf("\n\n************* The Misaligned Matrix in GRS **************\n");
489 Double_t ngA[3], ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
490 ng3->LocalToMaster(lA,ngA);
491 ng3->LocalToMaster(lB,ngB);
492 ng3->LocalToMaster(lC,ngC);
493 ng3->LocalToMaster(lD,ngD);
495 for(Int_t coord=0;coord<3;coord++){
496 fCombFMData[iSM*4][2*coord]=ngA[coord];
497 fCombFMData[iSM*4][2*coord+1]=1;
498 fCombFMData[iSM*4+1][2*coord]=ngB[coord];
499 fCombFMData[iSM*4+1][2*coord+1]=1;
500 fCombFMData[iSM*4+2][2*coord]=ngC[coord];
501 fCombFMData[iSM*4+2][2*coord+1]=1;
502 fCombFMData[iSM*4+3][2*coord]=ngD[coord];
503 fCombFMData[iSM*4+3][2*coord+1]=1;
509 //____________________________________________________________________________
510 void AliTOFAlignment::WriteCombData(const Char_t *nomefile, Int_t option)
512 // 1 for simulated data; 0 for data from survey file
513 // write combined data on a file
517 /* Open file in text mode: */
518 if( (data = fopen( nomefile, "w+t" )) != NULL ){
520 fprintf( data, "simulated data\n" );} else {
521 fprintf( data, "survey data\n" );}
523 fprintf( data, "data from InsertMisAlignmentBTOF method\n");}
524 else {fprintf( data, "real survey data from text file (coordinate in global RS)\n");}
525 fprintf( data, "Point Name,XPH,YPH,ZPH,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
526 fprintf( data, "> Data:\n");
527 for(Int_t i=0;i<72;i++){
528 if (fCombFMData[i][0]!=0){
529 fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1]*10,fCombFMData[i][3]*10,fCombFMData[i][5]*10);
535 printf( "Problem opening the file\n" );
541 //____________________________________________________________________________
542 void AliTOFAlignment::WriteSimSurveyData(const Char_t *nomefile)
544 // write sim data in standard format
549 /* Open file in text mode: */
550 if( (data = fopen( nomefile, "w+t" )) != NULL )
552 fprintf( data, "> Title:\n" );
553 fprintf( data, "simulated data\n" );
554 fprintf( data, "> Date:\n" );
555 fprintf( data, "24.09.2007\n" );
556 fprintf( data, "> Subdetector:\n" );
557 fprintf( data, "TOF\n" );
558 fprintf( data, "> Report URL:\n" );
559 fprintf( data, "https://edms.cern.ch/document/835615\n" );
560 fprintf( data, "> Version:\n" );
561 fprintf( data, "1\n");
562 fprintf( data, "> General Observations:\n");
563 fprintf( data, "data from InsertMisAlignmentBTOF method\n");
564 fprintf( data, "> Coordinate System:\n");
565 fprintf( data, "\\ALICEPH\n");
566 fprintf( data, "> Units:\n");
567 fprintf( data, "cm\n");
568 fprintf( data, "> Nr Columns:\n");
569 fprintf( data, "9\n");
570 fprintf( data, "> Column Names:\n");
571 fprintf( data, "Point Name,XPH,YPH,ZPH,Point Type,Target Used,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
572 fprintf( data, "> Data:\n");
573 for(Int_t i=0;i<72;i++)
574 if (fCombFMData[i][0]!=0)
575 fprintf( data, "SM%02iFM%i %f %f %f M Y %f %f %f\n", (i-i%4)/4, i%4, fCombFMData[i][0],fCombFMData[i][2],fCombFMData[i][4],fCombFMData[i][1],fCombFMData[i][3],fCombFMData[i][5]);
580 printf( "Problem opening the file\n" );
583 //____________________________________________________________________________
584 void AliTOFAlignment::MakeDefData(const Int_t nf,TString namefiles[])
586 //this method combines survey data from different files (namefiles[])
590 Float_t data[72][6][100];
591 for (Int_t i=0;i<72;i++)
592 for (Int_t j=0; j<6; j++)
593 for(Int_t k=0; k<100; k++)
597 Long64_t totdata[72]={0};
599 for (Int_t ii=0;ii<nf; ii++)
601 AliSurveyObj *so = new AliSurveyObj();
602 const Char_t *nome=namefiles[ii];
603 so->FillFromLocalFile(nome);
604 TObjArray *points = so->GetData();
605 Int_t nSurveyPoint=points->GetEntries();
606 for(Int_t jj=0;jj<nSurveyPoint;jj++){
607 const char* pointName= ((AliSurveyPoint *) points->At(jj))->GetPointName().Data();
608 nfm=atoi(&pointName[6]);
609 nsm=atoi(&pointName[2]);
610 data[nsm*4+nfm][0][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetX();
611 data[nsm*4+nfm][2][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetY();
612 data[nsm*4+nfm][4][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetZ();
613 data[nsm*4+nfm][1][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionX();
614 data[nsm*4+nfm][3][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionY();
615 data[nsm*4+nfm][5][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionZ();
616 totdata[nsm*4+nfm]=totdata[nsm*4+nfm]+1;
622 for(Int_t i=0; i<72 ;i++){
623 Float_t numx=0, numy=0,numz=0, comodox=0, comodoy=0, comodoz=0,denx=0, deny=0, denz=0;
625 for(Int_t j=0; j<totdata[i]; j++){
626 comodox=1/(data[i][1][j]/10*data[i][1][j]/10);//precision in mm, position in cm
627 numx=numx+data[i][0][j]*comodox;
629 comodoy=1/(data[i][3][j]/10*data[i][3][j]/10);
630 numy=numy+data[i][2][j]*comodoy;
632 comodoz=1/(data[i][5][j]/10*data[i][5][j]/10);
633 numz=numz+data[i][4][j]*comodoz;
636 fCombFMData[i][1]=TMath::Sqrt(1/denx); //error for x position
637 fCombFMData[i][3]=TMath::Sqrt(1/deny); //error for y position
638 fCombFMData[i][5]=TMath::Sqrt(1/denz); //error for z position
639 fCombFMData[i][0]=numx/denx; //combined survey data for x position of FM
640 fCombFMData[i][2]=numy/deny; //combined survey data for y position of FM
641 fCombFMData[i][4]=numz/denz; //combined survey data for z position of FM
645 for(Int_t i=0;i<72;i++)
646 if (fCombFMData[i][0]!=0){
647 fNFMforSM[(i-i%4)/4][i%4]=1;
648 fNFMforSM[(i-i%4)/4][4]=fNFMforSM[(i-i%4)/4][4]+1;
652 //_____________________________________________________________________________
653 void AliTOFAlignment::ReadSurveyDataAndAlign(){
655 // read the survey data and, if we know the positions of at least 3 FM
656 //for a SM, call the right Alignement procedure
658 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
660 Float_t deltaFM0=0, deltaFM1=0, deltaFM2=0, deltaFM3=0;
662 for(Int_t i=0; i<18; i++){
663 switch(fNFMforSM[i][4]){
665 printf("we don't know the position of any FM of SM %i\n",i);
668 printf("we know the position of only one FM for SM %i\n",i);
672 printf("we know the position of only 2 FM for SM %i\n",i);
676 if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1){
677 printf("we know the position of FM A B C for SM %i\n",i);
678 AliTOFAlignment::AlignFromSurveyABC(i);};
681 if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][3]==1){
682 printf("we know the position of FM A B D for SM %i\n",i);
683 AliTOFAlignment::AlignFromSurveyABD(i);};
686 if (fNFMforSM[i][0]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
687 printf("we know the position of FM A C D for SM %i\n",i);
688 AliTOFAlignment::AlignFromSurveyACD(i);};
691 if (fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
692 printf("we know the position of FM B C D for SM %i\n",i);
693 AliTOFAlignment::AlignFromSurveyBCD(i);};
698 printf("we know the position of all the 4 FM for SM %i\n",i);
699 //check the precision of the measurement
701 deltaFM0=fCombFMData[i*4][1]/TMath::Abs(fCombFMData[i*4][0])+fCombFMData[i*4][3]/TMath::Abs(fCombFMData[i*4][2])+fCombFMData[i*4][5]/TMath::Abs(fCombFMData[i*4][4]);
702 deltaFM1=fCombFMData[i*4+1][1]/TMath::Abs(fCombFMData[i*4+1][0])+fCombFMData[i*4+1][3]/TMath::Abs(fCombFMData[i*4+1][2])+fCombFMData[i*4+1][5]/TMath::Abs(fCombFMData[i*4+1][4]);
703 deltaFM2=fCombFMData[i*4+2][1]/TMath::Abs(fCombFMData[i*4+2][0])+fCombFMData[i*4+2][3]/TMath::Abs(fCombFMData[i*4+2][2])+fCombFMData[i*4+2][5]/TMath::Abs(fCombFMData[i*4+2][4]);
704 deltaFM3=fCombFMData[i*4+3][1]/TMath::Abs(fCombFMData[i*4+3][0])+fCombFMData[i*4+3][3]/TMath::Abs(fCombFMData[i*4+3][2])+fCombFMData[i*4+3][5]/TMath::Abs(fCombFMData[i*4+3][4]);
706 //to AlignFromSurvey we use the 3 FM whose positions are known with greatest precision
707 if(deltaFM0>=deltaFM1 && deltaFM0>=deltaFM2 && deltaFM0>=deltaFM3){
708 printf("to Align we use FM B,C,D");
709 AliTOFAlignment::AlignFromSurveyBCD(i);} else
710 if(deltaFM1>=deltaFM0 && deltaFM1>=deltaFM2 && deltaFM1>=deltaFM3){
711 printf("to Align we use FM A,C,D");
712 AliTOFAlignment::AlignFromSurveyACD(i);} else
713 if(deltaFM2>=deltaFM0 && deltaFM2>=deltaFM1 && deltaFM2>=deltaFM3){
714 printf("to Align we use FM A,B,D");
715 AliTOFAlignment::AlignFromSurveyABD(i);} else{
716 printf("to Align we use FM A,B,C");
717 AliTOFAlignment::AlignFromSurveyABC(i);}
724 // saving TOF AligObjs from survey on a file, for the moment..
725 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
726 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
727 TFile f("TOFAlignFromSurvey.root","RECREATE");
729 f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
735 //_____________________________________________________________________________
736 void AliTOFAlignment::AlignFromSurveyABC(Int_t iSM)
739 //From Survey data, derive the needed transformations to get the
741 //Again, highly "inspired" to Raffaele's example...
744 Double_t ngA[3], ngB[3], ngC[3]; // real FM point coord., global RS
745 // Get the 'realistic' input from the Survey Matrix
746 for(Int_t coord=0;coord<3;coord++){
747 ngA[coord]= fCombFMData[iSM*4][coord*2];
748 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
749 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
752 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
754 // From the real fiducial marks coordinates derive back the
755 // new global position of the surveyed volume
756 //*** What follows is the actual survey-to-alignment procedure
758 Double_t ab[3], bc[3], n[3];
759 Double_t plane[4], s=1.;
761 // first vector on the plane of the fiducial marks
762 for(Int_t i=0;i<3;i++){
763 ab[i] = (ngB[i] - ngA[i]);
766 // second vector on the plane of the fiducial marks
767 for(Int_t i=0;i<3;i++){
768 bc[i] = (ngC[i] - ngB[i]);
771 // vector normal to the plane of the fiducial marks obtained
772 // as cross product of the two vectors on the plane d0^d1
773 n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
774 n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
775 n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
777 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
779 s = Double_t(1.)/sizen ; //normalization factor
781 AliInfo("Problem in normalizing the vector");
784 // plane expressed in the hessian normal form, see:
785 // http://mathworld.wolfram.com/HessianNormalForm.html
786 // the first three are the coordinates of the orthonormal vector
787 // the fourth coordinate is equal to the distance from the origin
789 for(Int_t i=0;i<3;i++){
792 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
794 // The center of the square with fiducial marks as corners
795 // as the middle point of one diagonal - md
796 // Used below to get the center - orig - of the surveyed box
798 Double_t orig[3], md[3];
799 for(Int_t i=0;i<3;i++){
800 md[i] = (ngA[i] + ngC[i]) * 0.5;
803 // The center of the box, gives the global translation
804 for(Int_t i=0;i<3;i++){
805 orig[i] = md[i] - plane[i]*fgkZFM;
808 // get local directions needed to write the global rotation matrix
809 // for the surveyed volume by normalising vectors ab and bc
810 Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
814 for(Int_t i=0;i<3;i++){
818 Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
820 for(Int_t i=0;i<3;i++){
824 Double_t rot[9] = {ab[0],bc[0],plane[0],ab[1],bc[1],plane[1],ab[2],bc[2],plane[2]}; // the rotation matrix
825 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey
827 ng.SetTranslation(orig);
829 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
832 // Calculate the delta transformation wrt Ideal geometry
833 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
835 printf("\n\n**** The ideal matrix ***\n");
836 fTOFMatrixId[iSM]->Print();
838 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
839 printf("\n\n**** The inverse of the ideal matrix ***\n");
842 gdelta.MultiplyLeft(&ng);
843 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
844 gdelta.Print(); //this is the global delta trasformation
846 // Now Write the Alignment Objects....
847 Int_t index=0; //let all SM modules have index=0
848 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
849 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
850 TString symname(Form("TOF/sm%02d",iSM));
851 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
852 fTOFAlignObjArray->Add(o);
857 //_____________________________________________________________________________
858 void AliTOFAlignment::AlignFromSurveyABD(Int_t iSM)
861 //From Survey data, derive the needed transformations to get the
863 //Again, highly "inspired" to Raffaele's example...
866 Double_t ngA[3], ngB[3], ngD[3];// real FM point coord., global RS
868 // Get the 'realistic' input from the Survey Matrix
869 for(Int_t coord=0;coord<3;coord++){
870 ngA[coord]= fCombFMData[iSM*4][coord*2];
871 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
872 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
875 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
877 // From the new fiducial marks coordinates derive back the
878 // new global position of the surveyed volume
879 //*** What follows is the actual survey-to-alignment procedure
881 Double_t ab[3], ad[3], n[3];
882 Double_t plane[4], s=1.;
884 // first vector on the plane of the fiducial marks
885 for(Int_t i=0;i<3;i++){
886 ab[i] = (ngB[i] - ngA[i]);
889 // second vector on the plane of the fiducial marks
890 for(Int_t i=0;i<3;i++){
891 ad[i] = (ngD[i] - ngA[i]);
894 // vector normal to the plane of the fiducial marks obtained
895 // as cross product of the two vectors on the plane d0^d1
896 n[0] = (ab[1] * ad[2] - ab[2] * ad[1]);
897 n[1] = (ab[2] * ad[0] - ab[0] * ad[2]);
898 n[2] = (ab[0] * ad[1] - ab[1] * ad[0]);
900 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
902 s = Double_t(1.)/sizen ; //normalization factor
904 AliInfo("Problem in normalizing the vector");
907 // plane expressed in the hessian normal form, see:
908 // http://mathworld.wolfram.com/HessianNormalForm.html
909 // the first three are the coordinates of the orthonormal vector
910 // the fourth coordinate is equal to the distance from the origin
912 for(Int_t i=0;i<3;i++){
915 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
917 // The center of the square with fiducial marks as corners
918 // as the middle point of one diagonal - md
919 // Used below to get the center - orig - of the surveyed box
921 Double_t orig[3], md[3];
922 for(Int_t i=0;i<3;i++){
923 md[i] = (ngB[i] + ngD[i]) * 0.5;
926 // The center of the box, gives the global translation
927 for(Int_t i=0;i<3;i++){
928 orig[i] = md[i] - plane[i]*fgkZFM;
931 // get local directions needed to write the global rotation matrix
932 // for the surveyed volume by normalising vectors ab and bc
933 Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
935 for(Int_t i=0;i<3;i++){
939 Double_t sy = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
941 for(Int_t i=0;i<3;i++){
945 Double_t rot[9] = {ab[0],ad[0],plane[0],ab[1],ad[1],plane[1],ab[2],ad[2],plane[2]};
946 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
948 ng.SetTranslation(orig);
950 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
953 // Calculate the delta transformation wrt Ideal geometry
954 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
956 printf("\n\n**** The ideal matrix ***\n");
957 fTOFMatrixId[iSM]->Print();
959 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
960 printf("\n\n**** The inverse of the ideal matrix ***\n");
963 gdelta.MultiplyLeft(&ng);
964 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
965 gdelta.Print(); //global delta trasformation
967 // Now Write the Alignment Objects....
968 Int_t index=0; //let all SM modules have index=0
969 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
970 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
971 TString symname(Form("TOF/sm%02d",iSM));
972 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
973 fTOFAlignObjArray->Add(o);
976 //_____________________________________________________________________________
977 void AliTOFAlignment::AlignFromSurveyACD(Int_t iSM)
979 //From Survey data, derive the needed transformations to get the
981 //Again, highly "inspired" to Raffaele's example...
985 Double_t ngA[3], ngC[3], ngD[3];// real FM point coord., global RS
987 // Get the 'realistic' input from the Survey Matrix
988 for(Int_t coord=0;coord<3;coord++){
989 ngA[coord]= fCombFMData[iSM*4][coord*2];
990 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
991 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
994 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
996 // From the new fiducial marks coordinates derive back the
997 // new global position of the surveyed volume
998 //*** What follows is the actual survey-to-alignment procedure
1000 Double_t cd[3], ad[3], n[3];
1001 Double_t plane[4], s=1.;
1003 // first vector on the plane of the fiducial marks
1004 for(Int_t i=0;i<3;i++){
1005 cd[i] = (ngC[i] - ngD[i]);
1008 // second vector on the plane of the fiducial marks
1009 for(Int_t i=0;i<3;i++){
1010 ad[i] = (ngD[i] - ngA[i]);
1013 // vector normal to the plane of the fiducial marks obtained
1014 // as cross product of the two vectors on the plane d0^d1
1015 n[0] = (ad[1] * cd[2] - ad[2] * cd[1]);
1016 n[1] = (ad[2] * cd[0] - ad[0] * cd[2]);
1017 n[2] = (ad[0] * cd[1] - ad[1] * cd[0]);
1019 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1021 s = Double_t(1.)/sizen ; //normalization factor
1023 AliInfo("Problem in normalizing the vector");
1026 // plane expressed in the hessian normal form, see:
1027 // http://mathworld.wolfram.com/HessianNormalForm.html
1028 // the first three are the coordinates of the orthonormal vector
1029 // the fourth coordinate is equal to the distance from the origin
1031 for(Int_t i=0;i<3;i++){
1032 plane[i] = n[i] * s;
1034 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
1036 // The center of the square with fiducial marks as corners
1037 // as the middle point of one diagonal - md
1038 // Used below to get the center - orig - of the surveyed box
1040 Double_t orig[3], md[3];
1041 for(Int_t i=0;i<3;i++){
1042 md[i] = (ngA[i] + ngC[i]) * 0.5;
1045 // The center of the box, gives the global translation
1046 for(Int_t i=0;i<3;i++){
1047 orig[i] = md[i] + plane[i]*fgkZFM;
1050 // get local directions needed to write the global rotation matrix
1051 // for the surveyed volume by normalising vectors ab and bc
1052 Double_t sx = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
1054 for(Int_t i=0;i<3;i++){
1058 Double_t sy = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1060 for(Int_t i=0;i<3;i++){
1064 Double_t rot[9] = {cd[0],ad[0],-plane[0],cd[1],ad[1],-plane[1],cd[2],ad[2],-plane[2]};
1065 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1067 ng.SetTranslation(orig);
1068 ng.SetRotation(rot);
1069 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1072 // Calculate the delta transformation wrt Ideal geometry
1073 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1075 printf("\n\n**** The ideal matrix ***\n");
1076 fTOFMatrixId[iSM]->Print();
1078 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1079 printf("\n\n**** The inverse of the ideal matrix ***\n");
1082 gdelta.MultiplyLeft(&ng);
1083 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1084 gdelta.Print(); //global delta trasformation
1086 // Now Write the Alignment Objects....
1087 Int_t index=0; //let all SM modules have index=0
1088 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1089 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
1090 TString symname(Form("TOF/sm%02d",iSM));
1091 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1092 fTOFAlignObjArray->Add(o);
1095 //___________________________________________________________________________
1096 void AliTOFAlignment::AlignFromSurveyBCD(Int_t iSM)
1098 //From Survey data, derive the needed transformations to get the
1099 //Alignment Objects.
1100 //Again, highly "inspired" to Raffaele's example...
1103 Double_t ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
1106 // Get the 'realistic' input from the Survey Matrix
1107 for(Int_t coord=0;coord<3;coord++){
1108 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
1109 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
1110 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
1113 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
1115 // From the new fiducial marks coordinates derive back the
1116 // new global position of the surveyed volume
1117 //*** What follows is the actual survey-to-alignment procedure
1119 Double_t cd[3], bc[3], n[3];
1120 Double_t plane[4], s=1.;
1122 // first vector on the plane of the fiducial marks
1123 for(Int_t i=0;i<3;i++){
1124 cd[i] = (ngC[i] - ngD[i]);
1127 // second vector on the plane of the fiducial marks
1128 for(Int_t i=0;i<3;i++){
1129 bc[i] = (ngC[i] - ngB[i]);
1132 // vector normal to the plane of the fiducial marks obtained
1133 // as cross product of the two vectors on the plane d0^d1
1134 n[0] = (bc[1] * cd[2] - bc[2] * cd[1]);
1135 n[1] = (bc[2] * cd[0] - bc[0] * cd[2]);
1136 n[2] = (bc[0] * cd[1] - bc[1] * cd[0]);
1138 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1140 s = Double_t(1.)/sizen ; //normalization factor
1142 AliInfo("Problem in normalizing the vector");
1145 // plane expressed in the hessian normal form, see:
1146 // http://mathworld.wolfram.com/HessianNormalForm.html
1147 // the first three are the coordinates of the orthonormal vector
1148 // the fourth coordinate is equal to the distance from the origin
1150 for(Int_t i=0;i<3;i++){
1151 plane[i] = n[i] * s;
1153 plane[3] = ( plane[0] * ngB[0] + plane[1] * ngB[1] + plane[2] * ngB[2] );
1155 // The center of the square with fiducial marks as corners
1156 // as the middle point of one diagonal - md
1157 // Used below to get the center - orig - of the surveyed box
1159 Double_t orig[3], md[3];
1160 for(Int_t i=0;i<3;i++){
1161 md[i] = (ngB[i] + ngD[i]) * 0.5;
1164 // The center of the box, gives the global translation
1165 for(Int_t i=0;i<3;i++){
1166 orig[i] = md[i] + plane[i]*fgkZFM;
1169 // get local directions needed to write the global rotation matrix
1170 // for the surveyed volume by normalising vectors ab and bc
1171 Double_t sx = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1173 for(Int_t i=0;i<3;i++){
1177 Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
1179 for(Int_t i=0;i<3;i++){
1183 Double_t rot[9] = {cd[0],bc[0],-plane[0],cd[1],bc[1],-plane[1],cd[2],bc[2],-plane[2]};
1184 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1186 ng.SetTranslation(orig);
1187 ng.SetRotation(rot);
1188 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1191 // Calculate the delta transformation wrt Ideal geometry
1192 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1194 printf("\n\n**** The ideal matrix ***\n");
1195 fTOFMatrixId[iSM]->Print();
1197 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1198 printf("\n\n**** The inverse of the ideal matrix ***\n");
1201 gdelta.MultiplyLeft(&ng);
1202 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1203 gdelta.Print(); //global delta trasformation
1205 // Now Write the Alignment Objects....
1206 Int_t index=0; //let all SM modules have index=0
1207 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1208 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
1209 TString symname(Form("TOF/sm%02d",iSM));
1210 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1211 fTOFAlignObjArray->Add(o);