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