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