Skip to content

Commit f696f3f

Browse files
additional features added to TRK geometry
1 parent 68f955b commit f696f3f

File tree

2 files changed

+150
-2
lines changed

2 files changed

+150
-2
lines changed

Detectors/Upgrades/ALICE3/TRK/base/include/TRKBase/GeometryTGeo.h

Lines changed: 22 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -15,6 +15,7 @@
1515
#include <memory>
1616
#include <DetectorsCommonDataFormats/DetMatrixCache.h>
1717
#include "DetectorsCommonDataFormats/DetID.h"
18+
#include <iostream>
1819

1920
namespace o2
2021
{
@@ -58,7 +59,7 @@ class GeometryTGeo : public o2::detectors::DetMatrixCache
5859
int extractNumberOfLayersMLOT();
5960
int extractNumberOfLayersVD() const;
6061
int extractNumberOfPetalsVD() const;
61-
int extractNumberOfActivePartsVD() const;
62+
int extractNumberOfActivePartsVD() const;
6263
int extractNumberOfDisksVD() const;
6364
int extractNumberOfChipsPerPetalVD() const;
6465
int extractNumberOfStavesMLOT(int lay) const;
@@ -84,6 +85,22 @@ class GeometryTGeo : public o2::detectors::DetMatrixCache
8485
int getPetalCase(int index) const;
8586
int getDisk(int index) const;
8687

88+
void defineMLOTSensors();
89+
int getBarrelLayer(int) const;
90+
91+
//sensor ref X and alpha for ML & OT
92+
void extractSensorXAlphaMLOT(int, float&, float&);
93+
94+
//cache for tracking frames (ML & OT)
95+
bool isTrackingFrameCachedMLOT() const {return !mCacheRefXMLOT.empty();}
96+
void fillTrackingFramesCacheMLOT();
97+
98+
float getSensorRefAlphaMLOT(int index) const {return mCacheRefAlphaMLOT[index];}
99+
float getSensorXMLOT(int index) const {return mCacheRefXMLOT[index];}
100+
101+
//create matrix for tracking to local frame for MLOT
102+
TGeoHMatrix& createT2LMatrixMLOT(int);
103+
87104
/// This routine computes the chip index number from the subDetID, petal, disk, layer, stave /// TODO: retrieve also from chip when chips will be available
88105
/// \param int subDetID The subdetector ID, 0 for VD, 1 for MLOT
89106
/// \param int petalcase The petal case number for VD, from 0 to 3
@@ -173,6 +190,10 @@ class GeometryTGeo : public o2::detectors::DetMatrixCache
173190

174191
bool mOwner = true; //! is it owned by the singleton?
175192

193+
std::vector<int> sensorsMLOT;
194+
std::vector<float> mCacheRefXMLOT; ///cache for X of ML and OT
195+
std::vector<float> mCacheRefAlphaMLOT; ///cache for sensor ref alpha ML and OT
196+
176197
private:
177198
static std::unique_ptr<o2::trk::GeometryTGeo> sInstance;
178199
};

Detectors/Upgrades/ALICE3/TRK/base/src/GeometryTGeo.cxx

Lines changed: 128 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -12,6 +12,7 @@
1212
#include <TRKBase/GeometryTGeo.h>
1313
#include <TGeoManager.h>
1414
#include "TRKBase/SegmentationChip.h"
15+
#include <TMath.h>
1516

1617
using Segmentation = o2::trk::SegmentationChip;
1718

@@ -94,7 +95,7 @@ void GeometryTGeo::Build(int loadTrans)
9495
int numberOfChipsTotal = 0;
9596

9697
/// filling the information for the VD
97-
for (int i = 0; i < mNumberOfPetalsVD; i++) {
98+
for (int i = 0; i < mNumberOfPetalsVD; i++) {
9899
mNumberOfChipsPerPetalVD[i] = extractNumberOfChipsPerPetalVD();
99100
numberOfChipsTotal += mNumberOfChipsPerPetalVD[i];
100101
mLastChipIndex[i] = numberOfChipsTotal - 1;
@@ -111,6 +112,8 @@ void GeometryTGeo::Build(int loadTrans)
111112

112113
setSize(numberOfChipsTotal); /// temporary, number of chips = number of staves and active parts
113114
fillMatrixCache(loadTrans);
115+
defineMLOTSensors();
116+
fillTrackingFramesCacheMLOT();
114117
}
115118

116119
//__________________________________________________________________________
@@ -342,6 +345,32 @@ TGeoHMatrix* GeometryTGeo::extractMatrixSensor(int index) const
342345
return &matTmp;
343346
}
344347

348+
//__________________________________________________________________________
349+
void GeometryTGeo::defineMLOTSensors()
350+
{
351+
for(int i=0; i<mSize; i++){
352+
if (getSubDetID(i) == 0) {
353+
continue;
354+
}
355+
sensorsMLOT.push_back(i);
356+
}
357+
}
358+
359+
//__________________________________________________________________________
360+
void GeometryTGeo::fillTrackingFramesCacheMLOT()
361+
{
362+
//fill for every sensor of ML & OT its tracking frame parameters
363+
if(!isTrackingFrameCachedMLOT() && !sensorsMLOT.empty()){
364+
size_t newSize = sensorsMLOT.size();
365+
mCacheRefXMLOT.resize(newSize);
366+
mCacheRefAlphaMLOT.resize(newSize);
367+
for(int i=0; i< newSize; i++){
368+
int sensorId = sensorsMLOT[i];
369+
extractSensorXAlphaMLOT(sensorId, mCacheRefXMLOT[i], mCacheRefAlphaMLOT[i]);
370+
}
371+
}
372+
}
373+
345374
//__________________________________________________________________________
346375
void GeometryTGeo::fillMatrixCache(int mask)
347376
{
@@ -364,6 +393,21 @@ void GeometryTGeo::fillMatrixCache(int mask)
364393
}
365394
}
366395

396+
//build T2L matrices for ML & OT !! VD is yet to be implemented once its geometry will be more refined
397+
if ((mask & o2::math_utils::bit2Mask(o2::math_utils::TransformType::T2L)) && !getCacheT2L().isFilled()) {
398+
LOGP(info, "Loading {} T2L matrices from TGeo for ML & OT", getName());
399+
if(sensorsMLOT.size()){
400+
int m_Size = sensorsMLOT.size();
401+
auto& cacheT2L = getCacheT2L();
402+
cacheT2L.setSize(m_Size);
403+
for(int i=0; i<m_Size; i++) {
404+
int sensorID = sensorsMLOT[i];
405+
TGeoHMatrix& hm = createT2LMatrixMLOT(sensorID);
406+
cacheT2L.setMatrix(Mat3D(hm), i); //here, sensorIDs from 0 to 374, sensorIDs shifted to 36 !
407+
}
408+
}
409+
}
410+
367411
// TODO: build matrices for the cases T2L, T2G and T2GRot when needed
368412
}
369413

@@ -910,5 +954,88 @@ void GeometryTGeo::Print(Option_t*) const
910954
std::cout << "]" << std::endl;
911955
}
912956

957+
//__________________________________________________________________________
958+
int GeometryTGeo::getBarrelLayer(int chipID) const
959+
{
960+
//for barrel layers only,
961+
//so it would be consistent with number of layers i.e. from 0 to 10,
962+
//starting from VD0 to OT10;
963+
//skip the disks;
964+
965+
int subDetID = getSubDetID(chipID);
966+
int subLayerID = getLayer(chipID);
967+
968+
if (subDetID < 0 || subDetID > 1) {
969+
LOG(error) << "getBarrelLayer(): Invalid subDetID for barrel: " << subDetID
970+
<< ". Expected values are 0 or 1.";
971+
return -1;
972+
}
973+
974+
if (subLayerID < 0 || subLayerID > 7) {
975+
LOG(error) << "getBarrelLayer(): Invalid subLayerID for barrel: " << subDetID
976+
<< ". Expected values are between 0 and 7.";
977+
return -1;
978+
}
979+
980+
const int baseOffsets[] = {0, 3};
981+
982+
return baseOffsets[subDetID] + subLayerID;
983+
984+
}
985+
986+
//__________________________________________________________________________
987+
void GeometryTGeo::extractSensorXAlphaMLOT(int chipID, float& x, float& alp)
988+
{
989+
//works for ML and OT only, a.k.a flat sensors !!!
990+
double locA[3] = {-100., 0., 0.}, locB[3] = {100., 0., 0.}, gloA[3], gloB[3];
991+
double xp{0}, yp{0};
992+
993+
if(getSubDetID(chipID) == 0){
994+
995+
LOG(error) << "extractSensorXAlphaMLOT(): VD layers are not supported yet! chipID = " << chipID;
996+
return;
997+
998+
} else { //flat sensors, ML and OT
999+
const TGeoHMatrix* matL2G = extractMatrixSensor(chipID);
1000+
matL2G->LocalToMaster(locA, gloA);
1001+
matL2G->LocalToMaster(locB, gloB);
1002+
double dx = gloB[0] - gloA[0], dy = gloB[1] - gloA[1];
1003+
double t = (gloB[0] * dx + gloB[1] * dy) / (dx * dx + dy * dy);
1004+
xp = gloB[0] - dx * t;
1005+
yp = gloB[1] - dy * t;
1006+
}
1007+
1008+
alp = std::atan2(yp, xp);
1009+
x = std::hypot(xp, yp);
1010+
o2::math_utils::bringTo02Pi(alp);
1011+
1012+
///TODO:
1013+
//once the VD segmentation is done, VD should be added
1014+
}
1015+
1016+
//__________________________________________________________________________
1017+
TGeoHMatrix& GeometryTGeo::createT2LMatrixMLOT(int chipID)
1018+
{
1019+
//works only for ML & OT
1020+
//for VD is yet to be implemented once we have more refined geometry
1021+
if(getSubDetID(chipID) == 0){
1022+
1023+
LOG(error) << "createT2LMatrixMLOT(): VD layers are not supported yet! chipID = " << chipID
1024+
<< "returning dummy values! ";
1025+
static TGeoHMatrix dummy;
1026+
return dummy;
1027+
1028+
} else {
1029+
static TGeoHMatrix t2l;
1030+
t2l.Clear();
1031+
float alpha = getSensorRefAlphaMLOT(chipID);
1032+
t2l.RotateZ(alpha * TMath::RadToDeg());
1033+
const TGeoHMatrix* matL2G = extractMatrixSensor(chipID);
1034+
const TGeoHMatrix& matL2Gi = matL2G->Inverse();
1035+
t2l.MultiplyLeft(&matL2Gi);
1036+
return t2l;
1037+
}
1038+
}
1039+
9131040
} // namespace trk
9141041
} // namespace o2

0 commit comments

Comments
 (0)