]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSAlignMilleModule.cxx
LUTs mapping symnames and original global matrices removed from AliGeomManager, which...
[u/mrichter/AliRoot.git] / ITS / AliITSAlignMilleModule.cxx
1 /************************************************************************** \r
2  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. * \r
3  *                                                                        * \r
4  * Author: The ALICE Off-line Project.                                    * \r
5  * Contributors are mentioned in the code where appropriate.              * \r
6  *                                                                        * \r
7  * Permission to use, copy, modify and distribute this software and its   * \r
8  * documentation strictly for non-commercial purposes is hereby granted   * \r
9  * without fee, provided that the above copyright notice appears in all   * \r
10  * copies and that both the copyright notice and this permission notice   * \r
11  * appear in the supporting documentation. The authors make no claims     * \r
12  * about the suitability of this software for any purpose. It is          * \r
13  * provided "as is" without express or implied warranty.                  * \r
14  **************************************************************************/ \r
15  \r
16 /* $Id$    */ \r
17 //----------------------------------------------------------------------------- \r
18 /// \class AliITSAlignMilleModule\r
19 /// Alignment class for the ALICE ITS detector \r
20 /// \r
21 /// This class is used by AliITSAlignMille to build custom supermodules    \r
22 /// made of ITS sensitive modules. These supermodules are then aligned\r
23 /// \r
24 /// Custom supermodules must have VolumeID > 14335\r
25 ///\r
26 /// \author M. Lunardon  \r
27 //----------------------------------------------------------------------------- \r
28  \r
29 #include <TGeoManager.h> \r
30 #include <TGeoMatrix.h> \r
31  \r
32 #include "AliITSAlignMilleModule.h" \r
33 #include "AliITSgeomTGeo.h" \r
34 #include "AliGeomManager.h" \r
35 #include "AliAlignObjParams.h" \r
36 #include "AliLog.h" \r
37  \r
38 /// \cond CLASSIMP \r
39 ClassImp(AliITSAlignMilleModule) \r
40 /// \endcond \r
41     \r
42 //-------------------------------------------------------------\r
43 AliITSAlignMilleModule::AliITSAlignMilleModule() : TNamed(), \r
44   fNSensVol(0), \r
45   fIndex(-1),  \r
46   fVolumeID(0),  \r
47   fMatrix(NULL),\r
48   fSensVolMatrix(NULL),\r
49   fSensVolModifMatrix(NULL),\r
50   fTempAlignObj(NULL)\r
51\r
52   /// void constructor  \r
53   fMatrix = new TGeoHMatrix; \r
54   fSensVolMatrix = new TGeoHMatrix; \r
55   fSensVolModifMatrix = new TGeoHMatrix; \r
56   fTempAlignObj=new AliAlignObjParams;\r
57\r
58 //-------------------------------------------------------------\r
59 AliITSAlignMilleModule::AliITSAlignMilleModule(Int_t index, UShort_t volid, char* symname, TGeoHMatrix *m, Int_t nsv, UShort_t *volidsv) : TNamed(), \r
60   fNSensVol(0), \r
61   fIndex(-1),  \r
62   fVolumeID(0),  \r
63   fMatrix(NULL),\r
64   fSensVolMatrix(NULL),\r
65   fSensVolModifMatrix(NULL),\r
66   fTempAlignObj(NULL)\r
67\r
68   /// void constructor  \r
69   fMatrix = new TGeoHMatrix; \r
70   fSensVolMatrix = new TGeoHMatrix; \r
71   fSensVolModifMatrix = new TGeoHMatrix; \r
72   fTempAlignObj=new AliAlignObjParams;\r
73   if (Set(index,volid,symname,m,nsv,volidsv)) {\r
74     AliInfo("Error in AliITSAlignMilleModule::Set() - initializing void supermodule...");\r
75   }\r
76\r
77 //-------------------------------------------------------------\r
78 AliITSAlignMilleModule::AliITSAlignMilleModule(UShort_t volid) : TNamed(), \r
79   fNSensVol(0), \r
80   fIndex(-1),  \r
81   fVolumeID(0),  \r
82   fMatrix(NULL),\r
83   fSensVolMatrix(NULL),\r
84   fSensVolModifMatrix(NULL),\r
85   fTempAlignObj(NULL)\r
86\r
87   /// simple constructor building a supermodule from a single sensitive volume \r
88   fMatrix = new TGeoHMatrix; \r
89   fSensVolMatrix = new TGeoHMatrix; \r
90   fSensVolModifMatrix = new TGeoHMatrix;   \r
91   // temporary align object, just use the rotation...\r
92   fTempAlignObj=new AliAlignObjParams;\r
93 \r
94   fIndex = GetIndexFromVolumeID(volid);  \r
95   if (fIndex>=0 && gGeoManager) { // good sensitive module and geometry loaded\r
96     SetName(AliGeomManager::SymName(volid));\r
97     fVolumeID = volid;\r
98     AddSensitiveVolume(volid);\r
99     if (SensVolMatrix(volid, fMatrix))\r
100        AliInfo("Matrix not defined");\r
101   }\r
102   else {\r
103     AliInfo("Wrong VolumeID or Geometry not loaded - initializing void supermodule...");\r
104   }\r
105\r
106 //-------------------------------------------------------------\r
107 AliITSAlignMilleModule::~AliITSAlignMilleModule() { \r
108   /// Destructor \r
109   delete fMatrix; \r
110   delete fSensVolMatrix; \r
111   delete fSensVolModifMatrix; \r
112   delete fTempAlignObj;\r
113\r
114 //-------------------------------------------------------------\r
115 Int_t AliITSAlignMilleModule::Set(Int_t index, UShort_t volid, char* symname, TGeoHMatrix *m, Int_t nsv, UShort_t *volidsv) \r
116 {\r
117   // initialize a custom supermodule\r
118   // index, volid, symname and matrix must be given\r
119   // if (volidsv) add nsv sensitive volumes to the supermodules\r
120   // return 0 if success\r
121 \r
122   if (index<2198) {\r
123     AliInfo("Index must be >= 2198");\r
124     return -1;\r
125   }\r
126   if (volid<14336) {\r
127     AliInfo("VolumeID must be >= 14336");\r
128     return -2;\r
129   }\r
130   \r
131   if (!symname) return -3;\r
132   for (Int_t i=0; i<2198; i++) {\r
133     if (!strcmp(symname,AliITSgeomTGeo::GetSymName(i))) {\r
134       AliInfo("Symname already used by a Sensitive Volume");\r
135       return -3;\r
136     }\r
137   }\r
138   \r
139   if (!m) return -4;\r
140 \r
141   // can initialize needed stuffs\r
142   fIndex = index;\r
143   fVolumeID = volid;\r
144   SetName(symname);\r
145   (*fMatrix) = (*m);\r
146 \r
147   // add sensitive volumes\r
148   for (Int_t i=0; i<nsv; i++) AddSensitiveVolume(volidsv[i]);\r
149 \r
150   return 0;\r
151 }\r
152 //-------------------------------------------------------------\r
153 Int_t AliITSAlignMilleModule::GetIndexFromVolumeID(UShort_t voluid) {\r
154   /// index from volume ID\r
155   AliGeomManager::ELayerID lay = AliGeomManager::VolUIDToLayer(voluid);\r
156   if (lay<1|| lay>6) return -1;\r
157   Int_t idx=Int_t(voluid)-2048*lay;\r
158   if (idx>=AliGeomManager::LayerSize(lay)) return -1;\r
159   for (Int_t ilay=1; ilay<lay; ilay++) \r
160     idx += AliGeomManager::LayerSize(ilay);\r
161   return idx;\r
162 }\r
163 //-------------------------------------------------------------\r
164 void AliITSAlignMilleModule::AddSensitiveVolume(UShort_t voluid)\r
165 {\r
166   /// add a sensitive volume to this supermodule\r
167   if (GetIndexFromVolumeID(voluid)<0) return; // bad volid\r
168   fSensVolVolumeID[fNSensVol] = voluid;\r
169   fSensVolIndex[fNSensVol] = GetIndexFromVolumeID(voluid);\r
170   fNSensVol++;\r
171 }\r
172 //-------------------------------------------------------------\r
173 Bool_t AliITSAlignMilleModule::IsIn(UShort_t voluid) const \r
174 {\r
175   /// check if voluid is defined\r
176   if (!voluid) return kFALSE; // only positive voluid are accepted\r
177   for (Int_t i=0; i<fNSensVol; i++) {\r
178     if (fSensVolVolumeID[i]==voluid) return kTRUE;\r
179   }\r
180   return kFALSE;\r
181 }\r
182 //-------------------------------------------------------------\r
183 TGeoHMatrix *AliITSAlignMilleModule::GetSensitiveVolumeModifiedMatrix(UShort_t voluid, Double_t *deltalocal)\r
184 {\r
185   // modify the original TGeoHMatrix of the sensitive module 'voluid' according\r
186   // with a delta transform. applied to the supermodule matrix\r
187   // return NULL if error\r
188 \r
189   if (!IsIn(voluid)) return NULL;\r
190   if (!gGeoManager) return NULL;\r
191 \r
192   // prepare the TGeoHMatrix\r
193   Double_t tr[3],ang[3];\r
194   tr[0]=deltalocal[0]; // in centimeter\r
195   tr[1]=deltalocal[1]; \r
196   tr[2]=deltalocal[2];\r
197   ang[0]=deltalocal[3]; // psi   (X)  in deg\r
198   ang[1]=deltalocal[4]; // theta (Y)\r
199   ang[2]=deltalocal[5]; // phi   (Z)\r
200 \r
201   // reset align object (may not be needed...)\r
202   fTempAlignObj->SetTranslation(0,0,0);\r
203   fTempAlignObj->SetRotation(0,0,0);\r
204 \r
205   fTempAlignObj->SetRotation(ang[0],ang[1],ang[2]);\r
206   fTempAlignObj->SetTranslation(tr[0],tr[1],tr[2]);\r
207   AliDebug(3,Form("Delta angles: psi=%f  theta=%f   phi=%f",ang[0],ang[1],ang[2]));\r
208   TGeoHMatrix hm;\r
209   fTempAlignObj->GetMatrix(hm);\r
210   //printf("\n0: delta matrix\n");hm.Print();\r
211 \r
212   // 1) start setting fSensVolModif = fSensVol\r
213   if (SensVolMatrix(voluid, fSensVolModifMatrix)) return NULL;\r
214   //printf("\n1: modif=orig del sensvol\n");fSensVolModifMatrix->Print();\r
215 \r
216   // 2) set fSensVolModif = SensVolRel\r
217   fSensVolModifMatrix->MultiplyLeft( &fMatrix->Inverse() );\r
218   //printf("\n2: modif=relative del sensvol\n");fSensVolModifMatrix->Print();\r
219  \r
220   // 3) multiply left by delta\r
221   fSensVolModifMatrix->MultiplyLeft( &hm );\r
222   //printf("\n3: modif= delta*relative\n");fSensVolModifMatrix->Print();\r
223   \r
224   // 4) multiply left by fMatrix\r
225   fSensVolModifMatrix->MultiplyLeft( fMatrix );\r
226   //printf("\n4: modif=finale\n");fSensVolModifMatrix->Print();\r
227 \r
228   return fSensVolModifMatrix;\r
229 }\r
230 //-------------------------------------------------------------\r
231 AliAlignObjParams *AliITSAlignMilleModule::GetSensitiveVolumeMisalignment(UShort_t voluid, Double_t *deltalocal)\r
232 {\r
233   // calculate misalignment of sens.vol. 'voluid' according with a displacement 'deltalocal'\r
234   // of the mother volume. The misalignment is returned as AliAlignObjParams object\r
235 \r
236   if (!IsIn(voluid)) return NULL;\r
237   if (!gGeoManager) return NULL;\r
238   \r
239   // prepare the TGeoHMatrix\r
240   Double_t tr[3],ang[3];\r
241   tr[0]=deltalocal[0]; // in centimeter\r
242   tr[1]=deltalocal[1]; \r
243   tr[2]=deltalocal[2];\r
244   ang[0]=deltalocal[3]; // psi   (X)  in deg\r
245   ang[1]=deltalocal[4]; // theta (Y)\r
246   ang[2]=deltalocal[5]; // phi   (Z)\r
247 \r
248   // reset align object (may not be needed...)\r
249   fTempAlignObj->SetTranslation(0,0,0);\r
250   fTempAlignObj->SetRotation(0,0,0);\r
251 \r
252   fTempAlignObj->SetRotation(ang[0],ang[1],ang[2]);\r
253   fTempAlignObj->SetTranslation(tr[0],tr[1],tr[2]);\r
254   AliDebug(3,Form("Delta angles: psi=%f  theta=%f   phi=%f",ang[0],ang[1],ang[2]));\r
255   \r
256   return GetSensitiveVolumeMisalignment(voluid,fTempAlignObj);\r
257 }\r
258 //-------------------------------------------------------------\r
259 AliAlignObjParams *AliITSAlignMilleModule::GetSensitiveVolumeMisalignment(UShort_t voluid, AliAlignObjParams *a)\r
260 {\r
261   // return the misalignment of the sens. vol. 'voluid' corresponding with \r
262   // a misalignment 'a' in the mother volume\r
263   // return NULL if error\r
264 \r
265   // Gsv = Gg * Gg-1 * Gsv   -> Lsv,g = Gg-1 * Gsv\r
266   // G'sv = Gg * Dg * Lsv,g === Gsv * Dsv\r
267   // Gg * Dg * Gg-1 * Gsv = Gsv * Gsv-1 * Gg * Dg * Gg-1 * Gsv\r
268   //\r
269   // => Dsv = (Gsv-1 * Gg * Dg * Gg-1 * Gsv)\r
270   //\r
271 \r
272   if (!IsIn(voluid)) return NULL;\r
273   if (!gGeoManager) return NULL;\r
274 \r
275   //a->Print("");\r
276 \r
277   // prepare the Delta matrix Dg\r
278   TGeoHMatrix dg;\r
279   a->GetMatrix(dg);\r
280   //dg.Print();\r
281 \r
282   // 1) start setting fSensVolModif = Gsv\r
283   if (SensVolMatrix(voluid, fSensVolModifMatrix)) return NULL;\r
284   //printf("\n1: modif=orig del sensvol\n");fSensVolModifMatrix->Print();\r
285 \r
286   // 2) set fSensVolModif = Gg-1 * Gsv\r
287   fSensVolModifMatrix->MultiplyLeft( &fMatrix->Inverse() );\r
288   //printf("\n2: modif=relative del sensvol\n");fSensVolModifMatrix->Print();\r
289  \r
290   // 3) set fSensVolModif = Dg * Gg-1 * Gsv\r
291   fSensVolModifMatrix->MultiplyLeft( &dg );\r
292   //printf("\n3: modif= delta*relative\n");fSensVolModifMatrix->Print();\r
293   \r
294   // 4) set fSensVolModif = Gg * Dg * Gg-1 * Gsv\r
295   fSensVolModifMatrix->MultiplyLeft( fMatrix );\r
296   //printf("\n4: modif=quasi,manca il Gsv-1...\n");fSensVolModifMatrix->Print();\r
297 \r
298   // 5) set fSensVolModif = Gsv-1 * Gg * Dg * Gg-1 * Gsv\r
299   if (SensVolMatrix(voluid, &dg)) return NULL;\r
300   fSensVolModifMatrix->MultiplyLeft( &dg.Inverse() );\r
301   //printf("\n5: modif=finale\n");fSensVolModifMatrix->Print();\r
302 \r
303   // reset align object (may not be needed...)\r
304   fTempAlignObj->SetTranslation(0,0,0);\r
305   fTempAlignObj->SetRotation(0,0,0);\r
306 \r
307   if (!fTempAlignObj->SetMatrix(*fSensVolModifMatrix)) return NULL;\r
308   fTempAlignObj->SetVolUID(voluid);\r
309   fTempAlignObj->SetSymName(AliGeomManager::SymName(voluid));\r
310   \r
311   //fTempAlignObj->Print("");\r
312 \r
313   return fTempAlignObj;\r
314 }\r
315 //-------------------------------------------------------------\r
316 TGeoHMatrix *AliITSAlignMilleModule::GetSensitiveVolumeMatrix(UShort_t voluid)\r
317 {\r
318   // return TGeoHMatrix of the sens.vol. 'voluid' in the current geometry\r
319   if (SensVolMatrix(voluid,fSensVolMatrix)) return NULL;\r
320   return fSensVolMatrix;\r
321 }\r
322 //-------------------------------------------------------------\r
323 TGeoHMatrix *AliITSAlignMilleModule::GetSensitiveVolumeOrigGlobalMatrix(UShort_t voluid)\r
324 {\r
325   // return original ideal position (from AliGeomManager::GetOrigGlobalMatrix())\r
326   if (SensVolOrigGlobalMatrix(voluid,fSensVolMatrix)) return NULL;\r
327   return fSensVolMatrix;\r
328 }\r
329 //-------------------------------------------------------------\r
330 Int_t AliITSAlignMilleModule::SensVolMatrix(UShort_t volid, TGeoHMatrix *m) \r
331 {\r
332   // set matrix for sensitive modules (SPD corrected)\r
333   // return 0 if success\r
334   Double_t rot[9];\r
335   Int_t idx=GetIndexFromVolumeID(volid);\r
336   if (idx<0) return -1;\r
337   if (!AliITSgeomTGeo::GetRotation(idx,rot)) return -2;\r
338   m->SetRotation(rot);\r
339   Double_t oLoc[3]={0,0,0};\r
340   Double_t oGlo[3]={0,0,0};\r
341   if (!AliITSgeomTGeo::LocalToGlobal(idx,oLoc,oGlo)) return -3;\r
342   m->SetTranslation(oGlo);\r
343   return 0;\r
344 }\r
345 //-------------------------------------------------------------\r
346 Int_t AliITSAlignMilleModule::SensVolOrigGlobalMatrix(UShort_t volid, TGeoHMatrix *m) \r
347 {\r
348   // set original global matrix for sensitive modules (SPD corrected)\r
349   // return 0 if success\r
350   Int_t idx=GetIndexFromVolumeID(volid);\r
351   if (idx<0) return -1;\r
352   TGeoHMatrix mo;\r
353   if (!AliGeomManager::GetOrigGlobalMatrix(AliGeomManager::SymName(volid),mo));\r
354   (*m)=mo;\r
355 \r
356   // SPD y-shift by 81 mu\r
357   Double_t oLoc[3]={0.0,0.0081,0.0};\r
358   Double_t oGlo[3]={0,0,0};\r
359   m->LocalToMaster(oLoc,oGlo);\r
360   m->SetTranslation(oGlo);\r
361   return 0;\r
362 }\r
363 //-------------------------------------------------------------\r
364 UShort_t AliITSAlignMilleModule::GetVolumeIDFromSymname(const Char_t *symname) {\r
365   /// volume ID from symname\r
366   if (!symname) return 0;\r
367 \r
368   for (UShort_t voluid=2000; voluid<13300; voluid++) {\r
369     Int_t modId;\r
370     AliGeomManager::ELayerID layerId = AliGeomManager::VolUIDToLayer(voluid,modId);\r
371     if (layerId>0 && layerId<7 && modId>=0 && modId<AliGeomManager::LayerSize(layerId)) {\r
372       if (!strcmp(symname,AliGeomManager::SymName(layerId,modId))) return voluid;\r
373     }\r
374   }\r
375 \r
376   return 0;\r
377 }\r
378 \r
379 UShort_t AliITSAlignMilleModule::GetVolumeIDFromIndex(Int_t index) {\r
380   /// volume ID from index\r
381   if (index<0 || index>2197) return 0;\r
382   return GetVolumeIDFromSymname(AliITSgeomTGeo::GetSymName(index));\r
383 }\r
384 //-------------------------------------------------------------\r
385 void AliITSAlignMilleModule::Print(Option_t*) const \r
386 {\r
387   ///\r
388   printf("*** ITS SuperModule for AliITSAlignMille ***\n");\r
389   printf("symname  : %s\n",GetName());\r
390   printf("volumeID : %d\n",fVolumeID);\r
391   printf("index    : %d\n",fIndex);\r
392   fMatrix->Print();\r
393   printf("number of sensitive modules : %d\n",fNSensVol);\r
394   for (Int_t i=0; i<fNSensVol; i++) printf("   voluid[%d] = %d\n",i,fSensVolVolumeID[i]);\r
395 }\r
396 //_____________________________________________________________________________\r
397 AliITSAlignMilleModule::AliITSAlignMilleModule(const AliITSAlignMilleModule &m) :\r
398   TNamed(m),\r
399   fNSensVol(m.fNSensVol),\r
400   fIndex(m.fIndex),\r
401   fVolumeID(m.fVolumeID),\r
402   fMatrix(new TGeoHMatrix(*m.GetMatrix())),\r
403   fSensVolMatrix(new TGeoHMatrix),\r
404   fSensVolModifMatrix(new TGeoHMatrix),\r
405   fTempAlignObj(new AliAlignObjParams)\r
406 {\r
407   // Copy constructor\r
408   for (int i=0; i<fNSensVol; i++) {\r
409     fSensVolIndex[i]=m.fSensVolIndex[i];\r
410     fSensVolVolumeID[i]=m.fSensVolVolumeID[i];\r
411   }\r
412 }\r
413 //_____________________________________________________________________________\r
414 AliITSAlignMilleModule& AliITSAlignMilleModule::operator=(const AliITSAlignMilleModule &m)  \r
415 {\r
416   // operator =\r
417   //\r
418   if(this==&m) return *this;\r
419   ((TNamed *)this)->operator=(m);\r
420   \r
421   fNSensVol=m.fNSensVol;\r
422   fIndex=m.fIndex;\r
423   fVolumeID=m.fVolumeID;\r
424   delete fMatrix;\r
425   fMatrix=new TGeoHMatrix(*m.GetMatrix());\r
426   for (int i=0; i<fNSensVol; i++) {\r
427     fSensVolIndex[i]=m.fSensVolIndex[i];\r
428     fSensVolVolumeID[i]=m.fSensVolVolumeID[i];\r
429   }\r
430   return *this;\r
431 }\r
432 \r
433 //_____________________________________________________________________________\r
434 \r
435 \r