DD4hep  1.30.0
Detector Description Toolkit for High Energy Physics
Geant4ShapeConverter.cpp
Go to the documentation of this file.
1 //==========================================================================
2 // AIDA Detector description implementation
3 //--------------------------------------------------------------------------
4 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
5 // All rights reserved.
6 //
7 // For the licensing terms see $DD4hepINSTALL/LICENSE.
8 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
9 //
10 // Author : M.Frank
11 //
12 //==========================================================================
13 
14 // Framework include files
15 #include <DD4hep/Shapes.h>
16 #include <DD4hep/Printout.h>
18 
19 #include "Geant4ShapeConverter.h"
20 
21 // ROOT includes
22 #include <TClass.h>
23 #include <TGeoMatrix.h>
24 #include <TGeoBoolNode.h>
25 #include <TGeoScaledShape.h>
26 
27 // Geant4 include files
28 #include <G4Box.hh>
29 #include <G4Trd.hh>
30 #include <G4Tubs.hh>
31 #include <G4Trap.hh>
32 #include <G4Cons.hh>
33 #include <G4Hype.hh>
34 #include <G4Para.hh>
35 #include <G4Torus.hh>
36 #include <G4Sphere.hh>
37 #include <G4CutTubs.hh>
38 #include <G4Polycone.hh>
39 #include <G4Polyhedra.hh>
40 #include <G4Ellipsoid.hh>
41 #include <G4Paraboloid.hh>
42 #include <G4TwistedTubs.hh>
43 #include <G4GenericTrap.hh>
44 #include <G4ExtrudedSolid.hh>
45 #include <G4EllipticalTube.hh>
46 #include <G4TessellatedSolid.hh>
47 #include <G4TriangularFacet.hh>
48 #include <G4QuadrangularFacet.hh>
49 
50 // C/C++ include files
51 
52 using namespace dd4hep::detail;
53 
55 namespace dd4hep {
56 
58  namespace sim {
59 
60  static const double CM_2_MM = (CLHEP::centimeter/dd4hep::centimeter);
61 
63  template <typename T> G4VSolid* convertShape(const TGeoShape* shape) {
64  if ( shape ) {
65  except("convertShape","Unsupported shape: %s",shape->IsA()->GetName());
66  }
67  except("convertShape","Invalid shape conversion requested.");
68  return 0;
69  }
70 
71  template <> G4VSolid* convertShape<TGeoShapeAssembly>(const TGeoShape* /* shape */) {
72  return 0;
73  }
74 
75  template <> G4VSolid* convertShape<TGeoBBox>(const TGeoShape* shape) {
76  const TGeoBBox* sh = (const TGeoBBox*) shape;
77  return new G4Box(sh->GetName(), sh->GetDX() * CM_2_MM, sh->GetDY() * CM_2_MM, sh->GetDZ() * CM_2_MM);
78  }
79 
80  template <> G4VSolid* convertShape<TGeoTube>(const TGeoShape* shape) {
81  const TGeoTube* sh = (const TGeoTube*) shape;
82  return new G4Tubs(sh->GetName(), sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM, sh->GetDz() * CM_2_MM, 0, 2. * M_PI);
83  }
84 
85  template <> G4VSolid* convertShape<TGeoTubeSeg>(const TGeoShape* shape) {
86  const TGeoTubeSeg* sh = (const TGeoTubeSeg*) shape;
87  return new G4Tubs(sh->GetName(), sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM, sh->GetDz() * CM_2_MM,
88  sh->GetPhi1() * DEGREE_2_RAD, (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD);
89  }
90 
91  template <> G4VSolid* convertShape<TGeoCtub>(const TGeoShape* shape) {
92  const TGeoCtub* sh = (const TGeoCtub*) shape;
93  const Double_t* ln = sh->GetNlow();
94  const Double_t* hn = sh->GetNhigh();
95  G4ThreeVector lowNorm (ln[0], ln[1], ln[2]);
96  G4ThreeVector highNorm(hn[0], hn[1], hn[2]);
97  return new G4CutTubs(sh->GetName(),
98  sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM, sh->GetDz() * CM_2_MM,
99  sh->GetPhi1() * DEGREE_2_RAD, (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD,
100  std::move(lowNorm), std::move(highNorm));
101  }
102 
103  template <> G4VSolid* convertShape<TGeoEltu>(const TGeoShape* shape) {
104  const TGeoEltu* sh = (const TGeoEltu*) shape;
105  return new G4EllipticalTube(sh->GetName(),sh->GetA() * CM_2_MM, sh->GetB() * CM_2_MM, sh->GetDz() * CM_2_MM);
106  }
107 
108  template <> G4VSolid* convertShape<TwistedTubeObject>(const TGeoShape* shape) {
109  const TwistedTubeObject* sh = (const TwistedTubeObject*) shape;
110  if ( std::fabs(std::fabs(sh->GetNegativeEndZ()) - std::fabs(sh->GetPositiveEndZ())) < 1e-10 ) {
111  return new G4TwistedTubs(sh->GetName(),sh->GetPhiTwist() * DEGREE_2_RAD,
112  sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM,
113  sh->GetPositiveEndZ() * CM_2_MM,
114  sh->GetNsegments(),
115  (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD);
116  }
117  return new G4TwistedTubs(sh->GetName(),sh->GetPhiTwist() * DEGREE_2_RAD,
118  sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM,
120  sh->GetNsegments(),
121  (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD);
122  }
123 
124  template <> G4VSolid* convertShape<TGeoTrd1>(const TGeoShape* shape) {
125  const TGeoTrd1* sh = (const TGeoTrd1*) shape;
126  return new G4Trd(sh->GetName(),
127  sh->GetDx1() * CM_2_MM, sh->GetDx2() * CM_2_MM,
128  sh->GetDy() * CM_2_MM, sh->GetDy() * CM_2_MM,
129  sh->GetDz() * CM_2_MM);
130  }
131 
132  template <> G4VSolid* convertShape<TGeoTrd2>(const TGeoShape* shape) {
133  const TGeoTrd2* sh = (const TGeoTrd2*) shape;
134  return new G4Trd(sh->GetName(),
135  sh->GetDx1() * CM_2_MM, sh->GetDx2() * CM_2_MM,
136  sh->GetDy1() * CM_2_MM, sh->GetDy2() * CM_2_MM,
137  sh->GetDz() * CM_2_MM);
138  }
139 
140  template <> G4VSolid* convertShape<TGeoHype>(const TGeoShape* shape) {
141  const TGeoHype* sh = (const TGeoHype*) shape;
142  return new G4Hype(sh->GetName(), sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM,
143  sh->GetStIn() * DEGREE_2_RAD, sh->GetStOut() * DEGREE_2_RAD,
144  sh->GetDz() * CM_2_MM);
145  }
146 
147  template <> G4VSolid* convertShape<TGeoArb8>(const TGeoShape* shape) {
148  std::vector<G4TwoVector> vertices;
149  TGeoArb8* sh = (TGeoArb8*) shape;
150  Double_t* vtx_xy = sh->GetVertices();
151  for ( std::size_t i=0; i<8; ++i, vtx_xy +=2 )
152  vertices.emplace_back(vtx_xy[0] * CM_2_MM, vtx_xy[1] * CM_2_MM);
153  return new G4GenericTrap(sh->GetName(), sh->GetDz() * CM_2_MM, vertices);
154  }
155 
156  template <> G4VSolid* convertShape<TGeoPara>(const TGeoShape* shape) {
157  const auto* sh = static_cast<const TGeoPara*>(shape);
158  return new G4Para(sh->GetName(),
159  sh->GetX() * CM_2_MM, sh->GetY() * CM_2_MM, sh->GetZ() * CM_2_MM,
160  sh->GetAlpha() * DEGREE_2_RAD, sh->GetTheta() * DEGREE_2_RAD,
161  sh->GetPhi() * DEGREE_2_RAD);
162  }
163 
164  template <> G4VSolid* convertShape<TGeoXtru>(const TGeoShape* shape) {
165  const TGeoXtru* sh = (const TGeoXtru*) shape;
166  std::size_t nz = sh->GetNz();
167  std::vector<G4ExtrudedSolid::ZSection> z;
168  z.reserve(nz);
169  for(std::size_t i=0; i<nz; ++i) {
170  z.emplace_back(G4ExtrudedSolid::ZSection(sh->GetZ(i) * CM_2_MM, {sh->GetXOffset(i), sh->GetYOffset(i)}, sh->GetScale(i)));
171  }
172  std::size_t np = sh->GetNvert();
173  std::vector<G4TwoVector> polygon;
174  polygon.reserve(np);
175  for(std::size_t i=0; i<np; ++i) {
176  polygon.emplace_back(sh->GetX(i) * CM_2_MM,sh->GetY(i) * CM_2_MM);
177  }
178  return new G4ExtrudedSolid(sh->GetName(), polygon, z);
179  }
180 
181  template <> G4VSolid* convertShape<TGeoPgon>(const TGeoShape* shape) {
182  const TGeoPgon* sh = (const TGeoPgon*) shape;
183  std::vector<double> rmin, rmax, z;
184  for (Int_t i = 0; i < sh->GetNz(); ++i) {
185  rmin.emplace_back(sh->GetRmin(i) * CM_2_MM);
186  rmax.emplace_back(sh->GetRmax(i) * CM_2_MM);
187  z.emplace_back(sh->GetZ(i) * CM_2_MM);
188  }
189  return new G4Polyhedra(sh->GetName(), sh->GetPhi1() * DEGREE_2_RAD, sh->GetDphi() * DEGREE_2_RAD,
190  sh->GetNedges(), sh->GetNz(), &z[0], &rmin[0], &rmax[0]);
191  }
192 
193  template <> G4VSolid* convertShape<TGeoPcon>(const TGeoShape* shape) {
194  const TGeoPcon* sh = (const TGeoPcon*) shape;
195  std::vector<double> rmin, rmax, z;
196  for (Int_t i = 0; i < sh->GetNz(); ++i) {
197  rmin.emplace_back(sh->GetRmin(i) * CM_2_MM);
198  rmax.emplace_back(sh->GetRmax(i) * CM_2_MM);
199  z.emplace_back(sh->GetZ(i) * CM_2_MM);
200  }
201  return new G4Polycone(sh->GetName(), sh->GetPhi1() * DEGREE_2_RAD, sh->GetDphi() * DEGREE_2_RAD,
202  sh->GetNz(), &z[0], &rmin[0], &rmax[0]);
203  }
204 
205  template <> G4VSolid* convertShape<TGeoCone>(const TGeoShape* shape) {
206  const TGeoCone* sh = (const TGeoCone*) shape;
207  return new G4Cons(sh->GetName(), sh->GetRmin1() * CM_2_MM, sh->GetRmax1() * CM_2_MM, sh->GetRmin2() * CM_2_MM,
208  sh->GetRmax2() * CM_2_MM, sh->GetDz() * CM_2_MM, 0.0, 2.*M_PI);
209  }
210 
211  template <> G4VSolid* convertShape<TGeoConeSeg>(const TGeoShape* shape) {
212  const TGeoConeSeg* sh = (const TGeoConeSeg*) shape;
213  return new G4Cons(sh->GetName(), sh->GetRmin1() * CM_2_MM, sh->GetRmax1() * CM_2_MM,
214  sh->GetRmin2() * CM_2_MM, sh->GetRmax2() * CM_2_MM,
215  sh->GetDz() * CM_2_MM,
216  sh->GetPhi1() * DEGREE_2_RAD, (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD);
217  }
218 
219  template <> G4VSolid* convertShape<TGeoParaboloid>(const TGeoShape* shape) {
220  const TGeoParaboloid* sh = (const TGeoParaboloid*) shape;
221  return new G4Paraboloid(sh->GetName(), sh->GetDz() * CM_2_MM, sh->GetRlo() * CM_2_MM, sh->GetRhi() * CM_2_MM);
222  }
223 
224  template <> G4VSolid* convertShape<TGeoSphere>(const TGeoShape* shape) {
225  const TGeoSphere* sh = (const TGeoSphere*) shape;
226  return new G4Sphere(sh->GetName(), sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM,
227  sh->GetPhi1() * DEGREE_2_RAD, (sh->GetPhi2()-sh->GetPhi1()) * DEGREE_2_RAD,
228  sh->GetTheta1() * DEGREE_2_RAD, (sh->GetTheta2()- sh->GetTheta1()) * DEGREE_2_RAD);
229  }
230 
231  template <> G4VSolid* convertShape<TGeoTorus>(const TGeoShape* shape) {
232  const TGeoTorus* sh = (const TGeoTorus*) shape;
233  return new G4Torus(sh->GetName(), sh->GetRmin() * CM_2_MM, sh->GetRmax() * CM_2_MM, sh->GetR() * CM_2_MM,
234  sh->GetPhi1() * DEGREE_2_RAD, sh->GetDphi() * DEGREE_2_RAD);
235  }
236 
237  template <> G4VSolid* convertShape<TGeoTrap>(const TGeoShape* shape) {
238  const TGeoTrap* sh = (const TGeoTrap*) shape;
239  return new G4Trap(sh->GetName(), sh->GetDz() * CM_2_MM, sh->GetTheta() * DEGREE_2_RAD, sh->GetPhi() * DEGREE_2_RAD,
240  sh->GetH1() * CM_2_MM, sh->GetBl1() * CM_2_MM, sh->GetTl1() * CM_2_MM, sh->GetAlpha1() * DEGREE_2_RAD,
241  sh->GetH2() * CM_2_MM, sh->GetBl2() * CM_2_MM, sh->GetTl2() * CM_2_MM, sh->GetAlpha2() * DEGREE_2_RAD);
242  }
243 
244  template <> G4VSolid* convertShape<G4GenericTrap>(const TGeoShape* shape) {
245  std::vector<G4TwoVector> vertices;
246  TGeoTrap* sh = (TGeoTrap*) shape;
247  Double_t* vtx_xy = sh->GetVertices();
248  for ( std::size_t i=0; i<8; ++i, vtx_xy +=2 )
249  vertices.emplace_back(vtx_xy[0] * CM_2_MM, vtx_xy[1] * CM_2_MM);
250  return new G4GenericTrap(sh->GetName(), sh->GetDz() * CM_2_MM, vertices);
251  }
252 
253  template <> G4VSolid* convertShape<TGeoTessellated>(const TGeoShape* shape) {
254  TGeoTessellated* sh = (TGeoTessellated*) shape;
255  G4TessellatedSolid* g4 = new G4TessellatedSolid(sh->GetName());
256  int num_facet = sh->GetNfacets();
257 
258  printout(DEBUG,"TessellatedSolid","+++ %s> Converting %d facets", sh->GetName(), num_facet);
259  for(int i=0; i<num_facet; ++i) {
260  const TGeoFacet& facet = sh->GetFacet(i);
261  int nv = facet.GetNvert();
262 #if ROOT_VERSION_CODE >= ROOT_VERSION(6,31,1)
263  const auto& v0 = sh->GetVertex(facet[0]);
264  const auto& v1 = sh->GetVertex(facet[1]);
265  const auto& v2 = sh->GetVertex(facet[2]);
266 #else
267  const auto& v0 = sh->GetVertex(facet.GetVertexIndex(0));
268  const auto& v1 = sh->GetVertex(facet.GetVertexIndex(1));
269  const auto& v2 = sh->GetVertex(facet.GetVertexIndex(2));
270 #endif
271  G4VFacet* g4f = 0;
272  if ( nv == 3 ) {
273  g4f = new G4TriangularFacet(G4ThreeVector(v0.x() * CM_2_MM, v0.y() * CM_2_MM, v0.z() * CM_2_MM),
274  G4ThreeVector(v1.x() * CM_2_MM, v1.y() * CM_2_MM, v1.z() * CM_2_MM),
275  G4ThreeVector(v2.x() * CM_2_MM, v2.y() * CM_2_MM, v2.z() * CM_2_MM),
276  ABSOLUTE);
277  }
278  else if ( nv == 4 ) {
279 #if ROOT_VERSION_CODE >= ROOT_VERSION(6,31,1)
280  const auto& v3 = sh->GetVertex(facet[3]);
281 #else
282  const auto& v3 = sh->GetVertex(facet.GetVertexIndex(3));
283 #endif
284  g4f = new G4QuadrangularFacet(G4ThreeVector(v0.x() * CM_2_MM, v0.y() * CM_2_MM, v0.z() * CM_2_MM),
285  G4ThreeVector(v1.x() * CM_2_MM, v1.y() * CM_2_MM, v1.z() * CM_2_MM),
286  G4ThreeVector(v2.x() * CM_2_MM, v2.y() * CM_2_MM, v2.z() * CM_2_MM),
287  G4ThreeVector(v3.x() * CM_2_MM, v3.y() * CM_2_MM, v3.z() * CM_2_MM),
288  ABSOLUTE);
289  }
290  else {
291  except("TGeoTessellated", "Tessellated shape [%s] has facet with wrong number of vertices: %d",
292  sh->GetName(), nv);
293  }
294  g4->AddFacet(g4f);
295  }
296  g4->SetSolidClosed(sh->IsClosedBody());
297  return g4;
298  }
299 
300  } // End namespace sim
301 } // End namespace dd4hep
dd4hep::sim::convertShape< TGeoTrd2 >
G4VSolid * convertShape< TGeoTrd2 >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:132
dd4hep::sim::convertShape< TGeoPgon >
G4VSolid * convertShape< TGeoPgon >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:181
dd4hep::TwistedTubeObject::GetNsegments
int GetNsegments() const
Access the number of segments.
Definition: ShapesInterna.h:66
dd4hep::sim::convertShape< TGeoTrd1 >
G4VSolid * convertShape< TGeoTrd1 >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:124
dd4hep::sim::convertShape< TGeoPcon >
G4VSolid * convertShape< TGeoPcon >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:193
dd4hep::sim::convertShape< TGeoHype >
G4VSolid * convertShape< TGeoHype >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:140
M_PI
#define M_PI
Definition: Handle.h:33
dd4hep::sim::convertShape< TGeoBBox >
G4VSolid * convertShape< TGeoBBox >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:75
dd4hep::TwistedTubeObject::GetPhiTwist
double GetPhiTwist() const
Access twist angle.
Definition: ShapesInterna.h:60
dd4hep::sim::convertShape< TGeoTessellated >
G4VSolid * convertShape< TGeoTessellated >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:253
dd4hep::TwistedTubeObject::GetNegativeEndZ
double GetNegativeEndZ() const
Access the negative z.
Definition: ShapesInterna.h:62
CM_2_MM
#define CM_2_MM
Definition: ParticleActors.cpp:37
dd4hep::sim::convertShape< TGeoEltu >
G4VSolid * convertShape< TGeoEltu >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:103
Geant4ShapeConverter.h
dd4hep::TwistedTubeObject::GetPositiveEndZ
double GetPositiveEndZ() const
Access the positive z.
Definition: ShapesInterna.h:64
dd4hep::sim::convertShape< TGeoXtru >
G4VSolid * convertShape< TGeoXtru >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:164
dd4hep::sim::convertShape< TGeoCone >
G4VSolid * convertShape< TGeoCone >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:205
dd4hep::TwistedTubeObject
Concrete object implementation for the Header handle.
Definition: ShapesInterna.h:29
dd4hep::sim::convertShape< TwistedTubeObject >
G4VSolid * convertShape< TwistedTubeObject >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:108
dd4hep::sim::convertShape
G4VSolid * convertShape(const TGeoShape *shape)
Convert a specific TGeo shape into the geant4 equivalent.
Definition: Geant4ShapeConverter.cpp:63
dd4hep::detail
DD4hep internal namespace.
Definition: Alignments.h:32
dd4hep::sim::convertShape< G4GenericTrap >
G4VSolid * convertShape< G4GenericTrap >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:244
ShapesInterna.h
dd4hep::sim::convertShape< TGeoParaboloid >
G4VSolid * convertShape< TGeoParaboloid >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:219
dd4hep::sim::convertShape< TGeoPara >
G4VSolid * convertShape< TGeoPara >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:156
Shapes.h
dd4hep::sim::convertShape< TGeoTrap >
G4VSolid * convertShape< TGeoTrap >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:237
dd4hep::sim::convertShape< TGeoArb8 >
G4VSolid * convertShape< TGeoArb8 >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:147
dd4hep::sim::convertShape< TGeoShapeAssembly >
G4VSolid * convertShape< TGeoShapeAssembly >(const TGeoShape *)
Definition: Geant4ShapeConverter.cpp:71
dd4hep::sim::convertShape< TGeoConeSeg >
G4VSolid * convertShape< TGeoConeSeg >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:211
dd4hep::sim::convertShape< TGeoCtub >
G4VSolid * convertShape< TGeoCtub >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:91
dd4hep::sim::convertShape< TGeoTubeSeg >
G4VSolid * convertShape< TGeoTubeSeg >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:85
DEGREE_2_RAD
#define DEGREE_2_RAD
Definition: Handle.h:29
dd4hep::sim::convertShape< TGeoTorus >
G4VSolid * convertShape< TGeoTorus >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:231
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
TGeoConeSeg
Class of the ROOT toolkit. See http://root.cern.ch/root/htmldoc/ClassIndex.html.
Definition: ROOTClasses.h:17
Printout.h
dd4hep::sim::convertShape< TGeoSphere >
G4VSolid * convertShape< TGeoSphere >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:224
dd4hep::sim::convertShape< TGeoTube >
G4VSolid * convertShape< TGeoTube >(const TGeoShape *shape)
Definition: Geant4ShapeConverter.cpp:80