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 //_____________________________________________________________________________
136 AliTOFAlignment::AliTOFAlignment(const AliTOFAlignment &t):
138 fNTOFAlignObj(t.fNTOFAlignObj),
140 fTOFAlignObjArray(t.fTOFAlignObjArray)
142 //AliTOFAlignment copy Ctor
144 //AliTOFalignment main Ctor
145 for(Int_t i=0; i<18;i++)
146 for(Int_t j=0; j<5; j++)
147 fNFMforSM[i][j]=t.fNFMforSM[i][j];
148 for(Int_t i=0; i<72; i++)
149 for (Int_t j=0; j<6; j++)
150 fCombFMData[i][j]=t.fCombFMData[i][j];
152 //_____________________________________________________________________________
153 AliTOFAlignment& AliTOFAlignment::operator=(const AliTOFAlignment &t){
154 //AliTOFAlignment assignment operator
160 fNTOFAlignObj=t.fNTOFAlignObj;
162 fTOFAlignObjArray=t.fTOFAlignObjArray;
166 //_____________________________________________________________________________
167 AliTOFAlignment::~AliTOFAlignment() {
168 delete fTOFAlignObjArray;
172 //_____________________________________________________________________________
173 void AliTOFAlignment::Smear(Float_t * const tr, Float_t * const rot)
175 //Introduce Random Offset/Tilts
176 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
177 Float_t dx, dy, dz; // shifts
178 Float_t dpsi, dtheta, dphi; // angular displacements
179 TRandom *rnd = new TRandom(1567);
182 AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
183 UShort_t iIndex=0; //dummy volume index
184 // AliGeomManager::ELayerID iLayer = AliGeomManager::kTOF;
185 // Int_t iIndex=1; //dummy volume index
186 UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity
189 const Int_t kSize=100;
191 for (i = 0; i<nSMTOF ; i++) {
192 snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
194 dx = (rnd->Gaus(0.,1.))*tr[0];
195 dy = (rnd->Gaus(0.,1.))*tr[1];
196 dz = (rnd->Gaus(0.,1.))*tr[2];
200 AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
201 fTOFAlignObjArray->Add(o);
204 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
205 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
209 //_____________________________________________________________________________
210 void AliTOFAlignment::Align(Float_t * const tr, Float_t * const rot)
212 //Introduce Offset/Tilts
214 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
215 Float_t dx, dy, dz; // shifts
216 Float_t dpsi, dtheta, dphi; // angular displacements
220 AliGeomManager::ELayerID iLayer = AliGeomManager::kInvalidLayer;
221 UShort_t iIndex=0; //dummy volume index
222 UShort_t dvoluid = AliGeomManager::LayerToVolUID(iLayer,iIndex); //dummy volume identity
224 const Int_t kSize=100;
227 for (i = 0; i<nSMTOF ; i++) {
229 snprintf(path,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",i,i);
237 AliAlignObjParams *o =new AliAlignObjParams(path, dvoluid, dx, dy, dz, dpsi, dtheta, dphi, kTRUE);
238 fTOFAlignObjArray->Add(o);
240 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
241 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
243 //_____________________________________________________________________________
244 void AliTOFAlignment::WriteParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
246 //Write Align Par on CDB
247 AliCDBManager *man = AliCDBManager::Instance();
248 const Char_t *sel1 = "AlignPar" ;
249 const Int_t kSize=100;
251 snprintf(out,kSize,"%s/%s",sel,sel1);
252 AliCDBId idTOFAlign(out,minrun,maxrun);
253 AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
254 mdTOFAlign->SetResponsible("TOF");
255 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
256 man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
258 //_____________________________________________________________________________
259 void AliTOFAlignment::ReadParFromCDB(const Char_t *sel, Int_t nrun)
261 //Read Align Par from CDB
262 AliCDBManager *man = AliCDBManager::Instance();
263 const Char_t *sel1 = "AlignPar" ;
264 const Int_t kSize=100;
267 snprintf(out,kSize,"%s/%s",sel,sel1);
268 AliCDBEntry *entry = man->Get(out,nrun);
270 AliError(Form("Failed to get entry: %s",out));
273 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
274 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
275 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
278 //_____________________________________________________________________________
279 void AliTOFAlignment::WriteSimParOnCDB(const Char_t *sel, Int_t minrun, Int_t maxrun)
281 //Write Sim Align Par on CDB
282 AliCDBManager *man = AliCDBManager::Instance();
283 const Char_t *sel1 = "AlignSimPar" ;
284 const Int_t kSize=100;
286 snprintf(out,kSize,"%s/%s",sel,sel1);
287 AliCDBId idTOFAlign(out,minrun,maxrun);
288 AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
289 mdTOFAlign->SetResponsible("TOF");
290 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
291 man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
293 //_____________________________________________________________________________
294 void AliTOFAlignment::ReadSimParFromCDB(const Char_t *sel, Int_t nrun){
295 //Read Sim Align Par from CDB
296 AliCDBManager *man = AliCDBManager::Instance();
297 const Char_t *sel1 = "AlignSimPar" ;
298 const Int_t kSize=100;
300 snprintf(out,kSize,"%s/%s",sel,sel1);
301 AliCDBEntry *entry = man->Get(out,nrun);
303 AliError(Form("Failed to get entry: %s",out));
306 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
307 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
308 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
311 //_____________________________________________________________________________
312 void AliTOFAlignment::WriteOnCDBforDC()
314 //Write Align Par on CDB for DC06
315 AliCDBManager *man = AliCDBManager::Instance();
316 AliCDBId idTOFAlign("TOF/Align/Data",0,0);
317 AliCDBMetaData *mdTOFAlign = new AliCDBMetaData();
318 mdTOFAlign->SetComment("Alignment objects for ideal geometry, i.e. applying them to TGeo has to leave geometry unchanged");
319 mdTOFAlign->SetResponsible("TOF");
320 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
321 man->Put(fTOFAlignObjArray,idTOFAlign,mdTOFAlign);
323 //_____________________________________________________________________________
324 void AliTOFAlignment::ReadFromCDBforDC()
326 //Read Sim Align Par from CDB for DC06
327 AliCDBManager *man = AliCDBManager::Instance();
328 AliCDBEntry *entry = man->Get("TOF/Align/Data",0);
329 fTOFAlignObjArray=(TObjArray*)entry->GetObject();
330 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
331 AliInfo(Form("Number of Alignable Volumes from CDB: %d",fNTOFAlignObj));
335 //_____________________________________________________________________________
336 void AliTOFAlignment::BuildGeomForSurvey()
339 //Generates the ideal TOF structure with four Fiducial Marks in each
340 //supermodule (two on each z side) in their expected position.
343 fTOFmgr = new TGeoManager("Geom","survey to alignment for TOF");
344 TGeoMedium *medium = 0;
345 TGeoVolume *top = fTOFmgr->MakeBox("TOP",medium,1000,1000,1000);
346 fTOFmgr->SetTopVolume(top);
347 // make shape components:
348 // This is the BTOF containing the FTOA
349 TGeoTrd1 *strd1 = new TGeoTrd1(fgkX1BTOF*0.5,fgkX2BTOF*0.5, fgkYBTOF*0.5,fgkZBTOF*0.5);
350 TGeoVolume* trd1[18];
352 // Now four fiducial marks on SM, expressed in local coordinates
353 // They are positioned at x=+/- 38 cm, y=+/- 457.3 cm, z=11.2 cm
355 TGeoBBox *fmbox = new TGeoBBox(1,1,1);
356 TGeoVolume* fm = new TGeoVolume("FM",fmbox);
360 TGeoTranslation* mAtr = new TGeoTranslation("mAtr",-fgkXFM, -fgkYFM ,fgkZFM);
361 TGeoTranslation* mBtr = new TGeoTranslation("mBtr",fgkXFM, -fgkYFM ,fgkZFM );
362 TGeoTranslation* mCtr = new TGeoTranslation("mCtr",fgkXFM, fgkYFM ,fgkZFM );
363 TGeoTranslation* mDtr = new TGeoTranslation("mDtr",-fgkXFM, fgkYFM ,fgkZFM );
365 // position all this stuff in the global ALICE frame
367 const Int_t kSize=100;
372 Float_t smR = fgkRorigTOF;
373 for (Int_t iSM = 0; iSM < 18; iSM++) {
374 Int_t mod = iSM + 13;
375 if (mod > 17) mod -= 18;
376 snprintf(name,kSize, "BTOF%d",mod);
377 trd1[iSM] = new TGeoVolume(name,strd1);
378 Float_t phi = iSM * 20.;
379 Float_t phi2 = 270 + phi;
380 if (phi2 >= 360.) phi2 -= 360.;
381 smX = TMath::Sin(phi*TMath::Pi()/180.)*smR;
382 smY = -TMath::Cos(phi*TMath::Pi()/180.)*smR;
384 TGeoRotation* bTOFRot = new TGeoRotation("bTOFRot",phi,90,0.);
385 TGeoCombiTrans trans = *(new TGeoCombiTrans(smX,smY,smZ, bTOFRot));
386 TGeoMatrix* id = new TGeoHMatrix();
387 TGeoHMatrix transMat = *id * trans;
388 TGeoHMatrix *smTrans = new TGeoHMatrix(transMat);
390 trd1[iSM]->AddNode(fm,1,mAtr); //place FM in BTOF
391 trd1[iSM]->AddNode(fm,2,mBtr);
392 trd1[iSM]->AddNode(fm,3,mCtr);
393 trd1[iSM]->AddNode(fm,4,mDtr);
394 top->AddNode(trd1[iSM],1,smTrans); //place BTOF_iSM in ALICE
395 trd1[iSM]->SetVisDaughters();
396 trd1[iSM]->SetLineColor(iSM); //black
400 fTOFmgr->CloseGeometry();
401 fTOFmgr->GetTopVolume()->Draw();
402 fTOFmgr->SetVisOption(0);
403 fTOFmgr->SetVisLevel(6);
405 // Now Store the "Ideal" Global Matrices (local to global) for later use
407 for (Int_t iSM = 0; iSM < 18; iSM++) {
409 snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
410 printf("\n\n***************** TOF SuperModule: %s ****************** \n",name);
411 TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
412 fTOFMatrixId[iSM] = pn3->GetMatrix(); //save "ideal" global matrix
413 printf("\n\n*************** The Ideal Matrix in GRS *****************\n");
414 fTOFMatrixId[iSM]->Print();
419 //_____________________________________________________________________________
420 void AliTOFAlignment::InsertMisAlignment(Float_t * const mis)
422 // Now Apply the Displacements and store the misaligned FM positions...
426 Double_t lA[3]={-fgkXFM, -fgkYFM ,fgkZFM};
427 Double_t lB[3]={fgkXFM, -fgkYFM ,fgkZFM};
428 Double_t lC[3]={fgkXFM, fgkYFM ,fgkZFM};
429 Double_t lD[3]={-fgkXFM, fgkYFM ,fgkZFM};
431 const Int_t kSize=16;
434 for(Int_t iSM=0;iSM<18;iSM++){
435 snprintf(name,kSize, "TOP_1/BTOF%d_1", iSM);
437 printf("\n\n******Misaligning TOF SuperModule ************** %s \n",name);
439 // ************* get ideal global matrix *******************
440 TGeoHMatrix g3 = *fTOFmgr->GetCurrentMatrix();
441 AliInfo(Form("This is the ideal global trasformation of SM %i",iSM));
442 g3.Print(); // g3 is the local(BTOF) to global (ALICE) matrix and is the same of fTOFMatrixId
443 TGeoNode* n3 = fTOFmgr->GetCurrentNode();
444 TGeoMatrix* l3 = n3->GetMatrix();
446 Double_t gA[3], gB[3], gC[3], gD[3]; // ideal global FM point coord.
447 g3.LocalToMaster(lA,gA);
448 g3.LocalToMaster(lB,gB);
449 g3.LocalToMaster(lC,gC);
450 g3.LocalToMaster(lD,gD);
452 // We apply a delta transformation to the surveyed vol to represent
453 // its real position, given below by ng3 nl3, which differs from its
454 // ideal position saved above in g3 and l3
456 //we have to express the displacements as regards the old local RS (non misaligned BTOF)
457 Double_t dx = mis[0]; // shift along x
458 Double_t dy = mis[1]; // shift along y
459 Double_t dz = mis[2]; // shift along z
460 Double_t dphi = mis[3]; // rot around z
461 Double_t dtheta = mis[4]; // rot around x'
462 Double_t dpsi = mis[5]; // rot around z''
464 TGeoRotation* rrot = new TGeoRotation("rot",dphi,dtheta,dpsi);
465 TGeoCombiTrans localdelta = *(new TGeoCombiTrans(dx,dy,dz, rrot));
466 AliInfo(Form("This is the local delta trasformation for SM %i \n",iSM));
468 TGeoHMatrix nlocal = *l3 * localdelta;
469 TGeoHMatrix* nl3 = new TGeoHMatrix(nlocal); // new matrix, representing real position (from new local mis RS to the global one)
471 TGeoPhysicalNode* pn3 = fTOFmgr->MakePhysicalNode(name);
475 TGeoHMatrix* ng3 = pn3->GetMatrix(); //"real" global matrix, what survey sees
476 printf("\n\n************* The Misaligned Matrix in GRS **************\n");
478 Double_t ngA[3], ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
479 ng3->LocalToMaster(lA,ngA);
480 ng3->LocalToMaster(lB,ngB);
481 ng3->LocalToMaster(lC,ngC);
482 ng3->LocalToMaster(lD,ngD);
484 for(Int_t coord=0;coord<3;coord++){
485 fCombFMData[iSM*4][2*coord]=ngA[coord];
486 fCombFMData[iSM*4][2*coord+1]=1;
487 fCombFMData[iSM*4+1][2*coord]=ngB[coord];
488 fCombFMData[iSM*4+1][2*coord+1]=1;
489 fCombFMData[iSM*4+2][2*coord]=ngC[coord];
490 fCombFMData[iSM*4+2][2*coord+1]=1;
491 fCombFMData[iSM*4+3][2*coord]=ngD[coord];
492 fCombFMData[iSM*4+3][2*coord+1]=1;
498 //____________________________________________________________________________
499 void AliTOFAlignment::WriteCombData(const Char_t *nomefile, Int_t option)
501 // 1 for simulated data; 0 for data from survey file
502 // write combined data on a file
506 /* Open file in text mode: */
507 if( (data = fopen( nomefile, "w+t" )) != NULL ){
509 fprintf( data, "simulated data\n" );} else {
510 fprintf( data, "survey data\n" );}
512 fprintf( data, "data from InsertMisAlignmentBTOF method\n");}
513 else {fprintf( data, "real survey data from text file (coordinate in global RS)\n");}
514 fprintf( data, "Point Name,XPH,YPH,ZPH,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
515 fprintf( data, "> Data:\n");
516 for(Int_t i=0;i<72;i++){
517 if (fCombFMData[i][0]!=0){
518 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);
524 printf( "Problem opening the file\n" );
530 //____________________________________________________________________________
531 void AliTOFAlignment::WriteSimSurveyData(const Char_t *nomefile)
533 // write sim data in standard format
538 /* Open file in text mode: */
539 if( (data = fopen( nomefile, "w+t" )) != NULL )
541 fprintf( data, "> Title:\n" );
542 fprintf( data, "simulated data\n" );
543 fprintf( data, "> Date:\n" );
544 fprintf( data, "24.09.2007\n" );
545 fprintf( data, "> Subdetector:\n" );
546 fprintf( data, "TOF\n" );
547 fprintf( data, "> Report URL:\n" );
548 fprintf( data, "https://edms.cern.ch/document/835615\n" );
549 fprintf( data, "> Version:\n" );
550 fprintf( data, "1\n");
551 fprintf( data, "> General Observations:\n");
552 fprintf( data, "data from InsertMisAlignmentBTOF method\n");
553 fprintf( data, "> Coordinate System:\n");
554 fprintf( data, "\\ALICEPH\n");
555 fprintf( data, "> Units:\n");
556 fprintf( data, "cm\n");
557 fprintf( data, "> Nr Columns:\n");
558 fprintf( data, "9\n");
559 fprintf( data, "> Column Names:\n");
560 fprintf( data, "Point Name,XPH,YPH,ZPH,Point Type,Target Used,PrecisionX(mm),PrecisionY(mm),PrecisionZ(mm)\n");
561 fprintf( data, "> Data:\n");
562 for(Int_t i=0;i<72;i++)
563 if (fCombFMData[i][0]!=0)
564 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]);
569 printf( "Problem opening the file\n" );
572 //____________________________________________________________________________
573 void AliTOFAlignment::MakeDefData(const Int_t nf,TString namefiles[])
575 //this method combines survey data from different files (namefiles[])
579 Float_t data[72][6][100];
580 for (Int_t i=0;i<72;i++)
581 for (Int_t j=0; j<6; j++)
582 for(Int_t k=0; k<100; k++)
586 Long64_t totdata[72]={0};
588 for (Int_t ii=0;ii<nf; ii++)
590 AliSurveyObj *so = new AliSurveyObj();
591 const Char_t *nome=namefiles[ii];
592 so->FillFromLocalFile(nome);
593 TObjArray *points = so->GetData();
594 Int_t nSurveyPoint=points->GetEntries();
595 for(Int_t jj=0;jj<nSurveyPoint;jj++){
596 const char* pointName= ((AliSurveyPoint *) points->At(jj))->GetPointName().Data();
597 nfm=atoi(&pointName[6]);
598 nsm=atoi(&pointName[2]);
599 data[nsm*4+nfm][0][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetX();
600 data[nsm*4+nfm][2][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetY();
601 data[nsm*4+nfm][4][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetZ();
602 data[nsm*4+nfm][1][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionX();
603 data[nsm*4+nfm][3][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionY();
604 data[nsm*4+nfm][5][totdata[nsm*4+nfm]]=((AliSurveyPoint *) points->At(jj))->GetPrecisionZ();
605 totdata[nsm*4+nfm]=totdata[nsm*4+nfm]+1;
611 for(Int_t i=0; i<72 ;i++){
612 Float_t numx=0, numy=0,numz=0, comodox=0, comodoy=0, comodoz=0,denx=0, deny=0, denz=0;
614 for(Int_t j=0; j<totdata[i]; j++){
615 comodox=1/(data[i][1][j]/10*data[i][1][j]/10);//precision in mm, position in cm
616 numx=numx+data[i][0][j]*comodox;
618 comodoy=1/(data[i][3][j]/10*data[i][3][j]/10);
619 numy=numy+data[i][2][j]*comodoy;
621 comodoz=1/(data[i][5][j]/10*data[i][5][j]/10);
622 numz=numz+data[i][4][j]*comodoz;
625 fCombFMData[i][1]=TMath::Sqrt(1/denx); //error for x position
626 fCombFMData[i][3]=TMath::Sqrt(1/deny); //error for y position
627 fCombFMData[i][5]=TMath::Sqrt(1/denz); //error for z position
628 fCombFMData[i][0]=numx/denx; //combined survey data for x position of FM
629 fCombFMData[i][2]=numy/deny; //combined survey data for y position of FM
630 fCombFMData[i][4]=numz/denz; //combined survey data for z position of FM
634 for(Int_t i=0;i<72;i++)
635 if (fCombFMData[i][0]!=0){
636 fNFMforSM[(i-i%4)/4][i%4]=1;
637 fNFMforSM[(i-i%4)/4][4]=fNFMforSM[(i-i%4)/4][4]+1;
641 //_____________________________________________________________________________
642 void AliTOFAlignment::ReadSurveyDataAndAlign(){
644 // read the survey data and, if we know the positions of at least 3 FM
645 //for a SM, call the right Alignement procedure
647 fTOFAlignObjArray = new TObjArray(kMaxAlignObj);
649 Float_t deltaFM0=0, deltaFM1=0, deltaFM2=0, deltaFM3=0;
651 for(Int_t i=0; i<18; i++){
652 switch(fNFMforSM[i][4]){
654 printf("we don't know the position of any FM of SM %i\n",i);
657 printf("we know the position of only one FM for SM %i\n",i);
661 printf("we know the position of only 2 FM for SM %i\n",i);
665 if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1){
666 printf("we know the position of FM A B C for SM %i\n",i);
667 AliTOFAlignment::AlignFromSurveyABC(i);};
670 if (fNFMforSM[i][0]==1 && fNFMforSM[i][1]==1 && fNFMforSM[i][3]==1){
671 printf("we know the position of FM A B D for SM %i\n",i);
672 AliTOFAlignment::AlignFromSurveyABD(i);};
675 if (fNFMforSM[i][0]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
676 printf("we know the position of FM A C D for SM %i\n",i);
677 AliTOFAlignment::AlignFromSurveyACD(i);};
680 if (fNFMforSM[i][1]==1 && fNFMforSM[i][2]==1 && fNFMforSM[i][3]==1){
681 printf("we know the position of FM B C D for SM %i\n",i);
682 AliTOFAlignment::AlignFromSurveyBCD(i);};
687 printf("we know the position of all the 4 FM for SM %i\n",i);
688 //check the precision of the measurement
690 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]);
691 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]);
692 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]);
693 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]);
695 //to AlignFromSurvey we use the 3 FM whose positions are known with greatest precision
696 if(deltaFM0>=deltaFM1 && deltaFM0>=deltaFM2 && deltaFM0>=deltaFM3){
697 printf("to Align we use FM B,C,D");
698 AliTOFAlignment::AlignFromSurveyBCD(i);} else
699 if(deltaFM1>=deltaFM0 && deltaFM1>=deltaFM2 && deltaFM1>=deltaFM3){
700 printf("to Align we use FM A,C,D");
701 AliTOFAlignment::AlignFromSurveyACD(i);} else
702 if(deltaFM2>=deltaFM0 && deltaFM2>=deltaFM1 && deltaFM2>=deltaFM3){
703 printf("to Align we use FM A,B,D");
704 AliTOFAlignment::AlignFromSurveyABD(i);} else{
705 printf("to Align we use FM A,B,C");
706 AliTOFAlignment::AlignFromSurveyABC(i);}
713 // saving TOF AligObjs from survey on a file, for the moment..
714 fNTOFAlignObj=fTOFAlignObjArray->GetEntries();
715 AliInfo(Form("Number of Alignable Volumes: %d",fNTOFAlignObj));
716 TFile f("TOFAlignFromSurvey.root","RECREATE");
718 f.WriteObject(fTOFAlignObjArray,"TOFAlignObjs","kSingleKey");
724 //_____________________________________________________________________________
725 void AliTOFAlignment::AlignFromSurveyABC(Int_t iSM)
728 //From Survey data, derive the needed transformations to get the
730 //Again, highly "inspired" to Raffaele's example...
733 Double_t ngA[3], ngB[3], ngC[3]; // real FM point coord., global RS
734 // Get the 'realistic' input from the Survey Matrix
735 for(Int_t coord=0;coord<3;coord++){
736 ngA[coord]= fCombFMData[iSM*4][coord*2];
737 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
738 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
741 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
743 // From the real fiducial marks coordinates derive back the
744 // new global position of the surveyed volume
745 //*** What follows is the actual survey-to-alignment procedure
747 Double_t ab[3], bc[3], n[3];
748 Double_t plane[4], s=1.;
750 // first vector on the plane of the fiducial marks
751 for(Int_t i=0;i<3;i++){
752 ab[i] = (ngB[i] - ngA[i]);
755 // second vector on the plane of the fiducial marks
756 for(Int_t i=0;i<3;i++){
757 bc[i] = (ngC[i] - ngB[i]);
760 // vector normal to the plane of the fiducial marks obtained
761 // as cross product of the two vectors on the plane d0^d1
762 n[0] = (ab[1] * bc[2] - ab[2] * bc[1]);
763 n[1] = (ab[2] * bc[0] - ab[0] * bc[2]);
764 n[2] = (ab[0] * bc[1] - ab[1] * bc[0]);
766 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
768 s = Double_t(1.)/sizen ; //normalization factor
770 AliInfo("Problem in normalizing the vector");
773 // plane expressed in the hessian normal form, see:
774 // http://mathworld.wolfram.com/HessianNormalForm.html
775 // the first three are the coordinates of the orthonormal vector
776 // the fourth coordinate is equal to the distance from the origin
778 for(Int_t i=0;i<3;i++){
781 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
783 // The center of the square with fiducial marks as corners
784 // as the middle point of one diagonal - md
785 // Used below to get the center - orig - of the surveyed box
787 Double_t orig[3], md[3];
788 for(Int_t i=0;i<3;i++){
789 md[i] = (ngA[i] + ngC[i]) * 0.5;
792 // The center of the box, gives the global translation
793 for(Int_t i=0;i<3;i++){
794 orig[i] = md[i] - plane[i]*fgkZFM;
797 // get local directions needed to write the global rotation matrix
798 // for the surveyed volume by normalising vectors ab and bc
799 Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
803 for(Int_t i=0;i<3;i++){
807 Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
809 for(Int_t i=0;i<3;i++){
813 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
814 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey
816 ng.SetTranslation(orig);
818 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
821 // Calculate the delta transformation wrt Ideal geometry
822 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
824 printf("\n\n**** The ideal matrix ***\n");
825 fTOFMatrixId[iSM]->Print();
827 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
828 printf("\n\n**** The inverse of the ideal matrix ***\n");
831 gdelta.MultiplyLeft(&ng);
832 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
833 gdelta.Print(); //this is the global delta trasformation
835 // Now Write the Alignment Objects....
836 Int_t index=0; //let all SM modules have index=0
837 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
838 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
839 TString symname(Form("TOF/sm%02d",iSM));
840 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
841 fTOFAlignObjArray->Add(o);
846 //_____________________________________________________________________________
847 void AliTOFAlignment::AlignFromSurveyABD(Int_t iSM)
850 //From Survey data, derive the needed transformations to get the
852 //Again, highly "inspired" to Raffaele's example...
855 Double_t ngA[3], ngB[3], ngD[3];// real FM point coord., global RS
857 // Get the 'realistic' input from the Survey Matrix
858 for(Int_t coord=0;coord<3;coord++){
859 ngA[coord]= fCombFMData[iSM*4][coord*2];
860 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
861 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
864 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
866 // From the new fiducial marks coordinates derive back the
867 // new global position of the surveyed volume
868 //*** What follows is the actual survey-to-alignment procedure
870 Double_t ab[3], ad[3], n[3];
871 Double_t plane[4], s=1.;
873 // first vector on the plane of the fiducial marks
874 for(Int_t i=0;i<3;i++){
875 ab[i] = (ngB[i] - ngA[i]);
878 // second vector on the plane of the fiducial marks
879 for(Int_t i=0;i<3;i++){
880 ad[i] = (ngD[i] - ngA[i]);
883 // vector normal to the plane of the fiducial marks obtained
884 // as cross product of the two vectors on the plane d0^d1
885 n[0] = (ab[1] * ad[2] - ab[2] * ad[1]);
886 n[1] = (ab[2] * ad[0] - ab[0] * ad[2]);
887 n[2] = (ab[0] * ad[1] - ab[1] * ad[0]);
889 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
891 s = Double_t(1.)/sizen ; //normalization factor
893 AliInfo("Problem in normalizing the vector");
896 // plane expressed in the hessian normal form, see:
897 // http://mathworld.wolfram.com/HessianNormalForm.html
898 // the first three are the coordinates of the orthonormal vector
899 // the fourth coordinate is equal to the distance from the origin
901 for(Int_t i=0;i<3;i++){
904 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
906 // The center of the square with fiducial marks as corners
907 // as the middle point of one diagonal - md
908 // Used below to get the center - orig - of the surveyed box
910 Double_t orig[3], md[3];
911 for(Int_t i=0;i<3;i++){
912 md[i] = (ngB[i] + ngD[i]) * 0.5;
915 // The center of the box, gives the global translation
916 for(Int_t i=0;i<3;i++){
917 orig[i] = md[i] - plane[i]*fgkZFM;
920 // get local directions needed to write the global rotation matrix
921 // for the surveyed volume by normalising vectors ab and bc
922 Double_t sx = TMath::Sqrt(ab[0]*ab[0] + ab[1]*ab[1] + ab[2]*ab[2]);
924 for(Int_t i=0;i<3;i++){
928 Double_t sy = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
930 for(Int_t i=0;i<3;i++){
934 Double_t rot[9] = {ab[0],ad[0],plane[0],ab[1],ad[1],plane[1],ab[2],ad[2],plane[2]};
935 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
937 ng.SetTranslation(orig);
939 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
942 // Calculate the delta transformation wrt Ideal geometry
943 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
945 printf("\n\n**** The ideal matrix ***\n");
946 fTOFMatrixId[iSM]->Print();
948 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
949 printf("\n\n**** The inverse of the ideal matrix ***\n");
952 gdelta.MultiplyLeft(&ng);
953 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
954 gdelta.Print(); //global delta trasformation
956 // Now Write the Alignment Objects....
957 Int_t index=0; //let all SM modules have index=0
958 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
959 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
960 TString symname(Form("TOF/sm%02d",iSM));
961 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
962 fTOFAlignObjArray->Add(o);
965 //_____________________________________________________________________________
966 void AliTOFAlignment::AlignFromSurveyACD(Int_t iSM)
968 //From Survey data, derive the needed transformations to get the
970 //Again, highly "inspired" to Raffaele's example...
974 Double_t ngA[3], ngC[3], ngD[3];// real FM point coord., global RS
976 // Get the 'realistic' input from the Survey Matrix
977 for(Int_t coord=0;coord<3;coord++){
978 ngA[coord]= fCombFMData[iSM*4][coord*2];
979 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
980 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
983 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
985 // From the new fiducial marks coordinates derive back the
986 // new global position of the surveyed volume
987 //*** What follows is the actual survey-to-alignment procedure
989 Double_t cd[3], ad[3], n[3];
990 Double_t plane[4], s=1.;
992 // first vector on the plane of the fiducial marks
993 for(Int_t i=0;i<3;i++){
994 cd[i] = (ngC[i] - ngD[i]);
997 // second vector on the plane of the fiducial marks
998 for(Int_t i=0;i<3;i++){
999 ad[i] = (ngD[i] - ngA[i]);
1002 // vector normal to the plane of the fiducial marks obtained
1003 // as cross product of the two vectors on the plane d0^d1
1004 n[0] = (ad[1] * cd[2] - ad[2] * cd[1]);
1005 n[1] = (ad[2] * cd[0] - ad[0] * cd[2]);
1006 n[2] = (ad[0] * cd[1] - ad[1] * cd[0]);
1008 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1010 s = Double_t(1.)/sizen ; //normalization factor
1012 AliInfo("Problem in normalizing the vector");
1015 // plane expressed in the hessian normal form, see:
1016 // http://mathworld.wolfram.com/HessianNormalForm.html
1017 // the first three are the coordinates of the orthonormal vector
1018 // the fourth coordinate is equal to the distance from the origin
1020 for(Int_t i=0;i<3;i++){
1021 plane[i] = n[i] * s;
1023 plane[3] = ( plane[0] * ngA[0] + plane[1] * ngA[1] + plane[2] * ngA[2] );
1025 // The center of the square with fiducial marks as corners
1026 // as the middle point of one diagonal - md
1027 // Used below to get the center - orig - of the surveyed box
1029 Double_t orig[3], md[3];
1030 for(Int_t i=0;i<3;i++){
1031 md[i] = (ngA[i] + ngC[i]) * 0.5;
1034 // The center of the box, gives the global translation
1035 for(Int_t i=0;i<3;i++){
1036 orig[i] = md[i] + plane[i]*fgkZFM;
1039 // get local directions needed to write the global rotation matrix
1040 // for the surveyed volume by normalising vectors ab and bc
1041 Double_t sx = TMath::Sqrt(ad[0]*ad[0] + ad[1]*ad[1] + ad[2]*ad[2]);
1043 for(Int_t i=0;i<3;i++){
1047 Double_t sy = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1049 for(Int_t i=0;i<3;i++){
1053 Double_t rot[9] = {cd[0],ad[0],-plane[0],cd[1],ad[1],-plane[1],cd[2],ad[2],-plane[2]};
1054 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1056 ng.SetTranslation(orig);
1057 ng.SetRotation(rot);
1058 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1061 // Calculate the delta transformation wrt Ideal geometry
1062 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1064 printf("\n\n**** The ideal matrix ***\n");
1065 fTOFMatrixId[iSM]->Print();
1067 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1068 printf("\n\n**** The inverse of the ideal matrix ***\n");
1071 gdelta.MultiplyLeft(&ng);
1072 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1073 gdelta.Print(); //global delta trasformation
1075 // Now Write the Alignment Objects....
1076 Int_t index=0; //let all SM modules have index=0
1077 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1078 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
1079 TString symname(Form("TOF/sm%02d",iSM));
1080 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1081 fTOFAlignObjArray->Add(o);
1084 //___________________________________________________________________________
1085 void AliTOFAlignment::AlignFromSurveyBCD(Int_t iSM)
1087 //From Survey data, derive the needed transformations to get the
1088 //Alignment Objects.
1089 //Again, highly "inspired" to Raffaele's example...
1092 Double_t ngB[3], ngC[3], ngD[3];// real FM point coord., global RS
1095 // Get the 'realistic' input from the Survey Matrix
1096 for(Int_t coord=0;coord<3;coord++){
1097 ngB[coord]= fCombFMData[iSM*4+1][coord*2];
1098 ngC[coord]= fCombFMData[iSM*4+2][coord*2];
1099 ngD[coord]= fCombFMData[iSM*4+3][coord*2];
1102 printf("\n\n******Survey analysis for TOF SuperModule ************** %i \n",iSM);
1104 // From the new fiducial marks coordinates derive back the
1105 // new global position of the surveyed volume
1106 //*** What follows is the actual survey-to-alignment procedure
1108 Double_t cd[3], bc[3], n[3];
1109 Double_t plane[4], s=1.;
1111 // first vector on the plane of the fiducial marks
1112 for(Int_t i=0;i<3;i++){
1113 cd[i] = (ngC[i] - ngD[i]);
1116 // second vector on the plane of the fiducial marks
1117 for(Int_t i=0;i<3;i++){
1118 bc[i] = (ngC[i] - ngB[i]);
1121 // vector normal to the plane of the fiducial marks obtained
1122 // as cross product of the two vectors on the plane d0^d1
1123 n[0] = (bc[1] * cd[2] - bc[2] * cd[1]);
1124 n[1] = (bc[2] * cd[0] - bc[0] * cd[2]);
1125 n[2] = (bc[0] * cd[1] - bc[1] * cd[0]);
1127 Double_t sizen = TMath::Sqrt( n[0]*n[0] + n[1]*n[1] + n[2]*n[2] );
1129 s = Double_t(1.)/sizen ; //normalization factor
1131 AliInfo("Problem in normalizing the vector");
1134 // plane expressed in the hessian normal form, see:
1135 // http://mathworld.wolfram.com/HessianNormalForm.html
1136 // the first three are the coordinates of the orthonormal vector
1137 // the fourth coordinate is equal to the distance from the origin
1139 for(Int_t i=0;i<3;i++){
1140 plane[i] = n[i] * s;
1142 plane[3] = ( plane[0] * ngB[0] + plane[1] * ngB[1] + plane[2] * ngB[2] );
1144 // The center of the square with fiducial marks as corners
1145 // as the middle point of one diagonal - md
1146 // Used below to get the center - orig - of the surveyed box
1148 Double_t orig[3], md[3];
1149 for(Int_t i=0;i<3;i++){
1150 md[i] = (ngB[i] + ngD[i]) * 0.5;
1153 // The center of the box, gives the global translation
1154 for(Int_t i=0;i<3;i++){
1155 orig[i] = md[i] + plane[i]*fgkZFM;
1158 // get local directions needed to write the global rotation matrix
1159 // for the surveyed volume by normalising vectors ab and bc
1160 Double_t sx = TMath::Sqrt(cd[0]*cd[0] + cd[1]*cd[1] + cd[2]*cd[2]);
1162 for(Int_t i=0;i<3;i++){
1166 Double_t sy = TMath::Sqrt(bc[0]*bc[0] + bc[1]*bc[1] + bc[2]*bc[2]);
1168 for(Int_t i=0;i<3;i++){
1172 Double_t rot[9] = {cd[0],bc[0],-plane[0],cd[1],bc[1],-plane[1],cd[2],bc[2],-plane[2]};
1173 // the Aligned matrix for the current TOF SM in the Global RS, as derived from Survey:
1175 ng.SetTranslation(orig);
1176 ng.SetRotation(rot);
1177 printf("\n\n**** The Misaligned Matrix in GRS, as from Survey data ***\n");
1180 // Calculate the delta transformation wrt Ideal geometry
1181 // (Should be gdelta.rot ==I and gdelta.tr=0 if no misalignment is applied.)
1183 printf("\n\n**** The ideal matrix ***\n");
1184 fTOFMatrixId[iSM]->Print();
1186 TGeoHMatrix gdelta =fTOFMatrixId[iSM]->Inverse();
1187 printf("\n\n**** The inverse of the ideal matrix ***\n");
1190 gdelta.MultiplyLeft(&ng);
1191 printf("\n\n**** The Delta Matrix in GRS, as from Survey data ***\n");
1192 gdelta.Print(); //global delta trasformation
1194 // Now Write the Alignment Objects....
1195 Int_t index=0; //let all SM modules have index=0
1196 AliGeomManager::ELayerID layer = AliGeomManager::kInvalidLayer;
1197 UShort_t dvoluid = AliGeomManager::LayerToVolUID(layer,index); //dummy vol id
1198 TString symname(Form("TOF/sm%02d",iSM));
1199 AliAlignObjMatrix* o = new AliAlignObjMatrix(symname.Data(),dvoluid,gdelta,kTRUE);
1200 fTOFAlignObjArray->Add(o);