Logo ROOT   6.12/07
Reference Guide
TGDMLWrite.cxx
Go to the documentation of this file.
1 // @(#)root/gdml:$Id$
2 // Author: Anton Pytel 15/9/2011
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2011, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TGDMLWrite
13 \ingroup Geometry_gdml
14 
15  This class contains implementation of converting ROOT's gGeoManager
16 geometry to GDML file. gGeoManager is the instance of TGeoManager class
17 containing tree of geometries creating resulting geometry. GDML is xml
18 based format of file mirroring the tree of geometries according to GDML
19 schema rules. For more information about GDML see http://gdml.web.cern.ch.
20 Each object in ROOT is represented by xml tag (=xml node/element) in GDML.
21 
22  This class is not needed to be instanciated. It should always be called
23 by gGeoManager->Export("xyz.gdml") method. Export is driven by extenstion
24 that is why ".gdml" is important in resulting name.
25 
26  Whenever a new ROOT geometry object is implemented or there is a change
27 in GDML schema this class is needed to be updated to ensure proper mapping
28 between ROOT objects and GDML elements.
29 
30  Current status of mapping ROOT -> GDML is implemented in method called
31 TGDMLWrite::ChooseObject and it contains following "map":
32 
33 #### Solids:
34 
35 ~~~
36 TGeoBBox -> <box ... >
37 TGeoParaboloid -> <paraboloid ...>
38 TGeoSphere -> <sphere ...>
39 TGeoArb8 -> <arb8 ...>
40 TGeoConeSeg -> <cone ...>
41 TGeoCone -> <cone ...>
42 TGeoPara -> <para ...>
43 TGeoTrap -> <trap ...> or
44 - -> <arb8 ...>
45 TGeoGtra -> <twistedtrap ...> or
46 - -> <trap ...> or
47 - -> <arb8 ...>
48 TGeoTrd1 -> <trd ...>
49 TGeoTrd2 -> <trd ...>
50 TGeoTubeSeg -> <tube ...>
51 TGeoCtub -> <cutTube ...>
52 TGeoTube -> <tube ...>
53 TGeoPcon -> <polycone ...>
54 TGeoTorus -> <torus ...>
55 TGeoPgon -> <polyhedra ...>
56 TGeoEltu -> <eltube ...>
57 TGeoHype -> <hype ...>
58 TGeoXtru -> <xtru ...>
59 TGeoCompositeShape -> <union ...> or
60 - -> <subtraction ...> or
61 - -> <intersection ...>
62 
63 Special cases of solids:
64 TGeoScaledShape -> <elcone ...> if scaled TGeoCone or
65 - -> element without scale
66 TGeoCompositeShape -> <ellipsoid ...>
67 - intersection of:
68 - scaled TGeoSphere and TGeoBBox
69 ~~~
70 
71 #### Materials:
72 
73 ~~~
74 TGeoIsotope -> <isotope ...>
75 TGeoElement -> <element ...>
76 TGeoMaterial -> <material ...>
77 TGeoMixture -> <material ...>
78 ~~~
79 
80 #### Structure
81 
82 ~~~
83 TGeoVolume -> <volume ...> or
84 - -> <assembly ...>
85 TGeoNode -> <physvol ...>
86 TGeoPatternFinder -> <divisionvol ...>
87 ~~~
88 
89 There are options that can be set to change resulting document
90 
91 ##### Options:
92 
93 ~~~
94 g - is set by default in gGeoManager, this option ensures compatibility
95 - with Geant4. It means:
96 - -> atomic number of material will be changed if <1 to 1
97 - -> if polycone is set badly it will try to export it correctly
98 - -> if widht * ndiv + offset is more then width of object being divided
99 - (in divisions) then it will be rounded so it will not exceed or
100 - if kPhi divsion then it will keep range of offset in -360 -> 0
101 f - if this option is set then names of volumes and solids will have
102 - pointer as a suffix to ensure uniqness of names
103 n - if this option is set then names will not have suffix, but uniqness is
104 - of names is not secured
105 - - if none of this two options (f,n) is set then default behaviour is so
106 - that incremental suffix is added to the names.
107 - (eg. TGeoBBox_0x1, TGeoBBox_0x2 ...)
108 ~~~
109 
110 #### USAGE:
111 
112 ~~~
113 gGeoManager->Export("output.gdml");
114 gGeoManager->Export("output.gdml","","vg"); //the same as previous just
115  //options are set explicitly
116 gGeoManager->Export("output.gdml","","vgf");
117 gGeoManager->Export("output.gdml","","gn");
118 gGeoManager->Export("output.gdml","","f");
119 ...
120 ~~~
121 
122 #### Note:
123  Options discussed above are used only for TGDMLWrite class. There are
124 other options in the TGeoManager::Export(...) method that can be used.
125 See that function for details.
126 
127 */
128 
129 #include "TGeoManager.h"
130 #include "TGeoMaterial.h"
131 #include "TGeoMatrix.h"
132 #include "TXMLEngine.h"
133 #include "TGeoVolume.h"
134 #include "TGeoBBox.h"
135 #include "TGeoParaboloid.h"
136 #include "TGeoArb8.h"
137 #include "TGeoTube.h"
138 #include "TGeoCone.h"
139 #include "TGeoTrd1.h"
140 #include "TGeoTrd2.h"
141 #include "TGeoPcon.h"
142 #include "TGeoPgon.h"
143 #include "TGeoSphere.h"
144 #include "TGeoTorus.h"
145 #include "TGeoPara.h"
146 #include "TGeoHype.h"
147 #include "TGeoEltu.h"
148 #include "TGeoXtru.h"
149 #include "TGeoScaledShape.h"
150 #include "TGeoVolume.h"
151 #include "TROOT.h"
152 #include "TMath.h"
153 #include "TGeoBoolNode.h"
154 #include "TGeoMedium.h"
155 #include "TGeoElement.h"
156 #include "TGeoShape.h"
157 #include "TGeoCompositeShape.h"
158 #include "TGDMLWrite.h"
159 #include <stdlib.h>
160 #include <string>
161 #include <map>
162 #include <ctime>
163 
165 
167 
168 ////////////////////////////////////////////////////////////////////////////////
169 /// Default constructor.
170 
172  : TObject(),
173  fIsotopeList(0),
174  fElementList(0),
175  fAccPatt(0),
176  fRejShape(0),
177  fNameList(0),
178  fgNamingSpeed(0),
179  fgG4Compatibility(0),
180  fGdmlFile(0),
181  fTopVolumeName(0),
182  fGdmlE(0),
183  fDefineNode(0),
184  fMaterialsNode(0),
185  fSolidsNode(0),
186  fStructureNode(0),
187  fVolCnt(0),
188  fPhysVolCnt(0),
189  fActNameErr(0),
190  fSolCnt(0)
191 {
192  if (fgGDMLWrite) delete fgGDMLWrite;
193  fgGDMLWrite = this;
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// Destructor.
198 
200 {
201  delete fIsotopeList;
202  delete fElementList;
203  delete fAccPatt;
204  delete fRejShape;
205  delete fNameList;
206 
207  fgGDMLWrite = 0;
208 }
209 
210 ////////////////////////////////////////////////////////////////////////////////
211 /// Set convention of naming solids and volumes
212 
214 {
215  fgNamingSpeed = naming;
216 }
217 
218 ////////////////////////////////////////////////////////////////////////////////
219 /// Wrapper of all exporting methods
220 /// Creates blank GDML file and fills it with gGeoManager structure converted
221 /// to GDML structure of xml nodes
222 
223 void TGDMLWrite::WriteGDMLfile(TGeoManager * geomanager, const char* filename, TString option)
224 {
225  //option processing
226  option.ToLower();
227  if (option.Contains("g")) {
229  Info("WriteGDMLfile", "Geant4 compatibility mode set");
230  } else {
232  }
233  if (option.Contains("f")) {
235  Info("WriteGDMLfile", "Fast naming convention with pointer suffix set");
236  } else if (option.Contains("n")) {
238  Info("WriteGDMLfile", "Naming without prefix set - be careful uniqness of name is not ensured");
239  } else {
241  Info("WriteGDMLfile", "Potentially slow with incremental suffix naming convention set");
242  }
243 
244  //local variables
245  Int_t outputLayout = 1;
246  const char * krootNodeName = "gdml";
247  const char * knsRefGeneral = "http://www.w3.org/2001/XMLSchema-instance";
248  const char * knsNameGeneral = "xsi";
249  const char * knsRefGdml = "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd";
250  const char * knsNameGdml = "xsi:noNamespaceSchemaLocation";
251 
252  // First create engine
253  fGdmlE = new TXMLEngine;
255 
256  //create blank GDML file
257  fGdmlFile = fGdmlE->NewDoc();
258 
259  //create root node and add it to blank GDML file
260  XMLNodePointer_t rootNode = fGdmlE->NewChild(0, 0, krootNodeName, 0);
261  fGdmlE->DocSetRootElement(fGdmlFile, rootNode);
262 
263  //add namespaces to root node
264  fGdmlE->NewNS(rootNode, knsRefGeneral, knsNameGeneral);
265  fGdmlE->NewAttr(rootNode, 0, knsNameGdml, knsRefGdml);
266 
267  //initialize general lists and <define>, <solids>, <structure> nodes
268  fIsotopeList = new StructLst;
269  fElementList = new StructLst;
270 
271  fNameList = new NameLst;
272 
273  fDefineNode = fGdmlE->NewChild(0, 0, "define", 0);
274  fSolidsNode = fGdmlE->NewChild(0, 0, "solids", 0);
275  fStructureNode = fGdmlE->NewChild(0, 0, "structure", 0);
276  //========================
277 
278  //initialize list of accepted patterns for divisions (in ExtractVolumes)
279  fAccPatt = new StructLst;
280  fAccPatt->fLst["TGeoPatternX"] = kTRUE;
281  fAccPatt->fLst["TGeoPatternY"] = kTRUE;
282  fAccPatt->fLst["TGeoPatternZ"] = kTRUE;
283  fAccPatt->fLst["TGeoPatternCylR"] = kTRUE;
284  fAccPatt->fLst["TGeoPatternCylPhi"] = kTRUE;
285  //========================
286 
287  //initialize list of rejected shapes for divisions (in ExtractVolumes)
288  fRejShape = new StructLst;
289  //this shapes are rejected because, it is not possible to divide trd2
290  //in Y axis and while only trd2 object is imported from GDML
291  //it causes a problem when TGeoTrd1 is divided in Y axis
292  fRejShape->fLst["TGeoTrd1"] = kTRUE;
293  fRejShape->fLst["TGeoTrd2"] = kTRUE;
294  //=========================
295 
296  //Initialize global counters
297  fActNameErr = 0;
298  fVolCnt = 0;
299  fPhysVolCnt = 0;
300  fSolCnt = 0;
301  fTopVolumeName = "";
302 
303  //calling main extraction functions (with measuring time)
304  time_t startT, endT;
305  startT = time(NULL);
307 
308  Info("WriteGDMLfile", "Extracting volumes");
309  if (geomanager->GetTopVolume()) {
310  ExtractVolumes(geomanager->GetTopVolume());
311  } else {
312  Info("WriteGDMLfile", "Top volume does not exist!");
313  return;
314  }
315  Info("WriteGDMLfile", "%i solids added", fSolCnt);
316  Info("WriteGDMLfile", "%i volumes added", fVolCnt);
317  Info("WriteGDMLfile", "%i physvolumes added", fPhysVolCnt);
318  endT = time(NULL);
319  //<gdml>
320  fGdmlE->AddChild(rootNode, fDefineNode); // <define>...</define>
321  fGdmlE->AddChild(rootNode, fMaterialsNode); // <materials>...</materials>
322  fGdmlE->AddChild(rootNode, fSolidsNode); // <solids>...</solids>
323  fGdmlE->AddChild(rootNode, fStructureNode); // <structure>...</structure>
324  fGdmlE->AddChild(rootNode, CreateSetupN(fTopVolumeName.Data())); // <setup>...</setup>
325  //</gdml>
326  Double_t tdiffI = difftime(endT, startT);
327  TString tdiffS = (tdiffI == 0 ? TString("< 1 s") : TString::Format("%.0lf s", tdiffI));
328  Info("WriteGDMLfile", "Exporting time: %s", tdiffS.Data());
329  //=========================
330 
331  //Saving document
332  fGdmlE->SaveDoc(fGdmlFile, filename, outputLayout);
333  Info("WriteGDMLfile", "File %s saved", filename);
334  //cleaning
336  //unset processing bits:
337  UnsetTemporaryBits(geomanager);
338  delete fGdmlE;
339 }
340 
341 
342 ////////////////////////////////////////////////////////////////////////////////
343 /// Method exporting materials
344 
346 {
347  Info("ExtractMaterials", "Extracting materials");
348  //crate main <materials> node
349  XMLNodePointer_t materialsN = fGdmlE->NewChild(0, 0, "materials", 0);
350  Int_t matcnt = 0;
351 
352  //go through materials - iterator and object declaration
353  TIter next(materialsLst);
354  TGeoMaterial *lmaterial;
355 
356  while ((lmaterial = (TGeoMaterial *)next())) {
357  //generate uniq name
358  TString lname = GenName(lmaterial->GetName(), TString::Format("%p", lmaterial));
359 
360  if (lmaterial->IsMixture()) {
361  TGeoMixture *lmixture = (TGeoMixture *)lmaterial;
362  XMLNodePointer_t mixtureN = CreateMixtureN(lmixture, materialsN, lname);
363  fGdmlE->AddChild(materialsN, mixtureN);
364  } else {
365  XMLNodePointer_t materialN = CreateMaterialN(lmaterial, lname);
366  fGdmlE->AddChild(materialsN, materialN);
367  }
368  matcnt++;
369  }
370  Info("ExtractMaterials", "%i materials added", matcnt);
371  return materialsN;
372 }
373 
374 ////////////////////////////////////////////////////////////////////////////////
375 /// Method creating solid to xml file and returning its name
376 
378 {
379  XMLNodePointer_t solidN;
380  TString solname = "";
381  solidN = ChooseObject(volShape); //volume->GetShape()
382  fGdmlE->AddChild(fSolidsNode, solidN);
383  if (solidN != NULL) fSolCnt++;
384  solname = fNameList->fLst[TString::Format("%p", volShape)];
385  if (solname.Contains("missing_")) {
386  solname = "-1";
387  }
388  return solname;
389 }
390 
391 
392 ////////////////////////////////////////////////////////////////////////////////
393 /// Method extracting geometry structure recursively
394 
396 {
397  XMLNodePointer_t volumeN, childN;
398  TString volname, matname, solname, pattClsName, nodeVolNameBak;
399  TGeoPatternFinder *pattFinder = 0;
400  Bool_t isPattern = kFALSE;
401 
402  //create the name for volume/assembly
403  if (volume->IsTopVolume()) {
404  //not needed a special function for generating name
405  volname = volume->GetName();
406  fTopVolumeName = volname;
407  //register name to the pointer
408  fNameList->fLst[TString::Format("%p", volume)] = volname;
409  } else {
410  volname = GenName(volume->GetName(), TString::Format("%p", volume));
411  }
412 
413  //start to create main volume/assembly node
414  if (volume->IsAssembly()) {
415  volumeN = StartAssemblyN(volname);
416  } else {
417  //get reference material and add solid to <solids> + get name
418  matname = fNameList->fLst[TString::Format("%p", volume->GetMaterial())];
419  solname = ExtractSolid(volume->GetShape());
420  //If solid is not supported or corrupted
421  if (solname == "-1") {
422  Info("ExtractVolumes", "ERROR! %s volume was not added, because solid is either not supported or corrupted",
423  volname.Data());
424  //set volume as missing volume
425  fNameList->fLst[TString::Format("%p", volume)] = "missing_" + volname;
426  return;
427  }
428  volumeN = StartVolumeN(volname, solname, matname);
429 
430  //divisionvol can't be in assembly
431  pattFinder = volume->GetFinder();
432  //if found pattern
433  if (pattFinder) {
434  pattClsName = TString::Format("%s", pattFinder->ClassName());
435  TString shapeCls = TString::Format("%s", volume->GetShape()->ClassName());
436  //if pattern in accepted pattern list and not in shape rejected list
437  if ((fAccPatt->fLst[pattClsName] == kTRUE) &&
438  (fRejShape->fLst[shapeCls] != kTRUE)) {
439  isPattern = kTRUE;
440  }
441  }
442  }
443  //get all nodes in volume
444  TObjArray *nodeLst = volume->GetNodes();
445  TIter next(nodeLst);
446  TGeoNode *geoNode;
447  Int_t nCnt = 0;
448  //loop through all nodes
449  while ((geoNode = (TGeoNode *) next())) {
450  //get volume of current node and if not processed then process it
451  TGeoVolume * subvol = geoNode->GetVolume();
452  if (subvol->TestAttBit(fgkProcBitVol) == kFALSE) {
453  subvol->SetAttBit(fgkProcBitVol);
454  ExtractVolumes(subvol);
455  }
456 
457  //volume of this node has to exist because it was processed recursively
458  TString nodevolname = fNameList->fLst[TString::Format("%p", geoNode->GetVolume())];
459  if (nodevolname.Contains("missing_")) {
460  continue;
461  }
462  if (nCnt == 0) { //save name of the first node for divisionvol
463  nodeVolNameBak = nodevolname;
464  }
465 
466  if (isPattern == kFALSE) {
467  //create name for node
468  TString nodename, posname, rotname;
469  nodename = GenName(geoNode->GetName(), TString::Format("%p", geoNode));
470  nodename = nodename + "in" + volname;
471 
472  //create name for position and clear rotation
473  posname = nodename + "pos";
474  rotname = "";
475 
476  //position
477  const Double_t * pos = geoNode->GetMatrix()->GetTranslation();
478  Xyz nodPos;
479  nodPos.x = pos[0];
480  nodPos.y = pos[1];
481  nodPos.z = pos[2];
482  childN = CreatePositionN(posname.Data(), nodPos);
483  fGdmlE->AddChild(fDefineNode, childN); //adding node to <define> node
484  //Deal with reflection
485  XMLNodePointer_t scaleN = NULL;
486  Double_t lx, ly, lz;
487  Double_t xangle = 0;
488  Double_t zangle = 0;
489  lx = geoNode->GetMatrix()->GetRotationMatrix()[0];
490  ly = geoNode->GetMatrix()->GetRotationMatrix()[4];
491  lz = geoNode->GetMatrix()->GetRotationMatrix()[8];
492  if (geoNode->GetMatrix()->IsReflection()
493  && TMath::Abs(lx) == 1 && TMath::Abs(ly) == 1 && TMath::Abs(lz) == 1) {
494  scaleN = fGdmlE->NewChild(0, 0, "scale", 0);
495  fGdmlE->NewAttr(scaleN, 0, "name", (nodename + "scl").Data());
496  fGdmlE->NewAttr(scaleN, 0, "x", TString::Format("%.12g", lx));
497  fGdmlE->NewAttr(scaleN, 0, "y", TString::Format("%.12g", ly));
498  fGdmlE->NewAttr(scaleN, 0, "z", TString::Format("%.12g", lz));
499  //experimentally found out, that rotation should be updated like this
500  if (lx == -1) {
501  zangle = 180;
502  }
503  if (lz == -1) {
504  xangle = 180;
505  }
506  }
507 
508  //rotation
510  lxyz.x -= xangle;
511  lxyz.z -= zangle;
512  if ((lxyz.x != 0.0) || (lxyz.y != 0.0) || (lxyz.z != 0.0)) {
513  rotname = nodename + "rot";
514  childN = CreateRotationN(rotname.Data(), lxyz);
515  fGdmlE->AddChild(fDefineNode, childN); //adding node to <define> node
516  }
517 
518  //create physvol for main volume/assembly node
519  childN = CreatePhysVolN(geoNode->GetName(), geoNode->GetNumber(), nodevolname.Data(), posname.Data(), rotname.Data(), scaleN);
520  fGdmlE->AddChild(volumeN, childN);
521  }
522  nCnt++;
523  }
524  //create only one divisionvol node
525  if (isPattern && pattFinder) {
526  //retrieve attributes of division
527  Int_t ndiv, divaxis;
528  Double_t offset, width, xlo, xhi;
529  TString axis, unit;
530 
531  ndiv = pattFinder->GetNdiv();
532  width = pattFinder->GetStep();
533 
534  divaxis = pattFinder->GetDivAxis();
535  volume->GetShape()->GetAxisRange(divaxis, xlo, xhi);
536 
537  //compute relative start (not positional)
538  offset = pattFinder->GetStart() - xlo;
539  axis = GetPattAxis(divaxis, pattClsName, unit);
540 
541  //create division node
542  childN = CreateDivisionN(offset, width, ndiv, axis.Data(), unit.Data(), nodeVolNameBak.Data());
543  fGdmlE->AddChild(volumeN, childN);
544  }
545 
546  fVolCnt++;
547  //add volume/assembly node into the <structure> node
548  fGdmlE->AddChild(fStructureNode, volumeN);
549 
550 }
551 
552 ////////////////////////////////////////////////////////////////////////////////
553 /// Creates "atom" node for GDML
554 
556 {
557  XMLNodePointer_t atomN = fGdmlE->NewChild(0, 0, "atom", 0);
558  fGdmlE->NewAttr(atomN, 0, "unit", unit);
559  fGdmlE->NewAttr(atomN, 0, "value", TString::Format("%.12g", atom));
560  return atomN;
561 }
562 
563 ////////////////////////////////////////////////////////////////////////////////
564 /// Creates "D" density node for GDML
565 
566 XMLNodePointer_t TGDMLWrite::CreateDN(Double_t density, const char * unit)
567 {
568  XMLNodePointer_t densN = fGdmlE->NewChild(0, 0, "D", 0);
569  fGdmlE->NewAttr(densN, 0, "unit", unit);
570  fGdmlE->NewAttr(densN, 0, "value", TString::Format("%.12g", density));
571  return densN;
572 }
573 
574 ////////////////////////////////////////////////////////////////////////////////
575 /// Creates "fraction" node for GDML
576 
577 XMLNodePointer_t TGDMLWrite::CreateFractionN(Double_t percentage, const char * refName)
578 {
579  XMLNodePointer_t fractN = fGdmlE->NewChild(0, 0, "fraction", 0);
580  fGdmlE->NewAttr(fractN, 0, "n", TString::Format("%.12g", percentage));
581  fGdmlE->NewAttr(fractN, 0, "ref", refName);
582  return fractN;
583 }
584 
585 ////////////////////////////////////////////////////////////////////////////////
586 /// Creates "isotope" node for GDML
587 
589 {
590  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "isotope", 0);
591  fGdmlE->NewAttr(mainN, 0, "name", name);
592  fGdmlE->NewAttr(mainN, 0, "N", TString::Format("%i", isotope->GetN()));
593  fGdmlE->NewAttr(mainN, 0, "Z", TString::Format("%i", isotope->GetZ()));
594  fGdmlE->AddChild(mainN, CreateAtomN(isotope->GetA()));
595  return mainN;
596 }
597 
598 ////////////////////////////////////////////////////////////////////////////////
599 /// Creates "element" node for GDML
600 ///element node and attribute
601 
603 {
604  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "element", 0);
605  fGdmlE->NewAttr(mainN, 0, "name", name);
606  //local associative arrays for saving isotopes and their weight
607  //inside element
608  NameListF wPercentage;
609  NameListI wCounter;
610 
611  if (element->HasIsotopes()) {
612  Int_t nOfIso = element->GetNisotopes();
613  //go through isotopes
614  for (Int_t idx = 0; idx < nOfIso; idx++) {
615  TGeoIsotope *myIsotope = element->GetIsotope(idx);
616  if (!myIsotope) {
617  Fatal("CreateElementN", "Missing isotopes for element %s", element->GetName());
618  return mainN;
619  }
620 
621  //Get name of the Isotope (
622  TString lname = myIsotope->GetName();
623  //_iso suffix is added to avoid problems with same names
624  //for material, element and isotopes
625  lname = TString::Format("%s_iso", lname.Data());
626 
627  //cumulates abundance, in case 2 isotopes with same names
628  //within one element
629  wPercentage[lname] += element->GetRelativeAbundance(idx);
630  wCounter[lname]++;
631 
632  //check whether isotope name is not in list of isotopes
633  if (IsInList(fIsotopeList->fLst, lname)) {
634  continue;
635  }
636  //add isotope to list of isotopes and to main <materials> node
637  fIsotopeList->fLst[lname] = kTRUE;
638  XMLNodePointer_t isoNode = CreateIsotopN(myIsotope, lname);
639  fGdmlE->AddChild(materials, isoNode);
640  }
641  //loop through asoc array of isotopes
642  for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
643  if (itr->second > 1) {
644  Info("CreateMixtureN", "WARNING! 2 equal isotopes in one element. Check: %s isotope of %s element",
645  itr->first.Data(), name);
646  }
647  //add fraction child to element with reference to isotope
648  fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
649  }
650  } else {
651  fGdmlE->NewAttr(mainN, 0, "formula", element->GetName());
652  Int_t valZ = element->Z();
653  // Z can't be <1 in Geant4 and Z is optional parameter
654  if (valZ >= 1) {
655  fGdmlE->NewAttr(mainN, 0, "Z", TString::Format("%i", valZ));
656  }
657  fGdmlE->AddChild(mainN, CreateAtomN(element->A()));
658  }
659  return mainN;
660 }
661 
662 ////////////////////////////////////////////////////////////////////////////////
663 /// Creates "material" node for GDML with references to other sub elements
664 
666 {
667  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "material", 0);
668  fGdmlE->NewAttr(mainN, 0, "name", mname);
669  fGdmlE->AddChild(mainN, CreateDN(mixture->GetDensity()));
670  //local associative arrays for saving elements and their weight
671  //inside mixture
672  NameListF wPercentage;
673  NameListI wCounter;
674 
675  Int_t nOfElm = mixture->GetNelements();
676  //go through elements
677  for (Int_t idx = 0; idx < nOfElm; idx++) {
678  TGeoElement *myElement = mixture->GetElement(idx);
679 
680  //Get name of the element
681  //NOTE: that for element - GetTitle() returns the "name" tag
682  //and GetName() returns "formula" tag (see createElementN)
683  TString lname = myElement->GetTitle();
684  //_elm suffix is added to avoid problems with same names
685  //for material and element
686  lname = TString::Format("%s_elm", lname.Data());
687 
688  //cumulates percentage, in case 2 elements with same names within one mixture
689  wPercentage[lname] += mixture->GetWmixt()[idx];
690  wCounter[lname]++;
691 
692  //check whether element name is not in list of elements already created
693  if (IsInList(fElementList->fLst, lname)) {
694  continue;
695  }
696 
697  //add element to list of elements and to main <materials> node
698  fElementList->fLst[lname] = kTRUE;
699  XMLNodePointer_t elmNode = CreateElementN(myElement, materials, lname);
700  fGdmlE->AddChild(materials, elmNode);
701  }
702  //loop through asoc array
703  for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
704  if (itr->second > 1) {
705  Info("CreateMixtureN", "WARNING! 2 equal elements in one material. Check: %s element of %s material",
706  itr->first.Data(), mname.Data());
707  }
708  //add fraction child to material with reference to element
709  fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
710  }
711 
712  return mainN;
713 }
714 
715 ////////////////////////////////////////////////////////////////////////////////
716 /// Creates "material" node for GDML
717 
719 {
720  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "material", 0);
721  fGdmlE->NewAttr(mainN, 0, "name", mname);
722  Double_t valZ = material->GetZ();
723  //Z can't be zero in Geant4 so this is workaround for vacuum
724  TString tmpname = mname;
725  tmpname.ToLower();
726  if (valZ < 1) {
727  if (tmpname == "vacuum") {
728  valZ = 1;
729  } else {
730  if (fgG4Compatibility == kTRUE) {
731  Info("CreateMaterialN", "WARNING! value of Z in %s material can't be < 1 in Geant4, that is why it was changed to 1, please check it manually! ",
732  mname.Data());
733  valZ = 1;
734  } else {
735  Info("CreateMaterialN", "WARNING! value of Z in %s material can't be < 1 in Geant4", mname.Data());
736  }
737  }
738  }
739  fGdmlE->NewAttr(mainN, 0, "Z", TString::Format("%.12g", valZ)); //material->GetZ()));
740  fGdmlE->AddChild(mainN, CreateDN(material->GetDensity()));
741  fGdmlE->AddChild(mainN, CreateAtomN(material->GetA()));
742  return mainN;
743 }
744 
745 ////////////////////////////////////////////////////////////////////////////////
746 /// Creates "box" node for GDML
747 
749 {
750  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "box", 0);
751  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
752  fGdmlE->NewAttr(mainN, 0, "name", lname);
753  if (IsNullParam(geoShape->GetDX(), "DX", lname) ||
754  IsNullParam(geoShape->GetDY(), "DY", lname) ||
755  IsNullParam(geoShape->GetDZ(), "DZ", lname)) {
756  return NULL;
757  }
758  fGdmlE->NewAttr(mainN, 0, "x", TString::Format("%.12g", 2 * geoShape->GetDX()));
759  fGdmlE->NewAttr(mainN, 0, "y", TString::Format("%.12g", 2 * geoShape->GetDY()));
760  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDZ()));
761 
762  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
763  return mainN;
764 }
765 
766 ////////////////////////////////////////////////////////////////////////////////
767 /// Creates "paraboloid" node for GDML
768 
770 {
771  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "paraboloid", 0);
772  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
773  fGdmlE->NewAttr(mainN, 0, "name", lname);
774  if (IsNullParam(geoShape->GetRhi(), "Rhi", lname) ||
775  IsNullParam(geoShape->GetDz(), "Dz", lname)) {
776  return NULL;
777  }
778  fGdmlE->NewAttr(mainN, 0, "rlo", TString::Format("%.12g", geoShape->GetRlo()));
779  fGdmlE->NewAttr(mainN, 0, "rhi", TString::Format("%.12g", geoShape->GetRhi()));
780  fGdmlE->NewAttr(mainN, 0, "dz", TString::Format("%.12g", geoShape->GetDz()));
781 
782  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
783  return mainN;
784 }
785 
786 ////////////////////////////////////////////////////////////////////////////////
787 /// Creates "sphere" node for GDML
788 
790 {
791  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "sphere", 0);
792  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
793  fGdmlE->NewAttr(mainN, 0, "name", lname);
794  if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
795  return NULL;
796  }
797 
798  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
799  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
800  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
801  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetPhi2() - geoShape->GetPhi1()));
802  fGdmlE->NewAttr(mainN, 0, "starttheta", TString::Format("%.12g", geoShape->GetTheta1()));
803  fGdmlE->NewAttr(mainN, 0, "deltatheta", TString::Format("%.12g", geoShape->GetTheta2() - geoShape->GetTheta1()));
804 
805  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
806  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
807  return mainN;
808 }
809 
810 ////////////////////////////////////////////////////////////////////////////////
811 /// Creates "arb8" node for GDML
812 
814 {
815  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "arb8", 0);
816  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
817  fGdmlE->NewAttr(mainN, 0, "name", lname);
818  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
819  return NULL;
820  }
821 
822  fGdmlE->NewAttr(mainN, 0, "v1x", TString::Format("%.12g", geoShape->GetVertices()[0]));
823  fGdmlE->NewAttr(mainN, 0, "v1y", TString::Format("%.12g", geoShape->GetVertices()[1]));
824  fGdmlE->NewAttr(mainN, 0, "v2x", TString::Format("%.12g", geoShape->GetVertices()[2]));
825  fGdmlE->NewAttr(mainN, 0, "v2y", TString::Format("%.12g", geoShape->GetVertices()[3]));
826  fGdmlE->NewAttr(mainN, 0, "v3x", TString::Format("%.12g", geoShape->GetVertices()[4]));
827  fGdmlE->NewAttr(mainN, 0, "v3y", TString::Format("%.12g", geoShape->GetVertices()[5]));
828  fGdmlE->NewAttr(mainN, 0, "v4x", TString::Format("%.12g", geoShape->GetVertices()[6]));
829  fGdmlE->NewAttr(mainN, 0, "v4y", TString::Format("%.12g", geoShape->GetVertices()[7]));
830  fGdmlE->NewAttr(mainN, 0, "v5x", TString::Format("%.12g", geoShape->GetVertices()[8]));
831  fGdmlE->NewAttr(mainN, 0, "v5y", TString::Format("%.12g", geoShape->GetVertices()[9]));
832  fGdmlE->NewAttr(mainN, 0, "v6x", TString::Format("%.12g", geoShape->GetVertices()[10]));
833  fGdmlE->NewAttr(mainN, 0, "v6y", TString::Format("%.12g", geoShape->GetVertices()[11]));
834  fGdmlE->NewAttr(mainN, 0, "v7x", TString::Format("%.12g", geoShape->GetVertices()[12]));
835  fGdmlE->NewAttr(mainN, 0, "v7y", TString::Format("%.12g", geoShape->GetVertices()[13]));
836  fGdmlE->NewAttr(mainN, 0, "v8x", TString::Format("%.12g", geoShape->GetVertices()[14]));
837  fGdmlE->NewAttr(mainN, 0, "v8y", TString::Format("%.12g", geoShape->GetVertices()[15]));
838  fGdmlE->NewAttr(mainN, 0, "dz", TString::Format("%.12g", geoShape->GetDz()));
839 
840  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
841  return mainN;
842 }
843 
844 ////////////////////////////////////////////////////////////////////////////////
845 /// Creates "cone" node for GDML from TGeoConeSeg object
846 
848 {
849  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "cone", 0);
850  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
851  fGdmlE->NewAttr(mainN, 0, "name", lname);
852  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
853  return NULL;
854  }
855 
856  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
857  fGdmlE->NewAttr(mainN, 0, "rmin1", TString::Format("%.12g", geoShape->GetRmin1()));
858  fGdmlE->NewAttr(mainN, 0, "rmin2", TString::Format("%.12g", geoShape->GetRmin2()));
859  fGdmlE->NewAttr(mainN, 0, "rmax1", TString::Format("%.12g", geoShape->GetRmax1()));
860  fGdmlE->NewAttr(mainN, 0, "rmax2", TString::Format("%.12g", geoShape->GetRmax2()));
861  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
862  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetPhi2() - geoShape->GetPhi1()));
863 
864  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
865  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
866  return mainN;
867 }
868 
869 ////////////////////////////////////////////////////////////////////////////////
870 /// Creates "cone" node for GDML from TGeoCone object
871 
873 {
874  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "cone", 0);
875  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
876  fGdmlE->NewAttr(mainN, 0, "name", lname);
877  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
878  return NULL;
879  }
880 
881  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
882  fGdmlE->NewAttr(mainN, 0, "rmin1", TString::Format("%.12g", geoShape->GetRmin1()));
883  fGdmlE->NewAttr(mainN, 0, "rmin2", TString::Format("%.12g", geoShape->GetRmin2()));
884  fGdmlE->NewAttr(mainN, 0, "rmax1", TString::Format("%.12g", geoShape->GetRmax1()));
885  fGdmlE->NewAttr(mainN, 0, "rmax2", TString::Format("%.12g", geoShape->GetRmax2()));
886  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%i", 0));
887  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%i", 360));
888 
889  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
890  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
891  return mainN;
892 }
893 
894 ////////////////////////////////////////////////////////////////////////////////
895 /// Creates "para" node for GDML
896 
898 {
899  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "para", 0);
900  fGdmlE->NewAttr(mainN, 0, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
901 
902  fGdmlE->NewAttr(mainN, 0, "x", TString::Format("%.12g", 2 * geoShape->GetX()));
903  fGdmlE->NewAttr(mainN, 0, "y", TString::Format("%.12g", 2 * geoShape->GetY()));
904  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetZ()));
905  fGdmlE->NewAttr(mainN, 0, "alpha", TString::Format("%.12g", geoShape->GetAlpha()));
906  fGdmlE->NewAttr(mainN, 0, "theta", TString::Format("%.12g", geoShape->GetTheta()));
907  fGdmlE->NewAttr(mainN, 0, "phi", TString::Format("%.12g", geoShape->GetPhi()));
908 
909  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
910  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
911  return mainN;
912 }
913 
914 ////////////////////////////////////////////////////////////////////////////////
915 /// Creates "trap" node for GDML
916 
918 {
919  XMLNodePointer_t mainN;
920 
921  //if one base equals 0 create Arb8 instead of trap
922  if ((geoShape->GetBl1() == 0 || geoShape->GetTl1() == 0 || geoShape->GetH1() == 0) ||
923  (geoShape->GetBl2() == 0 || geoShape->GetTl2() == 0 || geoShape->GetH2() == 0)) {
924  mainN = CreateArb8N(geoShape);
925  return mainN;
926  }
927 
928  //if is twisted then create arb8
929  if (geoShape->IsTwisted()) {
930  mainN = CreateArb8N((TGeoArb8 *) geoShape);
931  return mainN;
932  }
933 
934  mainN = fGdmlE->NewChild(0, 0, "trap", 0);
935  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
936  fGdmlE->NewAttr(mainN, 0, "name", lname);
937  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
938  return NULL;
939  }
940 
941  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
942  fGdmlE->NewAttr(mainN, 0, "theta", TString::Format("%.12g", geoShape->GetTheta()));
943  fGdmlE->NewAttr(mainN, 0, "phi", TString::Format("%.12g", geoShape->GetPhi()));
944  fGdmlE->NewAttr(mainN, 0, "x1", TString::Format("%.12g", 2 * geoShape->GetBl1()));
945  fGdmlE->NewAttr(mainN, 0, "x2", TString::Format("%.12g", 2 * geoShape->GetTl1()));
946  fGdmlE->NewAttr(mainN, 0, "x3", TString::Format("%.12g", 2 * geoShape->GetBl2()));
947  fGdmlE->NewAttr(mainN, 0, "x4", TString::Format("%.12g", 2 * geoShape->GetTl2()));
948  fGdmlE->NewAttr(mainN, 0, "y1", TString::Format("%.12g", 2 * geoShape->GetH1()));
949  fGdmlE->NewAttr(mainN, 0, "y2", TString::Format("%.12g", 2 * geoShape->GetH2()));
950 
951  fGdmlE->NewAttr(mainN, 0, "alpha1", TString::Format("%.12g", geoShape->GetAlpha1()));
952  fGdmlE->NewAttr(mainN, 0, "alpha2", TString::Format("%.12g", geoShape->GetAlpha2()));
953 
954  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
955  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
956  return mainN;
957 }
958 
959 ////////////////////////////////////////////////////////////////////////////////
960 /// Creates "twistedtrap" node for GDML
961 
963 {
964  XMLNodePointer_t mainN;
965 
966  //if one base equals 0 create Arb8 instead of twisted trap
967  if ((geoShape->GetBl1() == 0 && geoShape->GetTl1() == 0 && geoShape->GetH1() == 0) ||
968  (geoShape->GetBl2() == 0 && geoShape->GetTl2() == 0 && geoShape->GetH2() == 0)) {
969  mainN = CreateArb8N((TGeoArb8 *) geoShape);
970  return mainN;
971  }
972 
973  //if is twisted then create arb8
974  if (geoShape->IsTwisted()) {
975  mainN = CreateArb8N((TGeoArb8 *) geoShape);
976  return mainN;
977  }
978 
979  //if parameter twistAngle (PhiTwist) equals zero create trap node
980  if (geoShape->GetTwistAngle() == 0) {
981  mainN = CreateTrapN((TGeoTrap *) geoShape);
982  return mainN;
983  }
984 
985  mainN = fGdmlE->NewChild(0, 0, "twistedtrap", 0);
986  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
987  fGdmlE->NewAttr(mainN, 0, "name", lname);
988  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
989  return NULL;
990  }
991 
992  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
993  fGdmlE->NewAttr(mainN, 0, "Theta", TString::Format("%.12g", geoShape->GetTheta()));
994  fGdmlE->NewAttr(mainN, 0, "Phi", TString::Format("%.12g", geoShape->GetPhi()));
995  fGdmlE->NewAttr(mainN, 0, "x1", TString::Format("%.12g", 2 * geoShape->GetBl1()));
996  fGdmlE->NewAttr(mainN, 0, "x2", TString::Format("%.12g", 2 * geoShape->GetTl1()));
997  fGdmlE->NewAttr(mainN, 0, "x3", TString::Format("%.12g", 2 * geoShape->GetBl2()));
998  fGdmlE->NewAttr(mainN, 0, "x4", TString::Format("%.12g", 2 * geoShape->GetTl2()));
999  fGdmlE->NewAttr(mainN, 0, "y1", TString::Format("%.12g", 2 * geoShape->GetH1()));
1000  fGdmlE->NewAttr(mainN, 0, "y2", TString::Format("%.12g", 2 * geoShape->GetH2()));
1001 
1002  fGdmlE->NewAttr(mainN, 0, "Alph", TString::Format("%.12g", geoShape->GetAlpha1()));
1003 
1004  //check if alpha1 equals to alpha2 (converting to string - to avoid problems with floats)
1005  if (TString::Format("%.12g", geoShape->GetAlpha1()) != TString::Format("%.12g", geoShape->GetAlpha2())) {
1006  Info("CreateTwistedTrapN",
1007  "ERROR! Object %s is not exported correctly because parameter Alpha2 is not declared in GDML schema",
1008  lname.Data());
1009  }
1010  //fGdmlE->NewAttr(mainN,0, "alpha2", TString::Format("%.12g", geoShape->GetAlpha2()));
1011  fGdmlE->NewAttr(mainN, 0, "PhiTwist", TString::Format("%.12g", geoShape->GetTwistAngle()));
1012 
1013  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1014  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1015  return mainN;
1016 }
1017 
1018 ////////////////////////////////////////////////////////////////////////////////
1019 /// Creates "trd" node for GDML from object TGeoTrd1
1020 
1022 {
1023  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "trd", 0);
1024  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1025  fGdmlE->NewAttr(mainN, 0, "name", lname);
1026  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1027  return NULL;
1028  }
1029 
1030  fGdmlE->NewAttr(mainN, 0, "x1", TString::Format("%.12g", 2 * geoShape->GetDx1()));
1031  fGdmlE->NewAttr(mainN, 0, "x2", TString::Format("%.12g", 2 * geoShape->GetDx2()));
1032  fGdmlE->NewAttr(mainN, 0, "y1", TString::Format("%.12g", 2 * geoShape->GetDy()));
1033  fGdmlE->NewAttr(mainN, 0, "y2", TString::Format("%.12g", 2 * geoShape->GetDy()));
1034  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1035 
1036  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1037  return mainN;
1038 }
1039 
1040 ////////////////////////////////////////////////////////////////////////////////
1041 /// Creates "trd" node for GDML from object TGeoTrd2
1042 
1044 {
1045  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "trd", 0);
1046  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1047  fGdmlE->NewAttr(mainN, 0, "name", lname);
1048  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1049  return NULL;
1050  }
1051 
1052  fGdmlE->NewAttr(mainN, 0, "x1", TString::Format("%.12g", 2 * geoShape->GetDx1()));
1053  fGdmlE->NewAttr(mainN, 0, "x2", TString::Format("%.12g", 2 * geoShape->GetDx2()));
1054  fGdmlE->NewAttr(mainN, 0, "y1", TString::Format("%.12g", 2 * geoShape->GetDy1()));
1055  fGdmlE->NewAttr(mainN, 0, "y2", TString::Format("%.12g", 2 * geoShape->GetDy2()));
1056  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1057 
1058  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1059  return mainN;
1060 }
1061 
1062 ////////////////////////////////////////////////////////////////////////////////
1063 /// Creates "tube" node for GDML from object TGeoTubeSeg
1064 
1066 {
1067  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "tube", 0);
1068  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1069  fGdmlE->NewAttr(mainN, 0, "name", lname);
1070  if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) ||
1071  IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1072  return NULL;
1073  }
1074 
1075  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
1076  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
1077  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1078  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
1079  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetPhi2() - geoShape->GetPhi1()));
1080 
1081  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1082  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1083  return mainN;
1084 }
1085 
1086 ////////////////////////////////////////////////////////////////////////////////
1087 /// Creates "cutTube" node for GDML
1088 
1090 {
1091  XMLNodePointer_t mainN;
1092 
1093  mainN = fGdmlE->NewChild(0, 0, "cutTube", 0);
1094  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1095  fGdmlE->NewAttr(mainN, 0, "name", lname);
1096  if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) ||
1097  IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1098  return NULL;
1099  }
1100  //This is not needed, because cutTube is already supported by Geant4 9.5
1101  if (fgG4Compatibility == kTRUE && kFALSE) {
1102  TGeoShape * fakeCtub = CreateFakeCtub(geoShape);
1103  mainN = ChooseObject(fakeCtub);
1104 
1105  //register name for cuttube shape (so it will be found during volume export)
1106  lname = fNameList->fLst[TString::Format("%p", fakeCtub)];
1107  fNameList->fLst[TString::Format("%p", geoShape)] = lname;
1108  Info("CreateCutTubeN", "WARNING! %s - CutTube was replaced by intersection of TGeoTubSeg and two TGeoBBoxes",
1109  lname.Data());
1110  return mainN;
1111  }
1112  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
1113  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
1114  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1115  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
1116  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetPhi2() - geoShape->GetPhi1()));
1117  fGdmlE->NewAttr(mainN, 0, "lowX", TString::Format("%.12g", geoShape->GetNlow()[0]));
1118  fGdmlE->NewAttr(mainN, 0, "lowY", TString::Format("%.12g", geoShape->GetNlow()[1]));
1119  fGdmlE->NewAttr(mainN, 0, "lowZ", TString::Format("%.12g", geoShape->GetNlow()[2]));
1120  fGdmlE->NewAttr(mainN, 0, "highX", TString::Format("%.12g", geoShape->GetNhigh()[0]));
1121  fGdmlE->NewAttr(mainN, 0, "highY", TString::Format("%.12g", geoShape->GetNhigh()[1]));
1122  fGdmlE->NewAttr(mainN, 0, "highZ", TString::Format("%.12g", geoShape->GetNhigh()[2]));
1123 
1124  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1125  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1126 
1127  return mainN;
1128 }
1129 
1130 ////////////////////////////////////////////////////////////////////////////////
1131 /// Creates "tube" node for GDML from object TGeoTube
1132 
1134 {
1135  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "tube", 0);
1136  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1137  fGdmlE->NewAttr(mainN, 0, "name", lname);
1138  if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) ||
1139  IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1140  return NULL;
1141  }
1142 
1143  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
1144  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
1145  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1146  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%i", 0));
1147  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%i", 360));
1148 
1149  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1150  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1151  return mainN;
1152 }
1153 
1154 ////////////////////////////////////////////////////////////////////////////////
1155 /// Creates "zplane" node for GDML
1156 
1158 {
1159  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "zplane", 0);
1160 
1161  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", z));
1162  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", rmin));
1163  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", rmax));
1164 
1165  return mainN;
1166 }
1167 
1168 ////////////////////////////////////////////////////////////////////////////////
1169 /// Creates "polycone" node for GDML
1170 
1172 {
1173  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "polycone", 0);
1174  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1175  fGdmlE->NewAttr(mainN, 0, "name", lname);
1176 
1177  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
1178  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetDphi()));
1179 
1180  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1181  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1182  Int_t nZPlns = geoShape->GetNz();
1183  for (Int_t it = 0; it < nZPlns; it++) {
1184  //add zplane child node
1185  fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1186  //compare actual plane and next plane
1187  if ((it < nZPlns - 1) && (geoShape->GetZ(it) == geoShape->GetZ(it + 1))) {
1188  //rmin of actual is greater then rmax of next one
1189  // | |rmax next
1190  // __ ...| |... __ < rmin actual
1191  // | | | |
1192  if (geoShape->GetRmin(it) > geoShape->GetRmax(it + 1)) {
1193  //adding plane from rmax next to rmin actual at the same z position
1194  if (fgG4Compatibility == kTRUE) {
1195  fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it + 1), geoShape->GetRmin(it)));
1196  Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4", lname.Data());
1197  } else {
1198  Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4", lname.Data());
1199  }
1200 
1201  }
1202  //rmin of next is greater then rmax of actual
1203  // | | | |
1204  // | |...___...| | rmin next
1205  // | | > rmax act
1206  if (geoShape->GetRmin(it + 1) > geoShape->GetRmax(it)) {
1207  //adding plane from rmax act to rmin next at the same z position
1208  if (fgG4Compatibility == kTRUE) {
1209  fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it), geoShape->GetRmin(it + 1)));
1210  Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4", lname.Data());
1211  } else {
1212  Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4", lname.Data());
1213  }
1214  }
1215  }
1216  }
1217  return mainN;
1218 }
1219 
1220 ////////////////////////////////////////////////////////////////////////////////
1221 /// Creates "torus" node for GDML
1222 
1224 {
1225  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "torus", 0);
1226  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1227  fGdmlE->NewAttr(mainN, 0, "name", lname);
1228  if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
1229  return NULL;
1230  }
1231 
1232  fGdmlE->NewAttr(mainN, 0, "rtor", TString::Format("%.12g", geoShape->GetR()));
1233  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
1234  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
1235  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
1236  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetDphi()));
1237 
1238  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1239  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1240 
1241  return mainN;
1242 }
1243 
1244 ////////////////////////////////////////////////////////////////////////////////
1245 /// Creates "polyhedra" node for GDML
1246 
1248 {
1249  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "polyhedra", 0);
1250  fGdmlE->NewAttr(mainN, 0, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1251 
1252  fGdmlE->NewAttr(mainN, 0, "startphi", TString::Format("%.12g", geoShape->GetPhi1()));
1253  fGdmlE->NewAttr(mainN, 0, "deltaphi", TString::Format("%.12g", geoShape->GetDphi()));
1254  fGdmlE->NewAttr(mainN, 0, "numsides", TString::Format("%i", geoShape->GetNedges()));
1255 
1256  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1257  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1258  for (Int_t it = 0; it < geoShape->GetNz(); it++) {
1259  //add zplane child node
1260  fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1261  }
1262  return mainN;
1263 }
1264 
1265 ////////////////////////////////////////////////////////////////////////////////
1266 /// Creates "eltube" node for GDML
1267 
1269 {
1270  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "eltube", 0);
1271  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1272  fGdmlE->NewAttr(mainN, 0, "name", lname);
1273  if (IsNullParam(geoShape->GetA(), "A", lname) ||
1274  IsNullParam(geoShape->GetB(), "B", lname) ||
1275  IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1276  return NULL;
1277  }
1278 
1279  fGdmlE->NewAttr(mainN, 0, "dx", TString::Format("%.12g", geoShape->GetA()));
1280  fGdmlE->NewAttr(mainN, 0, "dy", TString::Format("%.12g", geoShape->GetB()));
1281  fGdmlE->NewAttr(mainN, 0, "dz", TString::Format("%.12g", geoShape->GetDz()));
1282 
1283  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1284 
1285  return mainN;
1286 }
1287 
1288 ////////////////////////////////////////////////////////////////////////////////
1289 /// Creates "hype" node for GDML
1290 
1292 {
1293  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "hype", 0);
1294  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1295  fGdmlE->NewAttr(mainN, 0, "name", lname);
1296  if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1297  return NULL;
1298  }
1299 
1300 
1301  fGdmlE->NewAttr(mainN, 0, "rmin", TString::Format("%.12g", geoShape->GetRmin()));
1302  fGdmlE->NewAttr(mainN, 0, "rmax", TString::Format("%.12g", geoShape->GetRmax()));
1303  fGdmlE->NewAttr(mainN, 0, "inst", TString::Format("%.12g", geoShape->GetStIn()));
1304  fGdmlE->NewAttr(mainN, 0, "outst", TString::Format("%.12g", geoShape->GetStOut()));
1305  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", 2 * geoShape->GetDz()));
1306 
1307  fGdmlE->NewAttr(mainN, 0, "aunit", "deg");
1308  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1309 
1310  return mainN;
1311 }
1312 
1313 ////////////////////////////////////////////////////////////////////////////////
1314 /// Creates "xtru" node for GDML
1315 
1317 {
1318  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "xtru", 0);
1319  TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1320  fGdmlE->NewAttr(mainN, 0, "name", lname);
1321 
1322  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1323  XMLNodePointer_t childN;
1324  Int_t vertNum = geoShape->GetNvert();
1325  Int_t secNum = geoShape->GetNz();
1326  if (vertNum < 3 || secNum < 2) {
1327  Info("CreateXtrusionN", "ERROR! TGeoXtru %s has only %i vertices and %i sections. It was not exported",
1328  lname.Data(), vertNum, secNum);
1329  mainN = NULL;
1330  return mainN;
1331  }
1332  for (Int_t it = 0; it < vertNum; it++) {
1333  //add twoDimVertex child node
1334  childN = fGdmlE->NewChild(0, 0, "twoDimVertex", 0);
1335  fGdmlE->NewAttr(childN, 0, "x", TString::Format("%.12g", geoShape->GetX(it)));
1336  fGdmlE->NewAttr(childN, 0, "y", TString::Format("%.12g", geoShape->GetY(it)));
1337  fGdmlE->AddChild(mainN, childN);
1338  }
1339  for (Int_t it = 0; it < secNum; it++) {
1340  //add section child node
1341  childN = fGdmlE->NewChild(0, 0, "section", 0);
1342  fGdmlE->NewAttr(childN, 0, "zOrder", TString::Format("%i", it));
1343  fGdmlE->NewAttr(childN, 0, "zPosition", TString::Format("%.12g", geoShape->GetZ(it)));
1344  fGdmlE->NewAttr(childN, 0, "xOffset", TString::Format("%.12g", geoShape->GetXOffset(it)));
1345  fGdmlE->NewAttr(childN, 0, "yOffset", TString::Format("%.12g", geoShape->GetYOffset(it)));
1346  fGdmlE->NewAttr(childN, 0, "scalingFactor", TString::Format("%.12g", geoShape->GetScale(it)));
1347  fGdmlE->AddChild(mainN, childN);
1348  }
1349  return mainN;
1350 }
1351 
1352 ////////////////////////////////////////////////////////////////////////////////
1353 /// Creates "ellipsoid" node for GDML
1354 /// this is a special case, because ellipsoid is not defined in ROOT
1355 /// so when intersection of scaled sphere and TGeoBBox is found,
1356 /// it is considered as an ellipsoid
1357 
1359 {
1360  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "ellipsoid", 0);
1361  TGeoScaledShape *leftS = (TGeoScaledShape *)geoShape->GetBoolNode()->GetLeftShape(); //ScaledShape
1362  TGeoBBox *rightS = (TGeoBBox *)geoShape->GetBoolNode()->GetRightShape(); //BBox
1363 
1364 
1365  fGdmlE->NewAttr(mainN, 0, "name", elName.Data());
1366  Double_t sx = leftS->GetScale()->GetScale()[0];
1367  Double_t sy = leftS->GetScale()->GetScale()[1];
1368  Double_t radius = ((TGeoSphere *) leftS->GetShape())->GetRmax();
1369 
1370  Double_t ax, by, cz;
1371  cz = radius;
1372  ax = sx * radius;
1373  by = sy * radius;
1374 
1375  Double_t dz = rightS->GetDZ();
1376  Double_t zorig = rightS->GetOrigin()[2];
1377  Double_t zcut2 = dz + zorig;
1378  Double_t zcut1 = 2 * zorig - zcut2;
1379 
1380 
1381  fGdmlE->NewAttr(mainN, 0, "ax", TString::Format("%.12g", ax));
1382  fGdmlE->NewAttr(mainN, 0, "by", TString::Format("%.12g", by));
1383  fGdmlE->NewAttr(mainN, 0, "cz", TString::Format("%.12g", cz));
1384  fGdmlE->NewAttr(mainN, 0, "zcut1", TString::Format("%.12g", zcut1));
1385  fGdmlE->NewAttr(mainN, 0, "zcut2", TString::Format("%.12g", zcut2));
1386  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1387 
1388  return mainN;
1389 }
1390 
1391 ////////////////////////////////////////////////////////////////////////////////
1392 /// Creates "elcone" (elliptical cone) node for GDML
1393 /// this is a special case, because elliptical cone is not defined in ROOT
1394 /// so when scaled cone is found, it is considered as a elliptical cone
1395 
1397 {
1398  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "elcone", 0);
1399  fGdmlE->NewAttr(mainN, 0, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1400  Double_t zcut = ((TGeoCone *) geoShape->GetShape())->GetDz();
1401  Double_t rx1 = ((TGeoCone *) geoShape->GetShape())->GetRmax1();
1402  Double_t rx2 = ((TGeoCone *) geoShape->GetShape())->GetRmax2();
1403  Double_t zmax = zcut * ((rx1 + rx2) / (rx1 - rx2));
1404  Double_t z = zcut + zmax;
1405 
1406  Double_t sy = geoShape->GetScale()->GetScale()[1];
1407  Double_t ry1 = sy * rx1;
1408 
1409  fGdmlE->NewAttr(mainN, 0, "dx", TString::Format("%.12g/%.12g", rx1, z));
1410  fGdmlE->NewAttr(mainN, 0, "dy", TString::Format("%.12g/%.12g", ry1, z));
1411  fGdmlE->NewAttr(mainN, 0, "zmax", TString::Format("%.12g", zmax));
1412  fGdmlE->NewAttr(mainN, 0, "zcut", TString::Format("%.12g", zcut));
1413  fGdmlE->NewAttr(mainN, 0, "lunit", "cm");
1414 
1415  return mainN;
1416 }
1417 
1418 ////////////////////////////////////////////////////////////////////////////////
1419 /// Creates common part of union intersection and subtraction nodes
1420 
1422 {
1423  XMLNodePointer_t mainN, ndR, ndL, childN;
1424  TString nodeName = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1425  TString lboolType;
1426  TGeoBoolNode::EGeoBoolType boolType = geoShape->GetBoolNode()->GetBooleanOperator();
1427  switch (boolType) {
1429  lboolType = "union";
1430  break;
1432  lboolType = "subtraction";
1433  break;
1435  lboolType = "intersection";
1436  break;
1437  }
1438 
1440  const Double_t *ltr = geoShape->GetBoolNode()->GetLeftMatrix()->GetTranslation();
1442  const Double_t *rtr = geoShape->GetBoolNode()->GetRightMatrix()->GetTranslation();
1443 
1444  //specific case!
1445  //Ellipsoid tag preparing
1446  //if left == TGeoScaledShape AND right == TGeoBBox
1447  // AND if TGeoScaledShape->GetShape == TGeoSphere
1448  TGeoShape *leftS = geoShape->GetBoolNode()->GetLeftShape();
1449  TGeoShape *rightS = geoShape->GetBoolNode()->GetRightShape();
1450  if (strcmp(leftS->ClassName(), "TGeoScaledShape") == 0 &&
1451  strcmp(rightS->ClassName(), "TGeoBBox") == 0) {
1452  if (strcmp(((TGeoScaledShape *)leftS)->GetShape()->ClassName(), "TGeoSphere") == 0) {
1453  if (lboolType == "intersection") {
1454  mainN = CreateEllipsoidN(geoShape, nodeName);
1455  return mainN;
1456  }
1457  }
1458  }
1459 
1460  Xyz translL, translR;
1461  //translation
1462  translL.x = ltr[0];
1463  translL.y = ltr[1];
1464  translL.z = ltr[2];
1465  translR.x = rtr[0];
1466  translR.y = rtr[1];
1467  translR.z = rtr[2];
1468 
1469  //left and right nodes are created here also their names are created
1470  ndL = ChooseObject(geoShape->GetBoolNode()->GetLeftShape());
1471  ndR = ChooseObject(geoShape->GetBoolNode()->GetRightShape());
1472 
1473  //retrieve left and right node names by their pointer to make reference
1474  TString lname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetLeftShape())];
1475  TString rname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetRightShape())];
1476 
1477  //left and right nodes appended to main structure of nodes (if they are not already there)
1478  if (ndL != NULL) {
1479  fGdmlE->AddChild(fSolidsNode, ndL);
1480  fSolCnt++;
1481  } else {
1482  if (lname.Contains("missing_") || lname == "") {
1483  Info("CreateCommonBoolN", "ERROR! Left node is NULL - Boolean Shape will be skipped");
1484  return NULL;
1485  }
1486  }
1487  if (ndR != NULL) {
1488  fGdmlE->AddChild(fSolidsNode, ndR);
1489  fSolCnt++;
1490  } else {
1491  if (rname.Contains("missing_") || rname == "") {
1492  Info("CreateCommonBoolN", "ERROR! Right node is NULL - Boolean Shape will be skipped");
1493  return NULL;
1494  }
1495  }
1496 
1497  //create union node and its child nodes (or intersection or subtraction)
1498  /* <union name="...">
1499  * <first ref="left name" />
1500  * <second ref="right name" />
1501  * <firstposition .../>
1502  * <firstrotation .../>
1503  * <position .../>
1504  * <rotation .../>
1505  * </union>
1506  */
1507  mainN = fGdmlE->NewChild(0, 0, lboolType.Data(), 0);
1508  fGdmlE->NewAttr(mainN, 0, "name", nodeName);
1509 
1510  //<first> (left)
1511  childN = fGdmlE->NewChild(0, 0, "first", 0);
1512  fGdmlE->NewAttr(childN, 0, "ref", lname);
1513  fGdmlE->AddChild(mainN, childN);
1514 
1515  //<second> (right)
1516  childN = fGdmlE->NewChild(0, 0, "second", 0);
1517  fGdmlE->NewAttr(childN, 0, "ref", rname);
1518  fGdmlE->AddChild(mainN, childN);
1519 
1520  //<firstposition> (left)
1521  if ((translL.x != 0.0) || (translL.y != 0.0) || (translL.z != 0.0)) {
1522  childN = CreatePositionN((nodeName + lname + "pos").Data(), translL, "firstposition");
1523  fGdmlE->AddChild(mainN, childN);
1524  }
1525  //<firstrotation> (left)
1526  if ((lrot.x != 0.0) || (lrot.y != 0.0) || (lrot.z != 0.0)) {
1527  childN = CreateRotationN((nodeName + lname + "rot").Data(), lrot, "firstrotation");
1528  fGdmlE->AddChild(mainN, childN);
1529  }
1530  //<position> (right)
1531  if ((translR.x != 0.0) || (translR.y != 0.0) || (translR.z != 0.0)) {
1532  childN = CreatePositionN((nodeName + rname + "pos").Data(), translR, "position");
1533  fGdmlE->AddChild(mainN, childN);
1534  }
1535  //<rotation> (right)
1536  if ((rrot.x != 0.0) || (rrot.y != 0.0) || (rrot.z != 0.0)) {
1537  childN = CreateRotationN((nodeName + rname + "rot").Data(), rrot, "rotation");
1538  fGdmlE->AddChild(mainN, childN);
1539  }
1540 
1541  return mainN;
1542 }
1543 
1544 
1545 ////////////////////////////////////////////////////////////////////////////////
1546 /// Creates "position" kind of node for GDML
1547 
1548 XMLNodePointer_t TGDMLWrite::CreatePositionN(const char * name, Xyz position, const char * type, const char * unit)
1549 {
1550  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, type, 0);
1551  fGdmlE->NewAttr(mainN, 0, "name", name);
1552  fGdmlE->NewAttr(mainN, 0, "x", TString::Format("%.12g", position.x));
1553  fGdmlE->NewAttr(mainN, 0, "y", TString::Format("%.12g", position.y));
1554  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", position.z));
1555  fGdmlE->NewAttr(mainN, 0, "unit", unit);
1556  return mainN;
1557 }
1558 
1559 ////////////////////////////////////////////////////////////////////////////////
1560 /// Creates "rotation" kind of node for GDML
1561 
1562 XMLNodePointer_t TGDMLWrite::CreateRotationN(const char * name, Xyz rotation, const char * type, const char * unit)
1563 {
1564  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, type, 0);
1565  fGdmlE->NewAttr(mainN, 0, "name", name);
1566  fGdmlE->NewAttr(mainN, 0, "x", TString::Format("%.12g", rotation.x));
1567  fGdmlE->NewAttr(mainN, 0, "y", TString::Format("%.12g", rotation.y));
1568  fGdmlE->NewAttr(mainN, 0, "z", TString::Format("%.12g", rotation.z));
1569  fGdmlE->NewAttr(mainN, 0, "unit", unit);
1570  return mainN;
1571 }
1572 
1573 ////////////////////////////////////////////////////////////////////////////////
1574 /// Creates "setup" node for GDML
1575 
1576 XMLNodePointer_t TGDMLWrite::CreateSetupN(const char * topVolName, const char * name, const char * version)
1577 {
1578  XMLNodePointer_t setupN = fGdmlE->NewChild(0, 0, "setup", 0);
1579  fGdmlE->NewAttr(setupN, 0, "name", name);
1580  fGdmlE->NewAttr(setupN, 0, "version", version);
1581  XMLNodePointer_t fworldN = fGdmlE->NewChild(setupN, 0, "world", 0);
1582  fGdmlE->NewAttr(fworldN, 0, "ref", topVolName);
1583  return setupN;
1584 }
1585 
1586 ////////////////////////////////////////////////////////////////////////////////
1587 /// Creates "volume" node for GDML
1588 
1589 XMLNodePointer_t TGDMLWrite::StartVolumeN(const char * name, const char * solid, const char * material)
1590 {
1591  XMLNodePointer_t childN;
1592  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "volume", 0);
1593  fGdmlE->NewAttr(mainN, 0, "name", name);
1594 
1595  childN = fGdmlE->NewChild(0, 0, "materialref", 0);
1596  fGdmlE->NewAttr(childN, 0, "ref", material);
1597  fGdmlE->AddChild(mainN, childN);
1598 
1599  childN = fGdmlE->NewChild(0, 0, "solidref", 0);
1600  fGdmlE->NewAttr(childN, 0, "ref", solid);
1601  fGdmlE->AddChild(mainN, childN);
1602 
1603  return mainN;
1604 }
1605 
1606 ////////////////////////////////////////////////////////////////////////////////
1607 /// Creates "assembly" node for GDML
1608 
1610 {
1611  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "assembly", 0);
1612  fGdmlE->NewAttr(mainN, 0, "name", name);
1613 
1614  return mainN;
1615 }
1616 
1617 ////////////////////////////////////////////////////////////////////////////////
1618 /// Creates "physvol" node for GDML
1619 
1620 XMLNodePointer_t TGDMLWrite::CreatePhysVolN(const char *name, Int_t copyno, const char * volref, const char * posref, const char * rotref, XMLNodePointer_t scaleN)
1621 {
1622  fPhysVolCnt++;
1623  XMLNodePointer_t childN;
1624  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "physvol", 0);
1625  fGdmlE->NewAttr(mainN, 0, "name", name);
1626  fGdmlE->NewAttr(mainN, 0, "copynumber", TString::Format("%d",copyno));
1627 
1628  childN = fGdmlE->NewChild(0, 0, "volumeref", 0);
1629  fGdmlE->NewAttr(childN, 0, "ref", volref);
1630  fGdmlE->AddChild(mainN, childN);
1631 
1632  childN = fGdmlE->NewChild(0, 0, "positionref", 0);
1633  fGdmlE->NewAttr(childN, 0, "ref", posref);
1634  fGdmlE->AddChild(mainN, childN);
1635 
1636  //if is not empty string add this node
1637  if (strcmp(rotref, "") != 0) {
1638  childN = fGdmlE->NewChild(0, 0, "rotationref", 0);
1639  fGdmlE->NewAttr(childN, 0, "ref", rotref);
1640  fGdmlE->AddChild(mainN, childN);
1641  }
1642  if (scaleN != NULL) {
1643  fGdmlE->AddChild(mainN, scaleN);
1644  }
1645 
1646  return mainN;
1647 }
1648 
1649 ////////////////////////////////////////////////////////////////////////////////
1650 /// Creates "divisionvol" node for GDML
1651 
1652 XMLNodePointer_t TGDMLWrite::CreateDivisionN(Double_t offset, Double_t width, Int_t number, const char * axis, const char * unit, const char * volref)
1653 {
1654  XMLNodePointer_t childN = 0;
1655  XMLNodePointer_t mainN = fGdmlE->NewChild(0, 0, "divisionvol", 0);
1656  fGdmlE->NewAttr(mainN, 0, "axis", axis);
1657  fGdmlE->NewAttr(mainN, 0, "number", TString::Format("%i", number));
1658  if (fgG4Compatibility == kTRUE) {
1659  //if eg. full length is 20 and width * number = 20,0001 problem in geant4
1660  //unit is either in cm or degrees nothing else
1661  width = (floor(width * 1E4)) * 1E-4;
1662  if ((offset >= 0.) && (strcmp(axis, "kPhi") == 0)) {
1663  Int_t offsetI = (Int_t) offset;
1664  Double_t decimals = offset - offsetI;
1665  //put to range from 0 to 360 add decimals and then put to range 0 -> -360
1666  offset = (offsetI % 360) + decimals - 360;
1667  }
1668  }
1669  fGdmlE->NewAttr(mainN, 0, "width", TString::Format("%.12g", width));
1670 
1671  fGdmlE->NewAttr(mainN, 0, "offset", TString::Format("%.12g", offset));
1672  fGdmlE->NewAttr(mainN, 0, "unit", unit);
1673  if (strcmp(volref, "") != 0) {
1674  childN = fGdmlE->NewChild(0, 0, "volumeref", 0);
1675  fGdmlE->NewAttr(childN, 0, "ref", volref);
1676  }
1677  fGdmlE->AddChild(mainN, childN);
1678 
1679 
1680  return mainN;
1681 }
1682 
1683 ////////////////////////////////////////////////////////////////////////////////
1684 /// Chooses the object and method that should be used for processing object
1685 
1687 {
1688  const char * clsname = geoShape->ClassName();
1689  XMLNodePointer_t solidN;
1690 
1691  if (CanProcess((TObject *)geoShape) == kFALSE) {
1692  return NULL;
1693  }
1694 
1695  //process different shapes
1696  if (strcmp(clsname, "TGeoBBox") == 0) {
1697  solidN = CreateBoxN((TGeoBBox*) geoShape);
1698  } else if (strcmp(clsname, "TGeoParaboloid") == 0) {
1699  solidN = CreateParaboloidN((TGeoParaboloid*) geoShape);
1700  } else if (strcmp(clsname, "TGeoSphere") == 0) {
1701  solidN = CreateSphereN((TGeoSphere*) geoShape);
1702  } else if (strcmp(clsname, "TGeoArb8") == 0) {
1703  solidN = CreateArb8N((TGeoArb8*) geoShape);
1704  } else if (strcmp(clsname, "TGeoConeSeg") == 0) {
1705  solidN = CreateConeN((TGeoConeSeg*) geoShape);
1706  } else if (strcmp(clsname, "TGeoCone") == 0) {
1707  solidN = CreateConeN((TGeoCone*) geoShape);
1708  } else if (strcmp(clsname, "TGeoPara") == 0) {
1709  solidN = CreateParaN((TGeoPara*) geoShape);
1710  } else if (strcmp(clsname, "TGeoTrap") == 0) {
1711  solidN = CreateTrapN((TGeoTrap*) geoShape);
1712  } else if (strcmp(clsname, "TGeoGtra") == 0) {
1713  solidN = CreateTwistedTrapN((TGeoGtra*) geoShape);
1714  } else if (strcmp(clsname, "TGeoTrd1") == 0) {
1715  solidN = CreateTrdN((TGeoTrd1*) geoShape);
1716  } else if (strcmp(clsname, "TGeoTrd2") == 0) {
1717  solidN = CreateTrdN((TGeoTrd2*) geoShape);
1718  } else if (strcmp(clsname, "TGeoTubeSeg") == 0) {
1719  solidN = CreateTubeN((TGeoTubeSeg*) geoShape);
1720  } else if (strcmp(clsname, "TGeoCtub") == 0) {
1721  solidN = CreateCutTubeN((TGeoCtub*) geoShape);
1722  } else if (strcmp(clsname, "TGeoTube") == 0) {
1723  solidN = CreateTubeN((TGeoTube*) geoShape);
1724  } else if (strcmp(clsname, "TGeoPcon") == 0) {
1725  solidN = CreatePolyconeN((TGeoPcon*) geoShape);
1726  } else if (strcmp(clsname, "TGeoTorus") == 0) {
1727  solidN = CreateTorusN((TGeoTorus*) geoShape);
1728  } else if (strcmp(clsname, "TGeoPgon") == 0) {
1729  solidN = CreatePolyhedraN((TGeoPgon*) geoShape);
1730  } else if (strcmp(clsname, "TGeoEltu") == 0) {
1731  solidN = CreateEltubeN((TGeoEltu*) geoShape);
1732  } else if (strcmp(clsname, "TGeoHype") == 0) {
1733  solidN = CreateHypeN((TGeoHype*) geoShape);
1734  } else if (strcmp(clsname, "TGeoXtru") == 0) {
1735  solidN = CreateXtrusionN((TGeoXtru*) geoShape);
1736  } else if (strcmp(clsname, "TGeoScaledShape") == 0) {
1737  TGeoScaledShape * geoscale = (TGeoScaledShape *) geoShape;
1738  TString scaleObjClsName = geoscale->GetShape()->ClassName();
1739  if (scaleObjClsName == "TGeoCone") {
1740  solidN = CreateElConeN((TGeoScaledShape*) geoShape);
1741  } else {
1742  Info("ChooseObject",
1743  "ERROR! TGeoScaledShape object is not possible to process correctly. %s object is processed without scale",
1744  scaleObjClsName.Data());
1745  solidN = ChooseObject(geoscale->GetShape());
1746  //Name has to be propagated to geoscale level pointer
1747  fNameList->fLst[TString::Format("%p", geoscale)] =
1748  fNameList->fLst[TString::Format("%p", geoscale->GetShape())];
1749  }
1750  } else if (strcmp(clsname, "TGeoCompositeShape") == 0) {
1751  solidN = CreateCommonBoolN((TGeoCompositeShape*) geoShape);
1752  } else if (strcmp(clsname, "TGeoUnion") == 0) {
1753  solidN = CreateCommonBoolN((TGeoCompositeShape*) geoShape);
1754  } else if (strcmp(clsname, "TGeoIntersection") == 0) {
1755  solidN = CreateCommonBoolN((TGeoCompositeShape*) geoShape);
1756  } else if (strcmp(clsname, "TGeoSubtraction") == 0) {
1757  solidN = CreateCommonBoolN((TGeoCompositeShape*) geoShape);
1758  } else {
1759  Info("ChooseObject", "ERROR! %s Solid CANNOT be processed, solid is NOT supported",
1760  clsname);
1761  solidN = NULL;
1762  }
1763  if (solidN == NULL) {
1764  if (fNameList->fLst[TString::Format("%p", geoShape)] == "") {
1765  TString missingName = geoShape->GetName();
1766  GenName("missing_" + missingName, TString::Format("%p", geoShape));
1767  } else {
1768  fNameList->fLst[TString::Format("%p", geoShape)] = "missing_" + fNameList->fLst[TString::Format("%p", geoShape)];
1769  }
1770  }
1771 
1772  return solidN;
1773 }
1774 
1775 ////////////////////////////////////////////////////////////////////////////////
1776 /// Retrieves X Y Z angles from rotation matrix
1777 
1779 {
1780  TGDMLWrite::Xyz lxyz;
1781  Double_t a, b, c;
1782  Double_t rad = 180.0 / TMath::ACos(-1.0);
1783  const Double_t *r = rotationMatrix;
1784  Double_t cosb = TMath::Sqrt(r[0] * r[0] + r[1] * r[1]);
1785  if (cosb > 0.00001) {
1786  a = TMath::ATan2(r[5], r[8]) * rad;
1787  b = TMath::ATan2(-r[2], cosb) * rad;
1788  c = TMath::ATan2(r[1], r[0]) * rad;
1789  } else {
1790  a = TMath::ATan2(-r[7], r[4]) * rad;
1791  b = TMath::ATan2(-r[2], cosb) * rad;
1792  c = 0;
1793  }
1794  lxyz.x = a;
1795  lxyz.y = b;
1796  lxyz.z = c;
1797  return lxyz;
1798 }
1799 
1800 ////////////////////////////////////////////////////////////////////////////////
1801 /// Method creating cutTube as an intersection of tube and two boxes
1802 /// - not used anymore because cutube is supported in Geant4 9.5
1803 
1805 {
1806  Double_t rmin = geoShape->GetRmin();
1807  Double_t rmax = geoShape->GetRmax();
1808  Double_t z = geoShape->GetDz();
1809  Double_t startphi = geoShape->GetPhi1();
1810  Double_t deltaphi = geoShape->GetPhi2();
1811  Double_t x1 = geoShape->GetNlow()[0];
1812  Double_t y1 = geoShape->GetNlow()[1];
1813  Double_t z1 = geoShape->GetNlow()[2];
1814  Double_t x2 = geoShape->GetNhigh()[0];
1815  Double_t y2 = geoShape->GetNhigh()[1];
1816  Double_t z2 = geoShape->GetNhigh()[2];
1817  TString xname = geoShape->GetName();
1818 
1819 
1820  Double_t h0 = 2.*((TGeoBBox*)geoShape)->GetDZ();
1821  Double_t h1 = 2 * z;
1822  Double_t h2 = 2 * z;
1823  Double_t boxdx = 1E8 * (2 * rmax) + (2 * z);
1824 
1825  TGeoTubeSeg *T = new TGeoTubeSeg((xname + "T").Data(), rmin, rmax, h0, startphi, deltaphi);
1826  TGeoBBox *B1 = new TGeoBBox((xname + "B1").Data(), boxdx, boxdx, h1);
1827  TGeoBBox *B2 = new TGeoBBox((xname + "B2").Data(), boxdx, boxdx, h2);
1828 
1829 
1830  //first box position parameters
1831  Double_t phi1 = 360 - TMath::ATan2(x1, y1) * TMath::RadToDeg();
1832  Double_t theta1 = 360 - TMath::ATan2(sqrt(x1 * x1 + y1 * y1), z1) * TMath::RadToDeg();
1833 
1834  Double_t phi11 = TMath::ATan2(y1, x1) * TMath::RadToDeg() ;
1835  Double_t theta11 = TMath::ATan2(z1, sqrt(x1 * x1 + y1 * y1)) * TMath::RadToDeg() ;
1836 
1837  Double_t xpos1 = h1 * TMath::Cos((theta11) * TMath::DegToRad()) * TMath::Cos((phi11) * TMath::DegToRad()) * (-1);
1838  Double_t ypos1 = h1 * TMath::Cos((theta11) * TMath::DegToRad()) * TMath::Sin((phi11) * TMath::DegToRad()) * (-1);
1839  Double_t zpos1 = h1 * TMath::Sin((theta11) * TMath::DegToRad()) * (-1);
1840 
1841  //second box position parameters
1842  Double_t phi2 = 360 - TMath::ATan2(x2, y2) * TMath::RadToDeg();
1843  Double_t theta2 = 360 - TMath::ATan2(sqrt(x2 * x2 + y2 * y2), z2) * TMath::RadToDeg();
1844 
1845  Double_t phi21 = TMath::ATan2(y2, x2) * TMath::RadToDeg() ;
1846  Double_t theta21 = TMath::ATan2(z2, sqrt(x2 * x2 + y2 * y2)) * TMath::RadToDeg() ;
1847 
1848  Double_t xpos2 = h2 * TMath::Cos((theta21) * TMath::DegToRad()) * TMath::Cos((phi21) * TMath::DegToRad()) * (-1);
1849  Double_t ypos2 = h2 * TMath::Cos((theta21) * TMath::DegToRad()) * TMath::Sin((phi21) * TMath::DegToRad()) * (-1);
1850  Double_t zpos2 = h2 * TMath::Sin((theta21) * TMath::DegToRad()) * (-1);
1851 
1852 
1853  //positioning
1854  TGeoTranslation *t0 = new TGeoTranslation(0, 0, 0);
1855  TGeoTranslation *t1 = new TGeoTranslation(0 + xpos1, 0 + ypos1, 0 + (zpos1 - z));
1856  TGeoTranslation *t2 = new TGeoTranslation(0 + xpos2, 0 + ypos2, 0 + (zpos2 + z));
1857  TGeoRotation *r0 = new TGeoRotation((xname + "r0").Data());
1858  TGeoRotation *r1 = new TGeoRotation((xname + "r1").Data());
1859  TGeoRotation *r2 = new TGeoRotation((xname + "r2").Data());
1860 
1861  r1->SetAngles(phi1, theta1, 0);
1862  r2->SetAngles(phi2, theta2, 0);
1863 
1864  TGeoMatrix* m0 = new TGeoCombiTrans(*t0, *r0);
1865  TGeoMatrix* m1 = new TGeoCombiTrans(*t1, *r1);
1866  TGeoMatrix* m2 = new TGeoCombiTrans(*t2, *r2);
1867 
1868  TGeoCompositeShape *CS1 = new TGeoCompositeShape((xname + "CS1").Data(), new TGeoIntersection(T, B1, m0, m1));
1869  TGeoCompositeShape *cs = new TGeoCompositeShape((xname + "CS").Data(), new TGeoIntersection(CS1, B2, m0, m2));
1870  delete t0;
1871  delete t1;
1872  delete t2;
1873  delete r0;
1874  delete r1;
1875  delete r2;
1876  return cs;
1877 }
1878 
1879 ////////////////////////////////////////////////////////////////////////////////
1880 /// Checks whether name2check is in (NameList) list
1881 
1883 {
1884  Bool_t isIN = list[name2check];
1885  return isIN;
1886 }
1887 
1888 ////////////////////////////////////////////////////////////////////////////////
1889 ///NCNAME basic restrictions
1890 ///Replace "$" character with empty character etc.
1891 
1893 {
1894  TString newname = oldname.ReplaceAll("$", "");
1895  newname = newname.ReplaceAll(" ", "_");
1896  // :, @, $, %, &, /, +, ,, ;, whitespace characters or different parenthesis
1897  newname = newname.ReplaceAll(":", "");
1898  newname = newname.ReplaceAll("@", "");
1899  newname = newname.ReplaceAll("%", "");
1900  newname = newname.ReplaceAll("&", "");
1901  newname = newname.ReplaceAll("/", "");
1902  newname = newname.ReplaceAll("+", "");
1903  newname = newname.ReplaceAll(";", "");
1904  newname = newname.ReplaceAll("{", "");
1905  newname = newname.ReplaceAll("}", "");
1906  newname = newname.ReplaceAll("(", "");
1907  newname = newname.ReplaceAll(")", "");
1908  newname = newname.ReplaceAll("[", "");
1909  newname = newname.ReplaceAll("]", "");
1910  newname = newname.ReplaceAll("_refl", "");
1911  //workaround if first letter is digit than replace it to "O" (ou character)
1912  TString fstLet = newname(0, 1);
1913  if (fstLet.IsDigit()) {
1914  newname = "O" + newname(1, newname.Length());
1915  }
1916  return newname;
1917 }
1918 
1919 ////////////////////////////////////////////////////////////////////////////////
1920 /// Important function which is responsible for naming volumes, solids and materials
1921 
1923 {
1924  TString newname = GenName(oldname);
1925  if (newname != oldname) {
1926  if (fgkMaxNameErr > fActNameErr) {
1927  Info("GenName",
1928  "WARNING! Name of the object was changed because it failed to comply with NCNAME xml datatype restrictions.");
1929  } else if ((fgkMaxNameErr == fActNameErr)) {
1930  Info("GenName",
1931  "WARNING! Probably more names are going to be changed to comply with NCNAME xml datatype restriction, but it will not be displayed on the screen.");
1932  }
1933  fActNameErr++;
1934  }
1935  TString nameIter;
1936  Int_t iter = 0;
1937  switch (fgNamingSpeed) {
1938  case kfastButUglySufix:
1939  newname = newname + "0x" + objPointer;
1940  break;
1941  case kelegantButSlow:
1942  //0 means not in the list
1943  iter = fNameList->fLstIter[newname];
1944  if (iter == 0) {
1945  nameIter = "";
1946  } else {
1947  nameIter = TString::Format("0x%i", iter);
1948  }
1949  fNameList->fLstIter[newname]++;
1950  newname = newname + nameIter;
1951  break;
1952  case kwithoutSufixNotUniq:
1953  //no change
1954  newname = newname;
1955  break;
1956  }
1957  //store the name (mapped to pointer)
1958  fNameList->fLst[objPointer] = newname;
1959  return newname;
1960 }
1961 
1962 
1963 ////////////////////////////////////////////////////////////////////////////////
1964 /// Method which tests whether solids can be processed
1965 
1967 {
1968  Bool_t isProcessed = kFALSE;
1969  isProcessed = pointer->TestBit(fgkProcBit);
1970  pointer->SetBit(fgkProcBit, kTRUE);
1971  return !(isProcessed);
1972 }
1973 
1974 ////////////////////////////////////////////////////////////////////////////////
1975 /// Method that retrieves axis and unit along which object is divided
1976 
1977 TString TGDMLWrite::GetPattAxis(Int_t divAxis, const char * pattName, TString& unit)
1978 {
1979  TString resaxis;
1980  unit = "cm";
1981  switch (divAxis) {
1982  case 1:
1983  if (strcmp(pattName, "TGeoPatternX") == 0) {
1984  return "kXAxis";
1985  } else if (strcmp(pattName, "TGeoPatternCylR") == 0) {
1986  return "kRho";
1987  }
1988  break;
1989  case 2:
1990  if (strcmp(pattName, "TGeoPatternY") == 0) {
1991  return "kYAxis";
1992  } else if (strcmp(pattName, "TGeoPatternCylPhi") == 0) {
1993  unit = "deg";
1994  return "kPhi";
1995  }
1996  break;
1997  case 3:
1998  if (strcmp(pattName, "TGeoPatternZ") == 0) {
1999  return "kZAxis";
2000  }
2001  break;
2002  default:
2003  return "kUndefined";
2004  break;
2005  }
2006  return "kUndefined";
2007 }
2008 
2009 ////////////////////////////////////////////////////////////////////////////////
2010 /// Check for null parameter to skip the NULL objects
2011 
2013 {
2014  if (parValue == 0.) {
2015  Info("IsNullParam", "ERROR! %s is NULL due to %s = %.12g, Volume based on this shape will be skipped",
2016  objName.Data(),
2017  parName.Data(),
2018  parValue);
2019  return kTRUE;
2020  }
2021  return kFALSE;
2022 }
2023 
2024 ////////////////////////////////////////////////////////////////////////////////
2025 /// Unsetting bits that were changed in gGeoManager during export so that export
2026 /// can be run more times with the same instance of gGeoManager.
2027 
2029 {
2030  TIter next(geoMng->GetListOfVolumes());
2031  TGeoVolume *vol;
2032  while ((vol = (TGeoVolume *)next())) {
2033  ((TObject *)vol->GetShape())->SetBit(fgkProcBit, kFALSE);
2034  vol->SetAttBit(fgkProcBitVol, kFALSE);
2035  }
2036 
2037 }
Double_t GetPhi2() const
Definition: TGeoTube.h:149
Int_t GetZ() const
Definition: TGeoElement.h:120
Double_t GetDy2() const
Definition: TGeoTrd2.h:59
XMLNodePointer_t CreateMixtureN(TGeoMixture *mixture, XMLNodePointer_t materials, TString mname)
Creates "material" node for GDML with references to other sub elements.
Definition: TGDMLWrite.cxx:665
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Double_t GetTheta() const
Definition: TGeoArb8.h:124
Mixtures of elements.
Definition: TGeoMaterial.h:134
Spherical shell class.
Definition: TGeoSphere.h:17
Int_t GetNz() const
Definition: TGeoXtru.h:93
Cylindrical tube class.
Definition: TGeoTube.h:17
An array of TObjects.
Definition: TObjArray.h:37
Double_t GetDphi() const
Definition: TGeoTorus.h:73
Double_t GetAlpha2() const
Definition: TGeoArb8.h:133
Double_t * GetZ() const
Definition: TGeoPcon.h:79
XMLNodePointer_t CreateHypeN(TGeoHype *geoShape)
Creates "hype" node for GDML.
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:854
Bool_t TestAttBit(UInt_t f) const
Definition: TGeoAtt.h:68
Double_t GetXOffset(Int_t i) const
Definition: TGeoXtru.h:97
The manager class for any TGeo geometry.
Definition: TGeoManager.h:38
This class contains implementation of converting ROOT&#39;s gGeoManager geometry to GDML file...
Definition: TGDMLWrite.h:69
Box class.
Definition: TGeoBBox.h:17
Double_t GetDphi() const
Definition: TGeoPcon.h:72
TGeoShape * GetLeftShape() const
Definition: TGeoBoolNode.h:80
XMLNodePointer_t CreateFractionN(Double_t percentage, const char *refName)
Creates "fraction" node for GDML.
Definition: TGDMLWrite.cxx:577
Double_t GetPhi1() const
Definition: TGeoTube.h:148
virtual Double_t GetDX() const
Definition: TGeoBBox.h:70
XMLNodePointer_t CreateParaN(TGeoPara *geoShape)
Creates "para" node for GDML.
Definition: TGDMLWrite.cxx:897
XMLNodePointer_t fStructureNode
Definition: TGDMLWrite.h:136
static TGDMLWrite * fgGDMLWrite
Definition: TGDMLWrite.h:126
virtual Double_t GetDensity() const
Definition: TGeoMaterial.h:87
Gtra is a twisted trapezoid.
Definition: TGeoArb8.h:143
Double_t GetDy() const
Definition: TGeoTrd1.h:57
XMLNodePointer_t StartVolumeN(const char *name, const char *solid, const char *material)
Creates "volume" node for GDML.
Geometrical transformation package.
Definition: TGeoMatrix.h:40
StructLst * fRejShape
Definition: TGDMLWrite.h:121
virtual const Double_t * GetRotationMatrix() const
Definition: TGeoMatrix.h:468
double T(double x)
Definition: ChebyshevPol.h:34
Int_t Z() const
Definition: TGeoElement.h:73
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:638
A polycone.
Definition: TGeoPcon.h:17
virtual const Double_t * GetRotationMatrix() const =0
A polygone.
Definition: TGeoPgon.h:19
Double_t GetR() const
Definition: TGeoTorus.h:69
Double_t GetYOffset(Int_t i) const
Definition: TGeoXtru.h:98
XMLDocPointer_t NewDoc(const char *version="1.0")
creates new xml document with provided version
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition: TGeoVolume.h:48
Double_t GetPhi() const
Definition: TGeoArb8.h:125
Double_t GetPhi() const
Definition: TGeoPara.h:66
TString ExtractSolid(TGeoShape *volShape)
Method creating solid to xml file and returning its name.
Definition: TGDMLWrite.cxx:377
Int_t GetNedges() const
Definition: TGeoPgon.h:81
virtual Double_t GetB() const
Definition: TGeoEltu.h:44
virtual Double_t GetRmax() const
Definition: TGeoSphere.h:68
Double_t * GetWmixt() const
Definition: TGeoMaterial.h:175
Double_t GetRmin() const
Definition: TGeoTorus.h:70
virtual TGeoElement * GetElement(Int_t i=0) const
Retrieve the pointer to the element corresponding to component I.
const Double_t * GetNlow() const
Definition: TGeoTube.h:207
Class describing translations.
Definition: TGeoMatrix.h:121
Torus segment class.
Definition: TGeoTorus.h:17
void UnsetTemporaryBits(TGeoManager *geoMng)
Unsetting bits that were changed in gGeoManager during export so that export can be run more times wi...
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
XMLNodePointer_t CreateMaterialN(TGeoMaterial *material, TString mname)
Creates "material" node for GDML.
Definition: TGDMLWrite.cxx:718
Bool_t CanProcess(TObject *pointer)
Method which tests whether solids can be processed.
Basic string class.
Definition: TString.h:125
Base class describing materials.
Definition: TGeoMaterial.h:29
Double_t GetA() const
Definition: TGeoElement.h:122
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
XMLNodePointer_t CreatePositionN(const char *name, Xyz position, const char *type="position", const char *unit="cm")
Creates "position" kind of node for GDML.
Double_t GetPhi1() const
Definition: TGeoCone.h:160
void WriteGDMLfile(TGeoManager *geomanager, const char *filename="test.gdml", TString option="")
Wrapper of all exporting methods Creates blank GDML file and fills it with gGeoManager structure conv...
Definition: TGDMLWrite.cxx:223
void SetSkipComments(Bool_t on=kTRUE)
Definition: TXMLEngine.h:48
std::map< TString, Float_t > NameListF
Definition: TGDMLWrite.h:108
Int_t GetNdiv() const
NameListI fLstIter
Definition: TGDMLWrite.h:114
virtual Double_t GetRmin2() const
Definition: TGeoCone.h:75
Double_t GetDx2() const
Definition: TGeoTrd2.h:57
Int_t GetNvert() const
Definition: TGeoXtru.h:94
virtual EGeoBoolType GetBooleanOperator() const =0
XMLNodePointer_t CreateXtrusionN(TGeoXtru *geoShape)
Creates "xtru" node for GDML.
TGeoIsotope * GetIsotope(Int_t i) const
Return i-th isotope in the element.
A shape scaled by a TGeoScale transformation.
static constexpr double rad
XMLNodePointer_t CreateEllipsoidN(TGeoCompositeShape *geoShape, TString elName)
Creates "ellipsoid" node for GDML this is a special case, because ellipsoid is not defined in ROOT so...
virtual const Double_t * GetScale() const
Definition: TGeoMatrix.h:279
XMLNodePointer_t StartAssemblyN(const char *name)
Creates "assembly" node for GDML.
A trapezoid with only x length varying with z.
Definition: TGeoTrd1.h:17
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
void FreeDoc(XMLDocPointer_t xmldoc)
frees allocated document data and deletes document itself
Paraboloid class.
XMLNodePointer_t CreateTwistedTrapN(TGeoGtra *geoShape)
Creates "twistedtrap" node for GDML.
Definition: TGDMLWrite.cxx:962
Double_t GetPhi2() const
Definition: TGeoSphere.h:72
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
TList * GetListOfMaterials() const
Definition: TGeoManager.h:463
Double_t GetDx2() const
Definition: TGeoTrd1.h:56
XMLDocPointer_t fGdmlFile
Definition: TGDMLWrite.h:129
TObjArray * GetNodes()
Definition: TGeoVolume.h:170
Double_t GetPhi2() const
Definition: TGeoCone.h:161
Double_t GetPhi1() const
Definition: TGeoSphere.h:71
XMLNodePointer_t CreateElConeN(TGeoScaledShape *geoShape)
Creates "elcone" (elliptical cone) node for GDML this is a special case, because elliptical cone is n...
virtual Double_t GetRmax() const
Definition: TGeoTube.h:67
Double_t GetStOut() const
Definition: TGeoHype.h:69
Bool_t HasIsotopes() const
Definition: TGeoElement.h:85
XMLNodePointer_t CreateArb8N(TGeoArb8 *geoShape)
Creates "arb8" node for GDML.
Definition: TGDMLWrite.cxx:813
double sqrt(double)
static const double x2[5]
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
XMLNodePointer_t CreateZplaneN(Double_t z, Double_t rmin, Double_t rmax)
Creates "zplane" node for GDML.
TGeoPatternFinder * GetFinder() const
Definition: TGeoVolume.h:178
XMLNodePointer_t CreateConeN(TGeoConeSeg *geoShape)
Creates "cone" node for GDML from TGeoConeSeg object.
Definition: TGDMLWrite.cxx:847
Bool_t IsTwisted() const
Definition: TGeoArb8.h:73
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString...
Definition: TString.cxx:2365
Int_t GetNisotopes() const
Definition: TGeoElement.h:78
A phi segment of a conical tube.
Definition: TGeoCone.h:98
XMLNodePointer_t CreateTrapN(TGeoTrap *geoShape)
Creates "trap" node for GDML.
Definition: TGDMLWrite.cxx:917
Int_t GetN() const
Definition: TGeoElement.h:121
TGeoMatrix * GetLeftMatrix() const
Definition: TGeoBoolNode.h:78
Bool_t IsInList(NameList list, TString name2check)
Checks whether name2check is in (NameList) list.
XMLNodePointer_t CreateIsotopN(TGeoIsotope *isotope, const char *name)
Creates "isotope" node for GDML.
Definition: TGDMLWrite.cxx:588
void DocSetRootElement(XMLDocPointer_t xmldoc, XMLNodePointer_t xmlnode)
set main (root) node for document
Int_t fPhysVolCnt
Definition: TGDMLWrite.h:138
XMLNodePointer_t ExtractMaterials(TList *materialsLst)
Method exporting materials.
Definition: TGDMLWrite.cxx:345
TXMLEngine * fGdmlE
Definition: TGDMLWrite.h:131
UInt_t fActNameErr
Definition: TGDMLWrite.h:139
virtual ~TGDMLWrite()
Destructor.
Definition: TGDMLWrite.cxx:199
Int_t GetNz() const
Definition: TGeoPcon.h:73
Base class for chemical elements.
Definition: TGeoElement.h:36
constexpr Double_t DegToRad()
Definition: TMath.h:64
virtual Double_t GetA() const
Definition: TGeoMaterial.h:84
XMLNsPointer_t NewNS(XMLNodePointer_t xmlnode, const char *reference, const char *name=0)
create namespace attribute for xmlnode.
Definition: TXMLEngine.cxx:733
TGDMLWrite()
Default constructor.
Definition: TGDMLWrite.cxx:171
TRAP is a general trapezoid, i.e.
Definition: TGeoArb8.h:89
virtual TGeoMatrix * GetMatrix() const =0
TGeoMaterial * GetMaterial() const
Definition: TGeoVolume.h:175
void SetG4Compatibility(Bool_t G4Compatible)
Definition: TGDMLWrite.h:94
Double_t GetY(Int_t i) const
Definition: TGeoXtru.h:96
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:580
TH1F * h1
Definition: legend1.C:5
Double_t GetH2() const
Definition: TGeoArb8.h:130
XMLNodePointer_t CreateEltubeN(TGeoEltu *geoShape)
Creates "eltube" node for GDML.
Double_t GetX(Int_t i) const
Definition: TGeoXtru.h:95
virtual Double_t GetZ() const
Definition: TGeoMaterial.h:85
Bool_t IsNullParam(Double_t parValue, TString parName, TString objName)
Check for null parameter to skip the NULL objects.
A doubly linked list.
Definition: TList.h:44
XMLNodePointer_t fDefineNode
Definition: TGDMLWrite.h:133
Base finder class for patterns.
XMLNodePointer_t fSolidsNode
Definition: TGDMLWrite.h:135
virtual Double_t GetDz() const
Definition: TGeoCone.h:68
Class handling Boolean composition of shapes.
Double_t GetStep() const
virtual const char * GetName() const
Get the shape name.
Definition: TGeoShape.cxx:248
A trapezoid with both x and y lengths varying with z.
Definition: TGeoTrd2.h:17
Double_t * GetRmax() const
Definition: TGeoPcon.h:77
Int_t GetNumber() const
Definition: TGeoNode.h:92
Parallelepiped class.
Definition: TGeoPara.h:17
void SaveDoc(XMLDocPointer_t xmldoc, const char *filename, Int_t layout=1)
store document content to file if layout<=0, no any spaces or newlines will be placed between xmlnode...
Double_t GetX() const
Definition: TGeoPara.h:61
Double_t * GetRmin() const
Definition: TGeoPcon.h:75
virtual Double_t GetRmin1() const
Definition: TGeoCone.h:73
Xyz GetXYZangles(const Double_t *rotationMatrix)
Retrieves X Y Z angles from rotation matrix.
XMLNodePointer_t CreatePhysVolN(const char *name, Int_t copyno, const char *volref, const char *posref, const char *rotref, XMLNodePointer_t scaleN)
Creates "physvol" node for GDML.
Double_t GetTheta1() const
Definition: TGeoSphere.h:69
Base abstract class for all shapes.
Definition: TGeoShape.h:25
void AddChild(XMLNodePointer_t parent, XMLNodePointer_t child)
add child element to xmlnode
Definition: TXMLEngine.cxx:791
ROOT::R::TRInterface & r
Definition: Object.C:4
Class describing rotation + translation.
Definition: TGeoMatrix.h:291
void SetNamingSpeed(ENamingType naming)
Set convention of naming solids and volumes.
Definition: TGDMLWrite.cxx:213
XMLNodePointer_t CreateBoxN(TGeoBBox *geoShape)
Creates "box" node for GDML.
Definition: TGDMLWrite.cxx:748
auto * a
Definition: textangle.C:12
Double_t GetAlpha() const
Definition: TGeoPara.h:64
Double_t GetY() const
Definition: TGeoPara.h:62
std::map< TString, Bool_t > NameList
Definition: TGDMLWrite.h:105
XMLNodePointer_t CreatePolyconeN(TGeoPcon *geoShape)
Creates "polycone" node for GDML.
XMLNodePointer_t CreateParaboloidN(TGeoParaboloid *geoShape)
Creates "paraboloid" node for GDML.
Definition: TGDMLWrite.cxx:769
XMLNodePointer_t CreateTorusN(TGeoTorus *geoShape)
Creates "torus" node for GDML.
StructLst * fIsotopeList
Definition: TGDMLWrite.h:118
Double_t GetStart() const
virtual Double_t GetRmax1() const
Definition: TGeoCone.h:74
virtual Int_t GetNelements() const
Definition: TGeoMaterial.h:172
Double_t GetDz() const
Definition: TGeoTrd2.h:60
Double_t GetRhi() const
Hyperboloid class defined by 5 parameters.
Definition: TGeoHype.h:17
Double_t GetBl1() const
Definition: TGeoArb8.h:127
Ssiz_t Length() const
Definition: TString.h:386
virtual const Double_t * GetOrigin() const
Definition: TGeoBBox.h:73
double floor(double)
Double_t GetBl2() const
Definition: TGeoArb8.h:131
Double_t ACos(Double_t)
Definition: TMath.h:571
static constexpr double m2
TGeoCompositeShape * CreateFakeCtub(TGeoCtub *geoShape)
Method creating cutTube as an intersection of tube and two boxes.
Double_t GetScale(Int_t i) const
Definition: TGeoXtru.h:99
TString fTopVolumeName
Definition: TGDMLWrite.h:130
void * XMLNodePointer_t
Definition: TXMLEngine.h:17
virtual Double_t GetDY() const
Definition: TGeoBBox.h:71
Double_t GetDz() const
Bool_t fgG4Compatibility
Definition: TGDMLWrite.h:128
constexpr Double_t E()
Definition: TMath.h:74
Double_t Cos(Double_t)
Definition: TMath.h:550
Class describing rotations.
Definition: TGeoMatrix.h:174
const Bool_t kFALSE
Definition: RtypesCore.h:88
Double_t GetRlo() const
TString GenName(TString oldname)
NCNAME basic restrictions Replace "$" character with empty character etc.
Double_t GetPhi1() const
Definition: TGeoPcon.h:71
Bool_t IsTopVolume() const
True if this is the top volume of the geometry.
Definition: TGeoVolume.cxx:864
A tube segment cut with 2 planes.
Definition: TGeoTube.h:168
XMLNodePointer_t CreateDivisionN(Double_t offset, Double_t width, Int_t number, const char *axis, const char *unit, const char *volref)
Creates "divisionvol" node for GDML.
Double_t GetRmax() const
Definition: TGeoTorus.h:71
An extrusion with fixed outline shape in x-y and a sequence of z extents (segments).
Definition: TGeoXtru.h:21
XMLAttrPointer_t NewAttr(XMLNodePointer_t xmlnode, XMLNsPointer_t, const char *name, const char *value)
creates new attribute for xmlnode, namespaces are not supported for attributes
Definition: TXMLEngine.cxx:578
static const double x1[5]
auto * t1
Definition: textangle.C:20
#define ClassImp(name)
Definition: Rtypes.h:359
virtual Bool_t IsMixture() const
Definition: TGeoMaterial.h:108
double Double_t
Definition: RtypesCore.h:55
Int_t fVolCnt
Definition: TGDMLWrite.h:137
int type
Definition: TGX11.cxx:120
const Double_t * GetNhigh() const
Definition: TGeoTube.h:208
Conical tube class.
Definition: TGeoCone.h:17
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:570
Double_t GetDy1() const
Definition: TGeoTrd2.h:58
An arbitrary trapezoid with less than 8 vertices standing on two parallel planes perpendicular to Z a...
Definition: TGeoArb8.h:17
XMLNodePointer_t CreateTrdN(TGeoTrd1 *geoShape)
Creates "trd" node for GDML from object TGeoTrd1.
static const UInt_t fgkProcBitVol
Definition: TGDMLWrite.h:143
Double_t GetH1() const
Definition: TGeoArb8.h:126
Bool_t IsReflection() const
Definition: TGeoMatrix.h:69
XMLNodePointer_t CreatePolyhedraN(TGeoPgon *geoShape)
Creates "polyhedra" node for GDML.
virtual Double_t GetRmax2() const
Definition: TGeoCone.h:76
TGeoMatrix * GetRightMatrix() const
Definition: TGeoBoolNode.h:79
XMLNodePointer_t CreateAtomN(Double_t atom, const char *unit="g/mole")
Creates "atom" node for GDML.
Definition: TGDMLWrite.cxx:555
Double_t GetTwistAngle() const
Definition: TGeoArb8.h:166
TGeoScale * GetScale() const
virtual Int_t GetDivAxis()
Double_t * GetZ() const
Definition: TGeoXtru.h:100
Mother of all ROOT objects.
Definition: TObject.h:37
Double_t GetPhi1() const
Definition: TGeoTorus.h:72
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Double_t GetDx1() const
Definition: TGeoTrd2.h:56
constexpr Double_t RadToDeg()
Definition: TMath.h:60
void SetAttBit(UInt_t f)
Definition: TGeoAtt.h:65
virtual TGeoHMatrix Inverse() const =0
XMLNodePointer_t ChooseObject(TGeoShape *geoShape)
Chooses the object and method that should be used for processing object.
XMLNodePointer_t CreateSphereN(TGeoSphere *geoShape)
Creates "sphere" node for GDML.
Definition: TGDMLWrite.cxx:789
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition: TGeoNode.h:39
void ExtractVolumes(TGeoVolume *volume)
Method extracting geometry structure recursively.
Definition: TGDMLWrite.cxx:395
Double_t GetDz() const
Definition: TGeoArb8.h:63
Double_t GetTl1() const
Definition: TGeoArb8.h:128
XMLNodePointer_t fMaterialsNode
Definition: TGDMLWrite.h:134
Elliptical tube class.
Definition: TGeoEltu.h:17
virtual Double_t GetRmin() const
Definition: TGeoTube.h:66
NameLst * fNameList
Definition: TGDMLWrite.h:123
XMLNodePointer_t NewChild(XMLNodePointer_t parent, XMLNsPointer_t ns, const char *name, const char *content=0)
create new child element for parent node
Definition: TXMLEngine.cxx:707
XMLNodePointer_t CreateTubeN(TGeoTubeSeg *geoShape)
Creates "tube" node for GDML from object TGeoTubeSeg.
Int_t fgNamingSpeed
Definition: TGDMLWrite.h:127
StructLst * fElementList
Definition: TGDMLWrite.h:119
Double_t Sin(Double_t)
Definition: TMath.h:547
TObjArray * GetListOfVolumes() const
Definition: TGeoManager.h:465
Double_t GetAlpha1() const
Definition: TGeoArb8.h:129
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
virtual Double_t GetA() const
Definition: TGeoEltu.h:43
XMLNodePointer_t CreateDN(Double_t density, const char *unit="g/cm3")
Creates "D" density node for GDML.
Definition: TGDMLWrite.cxx:566
static const UInt_t fgkMaxNameErr
Definition: TGDMLWrite.h:144
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
XMLNodePointer_t CreateElementN(TGeoElement *element, XMLNodePointer_t materials, const char *name)
Creates "element" node for GDML element node and attribute.
Definition: TGDMLWrite.cxx:602
TGeoShape * GetShape() const
XMLNodePointer_t CreateCommonBoolN(TGeoCompositeShape *geoShape)
Creates common part of union intersection and subtraction nodes.
UInt_t fSolCnt
Definition: TGDMLWrite.h:140
Bool_t IsDigit() const
Returns true if all characters in string are digits (0-9) or white spaces, i.e.
Definition: TString.cxx:1817
Double_t A() const
Definition: TGeoElement.h:76
Double_t GetDz() const
Definition: TGeoTrd1.h:58
void SetAngles(Double_t phi, Double_t theta, Double_t psi)
Set matrix elements according to Euler angles.
Double_t GetTheta() const
Definition: TGeoPara.h:65
Double_t Sqrt(Double_t x)
Definition: TMath.h:590
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:908
Double_t GetRelativeAbundance(Int_t i) const
Return relative abundance of i-th isotope in this element.
TGeoVolume * GetTopVolume() const
Definition: TGeoManager.h:500
TGeoShape * GetShape() const
Definition: TGeoVolume.h:191
Double_t GetZ() const
Definition: TGeoPara.h:63
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
Double_t GetDx1() const
Definition: TGeoTrd1.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
TGeoBoolNode * GetBoolNode() const
TGeoShape * GetRightShape() const
Definition: TGeoBoolNode.h:81
virtual const Double_t * GetTranslation() const =0
TGeoVolume * GetVolume() const
Definition: TGeoNode.h:94
std::map< TString, Int_t > NameListI
Definition: TGDMLWrite.h:107
XMLNodePointer_t CreateSetupN(const char *topVolName, const char *name="default", const char *version="1.0")
Creates "setup" node for GDML.
virtual Double_t GetRmin() const
Definition: TGeoSphere.h:67
TString GetPattAxis(Int_t divAxis, const char *pattName, TString &unit)
Method that retrieves axis and unit along which object is divided.
virtual Double_t GetDz() const
Definition: TGeoTube.h:68
char name[80]
Definition: TGX11.cxx:109
XMLNodePointer_t CreateRotationN(const char *name, Xyz rotation, const char *type="rotation", const char *unit="deg")
Creates "rotation" kind of node for GDML.
StructLst * fAccPatt
Definition: TGDMLWrite.h:120
Double_t GetTl2() const
Definition: TGeoArb8.h:132
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
static const UInt_t fgkProcBit
Definition: TGDMLWrite.h:142
virtual Double_t GetDZ() const
Definition: TGeoBBox.h:72
Double_t * GetVertices()
Definition: TGeoArb8.h:67
const char * Data() const
Definition: TString.h:345
Double_t GetStIn() const
Definition: TGeoHype.h:68
Double_t GetTheta2() const
Definition: TGeoSphere.h:70
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
XMLNodePointer_t CreateCutTubeN(TGeoCtub *geoShape)
Creates "cutTube" node for GDML.
A phi segment of a tube.
Definition: TGeoTube.h:88