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