PandaRoot
PndMdtIGeometry.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 #ifndef PNDMDTIGEOMETRY_H
14 #define PNDMDTIGEOMETRY_H 1
15 
16 #include "TVector3.h"
17 #include "TGeoVolume.h"
18 #include <vector>
19 #include "TF1.h"
20 #include <map>
21 #include <set>
22 #include "PndMdtID.h"
23 #include "TGeoMatrix.h"
24 #include <iostream>
25 
26 using std::ostream;
27 ostream &operator<<(ostream &os, const TVector3 &v3);
28 
29 class PndMdtIGeometry : public TNamed {
30 
31  public:
32  static PndMdtIGeometry *Instance();
37 
38  void SetVerbose(Int_t _v) { fVerbose = _v; }
39  void AddSensor(TString _v) { fSensorSet.insert(_v); }
40  // version and Init() should be invoked first before using access functions
41  // now only support ROOT geometry version 0, 1, 2
42  Bool_t Init();
43  //================================================================
44  // Access to Tube Geometery for given parameters
45  Bool_t GetTubeCenter(Int_t iDetId, TVector3 &pos) const
46  {
47  std::map<Int_t, InfoType>::const_iterator it = fGeoMap.find(iDetId);
48  if (it == fGeoMap.end())
49  return kFALSE;
50  pos = it->second.Position;
51  return kTRUE;
52  }
53  //
54  Bool_t GetTubeLength(Int_t iDetId, Double_t &len) const
55  {
56  std::map<Int_t, InfoType>::const_iterator it = fGeoMap.find(iDetId);
57  if (it == fGeoMap.end())
58  return kFALSE;
59  len = it->second.Length;
60  return kTRUE;
61  }
62  //
63  Bool_t GetStripLength(Int_t iDetId, Double_t &len) const
64  {
65  std::map<Int_t, InfoType>::const_iterator it = fGeoMap.find(iDetId);
66  if (it == fGeoMap.end())
67  return kFALSE;
68  len = it->second.Length;
69  return kTRUE;
70  }
71  Bool_t GetStripCenter(Int_t iDetId, TVector3 &pos) const
72  {
73  std::map<Int_t, InfoType>::const_iterator it = fGeoMap.find(iDetId);
74  if (it == fGeoMap.end())
75  return kFALSE;
76  pos = it->second.Position;
77  return kTRUE;
78  }
79  Bool_t GetLayerBoundary(Int_t iLayer, TVector3 &LBPos, TVector3 &RTPos) const
80  {
81  std::map<Int_t, LayerBoundary>::const_iterator it = fLayerInfoMap.find(iLayer);
82  if (it == fLayerInfoMap.end())
83  return kFALSE;
84  LBPos = TVector3(it->second.minX, it->second.minY, it->second.minZ);
85  RTPos = TVector3(it->second.maxX, it->second.maxY, it->second.maxZ);
86  return kTRUE;
87  }
88  Bool_t MasterToLocal(Int_t iDetId, TVector3 &master, TVector3 &local)
89  {
90  std::map<Int_t, InfoType>::const_iterator it = fGeoMap.find(iDetId);
91  if (it == fGeoMap.end())
92  return kFALSE;
93  TVector3 center = it->second.Position;
94  TGeoMatrix *fMtrx = it->second.Matrix;
95  if (!fMtrx)
96  return kFALSE;
97  // fMtrx->Print();
98  // fMtrx->MasterToLocal(&center[0], &local[0]);
99  // std::cout<<"local 1 "<<local<<std::endl;
100  fMtrx->MasterToLocal(&master[0], &local[0]);
101  // std::cout<<"local 2 "<<local<<std::endl;
102  return kTRUE;
103  }
104  Bool_t MapWireToStrip(Int_t iDetId, const TVector3 &fEntryPos, Int_t &fStripDetID)
105  {
106  static const Double_t cSTRIPWIDTH = 1.; // cent meter;
107  static const Double_t cGAPOFSTRIPS = 0.0; // cent meter;
108  Short_t iMod = PndMdtID::Module(iDetId);
109  Short_t iSec = PndMdtID::Sector(iDetId);
110  Short_t iLayer = PndMdtID::Layer(iDetId);
111  Int_t layerID = PndMdtID::LayerID(iMod, iSec, iLayer);
112  LayerInfoMapIter it = fLayerInfoMap.find(layerID);
113  if (it == fLayerInfoMap.end())
114  return kFALSE;
115  Int_t fStripNo;
116  if (iMod == 1)
117  fStripNo = (fEntryPos.Z() - it->second.minZ) / (cSTRIPWIDTH + cGAPOFSTRIPS);
118  else if (iMod == 2) // endcap
119  {
120  if (iLayer < 2)
121  fStripNo = (fEntryPos.Y() - it->second.minY) / (cSTRIPWIDTH + cGAPOFSTRIPS);
122  else
123  fStripNo = (fEntryPos.X() - it->second.minX) / (cSTRIPWIDTH + cGAPOFSTRIPS);
124  } else
125  fStripNo = (fEntryPos.X() - it->second.minX) / (cSTRIPWIDTH + cGAPOFSTRIPS);
126  fStripDetID = PndMdtID::Identifier(iMod, iSec, iLayer, fStripNo);
127  return kTRUE;
128  }
129  void Print() const;
130 
131  private:
132  std::set<TString> fSensorSet;
133  Bool_t FindSensorinPath(TString path)
134  {
135  Bool_t ret = kFALSE;
136  std::set<TString>::iterator it = fSensorSet.begin();
137  std::set<TString>::iterator end = fSensorSet.end();
138  while (it != end) {
139  if (path.Contains(*it)) {
140  ret = kTRUE;
141  break;
142  }
143  ++it;
144  }
145  return ret;
146  }
147 
148  private:
149  Int_t fVerbose;
150  Bool_t fGoodGeometry;
151  TString fStateTip;
152 
153  public:
154  struct InfoType {
155  TVector3 Position;
156  Double_t Length;
157  TGeoMatrix *Matrix;
158  };
159 
160  private:
161  std::map<Int_t, InfoType> fGeoMap;
162  std::map<Int_t, TGeoMatrix *> fMatrixMap;
163 
164  public:
165  struct LayerBoundary {
166  Double_t minX;
167  Double_t minY;
168  Double_t minZ;
169  Double_t maxX;
170  Double_t maxY;
171  Double_t maxZ;
172  };
173 
174  private:
175  std::map<Int_t, LayerBoundary> fLayerInfoMap;
176  typedef std::map<Int_t, LayerBoundary>::iterator LayerInfoMapIter;
177  typedef std::pair<Int_t, LayerBoundary> LayerInfoMapValue;
178 
179  Int_t fBarrelVersion;
180  Int_t fEndcapVersion;
181  Int_t fForwardVersion;
182  Int_t fMuonFilterVersion;
183  Int_t fNumofMdt;
184  Int_t fNumofStrip;
185  Bool_t fInitialized;
186 
187  void LoadGeometry(TGeoNode *node);
188  Bool_t GetStripInfo();
189  Bool_t GetGeometryInfoV2(const TString &fFullPath, const TString &Name);
190  void InsertLayerInfo(Int_t layerID, const TVector3 &pos, const TVector3 &lwh);
191 
192  Bool_t AddTubeInfo(Int_t detID, const TVector3 &center, Double_t wirelen, TGeoMatrix *matrix);
193  Bool_t AddStripInfo(Int_t detID, TVector3 stripPos, Double_t stripLen, TGeoMatrix *matrix = 0);
194  Bool_t AddMatrix(Int_t layerID, TGeoMatrix *matrix);
195 
196  static PndMdtIGeometry *fInstance;
197 
198  ClassDef(PndMdtIGeometry, 1);
199 };
200 
201 #endif
ostream & operator<<(ostream &os, const TVector3 &v3)
static Short_t Module(Int_t detID)
Definition: PndMdtID.h:28
Bool_t MasterToLocal(Int_t iDetId, TVector3 &master, TVector3 &local)
static Short_t Layer(Int_t detID)
Definition: PndMdtID.h:30
Bool_t GetTubeLength(Int_t iDetId, Double_t &len) const
Bool_t GetLayerBoundary(Int_t iLayer, TVector3 &LBPos, TVector3 &RTPos) const
static Int_t LayerID(Int_t iMod, Int_t iOct, Int_t iLayer)
Definition: PndMdtID.h:25
static Int_t Identifier(Int_t iMod, Int_t iOct, Int_t iLayer, Int_t iBox, Int_t iWire)
Definition: PndMdtID.h:21
Bool_t MapWireToStrip(Int_t iDetId, const TVector3 &fEntryPos, Int_t &fStripDetID)
Bool_t GetStripCenter(Int_t iDetId, TVector3 &pos) const
static Short_t Sector(Int_t detID)
Definition: PndMdtID.h:29
Bool_t GetTubeCenter(Int_t iDetId, TVector3 &pos) const
void AddSensor(TString _v)
void SetVerbose(Int_t _v)
Bool_t GetStripLength(Int_t iDetId, Double_t &len) const
void Print() const
static PndMdtIGeometry * Instance()
basic_ostream< char, char_traits< char > > ostream