]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/UPGRADE/AliITSUGeomTGeo.h
Tracking2Local matrix is not stored in the geometry but precalculated
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUGeomTGeo.h
CommitLineData
451f5018 1#ifndef ALIITSUGEOMTGEO_H
2#define ALIITSUGEOMTGEO_H
a11ef2e4 3
451f5018 4/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
6
7/////////////////////////////////////////////////////////////////////////
8// AliITSUGeomTGeo is a simple interface class to TGeoManager //
9// It is used in the simulation and reconstruction in order to //
10// query the TGeo ITS geometry //
11// //
12// author - cvetan.cheshkov@cern.ch //
13// 15/02/2007 //
14// adapted to ITSupg 18/07/2012 - ruben.shahoyan@cern.ch //
15// RS: in order to preserve the static character of the class but //
16// make it dynamically access geometry, we need to check in every //
17// method if the structures are initialized. To be converted to //
18// singleton at later stage. //
19// //
20// Note on the upgrade detector types: //
21// The coarse type defines detectors served by different classes, //
75875328 22// like Pix. Each such a detector type can have kMaxSegmPerDetType //
451f5018 23// segmentations (pitch etc.) whose parameteres are stored in the //
75875328 24// AliITSsegmentation derived class (like AliITSUSegmentationPix) //
451f5018 25// This allows to have in the setup modules served by the same //
26// classes but with different segmentations. //
27// The full detector type is composed as: //
28// CoarseType*kMaxSegmPerDetType + segmentationType //
29// The only requirement on the segmentationType that should be //
30// < kMaxSegmPerDetType. //
31// The methods like GetLayerDetTypeID return the full detector type //
32// //
33// //
34/////////////////////////////////////////////////////////////////////////
35
36#include <TObject.h>
37#include <TGeoMatrix.h>
38#include <TString.h>
a11ef2e4 39#include <TObjArray.h>
451f5018 40
41class TGeoPNEntry;
42class TDatime;
43
44class AliITSUGeomTGeo : public TObject {
45
46 public:
47 enum {kITSVNA, kITSVUpg}; // ITS version
02d6eccc 48 enum {kDetTypePix=0, kNDetTypes, kMaxSegmPerDetType=10}; // defined detector types (each one can have different segmentations)
451f5018 49 //
50 AliITSUGeomTGeo(Bool_t build = kFALSE);
51 virtual ~AliITSUGeomTGeo();
52 AliITSUGeomTGeo(const AliITSUGeomTGeo &src);
53 AliITSUGeomTGeo& operator=(const AliITSUGeomTGeo &geom);
54 //
02d6eccc 55 Int_t GetNModules() const {return fNModules;}
56 Int_t GetNDetectors(Int_t lay) const {return fNDetectors[lay];}
57 Int_t GetNLadders(Int_t lay) const {return fNLadders[lay];}
58 Int_t GetNLayers() const {return fNLayers;}
451f5018 59
02d6eccc 60 Int_t GetModuleIndex(Int_t lay,int detInLay) const {return GetFirstModIndex(lay)+detInLay;}
451f5018 61 Int_t GetModuleIndex(Int_t lay,Int_t lad,Int_t det) const;
62 Bool_t GetModuleId(Int_t index,Int_t &lay,Int_t &lad,Int_t &det) const;
63 Int_t GetLayer(Int_t index) const;
64 Int_t GetLadder(Int_t index) const;
65 Int_t GetModIdInLayer(Int_t index) const;
66 Int_t GetModIdInLadder(Int_t index) const;
67 //
68 Int_t GetLastModIndex(Int_t lay) const {return fLastModIndex[lay];}
69 Int_t GetFirstModIndex(Int_t lay) const {return (lay==0) ? 0:fLastModIndex[lay-1]+1;}
70 //
71 const char *GetSymName(Int_t index) const;
72 const char *GetSymName(Int_t lay,Int_t lad,Int_t det) const;
73 //
74 // Attention: these are the matrices for the alignable volumes of the modules, i.e. not necessarily the sensors
75 TGeoHMatrix* GetMatrix(Int_t index) const;
76 TGeoHMatrix* GetMatrix(Int_t lay,Int_t lad,Int_t det) const;
77 Bool_t GetTranslation(Int_t index, Double_t t[3]) const;
78 Bool_t GetTranslation(Int_t lay,Int_t lad,Int_t det, Double_t t[3]) const;
79 Bool_t GetRotation(Int_t index, Double_t r[9]) const;
80 Bool_t GetRotation(Int_t lay,Int_t lad,Int_t det, Double_t r[9]) const;
81 Bool_t GetOrigMatrix(Int_t index, TGeoHMatrix &m) const;
82 Bool_t GetOrigMatrix(Int_t lay,Int_t lad,Int_t det, TGeoHMatrix &m) const;
83 Bool_t GetOrigTranslation(Int_t index, Double_t t[3]) const;
84 Bool_t GetOrigTranslation(Int_t lay,Int_t lad,Int_t det, Double_t t[3]) const;
85 Bool_t GetOrigRotation(Int_t index, Double_t r[9]) const;
86 Bool_t GetOrigRotation(Int_t lay,Int_t lad,Int_t det, Double_t r[9]) const;
87 //
cf457606 88 const TGeoHMatrix* GetMatrixT2L(Int_t index);
89 const TGeoHMatrix* GetMatrixT2L(Int_t lay,Int_t lad,Int_t det) {return GetMatrixT2L( GetModuleIndex(lay,lad,det) );}
90 const TGeoHMatrix* GetMatrixSens(Int_t index);
91 const TGeoHMatrix* GetMatrixSens(Int_t lay,Int_t lad,Int_t det) {return GetMatrixSens( GetModuleIndex(lay,lad,det) );}
451f5018 92 //
cf457606 93 Bool_t GetTrackingMatrix(Int_t index, TGeoHMatrix &m);
94 Bool_t GetTrackingMatrix(Int_t lay,Int_t lad,Int_t det, TGeoHMatrix &m);
451f5018 95 //
96 // Attention: these are transformations wrt sensitive volume!
cf457606 97 void LocalToGlobal(Int_t index, const Double_t *loc, Double_t *glob);
98 void LocalToGlobal(Int_t lay, Int_t lad, Int_t det,const Double_t *loc, Double_t *glob);
451f5018 99 //
cf457606 100 void GlobalToLocal(Int_t index, const Double_t *glob, Double_t *loc);
101 void GlobalToLocal(Int_t lay, Int_t lad, Int_t det,const Double_t *glob, Double_t *loc);
451f5018 102 //
cf457606 103 void LocalToGlobalVect(Int_t index, const Double_t *loc, Double_t *glob);
104 void GlobalToLocalVect(Int_t index, const Double_t *glob, Double_t *loc);
451f5018 105 Int_t GetLayerDetTypeID(Int_t lr) const;
106 Int_t GetModuleDetTypeID(Int_t id) const;
107 //
02d6eccc 108 virtual void Print(Option_t *opt="") const;
109
cf457606 110 static const char* GetITSVolPattern() {return fgkITSVolName;}
111 static const char* GetITSLayerPattern() {return fgkITSLrName;}
112 static const char* GetITSLadderPattern() {return fgkITSLadName;}
113 static const char* GetITSModulePattern() {return fgkITSModName;}
114 static const char* GetITSSensorPattern() {return fgkITSSensName;}
cf457606 115 static const char* GetITSsegmentationFileName() {return fgITSsegmFileName.Data();}
75875328 116 static const char* GetDetTypeName(Int_t i);
cf457606 117 static void SetITSsegmentationFileName(const char* nm) {fgITSsegmFileName = nm;}
451f5018 118 static UInt_t ComposeDetTypeID(UInt_t segmId);
119 //
392efe73 120 static const char *ComposeSymNameITS();
121 static const char *ComposeSymNameLayer(Int_t lr);
122 static const char *ComposeSymNameLadder(Int_t lr, Int_t lad);
123 static const char *ComposeSymNameModule(Int_t lr, Int_t lad, int det);
124 //
02d6eccc 125 // hack to avoid using AliGeomManager
cf457606 126 Int_t LayerToVolUID(Int_t lay,int detInLay) const {return GetModuleIndex(lay,detInLay);}
127 static Int_t ModuleVolUID(Int_t mod) {return mod;}
128 //
129 protected:
130 void FetchMatrices();
af4b47c9 131 void CreateT2LMatrices();
cf457606 132 TGeoHMatrix* ExtractMatrixT2L(Int_t index) const;
133 TGeoHMatrix* ExtractMatrixSens(Int_t index) const;
134 Bool_t GetLayer(Int_t index,Int_t &lay,Int_t &index2) const;
135 TGeoPNEntry* GetPNEntry(Int_t index) const;
136 Int_t ExtractNumberOfDetectors(Int_t lay) const;
137 Int_t ExtractNumberOfLadders(Int_t lay) const;
138 Int_t ExtractLayerDetType(Int_t lay) const;
139 Int_t ExtractNumberOfLayers() const;
140 void BuildITS();
02d6eccc 141 //
451f5018 142 private:
af4b47c9 143 //
451f5018 144 //
145 Int_t fVersion; // ITS Version
146 Int_t fNLayers; // number of layers
af4b47c9 147 Int_t fNModules; // The total number of modules
451f5018 148 Int_t *fNLadders; //[fNLayers] Array of the number of ladders/layer(layer)
149 Int_t *fLrDetType; //[fNLayers] Array of layer detector types
150 Int_t *fNDetectors; //[fNLayers] Array of the number of detector/ladder(layer)
151 Int_t *fLastModIndex; //[fNLayers] max ID of the detctor in the layer
152 //
cf457606 153 TObjArray* fMatSens; // Sensor's matrices pointers in the geometry
154 TObjArray* fMatT2L; // Tracking to Local matrices pointers in the geometry
155 //
451f5018 156 static const char* fgkITSVolName; // ITS mother volume name
157 static const char* fgkITSLrName; // ITS Layer name
158 static const char* fgkITSLadName; // ITS Ladder name
159 static const char* fgkITSModName; // ITS Module name
160 static const char* fgkITSSensName; // ITS Sensor name
161 static const char* fgkITSDetTypeName[kNDetTypes]; // ITS upg detType Names
162 //
163 static TString fgITSsegmFileName; // file name for segmentations
164 //
165 ClassDef(AliITSUGeomTGeo, 1) // ITS geometry based on TGeo
166};
167
168//_____________________________________________________________________________________________
169inline const char *AliITSUGeomTGeo::GetSymName(Int_t lay,Int_t lad,Int_t det) const
170{
171 // sym name
172 return GetSymName(GetModuleIndex(lay,lad,det));
173}
174
175//_____________________________________________________________________________________________
176inline TGeoHMatrix* AliITSUGeomTGeo::GetMatrix(Int_t lay,Int_t lad,Int_t det) const
177{
178 // module current matrix
179 return GetMatrix(GetModuleIndex(lay,lad,det));
180}
181
182//_____________________________________________________________________________________________
183inline Bool_t AliITSUGeomTGeo::GetTranslation(Int_t lay,Int_t lad,Int_t det, Double_t t[3]) const
184{
185 // translation
186 return GetTranslation(GetModuleIndex(lay,lad,det),t);
187}
188
189//_____________________________________________________________________________________________
190inline Bool_t AliITSUGeomTGeo::GetRotation(Int_t lay,Int_t lad,Int_t det, Double_t r[9]) const
191{
192 // rot
193 return GetRotation(GetModuleIndex(lay,lad,det),r);
194}
195
196//_____________________________________________________________________________________________
197inline Bool_t AliITSUGeomTGeo::GetOrigMatrix(Int_t lay,Int_t lad,Int_t det, TGeoHMatrix &m) const
198{
199 // orig matrix
200 return GetOrigMatrix(GetModuleIndex(lay,lad,det),m);
201}
202
203//_____________________________________________________________________________________________
204inline Bool_t AliITSUGeomTGeo::GetOrigTranslation(Int_t lay,Int_t lad,Int_t det, Double_t t[3]) const
205{
206 // orig trans
207 return GetOrigTranslation(GetModuleIndex(lay,lad,det),t);
208}
209
210//_____________________________________________________________________________________________
211inline Bool_t AliITSUGeomTGeo::GetOrigRotation(Int_t lay,Int_t lad,Int_t det, Double_t r[9]) const
212{
213 // orig rot
214 return GetOrigRotation(GetModuleIndex(lay,lad,det),r);
215}
216
217//_____________________________________________________________________________________________
cf457606 218inline Bool_t AliITSUGeomTGeo::GetTrackingMatrix(Int_t lay,Int_t lad,Int_t det, TGeoHMatrix &m)
451f5018 219{
cf457606 220 // tracking mat
221 return GetTrackingMatrix(GetModuleIndex(lay,lad,det),m);
451f5018 222}
223
224//_____________________________________________________________________________________________
cf457606 225inline Int_t AliITSUGeomTGeo::GetLayerDetTypeID(Int_t lr) const
451f5018 226{
cf457606 227 // detector type ID of layer
228 return fLrDetType[lr];
451f5018 229}
230
231//_____________________________________________________________________________________________
cf457606 232inline Int_t AliITSUGeomTGeo::GetModuleDetTypeID(Int_t id) const
451f5018 233{
cf457606 234 // detector type ID of module
235 return GetLayerDetTypeID(GetLayer(id));
236}
451f5018 237
238//_____________________________________________________________________________________________
cf457606 239inline const TGeoHMatrix* AliITSUGeomTGeo::GetMatrixSens(Int_t index)
451f5018 240{
cf457606 241 // access global to sensor matrix
242 if (!fMatSens) FetchMatrices();
243 return (TGeoHMatrix*)fMatSens->At(index);
451f5018 244}
245
246//_____________________________________________________________________________________________
cf457606 247inline const TGeoHMatrix* AliITSUGeomTGeo::GetMatrixT2L(Int_t index)
451f5018 248{
cf457606 249 // access tracking to local matrix
250 if (!fMatT2L) FetchMatrices();
251 return (TGeoHMatrix*)fMatT2L->At(index);
252}
253
254//______________________________________________________________________
255inline void AliITSUGeomTGeo::LocalToGlobal(Int_t index,const Double_t *loc, Double_t *glob)
256{
257 // sensor local to global
258 GetMatrixSens(index)->LocalToMaster(loc,glob);
259}
260
261//______________________________________________________________________
262inline void AliITSUGeomTGeo::GlobalToLocal(Int_t index, const Double_t *glob, Double_t *loc)
263{
264 // global to sensor local
265 GetMatrixSens(index)->MasterToLocal(glob,loc);
266}
267
268//______________________________________________________________________
269inline void AliITSUGeomTGeo::LocalToGlobalVect(Int_t index, const Double_t *loc, Double_t *glob)
270{
271 // sensor local to global
272 GetMatrixSens(index)->LocalToMasterVect(loc,glob);
273}
274
275//______________________________________________________________________
276inline void AliITSUGeomTGeo::GlobalToLocalVect(Int_t index, const Double_t *glob, Double_t *loc)
277{
278 // global to sensor local
279 GetMatrixSens(index)->MasterToLocalVect(glob,loc);
451f5018 280}
281
282//_____________________________________________________________________________________________
cf457606 283inline void AliITSUGeomTGeo::LocalToGlobal(Int_t lay, Int_t lad, Int_t det,const Double_t *loc, Double_t *glob)
451f5018 284{
cf457606 285 // Local2Master (sensor)
286 LocalToGlobal(GetModuleIndex(lay,lad,det), loc, glob);
451f5018 287}
288
289//_____________________________________________________________________________________________
cf457606 290inline void AliITSUGeomTGeo::GlobalToLocal(Int_t lay, Int_t lad, Int_t det,const Double_t *glob, Double_t *loc)
451f5018 291{
cf457606 292 // master2local (sensor)
293 GlobalToLocal(GetModuleIndex(lay,lad,det), glob, loc);
294}
451f5018 295
75875328 296//_____________________________________________________________________________________________
297inline const char* AliITSUGeomTGeo::GetDetTypeName(Int_t i)
298{
299 if (i>=kNDetTypes) i/=kMaxSegmPerDetType; // full type is provided
300 return fgkITSDetTypeName[i];
301}
451f5018 302
303#endif