DD4hep  1.31.0
Detector Description Toolkit for High Energy Physics
LCDDConverter.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 includes
15 #include "LCDDConverter.h"
16 #include <DD4hep/Plugins.h>
17 #include <DD4hep/Printout.h>
18 #include <DD4hep/Volumes.h>
19 #include <DD4hep/FieldTypes.h>
20 #include <DD4hep/DD4hepUnits.h>
21 #include <DD4hep/Segmentations.h>
24 #include <XML/DocumentHandler.h>
25 
26 // ROOT includes
27 #include <TROOT.h>
28 #include <TColor.h>
29 #include <TGeoShape.h>
30 
31 #include <TGeoArb8.h>
32 #include <TGeoBoolNode.h>
33 #include <TGeoCompositeShape.h>
34 #include <TGeoCone.h>
35 #include <TGeoEltu.h>
36 #include <TGeoHype.h>
37 #include <TGeoMatrix.h>
38 #include <TGeoParaboloid.h>
39 #include <TGeoPara.h>
40 #include <TGeoPcon.h>
41 #include <TGeoPgon.h>
42 #include <TGeoShapeAssembly.h>
43 #include <TGeoSphere.h>
44 #include <TGeoTorus.h>
45 #include <TGeoTrd1.h>
46 #include <TGeoTrd2.h>
47 #include <TGeoTube.h>
48 #include <TGeoScaledShape.h>
49 
50 #include <TGeoNode.h>
51 #include <TClass.h>
52 #include <TMath.h>
53 
55 #include <fstream>
56 #include <iostream>
57 #include <sstream>
58 
59 using namespace dd4hep;
60 using namespace dd4hep::detail;
61 
62 namespace {
63  typedef Position XYZRotation;
64 
65  XYZRotation getXYZangles(const Double_t* r) {
66  Double_t cosb = std::sqrt(r[0]*r[0] + r[1]*r[1]);
67  if (cosb > 0.00001) {
68  return XYZRotation(atan2(r[5], r[8]), atan2(-r[2], cosb), atan2(r[1], r[0]));
69  }
70  return XYZRotation(atan2(-r[7], r[4]),atan2(-r[2], cosb),0);
71  }
72 
73 #if 0
74  XYZRotation getXYZangles(const Double_t* rotationMatrix) {
75  Double_t a, b, c;
76  Double_t rad = 1.0; // RAD by default! 180.0 / TMath::ACos(-1.0);
77  const Double_t *r = rotationMatrix;
78  Double_t cosb = TMath::Sqrt(r[0] * r[0] + r[1] * r[1]);
79  if (cosb > 0.00001) {
80  a = TMath::ATan2(r[5], r[8]) * rad;
81  b = TMath::ATan2(-r[2], cosb) * rad;
82  c = TMath::ATan2(r[1], r[0]) * rad;
83  }
84  else {
85  a = TMath::ATan2(-r[7], r[4]) * rad;
86  b = TMath::ATan2(-r[2], cosb) * rad;
87  c = 0;
88  }
89  XYZRotation rr(a, b, c);
90  std::cout << " X:" << a << " " << rr.X() << " Y:" << b << " " << rr.Y() << " Z:" << c << " " << rr.Z()
91  << " lx:" << r[0] << " ly:" << r[4] << " lz:" << r[8] << std::endl;
92  return XYZRotation(a, b, c);
93  }
94 #endif
95 
96  bool is_volume(const TGeoVolume* volume) {
97  Volume v(volume);
98  return v.data() != 0;
99  }
100  bool is_placement(PlacedVolume node) {
101  return node.data() != 0;
102  }
103 
104  std::string genName(const std::string& n) { return n; }
105  std::string genName(const std::string& n, const void* ptr) {
106  std::string nn = genName(n);
107  char text[32];
108  ::snprintf(text,sizeof(text),"%p",ptr);
109  nn += "_";
110  nn += text;
111  return nn;
112  }
113 }
114 
115 void LCDDConverter::GeometryInfo::check(const std::string& name, const TNamed* _n, std::map<std::string, const TNamed*>& _m) const {
116  std::map<std::string, const TNamed*>::const_iterator i = _m.find(name);
117  if (i != _m.end()) {
118  const char* isa = _n ? _n->IsA()->GetName() : (*i).second ? (*i).second->IsA()->GetName() : "Unknown";
119  std::cout << isa << "(position): duplicate entry with name:" << name << " " << (void*) _n << " " << (void*) (*i).second << std::endl;
120  }
121  _m.insert(make_pair(name, _n));
122 }
123 
126  : m_detDesc(description), m_dataPtr(0) {
127 }
128 
130  if (m_dataPtr)
131  delete m_dataPtr;
132  m_dataPtr = 0;
133 }
134 
136 xml_h LCDDConverter::handleElement(const std::string& /* name */, Atom element) const {
137  GeometryInfo& geo = data();
138  xml_h e = geo.xmlElements[element];
139  if (!e) {
140  int Z = element->Z();
141  double A = element->A();
142  xml_elt_t atom(geo.doc, _U(atom));
143  // If we got an unphysical material (Z<1 or A<1)
144  // We pretend it is hydrogen and force Z=1 or A=1.00794 g/mole
145  geo.doc_materials.append(e = xml_elt_t(geo.doc, _U(element)));
146  e.append(atom);
147  e.setAttr(_U(name), element->GetName());
148  e.setAttr(_U(formula), element->GetName());
149  e.setAttr(_U(Z), Z>0 ? Z : 1);
150  atom.setAttr(_U(type), "A");
151  atom.setAttr(_U(unit), "g/mol");
152  atom.setAttr(_U(value), A>0.99 ? A : 1.00794 /* *(g/mole) */);
153  geo.xmlElements[element] = e;
154  }
155  return e;
156 }
157 
159 xml_h LCDDConverter::handleMaterial(const std::string& name, Material medium) const {
160  GeometryInfo& geo = data();
161  xml_h mat = geo.xmlMaterials[medium];
162  if (!mat) {
163  xml_h obj;
164  TGeoMaterial* geo_mat = medium->GetMaterial();
165  double d = geo_mat->GetDensity(); //*(gram/cm3);
166  if (d < 1e-10) d = 1e-10;
167  mat = xml_elt_t(geo.doc, _U(material));
168  mat.setAttr(_U(name), medium->GetName());
169  mat.append(obj = xml_elt_t(geo.doc, _U(D)));
170  obj.setAttr(_U(value), d /* *(g/cm3) */);
171  obj.setAttr(_U(unit), "g/cm3");
172  obj.setAttr(_U(type), "density");
173 
174  geo.checkMaterial(name, medium);
175 
176  if (geo_mat->IsMixture()) {
177  TGeoMixture *mix = (TGeoMixture*)geo_mat;
178  const double *wmix = mix->GetWmixt();
179  const int *nmix = mix->GetNmixt();
180  double sum = 0e0;
181  for (int i = 0, n = mix->GetNelements(); i < n; i++) {
182  TGeoElement *elt = mix->GetElement(i);
183  handleElement(elt->GetName(), Atom(elt));
184  sum += wmix[i];
185  }
186  for (int i = 0, n = mix->GetNelements(); i < n; i++) {
187  TGeoElement *elt = mix->GetElement(i);
188  //std::string formula = elt->GetTitle() + std::string("_elm");
189  if (nmix) {
190  mat.append(obj = xml_elt_t(geo.doc, _U(composite)));
191  obj.setAttr(_U(n), nmix[i]);
192  }
193  else {
194  mat.append(obj = xml_elt_t(geo.doc, _U(fraction)));
195  obj.setAttr(_U(n), wmix[i] / sum);
196  }
197  obj.setAttr(_U(ref), elt->GetName());
198  }
199  }
200  else if ( name != "dummy" ) {
201  // Do not exactly know where dummy comes from,
202  // but it causes havoc in Geant4 later
203  TGeoElement *elt = geo_mat->GetElement(0);
204  printout(INFO,"++ Converting non mixing material: %s",name.c_str());
205  xml_elt_t atom(geo.doc, _U(atom));
206  handleElement(elt->GetName(), Atom(elt));
207  mat.append(atom);
208  mat.setAttr(_U(Z), geo_mat->GetZ());
209  atom.setAttr(_U(type), "A");
210  atom.setAttr(_U(unit), "g/mol");
211  atom.setAttr(_U(value), geo_mat->GetA() /* *(g/mole) */);
212  }
213  geo.doc_materials.append(mat);
214  geo.xmlMaterials[medium] = mat;
215  }
216  return mat;
217 }
218 
220 xml_h LCDDConverter::handleSolid(const std::string& name, const TGeoShape* shape) const {
221  GeometryInfo& geo = data();
222  SolidMap::iterator sit = geo.xmlSolids.find(shape);
223  if (!shape) {
224  // This is an invalid volume. Let's pray returning nothing will work,
225  // and the non-existing solid is also nowhere referenced in the GDML.
226  return xml_h(0);
227  }
228  else if (sit != geo.xmlSolids.end()) {
229  // The solidis already registered. Return the reference
230  return (*sit).second;
231  }
232  else if (shape->IsA() == TGeoShapeAssembly::Class()) {
233  // Assemblies have no shape in GDML. Hence, return nothing.
234  return xml_h(0);
235  }
236  else {
237  xml_h solid(0);
238  xml_h zplane(0);
239  TClass* isa = shape->IsA();
240  std::string shape_name = shape->GetName(); //genName(shape->GetName(),shape);
241  geo.checkShape(name, shape);
242  if (isa == TGeoBBox::Class()) {
243  const TGeoBBox* sh = (const TGeoBBox*) shape;
244  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(box)));
245  solid.setAttr(_U(name), Unicode(shape_name));
246  solid.setAttr(_U(x), 2 * sh->GetDX());
247  solid.setAttr(_U(y), 2 * sh->GetDY());
248  solid.setAttr(_U(z), 2 * sh->GetDZ());
249  solid.setAttr(_U(lunit), "cm");
250  }
251  else if (isa == TGeoTube::Class()) {
252  const TGeoTube* sh = (const TGeoTube*) shape;
253  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(tube)));
254  solid.setAttr(_U(name), Unicode(shape_name));
255  solid.setAttr(_U(rmin), sh->GetRmin());
256  solid.setAttr(_U(rmax), sh->GetRmax());
257  solid.setAttr(_U(z), 2 * sh->GetDz());
258  solid.setAttr(_U(startphi), 0e0);
259  solid.setAttr(_U(deltaphi), 360.0);
260  solid.setAttr(_U(aunit), "deg");
261  solid.setAttr(_U(lunit), "cm");
262  }
263  else if (isa == TGeoTubeSeg::Class()) {
264  const TGeoTubeSeg* sh = (const TGeoTubeSeg*) shape;
265  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(tube)));
266  solid.setAttr(_U(name), Unicode(shape_name));
267  solid.setAttr(_U(rmin), sh->GetRmin());
268  solid.setAttr(_U(rmax), sh->GetRmax());
269  solid.setAttr(_U(z), 2 * sh->GetDz()); // Full zlen in GDML, half zlen in TGeo
270  solid.setAttr(_U(startphi), sh->GetPhi1());
271  solid.setAttr(_U(deltaphi), sh->GetPhi2());
272  solid.setAttr(_U(aunit), "deg");
273  solid.setAttr(_U(lunit), "cm");
274  }
275  else if (isa == TGeoEltu::Class()) {
276  const TGeoEltu* sh = (const TGeoEltu*) shape;
277  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(eltube)));
278  solid.setAttr(_U(name), Unicode(shape_name));
279  solid.setAttr(_U(dx), sh->GetA());
280  solid.setAttr(_U(dy), sh->GetB());
281  solid.setAttr(_U(dz), sh->GetDz());
282  solid.setAttr(_U(lunit), "cm");
283  }
284  else if (isa == TGeoTrd1::Class()) {
285  const TGeoTrd1* sh = (const TGeoTrd1*) shape;
286  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(trd)));
287  solid.setAttr(_U(name), Unicode(shape_name));
288  solid.setAttr(_U(x1), 2 * sh->GetDx1());
289  solid.setAttr(_U(x2), 2 * sh->GetDx2());
290  solid.setAttr(_U(y1), 2 * sh->GetDy());
291  solid.setAttr(_U(y2), 2 * sh->GetDy());
292  solid.setAttr(_U(z), 2 * sh->GetDz()); // Full zlen in GDML, half zlen in TGeo
293  solid.setAttr(_U(lunit), "cm");
294  }
295  else if (isa == TGeoTrd2::Class()) {
296  const TGeoTrd2* sh = (const TGeoTrd2*) shape;
297  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(trd)));
298  solid.setAttr(_U(name), Unicode(shape_name));
299  solid.setAttr(_U(x1), 2 * sh->GetDx1());
300  solid.setAttr(_U(x2), 2 * sh->GetDx2());
301  solid.setAttr(_U(y1), 2 * sh->GetDy1());
302  solid.setAttr(_U(y2), 2 * sh->GetDy2());
303  solid.setAttr(_U(z), 2 * sh->GetDz()); // Full zlen in GDML, half zlen in TGeo
304  solid.setAttr(_U(lunit), "cm");
305  }
306  else if (isa == TGeoHype::Class()) {
307  const TGeoHype* sh = (const TGeoHype*) shape;
308  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(hype)));
309  solid.setAttr(_U(name), Unicode(shape_name));
310  solid.setAttr(_U(rmin), sh->GetRmin());
311  solid.setAttr(_U(rmax), sh->GetRmax());
312  solid.setAttr(Unicode("inst"), sh->GetStIn());
313  solid.setAttr(_U(outst), sh->GetStOut());
314  solid.setAttr(_U(z), sh->GetDz()); // Full zlen in GDML, half zlen in TGeo
315  solid.setAttr(_U(aunit), "deg");
316  solid.setAttr(_U(lunit), "cm");
317  }
318  else if (isa == TGeoPgon::Class()) {
319  const TGeoPgon* sh = (const TGeoPgon*) shape;
320  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(polyhedra)));
321  solid.setAttr(_U(name), Unicode(shape_name));
322  solid.setAttr(_U(startphi), sh->GetPhi1());
323  solid.setAttr(_U(deltaphi), sh->GetDphi());
324  solid.setAttr(_U(numsides), sh->GetNedges());
325  solid.setAttr(_U(aunit), "deg");
326  solid.setAttr(_U(lunit), "cm");
327  for (Int_t i = 0; i < sh->GetNz(); ++i) {
328  zplane = xml_elt_t(geo.doc, _U(zplane));
329  zplane.setAttr(_U(z), sh->GetZ(i));
330  zplane.setAttr(_U(rmin), sh->GetRmin(i));
331  zplane.setAttr(_U(rmax), sh->GetRmax(i));
332  solid.append(zplane);
333  }
334  }
335  else if (isa == TGeoPcon::Class()) {
336  const TGeoPcon* sh = (const TGeoPcon*) shape;
337  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(polycone)));
338  solid.setAttr(_U(name), Unicode(shape_name));
339  solid.setAttr(_U(startphi), sh->GetPhi1());
340  solid.setAttr(_U(deltaphi), sh->GetDphi());
341  solid.setAttr(_U(aunit), "deg");
342  solid.setAttr(_U(lunit), "cm");
343  for (Int_t i = 0; i < sh->GetNz(); ++i) {
344  zplane = xml_elt_t(geo.doc, _U(zplane));
345  zplane.setAttr(_U(z), sh->GetZ(i));
346  zplane.setAttr(_U(rmin), sh->GetRmin(i));
347  zplane.setAttr(_U(rmax), sh->GetRmax(i));
348  solid.append(zplane);
349  }
350  solid.setAttr(_U(lunit), "cm");
351  }
352  else if (isa == TGeoCone::Class()) {
353  const TGeoCone* sh = (const TGeoCone*) shape;
354  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(cone)));
355  solid.setAttr(_U(name), Unicode(shape_name));
356  solid.setAttr(_U(z), 2 * sh->GetDz());
357  solid.setAttr(_U(rmin1), sh->GetRmin1());
358  solid.setAttr(_U(rmax1), sh->GetRmax1());
359  solid.setAttr(_U(rmin2), sh->GetRmin2());
360  solid.setAttr(_U(rmax2), sh->GetRmax2());
361  solid.setAttr(_U(startphi), 0e0);
362  solid.setAttr(_U(deltaphi), 360.0);
363  solid.setAttr(_U(aunit), "deg");
364  solid.setAttr(_U(lunit), "cm");
365  }
366  else if (isa == TGeoConeSeg::Class()) {
367  const TGeoConeSeg* sh = (const TGeoConeSeg*) shape;
368  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(cone)));
369  solid.setAttr(_U(name), Unicode(shape_name));
370  solid.setAttr(_U(z), 2*sh->GetDz());
371  solid.setAttr(_U(rmin1), sh->GetRmin1());
372  solid.setAttr(_U(rmin2), sh->GetRmin2());
373  solid.setAttr(_U(rmax1), sh->GetRmax1());
374  solid.setAttr(_U(rmax2), sh->GetRmax2());
375  solid.setAttr(_U(startphi), sh->GetPhi1());
376  solid.setAttr(_U(deltaphi), sh->GetPhi2() - sh->GetPhi1());
377  solid.setAttr(_U(aunit), "deg");
378  solid.setAttr(_U(lunit), "cm");
379  }
380  else if (isa == TGeoParaboloid::Class()) {
381  const TGeoParaboloid* sh = (const TGeoParaboloid*) shape;
382  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(paraboloid)));
383  solid.setAttr(_U(name), Unicode(shape_name));
384  solid.setAttr(_U(rlo), sh->GetRlo());
385  solid.setAttr(_U(rhi), sh->GetRhi());
386  solid.setAttr(_U(dz), sh->GetDz());
387  solid.setAttr(_U(lunit), "cm");
388  }
389 #if 0
390  else if (isa == TGeoEllipsoid::Class()) {
391  const TGeoEllipsoid* sh = (const TGeoEllipsoid*) shape;
392  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(ellipsoid)));
393  solid.setAttr(_U(lunit), "cm");
394  }
395 #endif
396  else if (isa == TGeoSphere::Class()) {
397  const TGeoSphere* sh = (const TGeoSphere*) shape;
398  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(sphere)));
399  solid.setAttr(_U(name), Unicode(shape_name));
400  solid.setAttr(_U(rmin), sh->GetRmin());
401  solid.setAttr(_U(rmax), sh->GetRmax());
402  solid.setAttr(_U(startphi), sh->GetPhi1());
403  solid.setAttr(_U(deltaphi), (sh->GetPhi2() - sh->GetPhi1()));
404  solid.setAttr(_U(starttheta), sh->GetTheta1());
405  solid.setAttr(_U(deltatheta), (sh->GetTheta2() - sh->GetTheta1()));
406  solid.setAttr(_U(aunit), "deg");
407  solid.setAttr(_U(lunit), "cm");
408  }
409  else if (isa == TGeoTorus::Class()) {
410  const TGeoTorus* sh = (const TGeoTorus*) shape;
411  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(torus)));
412  solid.setAttr(_U(name), Unicode(shape_name));
413  solid.setAttr(_U(rtor), sh->GetR());
414  solid.setAttr(_U(rmin), sh->GetRmin());
415  solid.setAttr(_U(rmax), sh->GetRmax());
416  solid.setAttr(_U(startphi), sh->GetPhi1());
417  solid.setAttr(_U(deltaphi), sh->GetDphi());
418  solid.setAttr(_U(aunit), "deg");
419  solid.setAttr(_U(lunit), "cm");
420  }
421  else if (isa == TGeoTrap::Class()) {
422  const TGeoTrap* sh = (const TGeoTrap*) shape;
423  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(trap)));
424  solid.setAttr(_U(name), Unicode(shape_name));
425  solid.setAttr(_U(z), 2 * sh->GetDz()); // Full zlen in GDML, half zlen in TGeo
426  solid.setAttr(_U(x1), 2 * sh->GetBl1());
427  solid.setAttr(_U(x2), 2 * sh->GetTl1());
428  solid.setAttr(_U(x3), 2 * sh->GetBl2());
429  solid.setAttr(_U(x4), 2 * sh->GetTl2());
430  solid.setAttr(_U(y1), 2 * sh->GetH1());
431  solid.setAttr(_U(y2), 2 * sh->GetH2());
432  solid.setAttr(_U(alpha1), sh->GetAlpha1());
433  solid.setAttr(_U(alpha2), sh->GetAlpha2());
434  solid.setAttr(_U(theta), sh->GetTheta());
435  solid.setAttr(_U(phi), sh->GetPhi());
436  solid.setAttr(_U(aunit), "deg");
437  solid.setAttr(_U(lunit), "cm");
438  }
439  else if (isa == TGeoPara::Class()) {
440  const TGeoPara* sh = (const TGeoPara*) shape;
441  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(para)));
442  solid.setAttr(_U(name), Unicode(shape_name));
443  solid.setAttr(_U(x), sh->GetX());
444  solid.setAttr(_U(y), sh->GetY());
445  solid.setAttr(_U(z), sh->GetZ());
446  solid.setAttr(_U(alpha), sh->GetAlpha());
447  solid.setAttr(_U(theta), sh->GetTheta());
448  solid.setAttr(_U(phi), sh->GetPhi());
449  solid.setAttr(_U(aunit), "deg");
450  solid.setAttr(_U(lunit), "cm");
451  }
452  else if (isa == TGeoArb8::Class()) {
453  TGeoArb8* sh = (TGeoArb8*) shape;
454  const double* vtx = sh->GetVertices();
455  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(arb8)));
456  solid.setAttr(_U(name), Unicode(shape_name));
457  solid.setAttr(_U(v1x), vtx[0]);
458  solid.setAttr(_U(v1y), vtx[1]);
459  solid.setAttr(_U(v2x), vtx[2]);
460  solid.setAttr(_U(v2y), vtx[3]);
461  solid.setAttr(_U(v3x), vtx[4]);
462  solid.setAttr(_U(v3y), vtx[5]);
463  solid.setAttr(_U(v4x), vtx[6]);
464  solid.setAttr(_U(v4y), vtx[7]);
465  solid.setAttr(_U(v5x), vtx[8]);
466  solid.setAttr(_U(v5y), vtx[9]);
467  solid.setAttr(_U(v6x), vtx[10]);
468  solid.setAttr(_U(v6y), vtx[11]);
469  solid.setAttr(_U(v7x), vtx[12]);
470  solid.setAttr(_U(v7y), vtx[13]);
471  solid.setAttr(_U(v8x), vtx[14]);
472  solid.setAttr(_U(v8y), vtx[15]);
473  solid.setAttr(_U(dz), sh->GetDz());
474  solid.setAttr(_U(lunit), "cm");
475  }
476  else if (isa == TGeoScaledShape::Class()) {
477  TGeoScaledShape* sh = (TGeoScaledShape*) shape;
478  const double* vals = sh->GetScale()->GetScale();
479  Solid s_sh(sh->GetShape());
480  handleSolid(s_sh.name(), s_sh.ptr());
481  geo.doc_solids.append(solid = xml_elt_t(geo.doc, _U(scale)));
482  solid.setAttr(_U(name), Unicode(shape_name));
483  solid.setAttr(_U(shape), s_sh.name());
484  solid.setAttr(_U(x), vals[0]);
485  solid.setAttr(_U(y), vals[1]);
486  solid.setAttr(_U(z), vals[2]);
487  solid.setAttr(_U(aunit), "deg");
488  solid.setAttr(_U(lunit), "cm");
489  }
490  else if (isa == TGeoCompositeShape::Class() ||
491  isa == TGeoUnion::Class() ||
492  isa == TGeoIntersection::Class() ||
493  isa == TGeoSubtraction::Class() ) {
494  const TGeoCompositeShape* sh = (const TGeoCompositeShape*) shape;
495  const TGeoBoolNode* boolean = sh->GetBoolNode();
496  TGeoBoolNode::EGeoBoolType oper = boolean->GetBooleanOperator();
497  TGeoMatrix* rm = boolean->GetRightMatrix();
498  TGeoMatrix* lm = boolean->GetLeftMatrix();
499  TGeoShape* ls = boolean->GetLeftShape();
500  TGeoShape* rs = boolean->GetRightShape();
501  xml_h left = handleSolid(ls->GetName(), ls);
502  xml_h right = handleSolid(rs->GetName(), rs);
503  xml_h first_solid(0), second_solid(0);
504  if (!left) {
505  throw std::runtime_error("G4Converter: No left Detector Solid present for composite shape:" + name);
506  }
507  if (!right) {
508  throw std::runtime_error("G4Converter: No right Detector Solid present for composite shape:" + name);
509  }
510 
511  //specific case!
512  //Ellipsoid tag preparing
513  //if left == TGeoScaledShape AND right == TGeoBBox
514  // AND if TGeoScaledShape->GetShape == TGeoSphere
515  if (strcmp(ls->ClassName(), "TGeoScaledShape") == 0 &&
516  strcmp(rs->ClassName(), "TGeoBBox") == 0) {
517  if (strcmp(((TGeoScaledShape *)ls)->GetShape()->ClassName(), "TGeoSphere") == 0) {
518  if (oper == TGeoBoolNode::kGeoIntersection) {
519  TGeoScaledShape* lls = (TGeoScaledShape *)ls;
520  TGeoBBox* rrs = (TGeoBBox*)rs;
521  solid = xml_elt_t(geo.doc,Unicode("ellipsoid"));
522  solid.setAttr(_U(name), Unicode(shape_name));
523  double sx = lls->GetScale()->GetScale()[0];
524  double sy = lls->GetScale()->GetScale()[1];
525  double radius = ((TGeoSphere *)lls->GetShape())->GetRmax();
526  double dz = rrs->GetDZ();
527  double zorig = rrs->GetOrigin()[2];
528  double zcut2 = dz + zorig;
529  double zcut1 = 2 * zorig - zcut2;
530  solid.setAttr(Unicode("ax"),sx * radius);
531  solid.setAttr(Unicode("by"),sy * radius);
532  solid.setAttr(Unicode("cz"),radius);
533  solid.setAttr(Unicode("zcut1"),zcut1);
534  solid.setAttr(Unicode("zcut2"),zcut2);
535  solid.setAttr(_U(lunit), "cm");
536  return data().xmlSolids[shape] = solid;
537  }
538  }
539  }
540 
541  if ( oper == TGeoBoolNode::kGeoSubtraction )
542  solid = xml_elt_t(geo.doc,_U(subtraction));
543  else if ( oper == TGeoBoolNode::kGeoUnion )
544  solid = xml_elt_t(geo.doc,_U(union));
545  else if ( oper == TGeoBoolNode::kGeoIntersection )
546  solid = xml_elt_t(geo.doc,_U(intersection));
547 
548  xml_h obj;
549  std::string lnam = left.attr<std::string>(_U(name));
550  std::string rnam = right.attr<std::string>(_U(name));
551 
552  geo.doc_solids.append(solid);
553  solid.append(first_solid = xml_elt_t(geo.doc, _U(first)));
554  solid.setAttr(_U(name), Unicode(shape_name));
555  first_solid.setAttr(_U(ref), lnam);
556  const double *tr = lm->GetTranslation();
557 
558  if ((tr[0] != 0.0) || (tr[1] != 0.0) || (tr[2] != 0.0)) {
559  first_solid.append(obj = xml_elt_t(geo.doc, _U(firstposition)));
560  obj.setAttr(_U(name), name+"_"+lnam+"_pos");
561  obj.setAttr(_U(x), tr[0]);
562  obj.setAttr(_U(y), tr[1]);
563  obj.setAttr(_U(z), tr[2]);
564  obj.setAttr(_U(unit), "cm");
565  }
566  if (lm->IsRotation()) {
567  TGeoMatrix const & linv = lm->Inverse();
568  XYZRotation rot = getXYZangles(linv.GetRotationMatrix());
569  if ((rot.X() != 0.0) || (rot.Y() != 0.0) || (rot.Z() != 0.0)) {
570  first_solid.append(obj = xml_elt_t(geo.doc, _U(firstrotation)));
571  obj.setAttr(_U(name), name+"_"+lnam+"_rot");
572  obj.setAttr(_U(x), rot.X());
573  obj.setAttr(_U(y), rot.Y());
574  obj.setAttr(_U(z), rot.Z());
575  obj.setAttr(_U(unit), "rad");
576  }
577  }
578  tr = rm->GetTranslation();
579  solid.append(second_solid = xml_elt_t(geo.doc, _U(second)));
580  second_solid.setAttr(_U(ref), rnam);
581  if ((tr[0] != 0.0) || (tr[1] != 0.0) || (tr[2] != 0.0)) {
582  xml_ref_t pos = handlePosition(rnam+"_pos", rm);
583  solid.setRef(_U(positionref), pos.name());
584  }
585  if (rm->IsRotation()) {
586  TGeoMatrix const & rinv = rm->Inverse();
587  XYZRotation rot = getXYZangles(rinv.GetRotationMatrix());
588  if ((rot.X() != 0.0) || (rot.Y() != 0.0) || (rot.Z() != 0.0)) {
589  xml_ref_t xml_rot = handleRotation(rnam+"_rot", &rinv);
590  solid.setRef(_U(rotationref), xml_rot.name());
591  }
592  }
593  }
594  if (!solid) {
595  std::string err = "Failed to handle unknown solid shape:" + name + " of type " + std::string(shape->IsA()->GetName());
596  throw std::runtime_error(err);
597  }
598  return data().xmlSolids[shape] = solid;
599  }
600 }
601 
603 xml_h LCDDConverter::handlePosition(const std::string& name, const TGeoMatrix* trafo) const {
604  GeometryInfo& geo = data();
605  xml_h pos = geo.xmlPositions[trafo];
606  if (!pos) {
607  const double* tr = trafo->GetTranslation();
608  if (tr[0] != 0.0 || tr[1] != 0.0 || tr[2] != 0.0) {
609  std::string gen_name = genName(name,trafo);
610  geo.checkPosition(gen_name, trafo);
611  geo.doc_define.append(pos = xml_elt_t(geo.doc, _U(position)));
612  pos.setAttr(_U(name), gen_name);
613  pos.setAttr(_U(x), tr[0]);
614  pos.setAttr(_U(y), tr[1]);
615  pos.setAttr(_U(z), tr[2]);
616  pos.setAttr(_U(unit), "cm");
617  }
618  else if (geo.identity_pos) {
619  pos = geo.identity_pos;
620  }
621  else {
622  geo.doc_define.append(geo.identity_pos = xml_elt_t(geo.doc, _U(position)));
623  geo.identity_pos.setAttr(_U(name), "identity_pos");
624  geo.identity_pos.setAttr(_U(x), 0);
625  geo.identity_pos.setAttr(_U(y), 0);
626  geo.identity_pos.setAttr(_U(z), 0);
627  geo.identity_pos.setAttr(_U(unit), "cm");
628  pos = geo.identity_pos;
629  geo.checkPosition("identity_pos", 0);
630  }
631  geo.xmlPositions[trafo] = pos;
632  }
633  return pos;
634 }
635 
637 xml_h LCDDConverter::handleRotation(const std::string& name, const TGeoMatrix* trafo) const {
638  GeometryInfo& geo = data();
639  xml_h rot = geo.xmlRotations[trafo];
640  if (!rot) {
641  XYZRotation r = getXYZangles(trafo->GetRotationMatrix());
642  if (!(r.X() == 0.0 && r.Y() == 0.0 && r.Z() == 0.0)) {
643  std::string gen_name = genName(name,trafo);
644  geo.checkRotation(gen_name, trafo);
645  geo.doc_define.append(rot = xml_elt_t(geo.doc, _U(rotation)));
646  rot.setAttr(_U(name), gen_name);
647  rot.setAttr(_U(x), r.X());
648  rot.setAttr(_U(y), r.Y());
649  rot.setAttr(_U(z), r.Z());
650  rot.setAttr(_U(unit), "rad");
651  }
652  else if (geo.identity_rot) {
653  rot = geo.identity_rot;
654  }
655  else {
656  geo.doc_define.append(geo.identity_rot = xml_elt_t(geo.doc, _U(rotation)));
657  geo.identity_rot.setAttr(_U(name), "identity_rot");
658  geo.identity_rot.setAttr(_U(x), 0);
659  geo.identity_rot.setAttr(_U(y), 0);
660  geo.identity_rot.setAttr(_U(z), 0);
661  geo.identity_rot.setAttr(_U(unit), "rad");
662  rot = geo.identity_rot;
663  geo.checkRotation("identity_rot", 0);
664  }
665  geo.xmlRotations[trafo] = rot;
666  }
667  return rot;
668 }
669 
671 xml_h LCDDConverter::handleVolume(const std::string& /* name */, Volume volume) const {
672  GeometryInfo& geo = data();
673  xml_h vol = geo.xmlVolumes[volume];
674  if (!vol) {
675  const TGeoVolume* v = volume;
676  Volume _v(v);
677  std::string n = genName(v->GetName(),v);
678  TGeoMedium* medium = v->GetMedium();
679  TGeoShape* sh = v->GetShape();
680  xml_ref_t sol = handleSolid(sh->GetName(), sh);
681 
682  geo.checkVolume(n, volume);
683  if (v->IsAssembly()) {
684  vol = xml_elt_t(geo.doc, _U(assembly));
685  vol.setAttr(_U(name), n);
686  }
687  else {
688  if (!sol)
689  throw std::runtime_error("G4Converter: No Detector Solid present for volume:" + n);
690  else if (!m)
691  throw std::runtime_error("G4Converter: No Detector material present for volume:" + n);
692 
693  vol = xml_elt_t(geo.doc, _U(volume));
694  vol.setAttr(_U(name), n);
695  if (m) {
696  std::string mat_name = medium->GetName();
697  xml_ref_t med = handleMaterial(mat_name, Material(medium));
698  vol.setRef(_U(materialref), med.name());
699  }
700  vol.setRef(_U(solidref), sol.name());
701  }
702  geo.doc_structure.append(vol);
703  geo.xmlVolumes[v] = vol;
704  const TObjArray* dau = const_cast<TGeoVolume*>(v)->GetNodes();
705  if (dau && dau->GetEntries() > 0) {
706  for (Int_t i = 0, n_dau = dau->GetEntries(); i < n_dau; ++i) {
707  TGeoNode* node = reinterpret_cast<TGeoNode*>(dau->At(i));
708  handlePlacement(node->GetName(), node);
709  }
710  }
711  if (geo.doc_header && is_volume(volume)) {
712  Region reg = _v.region();
713  LimitSet lim = _v.limitSet();
714  VisAttr vis = _v.visAttributes();
716  if (det.isValid()) {
717  xml_ref_t xml_data = handleSensitive(det.name(), det);
718  vol.setRef(_U(sdref), xml_data.name());
719  }
720  if (reg.isValid()) {
721  xml_ref_t xml_data = handleRegion(reg.name(), reg);
722  vol.setRef(_U(regionref), xml_data.name());
723  }
724  if (lim.isValid()) {
725  xml_ref_t xml_data = handleLimitSet(lim.name(), lim);
726  vol.setRef(_U(limitsetref), xml_data.name());
727  }
728  if (vis.isValid()) {
729  xml_ref_t xml_data = handleVis(vis.name(), vis);
730  vol.setRef(_U(visref), xml_data.name());
731  }
732  }
733  }
734  return vol;
735 }
736 
738 xml_h LCDDConverter::handleVolumeVis(const std::string& /* name */, const TGeoVolume* volume) const {
739  GeometryInfo& geo = data();
740  xml_h vol = geo.xmlVolumes[volume];
741  if (!vol) {
742  const TGeoVolume* v = volume;
743  Volume _v(volume);
744  if (is_volume(volume)) {
745  VisAttr vis = _v.visAttributes();
746  if (vis.isValid()) {
747  geo.doc_structure.append(vol = xml_elt_t(geo.doc, _U(volume)));
748  vol.setAttr(_U(name), v->GetName());
749  xml_ref_t xml_data = handleVis(vis.name(), vis);
750  vol.setRef(_U(visref), xml_data.name());
751  geo.xmlVolumes[v] = vol;
752  }
753  }
754  }
755  return vol;
756 }
757 
759 void LCDDConverter::collectVolume(const std::string& /* name */, const TGeoVolume* volume) const {
760  Volume v(volume);
761  if ( is_volume(volume) ) {
762  GeometryInfo& geo = data();
763  Region reg = v.region();
764  LimitSet lim = v.limitSet();
765  SensitiveDetector det = v.sensitiveDetector();
766  if (lim.isValid())
767  geo.limits.insert(lim);
768  if (reg.isValid())
769  geo.regions.insert(reg);
770  if (det.isValid())
771  geo.sensitives.insert(det);
772  }
773  else {
774  printout(WARNING,"LCDDConverter","++ CollectVolume: Skip volume: %s",volume->GetName());
775  }
776 }
777 
778 void LCDDConverter::checkVolumes(const std::string& /* name */, Volume v) const {
779  std::string n = v.name()+_toString(v.ptr(),"_%p");
780  NameSet::const_iterator i = m_checkNames.find(n);
781  if (i != m_checkNames.end()) {
782  std::stringstream str;
783  str << "++ CheckVolumes: Volume " << n << " ";
784  if (is_volume(v.ptr())) {
785  SensitiveDetector sd = v.sensitiveDetector();
786  VisAttr vis = v.visAttributes();
787  if (sd.isValid()) {
788  str << "of " << sd.name() << " ";
789  }
790  else if (vis.isValid()) {
791  str << "with VisAttrs " << vis.name() << " ";
792  }
793  }
794  str << "has duplicate entries." << std::endl;
795  printout(ERROR,"LCDDConverter",str.str().c_str());
796  return;
797  }
798  m_checkNames.insert(n);
799 }
800 
802 xml_h LCDDConverter::handlePlacement(const std::string& name,PlacedVolume node) const {
803  GeometryInfo& geo = data();
804  xml_h place = geo.xmlPlacements[node];
805  if (!place) {
806  TGeoMatrix* matrix = node->GetMatrix();
807  TGeoVolume* volume = node->GetVolume();
808  xml_ref_t vol = xml_h(geo.xmlVolumes[volume]);
809  xml_h mot = geo.xmlVolumes[node->GetMotherVolume()];
810 
811  place = xml_elt_t(geo.doc, _U(physvol));
812  if (mot) { // Beware of top level volume!
813  mot.append(place);
814  }
815  place.setRef(_U(volumeref), vol.name());
816  if (m) {
817  xml_ref_t pos = handlePosition(name+"_pos", matrix);
818  place.setRef(_U(positionref), pos.name());
819  if ( matrix->IsRotation() ) {
820  xml_ref_t rot = handleRotation(name+"_rot", matrix);
821  place.setRef(_U(rotationref), rot.name());
822  }
823  }
824  if (geo.doc_root.tag() != "gdml") {
825  if (is_placement(node)) {
826  const PlacedVolume::VolIDs& ids = node.volIDs();
827  for (PlacedVolume::VolIDs::const_iterator i = ids.begin(); i != ids.end(); ++i) {
828  xml_h pvid = xml_elt_t(geo.doc, _U(physvolid));
829  pvid.setAttr(_U(field_name), (*i).first);
830  pvid.setAttr(_U(value), (*i).second);
831  place.append(pvid);
832  }
833  }
834  }
835  geo.xmlPlacements[node] = place;
836  }
837  else {
838  std::cout << "Attempt to DOUBLE-place physical volume:" << name << " No:" << node->GetNumber() << std::endl;
839  }
840  return place;
841 }
842 
844 xml_h LCDDConverter::handleRegion(const std::string& /* name */, Region region) const {
845  GeometryInfo& geo = data();
846  xml_h reg = geo.xmlRegions[region];
847  if (!reg) {
848  geo.doc_regions.append(reg = xml_elt_t(geo.doc, _U(region)));
849  reg.setAttr(_U(name), region.name());
850  reg.setAttr(_U(cut), region.cut());
851  reg.setAttr(_U(eunit), "GeV"); // TGeo has energy in GeV
852  reg.setAttr(_U(lunit), "cm"); // TGeo has lengths in cm
853  reg.setAttr(_U(store_secondaries), region.storeSecondaries());
854  geo.xmlRegions[region] = reg;
855  }
856  return reg;
857 }
858 
860 xml_h LCDDConverter::handleLimitSet(const std::string& /* name */, LimitSet lim) const {
861  GeometryInfo& geo = data();
862  xml_h xml = geo.xmlLimits[lim];
863  if (!xml) {
864  geo.doc_limits.append(xml = xml_elt_t(geo.doc, _U(limitset)));
865  xml.setAttr(_U(name), lim.name());
866  const std::set<Limit>& obj = lim.limits();
867  for (std::set<Limit>::const_iterator i = obj.begin(); i != obj.end(); ++i) {
868  xml_h x = xml_elt_t(geo.doc, _U(limit));
869  const Limit& l = *i;
870  xml.append(x);
871  x.setAttr(_U(name), l.name);
872  x.setAttr(_U(unit), l.unit);
873  x.setAttr(_U(value), l.value);
874  x.setAttr(_U(particles), l.particles);
875  }
876  geo.xmlLimits[lim] = xml;
877  }
878  return xml;
879 }
880 
883  xml_h xml;
884  if (seg.isValid()) {
885  typedef DDSegmentation::Parameters _P;
886  std::string typ = seg.type();
887  _P p = seg.parameters();
888  xml = xml_elt_t(data().doc, Unicode(typ));
889  for (_P::const_iterator i = p.begin(); i != p.end(); ++i) {
890  const _P::value_type& v = *i;
891  if (v->name() == "lunit") {
892  std::string val = v->value() == _toDouble("mm") ? "mm" : v->value() == _toDouble("cm") ? "cm" :
893  v->value() == _toDouble("m") ? "m" : v->value() == _toDouble("micron") ? "micron" :
894  v->value() == _toDouble("nanometer") ? "namometer" : "??";
895  xml.setAttr(Unicode(v->name()), Unicode(val));
896  continue;
897  }
898  // translate from TGeo units to Geant4 units if necessary
900  double value = _toDouble(v->value()) / dd4hep::mm;
901  xml.setAttr(Unicode(v->name()), value);
902  } else if (v->unitType() == DDSegmentation::SegmentationParameter::AngleUnit) {
903  double value = _toDouble(v->value()) * DEGREE_2_RAD;
904  xml.setAttr(Unicode(v->name()), value);
905  } else {
906  xml.setAttr(Unicode(v->name()), v->value());
907  }
908  }
909  }
910  return xml;
911 }
912 
914 xml_h LCDDConverter::handleSensitive(const std::string& /* name */, SensitiveDetector sd) const {
915  GeometryInfo& geo = data();
916  xml_h sensdet = geo.xmlSensDets[sd];
917  if (!sensdet) {
918  geo.doc_detectors.append(sensdet = xml_elt_t(geo.doc, Unicode(sd.type())));
919  sensdet.setAttr(_U(name), sd.name());
920  sensdet.setAttr(_U(ecut), sd.energyCutoff());
921  sensdet.setAttr(_U(eunit), "MeV");
922  sensdet.setAttr(_U(verbose), int(sd.verbose() ? 1 : 0));
923  sensdet.setAttr(_U(hits_collection), sd.hitsCollection());
924  if (sd.combineHits())
925  sensdet.setAttr(_U(combine_hits), sd.combineHits());
926  Readout ro = sd.readout();
927  if (ro.isValid()) {
928  xml_ref_t ref = handleIdSpec(ro.idSpec().name(), ro.idSpec());
929  sensdet.setRef(_U(idspecref), ref.name());
931  if (seg)
932  sensdet.append(seg);
933  }
934  geo.xmlSensDets[sd] = sensdet;
935  }
936  return sensdet;
937 }
938 
940 xml_h LCDDConverter::handleIdSpec(const std::string& name, IDDescriptor id_spec) const {
941  GeometryInfo& geo = data();
942  xml_h id = geo.xmlIdSpecs[id_spec];
943  if (!id) {
944  int length = 0, start = 0;
945  IDDescriptor desc = id_spec;
946  geo.doc_idDict.append(id = xml_elt_t(geo.doc, _U(idspec)));
947  id.setAttr(_U(name), name);
948  const IDDescriptor::FieldMap& fm = desc.fields();
949  for (const auto& i : fm ) {
950  xml_h idfield = xml_elt_t(geo.doc, _U(idfield));
951 #if 0
952  const BitFieldElement* f = i.second;
953  start = f.first;
954  length = f.second<0 ? -f.second : f.second;
955  idfield.setAttr(_U(signed),f.second<0 ? true : false);
956  idfield.setAttr(_U(label),(*i).first);
957  idfield.setAttr(_U(length),length);
958  idfield.setAttr(_U(start),start);
959 #else
960  const BitFieldElement* f = i.second;
961  idfield.setAttr(_U(signed),f->isSigned() ? true : false);
962  idfield.setAttr(_U(label), f->name());
963  idfield.setAttr(_U(length), (int) f->width());
964  idfield.setAttr(_U(start), (int) f->offset());
965 #endif
966  id.append(idfield);
967  }
968  id.setAttr(_U(length), length + start);
969  geo.xmlIdSpecs[id_spec] = id;
970  }
971  return id;
972 }
973 
975 xml_h LCDDConverter::handleVis(const std::string& /* name */, VisAttr attr) const {
976  GeometryInfo& geo = data();
977  xml_h vis = geo.xmlVis[attr];
978  if (!vis) {
979  float red = 0, green = 0, blue = 0;
980  int style = attr.lineStyle();
981  int draw = attr.drawingStyle();
982 
983  geo.doc_display.append(vis = xml_elt_t(geo.doc, _U(vis)));
984  vis.setAttr(_U(name), attr.name());
985  vis.setAttr(_U(visible), attr.visible());
986  vis.setAttr(_U(show_daughters), attr.showDaughters());
987  if (style == VisAttr::SOLID)
988  vis.setAttr(_U(line_style), "unbroken");
989  else if (style == VisAttr::DASHED)
990  vis.setAttr(_U(line_style), "broken");
991  if (draw == VisAttr::SOLID)
992  vis.setAttr(_U(drawing_style), "solid");
993  else if (draw == VisAttr::WIREFRAME)
994  vis.setAttr(_U(drawing_style), "wireframe");
995 
996  xml_h col = xml_elt_t(geo.doc, _U(color));
997  attr.rgb(red, green, blue);
998  col.setAttr(_U(alpha), attr.alpha());
999  col.setAttr(_U(R), red);
1000  col.setAttr(_U(B), blue);
1001  col.setAttr(_U(G), green);
1002  vis.append(col);
1003  geo.xmlVis[attr] = vis;
1004  }
1005  return vis;
1006 }
1007 
1009 xml_h LCDDConverter::handleField(const std::string& /* name */, OverlayedField f) const {
1010  GeometryInfo& geo = data();
1011  xml_h field = geo.xmlFields[f];
1012  if (!field) {
1013  Handle<NamedObject> fld(f);
1014  std::string type = f->GetTitle();
1015  field = xml_elt_t(geo.doc, Unicode(type));
1016  field.setAttr(_U(name), f->GetName());
1017  fld = PluginService::Create<NamedObject*>(type + "_Convert2Detector", &m_detDesc, &field, &fld);
1018  printout(ALWAYS,"LCDDConverter","++ %s electromagnetic field:%s of type %s",
1019  (fld.isValid() ? "Converted" : "FAILED to convert "), f->GetName(), type.c_str());
1020  if (!fld.isValid()) {
1021  PluginDebug dbg;
1022  PluginService::Create<NamedObject*>(type + "_Convert2Detector", &m_detDesc, &field, &fld);
1023  except("LCDDConverter", "Failed to locate plugin to convert electromagnetic field:"
1024  + std::string(f->GetName()) + " of type " + type + ". "
1025  + dbg.missingFactory(type));
1026  }
1027  geo.doc_fields.append(field);
1028  }
1029  return field;
1030 }
1031 
1034  std::map<std::string, std::string> processors;
1035  static int s_idd = 9999999;
1036  std::string id;
1037  for (Detector::Properties::const_iterator i = prp.begin(); i != prp.end(); ++i) {
1038  const std::string& nam = (*i).first;
1039  const Detector::PropertyValues& vals = (*i).second;
1040  if (nam.substr(0, 6) == "geant4") {
1041  Detector::PropertyValues::const_iterator id_it = vals.find("id");
1042  if (id_it != vals.end()) {
1043  id = (*id_it).second;
1044  }
1045  else {
1046  char text[32];
1047  ::snprintf(text, sizeof(text), "%d", ++s_idd);
1048  id = text;
1049  }
1050  processors.insert(make_pair(id, nam));
1051  }
1052  }
1053  for (std::map<std::string, std::string>::const_iterator i = processors.begin(); i != processors.end(); ++i) {
1054  const GeoHandler* ptr = this;
1055  std::string nam = (*i).second;
1056  const Detector::PropertyValues& vals = prp[nam];
1057  auto iter = vals.find("type");
1058  if ( iter != vals.end() ) {
1059  std::string type = iter->second;
1060  std::string tag = type + "_Geant4_action";
1061  long result = PluginService::Create<long>(tag, &m_detDesc, ptr, &vals);
1062  if (0 == result) {
1063  PluginDebug dbg;
1064  result = PluginService::Create<long>(tag, &m_detDesc, ptr, &vals);
1065  if (0 == result) {
1066  except("LCDDConverter", "Failed to locate plugin to interprete files of type"
1067  " \"" + tag + "\" - no factory:" + type + ". " +
1068  dbg.missingFactory(tag));
1069  }
1070  }
1071  result = *(long*) result;
1072  if (result != 1) {
1073  except("LCDDConverter", "Failed to invoke the plugin " + tag + " of type " + type);
1074  }
1075  printout(INFO,"","+++ Executed Successfully Detector setup module %s.", type.c_str());
1076  continue;
1077  }
1078  printout(INFO,"","+++ FAILED to execute Detector setup module %s.", nam.c_str());
1079  }
1080 }
1081 
1084  GeometryInfo& geo = data();
1085  Header hdr = m_detDesc.header();
1086  if ( hdr.isValid() ) {
1087  xml_h obj;
1088  geo.doc_header.append(obj = xml_elt_t(geo.doc, _U(detector)));
1089  obj.setAttr(_U(name), hdr.name());
1090  geo.doc_header.append(obj = xml_elt_t(geo.doc, _U(generator)));
1091  obj.setAttr(_U(name), "LCDDConverter");
1092  obj.setAttr(_U(version), hdr.version());
1093  obj.setAttr(_U(file), hdr.url());
1094  obj.setAttr(_U(checksum), Unicode(m_detDesc.constantAsString("compact_checksum")));
1095  geo.doc_header.append(obj = xml_elt_t(geo.doc, _U(author)));
1096  obj.setAttr(_U(name), hdr.author());
1097  geo.doc_header.append(obj = xml_elt_t(geo.doc, _U(comment)));
1098  obj.setText(hdr.comment());
1099  return;
1100  }
1101  printout(WARNING,"LCDDConverter","+++ No Detector header information availible from the geometry description.");
1102 }
1103 
1104 template <typename O, typename C, typename F> void handle(const O* o, const C& c, F pmf) {
1105  for (typename C::const_iterator i = c.begin(); i != c.end(); ++i) {
1106  std::string n = (*i)->GetName();
1107  (o->*pmf)(n, *i);
1108  }
1109 }
1110 
1111 template <typename O, typename C, typename F> void handleMap(const O* o, const C& c, F pmf) {
1112  for (typename C::const_iterator i = c.begin(); i != c.end(); ++i)
1113  (o->*pmf)((*i).first, (*i).second);
1114 }
1115 
1116 template <typename O, typename C, typename F> void handleRMap(const O* o, const C& c, F pmf) {
1117  for (typename C::const_reverse_iterator i = c.rbegin(); i != c.rend(); ++i)
1118  handle(o, (*i).second, pmf);
1119 }
1120 
1123  Detector& description = m_detDesc;
1124  if (!top.isValid()) {
1125  throw std::runtime_error("Attempt to call createGDML with an invalid geometry!");
1126  }
1127  GeometryInfo& geo = *(m_dataPtr = new GeometryInfo);
1128  m_data->clear();
1129  collect(top, geo);
1130 
1131  printout(ALWAYS,"LCDDConverter","++ ==> Converting in memory detector description to GDML format...");
1132  xml::DocumentHandler docH;
1133  geo.doc = docH.create("gdml", docH.defaultComment());
1134  geo.doc_root = geo.doc.root();
1135  geo.doc_root.setAttr(Unicode("xmlns:xs"), "http://www.w3.org/2001/XMLSchema-instance");
1136  geo.doc_root.setAttr(Unicode("xs:noNamespaceSchemaLocation"),
1137  "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd");
1138  // geo.doc = docH.create("gdml_simple_extension",comment);
1139  // geo.doc_root.setAttr(Unicode("xmlns:gdml_simple_extension"),"http://www.example.org");
1140  // geo.doc_root.setAttr(Unicode("xs:noNamespaceSchemaLocation"),
1141  // "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd");
1142  Volume world_vol = description.worldVolume();
1143  geo.doc_root.append(geo.doc_define = xml_elt_t(geo.doc, _U(define)));
1144  geo.doc_root.append(geo.doc_materials = xml_elt_t(geo.doc, _U(materials)));
1145  geo.doc_root.append(geo.doc_solids = xml_elt_t(geo.doc, _U(solids)));
1146  geo.doc_root.append(geo.doc_structure = xml_elt_t(geo.doc, _U(structure)));
1147  geo.doc_root.append(geo.doc_setup = xml_elt_t(geo.doc, _U(setup)));
1148  geo.doc_setup.setRef(_U(world), genName(world_vol.name(),world_vol.ptr()));
1149  geo.doc_setup.setAttr(_U(name), Unicode("default"));
1150  geo.doc_setup.setAttr(_U(version), Unicode("1.0"));
1151 
1152  // Ensure that all required materials are present in the Detector material table
1153 #if 0
1154  const Detector::HandleMap& mat = description.materials();
1155  for(Detector::HandleMap::const_iterator i=mat.begin(); i!=mat.end(); ++i)
1156  geo.materials.insert(dynamic_cast<TGeoMedium*>((*i).second.ptr()));
1157 #endif
1158 
1159  // Start creating the objects for materials, solids and log volumes.
1161  printout(ALWAYS,"LCDDConverter","++ Handled %ld materials.",geo.materials.size());
1162 
1164  printout(ALWAYS,"LCDDConverter","++ Collected %ld volumes.",geo.volumes.size());
1165 
1167  printout(ALWAYS,"LCDDConverter","++ Handled %ld solids.",geo.solids.size());
1168 
1170  printout(ALWAYS,"LCDDConverter","++ Handled %ld volumes.",geo.volumes.size());
1171 
1172  m_checkNames.clear();
1174  return geo.doc;
1175 }
1176 
1179  if (!top.isValid()) {
1180  throw std::runtime_error("Attempt to call createDetector with an invalid geometry!");
1181  }
1182 
1183  GeometryInfo& geo = *(m_dataPtr = new GeometryInfo);
1184  m_data->clear();
1185  collect(top, geo);
1186  printout(ALWAYS,"LCDDConverter","++ ==> Dump visualisation attributes "
1187  "from in memory detector description...");
1188  xml::DocumentHandler docH;
1189  xml_elt_t elt(0);
1190  geo.doc = docH.create("visualization", docH.defaultComment());
1191  geo.doc_root = geo.doc.root();
1192  geo.doc_root.append(geo.doc_display = xml_elt_t(geo.doc, _U(display)));
1193  geo.doc_root.append(geo.doc_structure = xml_elt_t(geo.doc, _U(structure)));
1194 
1197  printout(ALWAYS,"LCDDConverter","++ Handled %ld volumes.",geo.volumes.size());
1198  return geo.doc;
1199 }
1200 
1203  Detector& description = m_detDesc;
1204  if (!top.isValid()) {
1205  throw std::runtime_error("Attempt to call createDetector with an invalid geometry!");
1206  }
1207 
1208  GeometryInfo& geo = *(m_dataPtr = new GeometryInfo);
1209  m_data->clear();
1210  collect(top, geo);
1211  xml::DocumentHandler docH;
1212  xml_elt_t elt(0);
1213  Volume world_vol = description.worldVolume();
1214  geo.doc = docH.create("description", docH.defaultComment());
1215  geo.doc_root = geo.doc.root();
1216  geo.doc_root.setAttr(Unicode("xmlns:description"), "http://www.lcsim.org/schemas/description/1.0");
1217  geo.doc_root.setAttr(Unicode("xmlns:xs"), "http://www.w3.org/2001/XMLSchema-instance");
1218  geo.doc_root.setAttr(Unicode("xs:noNamespaceSchemaLocation"),
1219  "http://www.lcsim.org/schemas/description/1.0/description.xsd");
1220 
1221  geo.doc_root.append(geo.doc_header = xml_elt_t(geo.doc, _U(header)));
1222  geo.doc_root.append(geo.doc_idDict = xml_elt_t(geo.doc, _U(iddict)));
1223  geo.doc_root.append(geo.doc_detectors = xml_elt_t(geo.doc, _U(sensitive_detectors)));
1224  geo.doc_root.append(geo.doc_limits = xml_elt_t(geo.doc, _U(limits)));
1225  geo.doc_root.append(geo.doc_regions = xml_elt_t(geo.doc, _U(regions)));
1226  geo.doc_root.append(geo.doc_display = xml_elt_t(geo.doc, _U(display)));
1227  geo.doc_root.append(geo.doc_gdml = xml_elt_t(geo.doc, _U(gdml)));
1228  geo.doc_root.append(geo.doc_fields = xml_elt_t(geo.doc, _U(fields)));
1229 
1230  geo.doc_gdml.append(geo.doc_define = xml_elt_t(geo.doc, _U(define)));
1231  geo.doc_gdml.append(geo.doc_materials = xml_elt_t(geo.doc, _U(materials)));
1232  geo.doc_gdml.append(geo.doc_solids = xml_elt_t(geo.doc, _U(solids)));
1233  geo.doc_gdml.append(geo.doc_structure = xml_elt_t(geo.doc, _U(structure)));
1234  geo.doc_gdml.append(geo.doc_setup = xml_elt_t(geo.doc, _U(setup)));
1235  geo.doc_setup.setRef(_U(world), genName(world_vol.name(),world_vol.ptr()));
1236  geo.doc_setup.setAttr(_U(name), Unicode("default"));
1237  geo.doc_setup.setAttr(_U(version), Unicode("1.0"));
1238 
1239  // Ensure that all required materials are present in the Detector material table
1240  const Detector::HandleMap& fld = description.fields();
1241  for (Detector::HandleMap::const_iterator i = fld.begin(); i != fld.end(); ++i)
1242  geo.fields.insert((*i).second);
1243 
1244  printout(ALWAYS,"LCDDConverter","++ ==> Converting in memory detector description to Detector format...");
1245  handleHeader();
1246  // Start creating the objects for materials, solids and log volumes.
1248  printout(ALWAYS,"LCDDConverter","++ Handled %ld materials.",geo.materials.size());
1249 
1251  printout(ALWAYS,"LCDDConverter","++ Collected %ld volumes.",geo.volumes.size());
1252 
1254  printout(ALWAYS,"LCDDConverter","++ Handled %ld solids.",geo.solids.size());
1255 
1256  handle(this, geo.vis, &LCDDConverter::handleVis);
1257  printout(ALWAYS,"LCDDConverter","++ Handled %ld visualization attributes.",geo.vis.size());
1258 
1260  printout(ALWAYS,"LCDDConverter","++ Handled %ld sensitive detectors.",geo.sensitives.size());
1261 
1263  printout(ALWAYS,"LCDDConverter","++ Handled %ld limit sets.",geo.limits.size());
1264 
1266  printout(ALWAYS,"LCDDConverter","++ Handled %ld regions.",geo.regions.size());
1267 
1269  printout(ALWAYS,"LCDDConverter","++ Handled %ld volumes.",geo.volumes.size());
1270 
1272  printout(ALWAYS,"LCDDConverter","++ Handled %ld fields.",geo.fields.size());
1273 
1274  m_checkNames.clear();
1276 #if 0
1277  //==================== Fields
1279 #endif
1280  return geo.doc;
1281 }
1282 
1285  : doc(0), doc_root(0), doc_header(0), doc_idDict(0), doc_detectors(0), doc_limits(0),
1286  doc_regions(0), doc_display(0), doc_gdml(0), doc_fields(0), doc_define(0),
1287  doc_materials(0), doc_solids(0), doc_structure(0), doc_setup(0)
1288 {
1289 }
1290 
1291 static long dump_output(xml_doc_t doc, int argc, char** argv) {
1292  xml::DocumentHandler docH;
1293  return docH.output(doc, argc > 0 ? argv[0] : "");
1294 }
1295 
1296 long create_gdml_from_dd4hep(Detector& description, int argc, char** argv) {
1297  LCDDConverter wr(description);
1298  return dump_output(wr.createGDML(description.world()), argc, argv);
1299 }
1300 
1301 static long create_description(Detector& description, int argc, char** argv) {
1302  LCDDConverter wr(description);
1303  return dump_output(wr.createDetector(description.world()), argc, argv);
1304 }
1305 
1306 static long create_vis(Detector& description, int argc, char** argv) {
1307  LCDDConverter wr(description);
1308  return dump_output(wr.createVis(description.world()), argc, argv);
1309 }
1310 
1311 static long create_visASCII(Detector& description, int /* argc */, char** argv) {
1312  LCDDConverter wr(description);
1313  /* xml_doc_t doc = */ wr.createVis(description.world());
1314  LCDDConverter::GeometryInfo& geo = wr.data();
1315  std::map<std::string, xml_comp_t> vis_map;
1316  for (xml_coll_t c(geo.doc_display, _U(vis)); c; ++c)
1317  vis_map.insert(make_pair(xml_comp_t(c).nameStr(), xml_comp_t(c)));
1318 
1319  const char* sep = ";";
1320  std::ofstream os(argv[0]);
1321  for (xml_coll_t c(geo.doc_structure, _U(volume)); c; ++c) {
1322  xml_comp_t vol = c;
1323  xml_comp_t ref = c.child(_U(visref));
1324  auto iter = vis_map.find(ref.refStr());
1325  if ( iter != vis_map.end() ) {
1326  xml_comp_t vis = iter->second;
1327  xml_comp_t col = vis.child(_U(color));
1328  os << "vol:" << vol.nameStr() << sep << "vis:" << vis.nameStr() << sep
1329  << "visible:" << vis.visible() << sep << "r:"
1330  << col.R() << sep << "g:" << col.G() << sep << "b:" << col.B() << sep
1331  << "alpha:" << col.alpha() << sep << "line_style:"
1332  << vis.attr < std::string > (_U(line_style)) << sep
1333  << "drawing_style:" << vis.attr < std::string> (_U(drawing_style)) << sep
1334  << "show_daughters:" << vis.show_daughters() << sep << std::endl;
1335  }
1336  }
1337  os.close();
1338  return 1;
1339 }
1340 
1341 DECLARE_APPLY(DD4hepGeometry2VIS, create_vis)
1342 DECLARE_APPLY(DD4hepGeometry2VISASCII, create_visASCII)
1343 //DECLARE_APPLY(DD4hepGeometry2GDML, create_gdml_from_dd4hep)
1344 DECLARE_APPLY(DD4hepGeometry2Detector, create_description)
dd4hep::xml::DocumentHandler::create
Document create(const char *tag, const char *comment=0) const
Create new XML document by parsing empty xml buffer.
Definition: DocumentHandler.cpp:680
dd4hep::xml::Collection_t
Class to support the access to collections of XmlNodes (or XmlElements)
Definition: XMLElements.h:636
dd4hep::DDSegmentation::BitFieldElement::offset
unsigned offset() const
Definition: BitFieldCoder.h:63
dd4hep::Header::comment
const std::string & comment() const
Accessor to object comment.
Definition: Objects.cpp:127
dd4hep::Detector::world
virtual DetElement world() const =0
Return reference to the top-most (world) detector element.
dd4hep::detail::GeoHandlerTypes::GeometryInfo::solids
std::vector< TGeoShape * > solids
Definition: GeoHandler.h:62
dd4hep::Atom
Handle class describing an element in the periodic table.
Definition: Objects.h:241
Segmentations.h
dd4hep::detail::LCDDConverter::GeometryInfo::doc_regions
xml_elt_t doc_regions
Definition: LCDDConverter.h:113
dd4hep::DDSegmentation::BitFieldElement
Helper class for BitFieldCoder that corresponds to one field value.
Definition: BitFieldCoder.h:31
Volumes.h
dd4hep::IDDescriptor::FieldMap
std::vector< std::pair< std::string, const Field * > > FieldMap
Definition: IDDescriptor.h:39
dd4hep::xml::Element::append
void append(Handle_t handle) const
Append a new element to the existing tree.
Definition: XMLElements.h:839
v
View * v
Definition: MultiView.cpp:28
dd4hep::detail::LCDDConverter::m_detDesc
Detector & m_detDesc
Reference to detector description.
Definition: LCDDConverter.h:120
dd4hep::SensitiveDetector
Handle class to hold the information of a sensitive detector.
Definition: DetElement.h:43
dd4hep::DDSegmentation::BitFieldElement::isSigned
bool isSigned() const
Definition: BitFieldCoder.h:69
dd4hep::Limit::particles
std::string particles
Particle the limit should be applied to.
Definition: Objects.h:391
dd4hep::detail::LCDDConverter::m_checkNames
NameSet m_checkNames
Definition: LCDDConverter.h:121
dd4hep::Detector::Properties
std::map< std::string, PropertyValues > Properties
Definition: Detector.h:95
dd4hep::DDSegmentation::SegmentationParameter::LengthUnit
@ LengthUnit
Definition: SegmentationParameter.h:111
dd4hep::PlacedVolume
Handle class holding a placed volume (also called physical volume)
Definition: Volumes.h:164
dd4hep::VisAttr
Handle class describing visualization attributes.
Definition: Objects.h:323
dd4hep::Region::storeSecondaries
bool storeSecondaries() const
Access secondaries flag.
Definition: Objects.cpp:524
dd4hep::DDSegmentation::SegmentationParameter::AngleUnit
@ AngleUnit
Definition: SegmentationParameter.h:111
dd4hep::xml::Element::tag
std::string tag() const
Access the tag name of this DOM element.
Definition: XMLElements.h:823
dd4hep::detail::LCDDConverter::GeometryInfo::doc_gdml
xml_elt_t doc_gdml
Definition: LCDDConverter.h:113
dd4hep::xml::Handle_t::setRef
Handle_t setRef(const XmlChar *tag, const XmlChar *ref)
Add reference child as a new child node. The obj must have the "name" attribute!
Definition: XMLElements.cpp:940
dd4hep::xml::Handle_t::append
void append(Handle_t e) const
Append a DOM element to the current node.
Definition: XMLElements.cpp:730
dd4hep::LimitSet::limits
const std::set< Limit > & limits() const
Accessor to limits container.
Definition: Objects.cpp:463
dd4hep::detail::LCDDConverter::GeometryInfo::xmlSensDets
SensDetMap xmlSensDets
Definition: LCDDConverter.h:81
dd4hep::Limit
Small object describing a limit structure acting on a particle type.
Definition: Objects.h:388
dd4hep::IDDescriptor
Class implementing the ID encoding of the detector response.
Definition: IDDescriptor.h:36
DECLARE_APPLY
#define DECLARE_APPLY(name, func)
Definition: Factories.h:281
dd4hep::detail::LCDDConverter::createDetector
xml_doc_t createDetector(DetElement top)
Create geometry conversion in Detector format.
Definition: LCDDConverter.cpp:1202
dd4hep::detail::LCDDConverter::handleField
virtual xml_h handleField(const std::string &name, OverlayedField field) const
Convert the electric or magnetic fields into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:1009
dd4hep::detail::LCDDConverter::GeometryInfo::doc_solids
xml_elt_t doc_solids
Definition: LCDDConverter.h:114
dd4hep::Detector::PropertyValues
std::map< std::string, std::string > PropertyValues
Definition: Detector.h:94
dd4hep::Handle::isValid
bool isValid() const
Check the validity of the object held by the handle.
Definition: Handle.h:126
dd4hep::Handle< NamedObject >
DetectorInterna.h
dd4hep::_toString
std::string _toString(bool value)
String conversions: boolean value to string.
Definition: Handle.cpp:332
dd4hep::detail::LCDDConverter::GeometryInfo::limits
std::set< LimitSet > limits
Definition: LCDDConverter.h:87
dd4hep::detail::LCDDConverter::handleRotation
virtual xml_h handleRotation(const std::string &name, const TGeoMatrix *trafo) const
Convert the Rotation into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:637
dd4hep::Solid_type< TGeoShape >
dd4hep::Detector::fields
virtual const HandleMap & fields() const =0
Accessor to the map of field entries, which together form the global field.
dd4hep::xml::Handle_t
Class to easily access the properties of single XmlElements.
Definition: XMLElements.h:380
dd4hep::Header::url
const std::string & url() const
Accessor to object url.
Definition: Objects.cpp:87
dd4hep::detail::LCDDConverter::GeometryInfo::xmlFields
FieldMap xmlFields
Definition: LCDDConverter.h:84
dd4hep::Detector::worldVolume
virtual Volume worldVolume() const =0
Return handle to the world volume containing everything.
dd4hep::detail::LCDDConverter::handleElement
virtual xml_h handleElement(const std::string &name, Atom element) const
Convert the geometry type element into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:136
dd4hep::detail::GeoHandlerTypes::GeometryInfo::vis
std::set< VisAttr > vis
Definition: GeoHandler.h:66
dd4hep::Handle::name
const char * name() const
Access the object name (or "" if not supported by the object)
dd4hep::detail::LCDDConverter::data
GeometryInfo & data() const
Definition: LCDDConverter.h:124
dd4hep::detail::LCDDConverter
Geometry converter from dd4hep to Geant 4 in Detector format.
Definition: LCDDConverter.h:49
dd4hep::detail::LCDDConverter::GeometryInfo::doc_display
xml_elt_t doc_display
Definition: LCDDConverter.h:113
dd4hep::detail::LCDDConverter::handleSolid
virtual xml_h handleSolid(const std::string &name, const TGeoShape *volume) const
Convert the geometry type solid into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:220
dd4hep::VisAttr::alpha
float alpha() const
Get alpha value.
Definition: Objects.cpp:370
dd4hep::detail::LCDDConverter::GeometryInfo::xmlRegions
RegionMap xmlRegions
Definition: LCDDConverter.h:77
dd4hep::OverlayedField
Class describing a field overlay with several sources.
Definition: Fields.h:138
dd4hep::detail::LCDDConverter::handleMaterial
virtual xml_h handleMaterial(const std::string &name, Material medium) const
Convert the geometry type material into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:159
xml_comp_t
dd4hep::xml::Component xml_comp_t
Definition: XML.h:33
dd4hep::detail::LCDDConverter::GeometryInfo::doc_fields
xml_elt_t doc_fields
Definition: LCDDConverter.h:113
dd4hep::detail::LCDDConverter::handlePlacement
virtual xml_h handlePlacement(const std::string &name, PlacedVolume node) const
Convert the geometry type volume placement into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:802
dd4hep::detail::LCDDConverter::createVis
xml_doc_t createVis(DetElement top)
Create geometry conversion in Vis format.
Definition: LCDDConverter.cpp:1178
dd4hep::Header::name
const std::string name() const
Accessor to object name.
Definition: Objects.cpp:67
DocumentHandler.h
dd4hep::detail::LCDDConverter::GeometryInfo::doc_idDict
xml_elt_t doc_idDict
Definition: LCDDConverter.h:113
dd4hep::VisAttr::lineStyle
int lineStyle() const
Get line style.
Definition: Objects.cpp:350
dd4hep::detail::LCDDConverter::GeometryInfo::checkPosition
void checkPosition(const std::string &name, const TNamed *n) const
Definition: LCDDConverter.h:95
dd4hep::Material
Handle class describing a material.
Definition: Objects.h:271
dd4hep::VisAttr::drawingStyle
int drawingStyle() const
Get drawing style.
Definition: Objects.cpp:360
dd4hep::VisAttr::showDaughters
bool showDaughters() const
Get Flag to show/hide daughter elements.
Definition: Objects.cpp:330
dd4hep::detail::LCDDConverter::GeometryInfo::xmlSolids
SolidMap xmlSolids
Definition: LCDDConverter.h:74
dd4hep::detail::LCDDConverter::collectVolume
virtual void collectVolume(const std::string &name, const TGeoVolume *volume) const
Dump logical volume in GDML format to output stream.
Definition: LCDDConverter.cpp:759
dd4hep::detail::LCDDConverter::GeometryInfo::identity_pos
xml_h identity_pos
Definition: LCDDConverter.h:112
dd4hep::detail::LCDDConverter::GeometryInfo::xmlLimits
LimitMap xmlLimits
Definition: LCDDConverter.h:79
dd4hep::View::name
const std::string & name() const
Access to the view name/title.
Definition: View.h:81
dd4hep::detail::LCDDConverter::GeometryInfo::checkMaterial
void checkMaterial(const std::string &name, Material n) const
Definition: LCDDConverter.h:107
dd4hep::DetElement
Handle class describing a detector element.
Definition: DetElement.h:187
dd4hep::detail::LCDDConverter::GeometryInfo::doc_define
xml_elt_t doc_define
Definition: LCDDConverter.h:114
dd4hep::PlacedVolumeExtension::VolIDs
Volume ID container.
Definition: Volumes.h:89
dd4hep::detail::LCDDConverter::createGDML
xml_doc_t createGDML(DetElement top)
Create geometry conversion in GDML format.
Definition: LCDDConverter.cpp:1122
dd4hep::Volume
Handle class holding a placed volume (also called physical volume)
Definition: Volumes.h:371
dd4hep::Solid_type::name
const char * name() const
Access to shape name.
Definition: Shapes.cpp:54
dd4hep::detail::GeoHandler
The base class for all dd4hep geometry crawlers.
Definition: GeoHandler.h:87
dd4hep::detail::LCDDConverter::GeometryInfo::xmlVis
VisMap xmlVis
Definition: LCDDConverter.h:78
dd4hep::detail
DD4hep internal namespace.
Definition: Alignments.h:32
dd4hep::DDSegmentation::BitFieldElement::name
const std::string & name() const
Definition: BitFieldCoder.h:60
mix
#define mix(a, b, c)
Definition: Geant4EventSeed.h:113
dd4hep::LimitSet
Handle class describing a set of limits as they are used for simulation.
Definition: Objects.h:424
dd4hep::SensitiveDetector::hitsCollection
const std::string & hitsCollection() const
Access the hits collection name.
Definition: DetElement.cpp:447
dd4hep::Volume::limitSet
LimitSet limitSet() const
Access to the limit set.
Definition: Volumes.cpp:1304
dd4hep::Limit::name
std::string name
Limit name.
Definition: Objects.h:393
dd4hep::detail::LCDDConverter::GeometryInfo::doc_header
xml_elt_t doc_header
Definition: LCDDConverter.h:113
dd4hep::VisAttr::DASHED
@ DASHED
Definition: Objects.h:326
dd4hep::detail::LCDDConverter::GeometryInfo
Data structure of the geometry converter from dd4hep to Geant 4 in Detector format.
Definition: LCDDConverter.h:70
dd4hep::DDSegmentation::Parameters
std::vector< Parameter > Parameters
Definition: SegmentationsInterna.h:35
dd4hep::detail::LCDDConverter::GeometryInfo::checkVolume
void checkVolume(const std::string &name, const TNamed *n) const
Definition: LCDDConverter.h:101
dd4hep::detail::LCDDConverter::checkVolumes
void checkVolumes(const std::string &name, Volume volume) const
Data integrity checker.
Definition: LCDDConverter.cpp:778
dd4hep::Volume::region
Region region() const
Access to the handle to the region structure.
Definition: Volumes.cpp:1285
Unicode
dd4hep::xml::Strng_t Unicode
Definition: XML.h:24
dd4hep::xml
Namespace for the AIDA detector description toolkit supporting XML utilities.
Definition: ConditionsTags.h:27
dd4hep::Detector::constantAsString
virtual std::string constantAsString(const std::string &name) const =0
Typed access to constants: access string values.
dd4hep::detail::LCDDConverter::GeometryInfo::doc_detectors
xml_elt_t doc_detectors
Definition: LCDDConverter.h:113
dd4hep::detail::LCDDConverter::GeometryInfo::regions
std::set< Region > regions
Definition: LCDDConverter.h:86
dd4hep::xml::Element::setRef
Attribute setRef(const XmlChar *tag, const XmlChar *refname) const
Set the reference attribute to the node (adds attribute ref="ref-name")
Definition: XMLElements.cpp:1109
dd4hep::detail::LCDDConverter::handleVolume
virtual xml_h handleVolume(const std::string &name, Volume volume) const
Convert the geometry type logical volume into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:671
dd4hep::detail::LCDDConverter::GeometryInfo::identity_rot
xml_h identity_rot
Definition: LCDDConverter.h:112
dd4hep::detail::LCDDConverter::m_dataPtr
GeometryInfo * m_dataPtr
Definition: LCDDConverter.h:122
dd4hep::DDSegmentation::BitFieldElement::width
unsigned width() const
Definition: BitFieldCoder.h:66
_U
#define _U(a)
Definition: Tags.h:23
dd4hep::SensitiveDetector::combineHits
bool combineHits() const
Access flag to combine hist.
Definition: DetElement.cpp:471
TNamed
Class of the ROOT toolkit. See http://root.cern.ch/root/htmldoc/ClassIndex.html.
Definition: ROOTClasses.h:37
dd4hep::Header::author
const std::string & author() const
Accessor to object author.
Definition: Objects.cpp:97
dd4hep::xml::RefElement::name
const XmlChar * name() const
Access the object's name in unicode.
Definition: XMLElements.cpp:1180
Plugins.h
dd4hep::Region
Handle class describing a region as used in simulation.
Definition: Objects.h:461
dd4hep::PluginDebug
Helper to debug plugin manager calls.
Definition: Plugins.h:73
dd4hep::VisAttr::rgb
bool rgb(float &red, float &green, float &blue) const
Get RGB values of the color (if valid)
Definition: Objects.cpp:395
dd4hep::detail::GeoHandler::m_data
std::map< int, std::vector< const TGeoNode * > > * m_data
actual container with std::vector (preserves order)
Definition: GeoHandler.h:93
dd4hep::Header
Handle class describing the basic information about geometry objects as it is defined in Detector.
Definition: Objects.h:157
dd4hep::Readout::segmentation
Segmentation segmentation() const
Access segmentation structure.
Definition: Readout.cpp:134
dd4hep::detail::LCDDConverter::GeometryInfo::xmlVolumes
VolumeMap xmlVolumes
Definition: LCDDConverter.h:75
dd4hep::Limit::unit
std::string unit
Units.
Definition: Objects.h:395
dd4hep::detail::LCDDConverter::~LCDDConverter
virtual ~LCDDConverter()
Standard destructor.
Definition: LCDDConverter.cpp:129
dd4hep::xml::Strng_t
Helper class to encapsulate a unicode string.
Definition: XMLElements.h:170
dd4hep::PluginDebug::missingFactory
std::string missingFactory(const std::string &name) const
Helper to check factory existence.
Definition: Plugins.cpp:147
dd4hep::xml::DocumentHandler::defaultComment
static std::string defaultComment()
Default comment string.
Definition: DocumentHandler.cpp:629
dd4hep::detail::GeoHandlerTypes::GeometryInfo::fields
std::set< Ref_t > fields
Definition: GeoHandler.h:67
dd4hep::verbose
std::size_t verbose(const std::string &src, const std::string &msg)
Definition: RootDictionary.h:63
dd4hep::Limit::value
double value
Double value.
Definition: Objects.h:399
dd4hep::detail::LCDDConverter::GeometryInfo::doc
xml_doc_t doc
Definition: LCDDConverter.h:111
dd4hep::detail::LCDDConverter::handleProperties
void handleProperties(Detector::Properties &prp) const
Handle the geant 4 specific properties.
Definition: LCDDConverter.cpp:1033
dd4hep::xml::Handle_t::setText
void setText(const XmlChar *text) const
Set the element's text.
Definition: XMLElements.cpp:786
dd4hep::detail::LCDDConverter::GeometryInfo::xmlElements
ElementMap xmlElements
Definition: LCDDConverter.h:72
dd4hep::Region::cut
double cut() const
Access cut value.
Definition: Objects.cpp:514
dd4hep::detail::LCDDConverter::GeometryInfo::doc_limits
xml_elt_t doc_limits
Definition: LCDDConverter.h:113
dd4hep::SensitiveDetector::type
std::string type() const
Access the type of the sensitive detector.
Definition: DetElement.cpp:409
dd4hep::IDDescriptor::fields
const FieldMap & fields() const
Access the fieldmap container.
Definition: IDDescriptor.cpp:87
dd4hep::xml::Document
Class supporting the basic functionality of an XML document.
Definition: XMLElements.h:697
dd4hep::Detector::HandleMap
std::map< std::string, Handle< NamedObject > > HandleMap
Type definition of a map of named handles.
Definition: Detector.h:93
dd4hep::detail::LCDDConverter::handleRegion
virtual xml_h handleRegion(const std::string &name, Region region) const
Convert the geometry type region into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:844
dd4hep::SensitiveDetector::verbose
bool verbose() const
Access flag to combine hist.
Definition: DetElement.cpp:459
xml_h
dd4hep::xml::Handle_t xml_h
Definition: ConditionsRepository.cpp:32
dd4hep::detail::LCDDConverter::GeometryInfo::sensitives
std::set< SensitiveDetector > sensitives
Definition: LCDDConverter.h:85
dd4hep::detail::LCDDConverter::GeometryInfo::xmlPositions
TrafoMap xmlPositions
Definition: LCDDConverter.h:82
dd4hep::Position
ROOT::Math::XYZVector Position
Definition: Objects.h:80
dd4hep::xml::RefElement
User abstraction class to manipulate named XML elements (references) within a document.
Definition: XMLElements.h:953
dd4hep::Segmentation::type
std::string type() const
Accessor: Segmentation type.
Definition: Segmentations.cpp:48
dd4hep::detail::LCDDConverter::GeometryInfo::checkShape
void checkShape(const std::string &name, const TNamed *n) const
Definition: LCDDConverter.h:104
dd4hep::detail::LCDDConverter::handleIdSpec
virtual xml_h handleIdSpec(const std::string &name, IDDescriptor idspec) const
Convert the geometry id dictionary entry to the corresponding Xml object(s).
Definition: LCDDConverter.cpp:940
dd4hep::detail::LCDDConverter::GeometryInfo::doc_setup
xml_elt_t doc_setup
Definition: LCDDConverter.h:114
dd4hep::_toDouble
double _toDouble(const std::string &value)
String conversions: string to double value.
Definition: Handle.cpp:116
dd4hep::xml::DocumentHandler::output
virtual int output(Document doc, const std::string &fname) const
Write xml document to output file (stdout if file name empty)
Definition: DocumentHandler.cpp:395
DEGREE_2_RAD
#define DEGREE_2_RAD
Definition: Handle.h:27
dd4hep::detail::LCDDConverter::handleLimitSet
virtual xml_h handleLimitSet(const std::string &name, LimitSet limitset) const
Convert the geometry type LimitSet into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:860
dd4hep::VisAttr::visible
bool visible() const
Get visibility flag.
Definition: Objects.cpp:340
dd4hep::detail::LCDDConverter::handleVolumeVis
virtual xml_h handleVolumeVis(const std::string &name, const TGeoVolume *volume) const
Dump logical volume in GDML format to output stream.
Definition: LCDDConverter.cpp:738
dd4hep::xml::DocumentHandler
Class supporting to read and parse XML documents.
Definition: DocumentHandler.h:39
dd4hep::VisAttr::SOLID
@ SOLID
Definition: Objects.h:326
dd4hep::Handle::ptr
T * ptr() const
Access to the held object.
Definition: Handle.h:151
dd4hep::detail::LCDDConverter::handleHeader
virtual void handleHeader() const
Add header information in Detector format.
Definition: LCDDConverter.cpp:1083
ObjectsInterna.h
dd4hep::detail::LCDDConverter::GeometryInfo::xmlIdSpecs
IdSpecMap xmlIdSpecs
Definition: LCDDConverter.h:80
dd4hep::xml::Element
User abstraction class to manipulate XML elements within a document.
Definition: XMLElements.h:769
dd4hep::Segmentation::parameters
DDSegmentation::Parameters parameters() const
Access to the parameters.
Definition: Segmentations.cpp:57
dd4hep::SensitiveDetector::energyCutoff
double energyCutoff() const
Access energy cut off.
Definition: DetElement.cpp:436
dd4hep::detail::LCDDConverter::GeometryInfo::check
void check(const std::string &name, const TNamed *n, std::map< std::string, const TNamed * > &array) const
Definition: LCDDConverter.cpp:115
dd4hep
Namespace for the AIDA detector description toolkit.
Definition: AlignmentsCalib.h:28
dd4hep::VisAttr::WIREFRAME
@ WIREFRAME
Definition: Objects.h:326
dd4hep::xml::Element::setAttr
Attribute setAttr(const XmlChar *nam, const T &val) const
Set single attribute.
Definition: XMLElements.h:900
det
DetElement::Object * det
Definition: AlignmentsCalculator.cpp:66
dd4hep::SensitiveDetector::readout
Readout readout() const
Access readout structure of the sensitive detector.
Definition: DetElement.cpp:420
xml_elt_t
dd4hep::xml::Element xml_elt_t
Definition: ConditionsRepository.cpp:33
dd4hep::PlacedVolume::data
Object * data() const
Check if placement is properly instrumented.
Definition: Volumes.cpp:447
dd4hep::detail::LCDDConverter::GeometryInfo::xmlPlacements
PlacementMap xmlPlacements
Definition: LCDDConverter.h:76
dd4hep::detail::LCDDConverter::LCDDConverter
LCDDConverter(Detector &description)
Initializing Constructor.
Definition: LCDDConverter.cpp:125
LCDDConverter.h
dd4hep::Detector
The main interface to the dd4hep detector description package.
Definition: Detector.h:90
dd4hep::detail::LCDDConverter::handleSensitive
virtual xml_h handleSensitive(const std::string &name, SensitiveDetector sens_det) const
Convert the geometry type SensitiveDetector into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:914
TGeoConeSeg
Class of the ROOT toolkit. See http://root.cern.ch/root/htmldoc/ClassIndex.html.
Definition: ROOTClasses.h:17
dd4hep::Readout
Handle to the implementation of the readout structure of a subdetector.
Definition: Readout.h:38
dd4hep::Segmentation
Handle class supporting generic Segmentations of sensitive detectors.
Definition: Segmentations.h:41
dd4hep::detail::LCDDConverter::GeometryInfo::xmlMaterials
MaterialMap xmlMaterials
Definition: LCDDConverter.h:73
dd4hep::xml::Handle_t::setAttr
Attribute setAttr(const XmlChar *t, const XmlChar *v) const
Generic attribute setter with unicode value.
Definition: XMLElements.cpp:922
dd4hep::detail::LCDDConverter::handleSegmentation
virtual xml_h handleSegmentation(Segmentation seg) const
Convert the segmentation of a SensitiveDetector into the corresponding Detector object.
Definition: LCDDConverter.cpp:882
dd4hep::detail::LCDDConverter::GeometryInfo::doc_root
xml_elt_t doc_root
Definition: LCDDConverter.h:113
dd4hep::detail::LCDDConverter::GeometryInfo::doc_structure
xml_elt_t doc_structure
Definition: LCDDConverter.h:114
dd4hep::Readout::idSpec
IDDescriptor idSpec() const
Access IDDescription structure.
Definition: Readout.cpp:112
create_gdml_from_dd4hep
long create_gdml_from_dd4hep(Detector &description, int argc, char **argv)
Definition: LCDDConverter.cpp:1296
dd4hep::PlacedVolume::volIDs
const PlacedVolumeExtension::VolIDs & volIDs() const
Access to the volume IDs.
Definition: Volumes.cpp:496
dd4hep::Volume::visAttributes
VisAttr visAttributes() const
Access the visualisation attributes.
Definition: Volumes.cpp:1239
DD4hepUnits.h
dd4hep::detail::GeoHandlerTypes::GeometryInfo::volumes
std::vector< Volume > volumes
Definition: GeoHandler.h:65
FieldTypes.h
dd4hep::detail::LCDDConverter::GeometryInfo::GeometryInfo
GeometryInfo()
Helper constructor.
Definition: LCDDConverter.cpp:1284
dd4hep::xml::Handle_t::attr
T attr(const Attribute a) const
Access typed attribute value by the XmlAttr.
Definition: XMLElements.h:454
dd4hep::detail::GeoHandlerTypes::GeometryInfo::materials
std::set< Material > materials
Definition: GeoHandler.h:68
Printout.h
dd4hep::detail::LCDDConverter::GeometryInfo::xmlRotations
TrafoMap xmlRotations
Definition: LCDDConverter.h:83
dd4hep::detail::LCDDConverter::handlePosition
virtual xml_h handlePosition(const std::string &name, const TGeoMatrix *trafo) const
Convert the Position into the corresponding Xml object(s).
Definition: LCDDConverter.cpp:603
dd4hep::detail::LCDDConverter::GeometryInfo::doc_materials
xml_elt_t doc_materials
Definition: LCDDConverter.h:114
dd4hep::Header::version
const std::string & version() const
Accessor to object version.
Definition: Objects.cpp:117
dd4hep::Detector::header
virtual Header header() const =0
Accessor to the map of header entries.
handle
void handle(const O *o, const C &c, F pmf)
Definition: LCDDConverter.cpp:1104
dd4hep::detail::LCDDConverter::GeometryInfo::checkRotation
void checkRotation(const std::string &name, const TNamed *n) const
Definition: LCDDConverter.h:98
dd4hep::Volume::sensitiveDetector
Handle< NamedObject > sensitiveDetector() const
Access to the handle to the sensitive detector.
Definition: Volumes.cpp:1316
dd4hep::detail::LCDDConverter::handleVis
virtual xml_h handleVis(const std::string &name, VisAttr vis) const
Convert the geometry visualisation attributes to the corresponding Xml object(s).
Definition: LCDDConverter.cpp:975
handleRMap
void handleRMap(const O *o, const C &c, F pmf)
Definition: LCDDConverter.cpp:1116
dd4hep::Detector::properties
virtual Properties & properties() const =0
Access to properties map.
dd4hep::detail::GeoHandler::collect
GeoHandler & collect(DetElement top)
Collect geometry information from traversal.
Definition: GeoHandler.cpp:99
handleMap
void handleMap(const O *o, const C &c, F pmf)
Definition: LCDDConverter.cpp:1111
dd4hep::xml::Document::root
Handle_t root() const
Access the ROOT eleemnt of the DOM document.
Definition: XMLElements.cpp:1043