Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
16geometry to GDML file. gGeoManager is the instance of TGeoManager class
17containing tree of geometries creating resulting geometry. GDML is xml
18based format of file mirroring the tree of geometries according to GDML
19schema rules. For more information about GDML see http://gdml.web.cern.ch.
20Each 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
23by gGeoManager->Export("xyz.gdml") method. Export is driven by extenstion
24that is why ".gdml" is important in resulting name.
25
26 Whenever a new ROOT geometry object is implemented or there is a change
27in GDML schema this class is needed to be updated to ensure proper mapping
28between ROOT objects and GDML elements.
29
30 Current status of mapping ROOT -> GDML is implemented in method called
31TGDMLWrite::ChooseObject and it contains following "map":
32
33#### Solids:
34
35~~~
36TGeoBBox -> <box ... >
37TGeoParaboloid -> <paraboloid ...>
38TGeoSphere -> <sphere ...>
39TGeoArb8 -> <arb8 ...>
40TGeoConeSeg -> <cone ...>
41TGeoCone -> <cone ...>
42TGeoPara -> <para ...>
43TGeoTrap -> <trap ...> or
44- -> <arb8 ...>
45TGeoGtra -> <twistedtrap ...> or
46- -> <trap ...> or
47- -> <arb8 ...>
48TGeoTrd1 -> <trd ...>
49TGeoTrd2 -> <trd ...>
50TGeoTubeSeg -> <tube ...>
51TGeoCtub -> <cutTube ...>
52TGeoTube -> <tube ...>
53TGeoPcon -> <polycone ...>
54TGeoTorus -> <torus ...>
55TGeoPgon -> <polyhedra ...>
56TGeoEltu -> <eltube ...>
57TGeoHype -> <hype ...>
58TGeoXtru -> <xtru ...>
59TGeoTessellated -> <tessellated ...>
60TGeoCompositeShape -> <union ...> or
61- -> <subtraction ...> or
62- -> <intersection ...>
63
64Special cases of solids:
65TGeoScaledShape -> <elcone ...> if scaled TGeoCone or
66- -> element without scale
67TGeoCompositeShape -> <ellipsoid ...>
68- intersection of:
69- scaled TGeoSphere and TGeoBBox
70~~~
71
72#### Materials:
73
74~~~
75TGeoIsotope -> <isotope ...>
76TGeoElement -> <element ...>
77TGeoMaterial -> <material ...>
78TGeoMixture -> <material ...>
79~~~
80
81#### Structure
82
83~~~
84TGeoVolume -> <volume ...> or
85- -> <assembly ...>
86TGeoNode -> <physvol ...>
87TGeoPatternFinder -> <divisionvol ...>
88~~~
89
90There are options that can be set to change resulting document
91
92##### Options:
93
94~~~
95g - is set by default in gGeoManager, this option ensures compatibility
96- with Geant4. It means:
97- -> atomic number of material will be changed if <1 to 1
98- -> if polycone is set badly it will try to export it correctly
99- -> if widht * ndiv + offset is more then width of object being divided
100- (in divisions) then it will be rounded so it will not exceed or
101- if kPhi divsion then it will keep range of offset in -360 -> 0
102f - if this option is set then names of volumes and solids will have
103- pointer as a suffix to ensure uniqness of names
104n - if this option is set then names will not have suffix, but uniqness is
105- of names is not secured
106- - if none of this two options (f,n) is set then default behaviour is so
107- that incremental suffix is added to the names.
108- (eg. TGeoBBox_0x1, TGeoBBox_0x2 ...)
109~~~
110
111#### USAGE:
112
113~~~
114gGeoManager->Export("output.gdml");
115gGeoManager->Export("output.gdml","","vg"); //the same as previous just
116 //options are set explicitly
117gGeoManager->Export("output.gdml","","vgf");
118gGeoManager->Export("output.gdml","","gn");
119gGeoManager->Export("output.gdml","","f");
120...
121~~~
122
123#### Note:
124 Options discussed above are used only for TGDMLWrite class. There are
125other options in the TGeoManager::Export(...) method that can be used.
126See that function for details.
127
128*/
129
130#include "TGDMLWrite.h"
131
132#include "TGeoManager.h"
133#include "TGeoMaterial.h"
134#include "TGeoMatrix.h"
135#include "TXMLEngine.h"
136#include "TGeoVolume.h"
137#include "TGeoBBox.h"
138#include "TGeoParaboloid.h"
139#include "TGeoArb8.h"
140#include "TGeoTube.h"
141#include "TGeoCone.h"
142#include "TGeoTrd1.h"
143#include "TGeoTrd2.h"
144#include "TGeoPcon.h"
145#include "TGeoPgon.h"
146#include "TGeoSphere.h"
147#include "TGeoTorus.h"
148#include "TGeoPara.h"
149#include "TGeoHype.h"
150#include "TGeoEltu.h"
151#include "TGeoXtru.h"
152#include "TGeoScaledShape.h"
153#include "TMath.h"
154#include "TGeoBoolNode.h"
155#include "TGeoMedium.h"
156#include "TGeoElement.h"
157#include "TGeoShape.h"
158#include "TGeoCompositeShape.h"
159#include "TGeoOpticalSurface.h"
160#include <cstdlib>
161#include <string>
162#include <map>
163#include <set>
164#include <ctime>
165#include <sstream>
166
168
170
171namespace {
172
173// Helper to replace string patterns
174std::string str_replace(const std::string &str, const std::string &pattern, const std::string &replacement)
175{
176 std::string res = str;
177 for (size_t id = res.find(pattern); id != std::string::npos; id = res.find(pattern))
178 res.replace(id, pattern.length(), replacement);
179 return res;
180}
181
182// Create a NCN compliant name (no '/' and no '#')
183std::string make_NCName(const std::string &in)
184{
185 std::string res = str_replace(in, "/", "_");
186 res = str_replace(res, "#", "_");
187 return res;
188}
189
190// Materials extractor from a volume tree
191struct MaterialExtractor {
192 std::set<TGeoMaterial *> materials;
193 void operator()(const TGeoVolume *v)
194 {
195 materials.insert(v->GetMaterial());
196 for (Int_t i = 0; i < v->GetNdaughters(); ++i)
197 (*this)(v->GetNode(i)->GetVolume());
198 }
199};
200} // namespace
201
202////////////////////////////////////////////////////////////////////////////////
203/// Default constructor.
204
206 : TObject(),
207 fIsotopeList(nullptr),
208 fElementList(nullptr),
209 fAccPatt(nullptr),
210 fRejShape(nullptr),
211 fNameList(nullptr),
212 fgNamingSpeed(0),
213 fgG4Compatibility(false),
214 fGdmlFile(nullptr),
215 fTopVolumeName(0),
216 fGdmlE(nullptr),
217 fDefineNode(nullptr),
218 fMaterialsNode(nullptr),
219 fSolidsNode(nullptr),
220 fStructureNode(nullptr),
221 fVolCnt(0),
222 fPhysVolCnt(0),
223 fActNameErr(0),
224 fSolCnt(0),
225 fFltPrecision(17) // %.17g
226{
227 if (fgGDMLWrite)
228 delete fgGDMLWrite;
229 fgGDMLWrite = this;
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// Destructor.
234
236{
237 delete fIsotopeList;
238 delete fElementList;
239 delete fAccPatt;
240 delete fRejShape;
241 delete fNameList;
242
243 fgGDMLWrite = nullptr;
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Set convention of naming solids and volumes
248
253
254////////////////////////////////////////////////////////////////////////////////
255/// Ignore dummy material instance, which causes trouble reading GDML in Geant4
256
261
262////////////////////////////////////////////////////////////////////////////////
263// wrapper of all main methods for extraction
265{
266 TList *materials = geomanager->GetListOfMaterials();
267 TGeoNode *node = geomanager->GetTopNode();
268 if (!node) {
269 Info("WriteGDMLfile", "Top volume does not exist!");
270 return;
271 }
272 fTopVolume = node->GetVolume();
274 WriteGDMLfile(geomanager, node, materials, filename, option);
275}
276
277////////////////////////////////////////////////////////////////////////////////
278// Wrapper to only selectively write one branch of the volume hierarchy to file
280{
281 TGeoVolume *volume = node->GetVolume();
282 TList materials, volumes, nodes;
283 MaterialExtractor extract;
284 if (!volume) {
285 Info("WriteGDMLfile", "Invalid Volume reference to extract GDML information!");
286 return;
287 }
288 extract(volume);
289 for (TGeoMaterial *m : extract.materials)
290 materials.Add(m);
291 fTopVolumeName = volume->GetName();
292 fTopVolume = volume;
293 fSurfaceList.clear();
294 fVolumeList.clear();
295 fNodeList.clear();
296 WriteGDMLfile(geomanager, node, &materials, filename, option);
297 materials.Clear("nodelete");
298 volumes.Clear("nodelete");
299 nodes.Clear("nodelete");
300}
301
302////////////////////////////////////////////////////////////////////////////////
303/// Wrapper of all exporting methods
304/// Creates blank GDML file and fills it with gGeoManager structure converted
305/// to GDML structure of xml nodes
306
309{
310 // option processing
311 option.ToLower();
312 if (option.Contains("g")) {
314 Info("WriteGDMLfile", "Geant4 compatibility mode set");
315 } else {
317 }
318 if (option.Contains("f")) {
320 Info("WriteGDMLfile", "Fast naming convention with pointer suffix set");
321 } else if (option.Contains("n")) {
323 Info("WriteGDMLfile", "Naming without prefix set - be careful uniqness of name is not ensured");
324 } else {
326 Info("WriteGDMLfile", "Potentially slow with incremental suffix naming convention set");
327 }
328
330
332 switch (def_units) {
333 case TGeoManager::kG4Units: fDefault_lunit = "mm"; break;
334 case TGeoManager::kRootUnits: fDefault_lunit = "cm"; break;
335 default: // G4 units
336 fDefault_lunit = "mm";
337 break;
338 }
339
340 // local variables
342 const char *krootNodeName = "gdml";
343 const char *knsRefGeneral = "http://www.w3.org/2001/XMLSchema-instance";
344 const char *knsNameGeneral = "xsi";
345 const char *knsRefGdml = "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd";
346 const char *knsNameGdml = "xsi:noNamespaceSchemaLocation";
347
348 // First create engine
349 fGdmlE = new TXMLEngine;
351
352 // create blank GDML file
354
355 // create root node and add it to blank GDML file
356 XMLNodePointer_t rootNode = fGdmlE->NewChild(nullptr, nullptr, krootNodeName, nullptr);
358
359 // add namespaces to root node
362
363 // initialize general lists and <define>, <solids>, <structure> nodes
366
367 fNameList = new NameLst;
368
369 fDefineNode = fGdmlE->NewChild(nullptr, nullptr, "define", nullptr);
370 fSolidsNode = fGdmlE->NewChild(nullptr, nullptr, "solids", nullptr);
371 fStructureNode = fGdmlE->NewChild(nullptr, nullptr, "structure", nullptr);
372 //========================
373
374 // initialize list of accepted patterns for divisions (in ExtractVolumes)
375 fAccPatt = new StructLst;
376 fAccPatt->fLst["TGeoPatternX"] = kTRUE;
377 fAccPatt->fLst["TGeoPatternY"] = kTRUE;
378 fAccPatt->fLst["TGeoPatternZ"] = kTRUE;
379 fAccPatt->fLst["TGeoPatternCylR"] = kTRUE;
380 fAccPatt->fLst["TGeoPatternCylPhi"] = kTRUE;
381 //========================
382
383 // initialize list of rejected shapes for divisions (in ExtractVolumes)
384 fRejShape = new StructLst;
385 // this shapes are rejected because, it is not possible to divide trd2
386 // in Y axis and while only trd2 object is imported from GDML
387 // it causes a problem when TGeoTrd1 is divided in Y axis
388 fRejShape->fLst["TGeoTrd1"] = kTRUE;
389 fRejShape->fLst["TGeoTrd2"] = kTRUE;
390 //=========================
391
392 // Initialize global counters
393 fActNameErr = 0;
394 fVolCnt = 0;
395 fPhysVolCnt = 0;
396 fSolCnt = 0;
397
398 // calling main extraction functions (with measuring time)
399 time_t startT, endT;
400 startT = time(nullptr);
401 ExtractMatrices(geomanager->GetListOfGDMLMatrices());
404
405 Info("WriteGDMLfile", "Extracting volumes");
406 ExtractVolumes(node);
407 Info("WriteGDMLfile", "%i solids added", fSolCnt);
408 Info("WriteGDMLfile", "%i volumes added", fVolCnt);
409 Info("WriteGDMLfile", "%i physvolumes added", fPhysVolCnt);
410 ExtractSkinSurfaces(geomanager->GetListOfSkinSurfaces());
411 ExtractBorderSurfaces(geomanager->GetListOfBorderSurfaces());
412 ExtractOpticalSurfaces(geomanager->GetListOfOpticalSurfaces());
413 endT = time(nullptr);
414 //<gdml>
415 fGdmlE->AddChild(rootNode, fDefineNode); // <define>...</define>
416 fGdmlE->AddChild(rootNode, fMaterialsNode); // <materials>...</materials>
417 fGdmlE->AddChild(rootNode, fSolidsNode); // <solids>...</solids>
418 fGdmlE->AddChild(rootNode, fStructureNode); // <structure>...</structure>
419 fGdmlE->AddChild(rootNode, CreateSetupN(fTopVolumeName.Data())); // <setup>...</setup>
420 //</gdml>
422 TString tdiffS = (tdiffI == 0 ? TString("< 1 s") : TString::Format("%.0lf s", tdiffI));
423 Info("WriteGDMLfile", "Exporting time: %s", tdiffS.Data());
424 //=========================
425
426 // Saving document
428 Info("WriteGDMLfile", "File %s saved", filename);
429 // cleaning
431 // unset processing bits:
433 delete fGdmlE;
434}
435
436////////////////////////////////////////////////////////////////////////////////
437/// Method exporting GDML matrices
438
440{
441 if (!matrixList->GetEntriesFast())
442 return;
444 TIter next(matrixList);
446 while ((matrix = (TGDMLMatrix *)next())) {
449 }
450}
451
452////////////////////////////////////////////////////////////////////////////////
453/// Method exporting GDML matrices
454
456{
457 if (!geom->GetNproperties())
458 return;
462 for (Int_t i = 0; i < geom->GetNproperties(); ++i) {
463 value = geom->GetProperty(i, property);
466 }
467}
468
469////////////////////////////////////////////////////////////////////////////////
470/// Method exporting optical surfaces
471
surfaces)
473{
474 if (!surfaces->GetEntriesFast())
475 return;
477 surfaces);
479 while ((surf = (TGeoOpticalSurface *)next())) {
480 if (fSurfaceList.find(surf) == fSurfaceList.end())
481 continue;
484 // Info("ExtractSkinSurfaces", "Extracted optical surface: %s",surf->GetName());
485 }
486}
487
488////////////////////////////////////////////////////////////////////////////////
489/// Method exporting skin surfaces
490
surfaces)
492{
493 if (!surfaces->GetEntriesFast())
494 return;
496 surfaces);
498 while ((surf = (TGeoSkinSurface *)next())) {
499 if (fVolumeList.find(surf->GetVolume()) == fVolumeList.end())
500 continue;
503 fSurfaceList.insert(surf->GetSurface());
504 // Info("ExtractSkinSurfaces", "Extracted skin surface: %s",surf->GetName());
505 }
506}
507
508////////////////////////////////////////////////////////////////////////////////
509/// Method exporting border surfaces
510
surfaces)
512{
513 if (!surfaces->GetEntriesFast())
514 return;
516 surfaces);
518 while ((surf = (TGeoBorderSurface *)next())) {
519 auto ia = fNodeList.find(surf->GetNode1());
520 auto ib = fNodeList.find(surf->GetNode2());
521 if (ia == fNodeList.end() && ib == fNodeList.end()) {
522 continue;
523 } else if (ia == fNodeList.end() && ib != fNodeList.end()) {
524 Warning("ExtractBorderSurfaces",
525 "Inconsistent border surface extraction %s: Node %s"
526 " is not part of GDML!",
527 surf->GetName(), surf->GetNode1()->GetName());
528 continue;
529 } else if (ia != fNodeList.end() && ib == fNodeList.end()) {
530 Warning("ExtractBorderSurfaces",
531 "Inconsistent border surface extraction %s: Node %s"
532 " is not part of GDML!",
533 surf->GetName(), surf->GetNode2()->GetName());
534 continue;
535 }
538 fSurfaceList.insert(surf->GetSurface());
539 // Info("ExtractBorderSurfaces", "Extracted border surface: %s",surf->GetName());
540 }
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// Method exporting materials
545
547{
548 Info("ExtractMaterials", "Extracting materials");
549 // crate main <materials> node
550 XMLNodePointer_t materialsN = fGdmlE->NewChild(nullptr, nullptr, "materials", nullptr);
551 Int_t matcnt = 0;
552
553 // go through materials - iterator and object declaration
554 TIter next(materialsLst);
557 TGeoMaterial *dummy_mat = dummy_med ? dummy_med->GetMaterial() : nullptr;
558 std::string dummy_nam = dummy_mat ? dummy_mat->GetName() : "dummy";
559
560 while ((lmaterial = (TGeoMaterial *)next())) {
561 // check for dummy material: if requested, ignore it
562 std::string mname = lmaterial->GetName();
564 Info("ExtractMaterials", "Skip dummy material: %s", dummy_nam.c_str());
565 continue;
566 }
567 // generate uniq name
569
570 if (lmaterial->IsMixture()) {
574 } else {
577 }
578 matcnt++;
579 }
580 Info("ExtractMaterials", "%i materials added", matcnt);
581 return materialsN;
582}
583
584////////////////////////////////////////////////////////////////////////////////
585/// Method creating solid to xml file and returning its name
586
588{
590 TString solname = "";
591 solidN = ChooseObject(volShape); // volume->GetShape()
593 if (solidN != nullptr)
594 fSolCnt++;
596 if (solname.Contains("missing_")) {
597 solname = "-1";
598 }
599 return solname;
600}
601
602////////////////////////////////////////////////////////////////////////////////
603/// Method extracting geometry structure recursively
604
606{
608 TGeoVolume *volume = node->GetVolume();
610 TGeoPatternFinder *pattFinder = nullptr;
613
614 fNodeList.insert(node);
615 fVolumeList.insert(volume);
616 // create the name for volume/assembly
617 if (volume == fTopVolume) {
618 // not needed a special function for generating name
619 volname = volume->GetName();
621 // register name to the pointer
622 fNameList->fLst[TString::Format("%p", volume)] = volname;
623 } else {
624 volname = GenName(volume->GetName(), TString::Format("%p", volume));
625 }
626
627 // start to create main volume/assembly node
628 if (volume->IsAssembly()) {
630 } else {
631 // get reference material and add solid to <solids> + get name
632 matname = fNameList->fLst[TString::Format("%p", volume->GetMaterial())];
633 solname = ExtractSolid(volume->GetShape());
634 // If solid is not supported or corrupted
635 if (solname == "-1") {
636 Info("ExtractVolumes", "ERROR! %s volume was not added, because solid is either not supported or corrupted",
637 volname.Data());
638 // set volume as missing volume
639 fNameList->fLst[TString::Format("%p", volume)] = "missing_" + volname;
640 return;
641 }
643
644 // divisionvol can't be in assembly
645 pattFinder = volume->GetFinder();
646 // if found pattern
647 if (pattFinder) {
648 pattClsName = TString::Format("%s", pattFinder->ClassName());
649 TString shapeCls = TString::Format("%s", volume->GetShape()->ClassName());
650 // if pattern in accepted pattern list and not in shape rejected list
651 if ((fAccPatt->fLst[pattClsName] == kTRUE) && (fRejShape->fLst[shapeCls] != kTRUE)) {
653 }
654 }
655 }
656 // get all nodes in volume
657 TObjArray *nodeLst = volume->GetNodes();
658 TIter next(nodeLst);
661 Int_t nCnt = 0;
662 // loop through all nodes
663 while ((geoNode = (TGeoNode *)next())) {
664 // get volume of current node and if not processed then process it
665 TGeoVolume *subvol = geoNode->GetVolume();
666 fNodeList.insert(geoNode);
667 if (subvol->TestAttBit(fgkProcBitVol) == kFALSE) {
668 subvol->SetAttBit(fgkProcBitVol);
670 }
671
672 // volume of this node has to exist because it was processed recursively
674 if (nodevolname.Contains("missing_")) {
675 continue;
676 }
677 if (nCnt == 0) { // save name of the first node for divisionvol
679 }
680
681 if (isPattern == kFALSE) {
682 // create name for node
683 TString nodename, posname, rotname;
684 nodename = GenName(geoNode->GetName(), TString::Format("%p", geoNode));
685 nodename = nodename + "in" + volname;
686
687 // create name for position and clear rotation
688 posname = nodename + "pos";
689 rotname = "";
690
691 // position
692 const Double_t *pos = geoNode->GetMatrix()->GetTranslation();
693 Xyz nodPos;
694 nodPos.x = pos[0];
695 nodPos.y = pos[1];
696 nodPos.z = pos[2];
697 childN = CreatePositionN(posname.Data(), nodPos, "position", fDefault_lunit);
698 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
699 // Deal with reflection
700 XMLNodePointer_t scaleN = nullptr;
701 Double_t lx, ly, lz;
702 Double_t xangle = 0;
703 Double_t zangle = 0;
704 lx = geoNode->GetMatrix()->GetRotationMatrix()[0];
705 ly = geoNode->GetMatrix()->GetRotationMatrix()[4];
706 lz = geoNode->GetMatrix()->GetRotationMatrix()[8];
707 if (geoNode->GetMatrix()->IsReflection() && TMath::Abs(lx) == 1 && TMath::Abs(ly) == 1 &&
708 TMath::Abs(lz) == 1) {
709 scaleN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
710 fGdmlE->NewAttr(scaleN, nullptr, "name", (nodename + "scl").Data());
711 fGdmlE->NewAttr(scaleN, nullptr, "x", TString::Format(fltPrecision.Data(), lx));
712 fGdmlE->NewAttr(scaleN, nullptr, "y", TString::Format(fltPrecision.Data(), ly));
713 fGdmlE->NewAttr(scaleN, nullptr, "z", TString::Format(fltPrecision.Data(), lz));
714 // experimentally found out, that rotation should be updated like this
715 if (lx == -1) {
716 zangle = 180;
717 }
718 if (lz == -1) {
719 xangle = 180;
720 }
721 }
722
723 // rotation
724 TGDMLWrite::Xyz lxyz = GetXYZangles(geoNode->GetMatrix()->GetRotationMatrix());
725 lxyz.x -= xangle;
726 lxyz.z -= zangle;
727 if ((lxyz.x != 0.0) || (lxyz.y != 0.0) || (lxyz.z != 0.0)) {
728 rotname = nodename + "rot";
729 childN = CreateRotationN(rotname.Data(), lxyz);
730 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
731 }
732
733 // create physvol for main volume/assembly node
735 childN = CreatePhysVolN(physvolname, geoNode->GetNumber(), nodevolname.Data(), posname.Data(), rotname.Data(),
736 scaleN);
738 }
739 nCnt++;
740 }
741 // create only one divisionvol node
742 if (isPattern && pattFinder) {
743 // retrieve attributes of division
744 Int_t ndiv, divaxis;
745 Double_t offset, width, xlo, xhi;
746 TString axis, unit;
747
748 ndiv = pattFinder->GetNdiv();
749 width = pattFinder->GetStep();
750
751 divaxis = pattFinder->GetDivAxis();
752 volume->GetShape()->GetAxisRange(divaxis, xlo, xhi);
753
754 // compute relative start (not positional)
755 offset = pattFinder->GetStart() - xlo;
756 axis = GetPattAxis(divaxis, pattClsName, unit);
757
758 // create division node
759 childN = CreateDivisionN(offset, width, ndiv, axis.Data(), unit.Data(), nodeVolNameBak.Data());
761 }
762
763 fVolCnt++;
764 // add volume/assembly node into the <structure> node
766}
767
768////////////////////////////////////////////////////////////////////////////////
769/// Creates "atom" node for GDML
770
772{
774 XMLNodePointer_t atomN = fGdmlE->NewChild(nullptr, nullptr, "atom", nullptr);
775 fGdmlE->NewAttr(atomN, nullptr, "unit", unit);
776 fGdmlE->NewAttr(atomN, nullptr, "value", TString::Format(fltPrecision.Data(), atom));
777 return atomN;
778}
779
780////////////////////////////////////////////////////////////////////////////////
781/// Creates "D" density node for GDML
782
784{
786 XMLNodePointer_t densN = fGdmlE->NewChild(nullptr, nullptr, "D", nullptr);
787 fGdmlE->NewAttr(densN, nullptr, "unit", unit);
788 fGdmlE->NewAttr(densN, nullptr, "value", TString::Format(fltPrecision.Data(), density));
789 return densN;
790}
791
792////////////////////////////////////////////////////////////////////////////////
793/// Creates "fraction" node for GDML
794
796{
798 XMLNodePointer_t fractN = fGdmlE->NewChild(nullptr, nullptr, "fraction", nullptr);
800 fGdmlE->NewAttr(fractN, nullptr, "ref", refName);
801 return fractN;
802}
803
804////////////////////////////////////////////////////////////////////////////////
805/// Creates "property" node for GDML
806
808{
809 XMLNodePointer_t propertyN = fGdmlE->NewChild(nullptr, nullptr, "property", nullptr);
810 fGdmlE->NewAttr(propertyN, nullptr, "name", property.GetName());
811 fGdmlE->NewAttr(propertyN, nullptr, "ref", property.GetTitle());
812 return propertyN;
813}
814
815////////////////////////////////////////////////////////////////////////////////
816/// Creates "isotope" node for GDML
817
819{
820 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "isotope", nullptr);
821 fGdmlE->NewAttr(mainN, nullptr, "name", name);
822 fGdmlE->NewAttr(mainN, nullptr, "N", TString::Format("%i", isotope->GetN()));
823 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format("%i", isotope->GetZ()));
825 return mainN;
826}
827
828////////////////////////////////////////////////////////////////////////////////
829/// Creates "element" node for GDML
830/// element node and attribute
831
833{
834 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "element", nullptr);
835 fGdmlE->NewAttr(mainN, nullptr, "name", name);
836 // local associative arrays for saving isotopes and their weight
837 // inside element
840
841 if (element->HasIsotopes()) {
842 Int_t nOfIso = element->GetNisotopes();
843 // go through isotopes
844 for (Int_t idx = 0; idx < nOfIso; idx++) {
845 TGeoIsotope *myIsotope = element->GetIsotope(idx);
846 if (!myIsotope) {
847 Fatal("CreateElementN", "Missing isotopes for element %s", element->GetName());
848 return mainN;
849 }
850
851 // Get name of the Isotope (
852 TString lname = myIsotope->GetName();
853 //_iso suffix is added to avoid problems with same names
854 // for material, element and isotopes
855 lname = TString::Format("%s_iso", lname.Data());
856
857 // cumulates abundance, in case 2 isotopes with same names
858 // within one element
859 wPercentage[lname] += element->GetRelativeAbundance(idx);
860 wCounter[lname]++;
861
862 // check whether isotope name is not in list of isotopes
864 continue;
865 }
866 // add isotope to list of isotopes and to main <materials> node
869 fGdmlE->AddChild(materials, isoNode);
870 }
871 // loop through asoc array of isotopes
872 for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
873 if (itr->second > 1) {
874 Info("CreateMixtureN", "WARNING! 2 equal isotopes in one element. Check: %s isotope of %s element",
875 itr->first.Data(), name);
876 }
877 // add fraction child to element with reference to isotope
878 fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
879 }
880 } else {
881 fGdmlE->NewAttr(mainN, nullptr, "formula", element->GetName());
882 Int_t valZ = element->Z();
883 // Z can't be <1 in Geant4 and Z is optional parameter
884 if (valZ >= 1) {
885 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format("%i", valZ));
886 }
887 Int_t valN = element->N();
888 fGdmlE->NewAttr(mainN, nullptr, "N", TString::Format("%i", valN));
890 }
891 return mainN;
892}
893
894////////////////////////////////////////////////////////////////////////////////
895/// Creates "material" node for GDML with references to other sub elements
896
898{
899 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "material", nullptr);
900 fGdmlE->NewAttr(mainN, nullptr, "name", mname);
901
902 // Write properties
903 TList const &properties = mixture->GetProperties();
904 if (properties.GetSize()) {
905 TIter next(&properties);
907 while ((property = (TNamed *)next()))
909 }
910 // Write CONST properties
911 TList const &const_properties = mixture->GetConstProperties();
912 if (const_properties.GetSize()) {
913 TIter next(&const_properties);
915 while ((property = (TNamed *)next()))
917 }
918
919 fGdmlE->AddChild(mainN, CreateDN(mixture->GetDensity()));
920 // local associative arrays for saving elements and their weight
921 // inside mixture
924
925 Int_t nOfElm = mixture->GetNelements();
926 // go through elements
927 for (Int_t idx = 0; idx < nOfElm; idx++) {
928 TGeoElement *myElement = mixture->GetElement(idx);
929
930 // Get name of the element
931 // NOTE: that for element - GetTitle() returns the "name" tag
932 // and GetName() returns "formula" tag (see createElementN)
933 TString lname = myElement->GetTitle();
934 //_elm suffix is added to avoid problems with same names
935 // for material and element
936 lname = TString::Format("%s_elm", lname.Data());
937
938 // cumulates percentage, in case 2 elements with same names within one mixture
939 wPercentage[lname] += mixture->GetWmixt()[idx];
940 wCounter[lname]++;
941
942 // check whether element name is not in list of elements already created
944 continue;
945 }
946
947 // add element to list of elements and to main <materials> node
950 fGdmlE->AddChild(materials, elmNode);
951 }
952 // loop through asoc array
953 for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
954 if (itr->second > 1) {
955 Info("CreateMixtureN", "WARNING! 2 equal elements in one material. Check: %s element of %s material",
956 itr->first.Data(), mname.Data());
957 }
958 // add fraction child to material with reference to element
959 fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
960 }
961
962 return mainN;
963}
964
965////////////////////////////////////////////////////////////////////////////////
966/// Creates "material" node for GDML
967
969{
970 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "material", nullptr);
971 fGdmlE->NewAttr(mainN, nullptr, "name", mname);
972 Double_t valZ = material->GetZ();
973 // Z can't be zero in Geant4 so this is workaround for vacuum
976 tmpname.ToLower();
977 if (valZ < 1) {
978 if (tmpname == "vacuum") {
979 valZ = 1;
980 } else {
981 if (fgG4Compatibility == kTRUE) {
982 Info("CreateMaterialN",
983 "WARNING! value of Z in %s material can't be < 1 in Geant4, that is why it was changed to 1, please "
984 "check it manually! ",
985 mname.Data());
986 valZ = 1;
987 } else {
988 Info("CreateMaterialN", "WARNING! value of Z in %s material can't be < 1 in Geant4", mname.Data());
989 }
990 }
991 }
992 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format(fltPrecision.Data(), valZ)); // material->GetZ()));
993
994 // Create properties if any: Properties according to the GDML schema MUST come first!
995 TList const &properties = material->GetProperties();
996 if (properties.GetSize()) {
997 TIter next(&properties);
999 while ((property = (TNamed *)next()))
1001 }
1002 // Write CONST properties
1003 TList const &const_properties = material->GetConstProperties();
1004 if (const_properties.GetSize()) {
1005 TIter next(&const_properties);
1007 while ((property = (TNamed *)next()))
1009 }
1010
1011 // Now add the other children
1012 fGdmlE->AddChild(mainN, CreateDN(material->GetDensity()));
1013 fGdmlE->AddChild(mainN, CreateAtomN(material->GetA()));
1014 return mainN;
1015}
1016
1017////////////////////////////////////////////////////////////////////////////////
1018/// Creates "box" node for GDML
1019
1021{
1022 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "box", nullptr);
1024 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1025 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1026 if (IsNullParam(geoShape->GetDX(), "DX", lname) || IsNullParam(geoShape->GetDY(), "DY", lname) ||
1027 IsNullParam(geoShape->GetDZ(), "DZ", lname)) {
1028 return nullptr;
1029 }
1030 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDX()));
1031 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDY()));
1032 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDZ()));
1033
1034 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1035 return mainN;
1036}
1037
1038////////////////////////////////////////////////////////////////////////////////
1039/// Creates "paraboloid" node for GDML
1040
1042{
1043 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "paraboloid", nullptr);
1045 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1046 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1047 if (IsNullParam(geoShape->GetRhi(), "Rhi", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1048 return nullptr;
1049 }
1050 fGdmlE->NewAttr(mainN, nullptr, "rlo", TString::Format(fltPrecision.Data(), geoShape->GetRlo()));
1051 fGdmlE->NewAttr(mainN, nullptr, "rhi", TString::Format(fltPrecision.Data(), geoShape->GetRhi()));
1052 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1053
1054 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1055 return mainN;
1056}
1057
1058////////////////////////////////////////////////////////////////////////////////
1059/// Creates "sphere" node for GDML
1060
1062{
1063 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "sphere", nullptr);
1065 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1066 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1067 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
1068 return nullptr;
1069 }
1070
1071 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1072 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1073 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1074 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1075 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1076 fGdmlE->NewAttr(mainN, nullptr, "starttheta", TString::Format(fltPrecision.Data(), geoShape->GetTheta1()));
1077 fGdmlE->NewAttr(mainN, nullptr, "deltatheta",
1078 TString::Format(fltPrecision.Data(), geoShape->GetTheta2() - geoShape->GetTheta1()));
1079
1080 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1081 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1082 return mainN;
1083}
1084
1085////////////////////////////////////////////////////////////////////////////////
1086/// Creates "arb8" node for GDML
1087
1089{
1090 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "arb8", nullptr);
1092 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1093 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1094 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1095 return nullptr;
1096 }
1097
1098 fGdmlE->NewAttr(mainN, nullptr, "v1x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[0]));
1099 fGdmlE->NewAttr(mainN, nullptr, "v1y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[1]));
1100 fGdmlE->NewAttr(mainN, nullptr, "v2x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[2]));
1101 fGdmlE->NewAttr(mainN, nullptr, "v2y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[3]));
1102 fGdmlE->NewAttr(mainN, nullptr, "v3x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[4]));
1103 fGdmlE->NewAttr(mainN, nullptr, "v3y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[5]));
1104 fGdmlE->NewAttr(mainN, nullptr, "v4x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[6]));
1105 fGdmlE->NewAttr(mainN, nullptr, "v4y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[7]));
1106 fGdmlE->NewAttr(mainN, nullptr, "v5x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[8]));
1107 fGdmlE->NewAttr(mainN, nullptr, "v5y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[9]));
1108 fGdmlE->NewAttr(mainN, nullptr, "v6x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[10]));
1109 fGdmlE->NewAttr(mainN, nullptr, "v6y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[11]));
1110 fGdmlE->NewAttr(mainN, nullptr, "v7x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[12]));
1111 fGdmlE->NewAttr(mainN, nullptr, "v7y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[13]));
1112 fGdmlE->NewAttr(mainN, nullptr, "v8x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[14]));
1113 fGdmlE->NewAttr(mainN, nullptr, "v8y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[15]));
1114 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1115
1116 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1117 return mainN;
1118}
1119
1120////////////////////////////////////////////////////////////////////////////////
1121/// Creates "cone" node for GDML from TGeoConeSeg object
1122
1124{
1125 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "cone", nullptr);
1127 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1128 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1129 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1130 return nullptr;
1131 }
1132
1133 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1134 fGdmlE->NewAttr(mainN, nullptr, "rmin1", TString::Format(fltPrecision.Data(), geoShape->GetRmin1()));
1135 fGdmlE->NewAttr(mainN, nullptr, "rmin2", TString::Format(fltPrecision.Data(), geoShape->GetRmin2()));
1136 fGdmlE->NewAttr(mainN, nullptr, "rmax1", TString::Format(fltPrecision.Data(), geoShape->GetRmax1()));
1137 fGdmlE->NewAttr(mainN, nullptr, "rmax2", TString::Format(fltPrecision.Data(), geoShape->GetRmax2()));
1138 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1139 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1140 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1141
1142 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1143 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1144 return mainN;
1145}
1146
1147////////////////////////////////////////////////////////////////////////////////
1148/// Creates "cone" node for GDML from TGeoCone object
1149
1151{
1152 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "cone", nullptr);
1154 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1155 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1156 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1157 return nullptr;
1158 }
1159
1160 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1161 fGdmlE->NewAttr(mainN, nullptr, "rmin1", TString::Format(fltPrecision.Data(), geoShape->GetRmin1()));
1162 fGdmlE->NewAttr(mainN, nullptr, "rmin2", TString::Format(fltPrecision.Data(), geoShape->GetRmin2()));
1163 fGdmlE->NewAttr(mainN, nullptr, "rmax1", TString::Format(fltPrecision.Data(), geoShape->GetRmax1()));
1164 fGdmlE->NewAttr(mainN, nullptr, "rmax2", TString::Format(fltPrecision.Data(), geoShape->GetRmax2()));
1165 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format("%i", 0));
1166 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format("%i", 360));
1167
1168 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1169 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1170 return mainN;
1171}
1172
1173////////////////////////////////////////////////////////////////////////////////
1174/// Creates "para" node for GDML
1175
1177{
1178 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "para", nullptr);
1180 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1181
1182 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), 2 * geoShape->GetX()));
1183 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), 2 * geoShape->GetY()));
1184 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetZ()));
1185 fGdmlE->NewAttr(mainN, nullptr, "alpha", TString::Format(fltPrecision.Data(), geoShape->GetAlpha()));
1186 fGdmlE->NewAttr(mainN, nullptr, "theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1187 fGdmlE->NewAttr(mainN, nullptr, "phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1188
1189 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1190 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1191 return mainN;
1192}
1193
1194////////////////////////////////////////////////////////////////////////////////
1195/// Creates "trap" node for GDML
1196
1198{
1201
1202 // if one base equals 0 create Arb8 instead of trap
1203 if ((geoShape->GetBl1() == 0 || geoShape->GetTl1() == 0 || geoShape->GetH1() == 0) ||
1204 (geoShape->GetBl2() == 0 || geoShape->GetTl2() == 0 || geoShape->GetH2() == 0)) {
1206 return mainN;
1207 }
1208
1209 // if is twisted then create arb8
1210 if (geoShape->IsTwisted()) {
1212 return mainN;
1213 }
1214
1215 mainN = fGdmlE->NewChild(nullptr, nullptr, "trap", nullptr);
1216 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1217 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1218 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1219 return nullptr;
1220 }
1221
1222 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1223 fGdmlE->NewAttr(mainN, nullptr, "theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1224 fGdmlE->NewAttr(mainN, nullptr, "phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1225 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl1()));
1226 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl1()));
1227 fGdmlE->NewAttr(mainN, nullptr, "x3", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl2()));
1228 fGdmlE->NewAttr(mainN, nullptr, "x4", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl2()));
1229 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH1()));
1230 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH2()));
1231
1232 fGdmlE->NewAttr(mainN, nullptr, "alpha1", TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()));
1233 fGdmlE->NewAttr(mainN, nullptr, "alpha2", TString::Format(fltPrecision.Data(), geoShape->GetAlpha2()));
1234
1235 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1236 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1237 return mainN;
1238}
1239
1240////////////////////////////////////////////////////////////////////////////////
1241/// Creates "twistedtrap" node for GDML
1242
1244{
1247
1248 // if one base equals 0 create Arb8 instead of twisted trap
1249 if ((geoShape->GetBl1() == 0 && geoShape->GetTl1() == 0 && geoShape->GetH1() == 0) ||
1250 (geoShape->GetBl2() == 0 && geoShape->GetTl2() == 0 && geoShape->GetH2() == 0)) {
1252 return mainN;
1253 }
1254
1255 // if is twisted then create arb8
1256 if (geoShape->IsTwisted()) {
1258 return mainN;
1259 }
1260
1261 // if parameter twistAngle (PhiTwist) equals zero create trap node
1262 if (geoShape->GetTwistAngle() == 0) {
1264 return mainN;
1265 }
1266
1267 mainN = fGdmlE->NewChild(nullptr, nullptr, "twistedtrap", nullptr);
1268 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1269 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1270 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1271 return nullptr;
1272 }
1273
1274 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1275 fGdmlE->NewAttr(mainN, nullptr, "Theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1276 fGdmlE->NewAttr(mainN, nullptr, "Phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1277 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl1()));
1278 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl1()));
1279 fGdmlE->NewAttr(mainN, nullptr, "x3", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl2()));
1280 fGdmlE->NewAttr(mainN, nullptr, "x4", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl2()));
1281 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH1()));
1282 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH2()));
1283
1284 fGdmlE->NewAttr(mainN, nullptr, "Alph", TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()));
1285
1286 // check if alpha1 equals to alpha2 (converting to string - to avoid problems with floats)
1287 if (TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()) !=
1288 TString::Format(fltPrecision.Data(), geoShape->GetAlpha2())) {
1289 Info("CreateTwistedTrapN",
1290 "ERROR! Object %s is not exported correctly because parameter Alpha2 is not declared in GDML schema",
1291 lname.Data());
1292 }
1293 // fGdmlE->NewAttr(mainN,0, "alpha2", TString::Format(fltPrecision.Data(), geoShape->GetAlpha2()));
1294 fGdmlE->NewAttr(mainN, nullptr, "PhiTwist", TString::Format(fltPrecision.Data(), geoShape->GetTwistAngle()));
1295
1296 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1297 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1298 return mainN;
1299}
1300
1301////////////////////////////////////////////////////////////////////////////////
1302/// Creates "trd" node for GDML from object TGeoTrd1
1303
1305{
1306 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "trd", nullptr);
1308 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1309 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1310 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1311 return nullptr;
1312 }
1313
1314 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx1()));
1315 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx2()));
1316 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy()));
1317 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy()));
1318 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1319
1320 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1321 return mainN;
1322}
1323
1324////////////////////////////////////////////////////////////////////////////////
1325/// Creates "trd" node for GDML from object TGeoTrd2
1326
1328{
1329 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "trd", nullptr);
1331 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1332 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1333 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1334 return nullptr;
1335 }
1336
1337 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx1()));
1338 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx2()));
1339 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy1()));
1340 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy2()));
1341 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1342
1343 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1344 return mainN;
1345}
1346
1347////////////////////////////////////////////////////////////////////////////////
1348/// Creates "tube" node for GDML from object TGeoTubeSeg
1349
1351{
1352 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tube", nullptr);
1354 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1355 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1356 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1357 return nullptr;
1358 }
1359
1360 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1361 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1362 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1363 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1364 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1365 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1366
1367 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1368 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1369 return mainN;
1370}
1371
1372////////////////////////////////////////////////////////////////////////////////
1373/// Creates "cutTube" node for GDML
1374
1376{
1379
1380 mainN = fGdmlE->NewChild(nullptr, nullptr, "cutTube", nullptr);
1381 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1382 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1383 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1384 return nullptr;
1385 }
1386 // This is not needed, because cutTube is already supported by Geant4 9.5
1387 if (fgG4Compatibility == kTRUE && kFALSE) {
1390
1391 // register name for cuttube shape (so it will be found during volume export)
1394 Info("CreateCutTubeN", "WARNING! %s - CutTube was replaced by intersection of TGeoTubSeg and two TGeoBBoxes",
1395 lname.Data());
1396 return mainN;
1397 }
1398 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1399 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1400 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1401 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1402 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1403 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1404 fGdmlE->NewAttr(mainN, nullptr, "lowX", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[0]));
1405 fGdmlE->NewAttr(mainN, nullptr, "lowY", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[1]));
1406 fGdmlE->NewAttr(mainN, nullptr, "lowZ", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[2]));
1407 fGdmlE->NewAttr(mainN, nullptr, "highX", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[0]));
1408 fGdmlE->NewAttr(mainN, nullptr, "highY", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[1]));
1409 fGdmlE->NewAttr(mainN, nullptr, "highZ", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[2]));
1410
1411 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1412 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1413
1414 return mainN;
1415}
1416
1417////////////////////////////////////////////////////////////////////////////////
1418/// Creates "tube" node for GDML from object TGeoTube
1419
1421{
1422 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tube", nullptr);
1424 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1425 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1426 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1427 return nullptr;
1428 }
1429
1430 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1431 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1432 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1433 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format("%i", 0));
1434 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format("%i", 360));
1435
1436 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1437 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1438 return mainN;
1439}
1440
1441////////////////////////////////////////////////////////////////////////////////
1442/// Creates "zplane" node for GDML
1443
1445{
1446 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "zplane", nullptr);
1448
1449 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), z));
1450 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), rmin));
1451 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), rmax));
1452
1453 return mainN;
1454}
1455
1456////////////////////////////////////////////////////////////////////////////////
1457/// Creates "polycone" node for GDML
1458
1460{
1461 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "polycone", nullptr);
1463 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1464 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1465
1466 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1467 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1468
1469 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1470 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1471 Int_t nZPlns = geoShape->GetNz();
1472 for (Int_t it = 0; it < nZPlns; it++) {
1473 // add zplane child node
1474 fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1475 // compare actual plane and next plane
1476 if ((it < nZPlns - 1) && (geoShape->GetZ(it) == geoShape->GetZ(it + 1))) {
1477 // rmin of actual is greater then rmax of next one
1478 // | |rmax next
1479 // __ ...| |... __ < rmin actual
1480 // | | | |
1481 if (geoShape->GetRmin(it) > geoShape->GetRmax(it + 1)) {
1482 // adding plane from rmax next to rmin actual at the same z position
1483 if (fgG4Compatibility == kTRUE) {
1485 CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it + 1), geoShape->GetRmin(it)));
1486 Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4",
1487 lname.Data());
1488 } else {
1489 Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4",
1490 lname.Data());
1491 }
1492 }
1493 // rmin of next is greater then rmax of actual
1494 // | | | |
1495 // | |...___...| | rmin next
1496 // | | > rmax act
1497 if (geoShape->GetRmin(it + 1) > geoShape->GetRmax(it)) {
1498 // adding plane from rmax act to rmin next at the same z position
1499 if (fgG4Compatibility == kTRUE) {
1501 CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it), geoShape->GetRmin(it + 1)));
1502 Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4",
1503 lname.Data());
1504 } else {
1505 Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4",
1506 lname.Data());
1507 }
1508 }
1509 }
1510 }
1511 return mainN;
1512}
1513
1514////////////////////////////////////////////////////////////////////////////////
1515/// Creates "torus" node for GDML
1516
1518{
1519 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "torus", nullptr);
1521 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1522 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1523 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
1524 return nullptr;
1525 }
1526
1527 fGdmlE->NewAttr(mainN, nullptr, "rtor", TString::Format(fltPrecision.Data(), geoShape->GetR()));
1528 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1529 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1530 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1531 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1532
1533 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1534 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1535
1536 return mainN;
1537}
1538
1539////////////////////////////////////////////////////////////////////////////////
1540/// Creates "polyhedra" node for GDML
1541
1543{
1544 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "polyhedra", nullptr);
1546 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1547
1548 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1549 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1550 fGdmlE->NewAttr(mainN, nullptr, "numsides", TString::Format("%i", geoShape->GetNedges()));
1551
1552 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1553 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1554 for (Int_t it = 0; it < geoShape->GetNz(); it++) {
1555 // add zplane child node
1556 fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1557 }
1558 return mainN;
1559}
1560
1561////////////////////////////////////////////////////////////////////////////////
1562/// Creates "eltube" node for GDML
1563
1565{
1566 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "eltube", nullptr);
1568 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1569 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1570 if (IsNullParam(geoShape->GetA(), "A", lname) || IsNullParam(geoShape->GetB(), "B", lname) ||
1571 IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1572 return nullptr;
1573 }
1574
1575 fGdmlE->NewAttr(mainN, nullptr, "dx", TString::Format(fltPrecision.Data(), geoShape->GetA()));
1576 fGdmlE->NewAttr(mainN, nullptr, "dy", TString::Format(fltPrecision.Data(), geoShape->GetB()));
1577 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1578
1579 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1580
1581 return mainN;
1582}
1583
1584////////////////////////////////////////////////////////////////////////////////
1585/// Creates "hype" node for GDML
1586
1588{
1589 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "hype", nullptr);
1591 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1592 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1593 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1594 return nullptr;
1595 }
1596
1597 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1598 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1599 fGdmlE->NewAttr(mainN, nullptr, "inst", TString::Format(fltPrecision.Data(), geoShape->GetStIn()));
1600 fGdmlE->NewAttr(mainN, nullptr, "outst", TString::Format(fltPrecision.Data(), geoShape->GetStOut()));
1601 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1602
1603 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1604 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1605
1606 return mainN;
1607}
1608
1609////////////////////////////////////////////////////////////////////////////////
1610/// Creates "xtru" node for GDML
1611
1613{
1614 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "xtru", nullptr);
1616 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1617 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1618
1619 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1621 Int_t vertNum = geoShape->GetNvert();
1622 Int_t secNum = geoShape->GetNz();
1623 if (vertNum < 3 || secNum < 2) {
1624 Info("CreateXtrusionN", "ERROR! TGeoXtru %s has only %i vertices and %i sections. It was not exported",
1625 lname.Data(), vertNum, secNum);
1626 mainN = nullptr;
1627 return mainN;
1628 }
1629 for (Int_t it = 0; it < vertNum; it++) {
1630 // add twoDimVertex child node
1631 childN = fGdmlE->NewChild(nullptr, nullptr, "twoDimVertex", nullptr);
1632 fGdmlE->NewAttr(childN, nullptr, "x", TString::Format(fltPrecision.Data(), geoShape->GetX(it)));
1633 fGdmlE->NewAttr(childN, nullptr, "y", TString::Format(fltPrecision.Data(), geoShape->GetY(it)));
1635 }
1636 for (Int_t it = 0; it < secNum; it++) {
1637 // add section child node
1638 childN = fGdmlE->NewChild(nullptr, nullptr, "section", nullptr);
1639 fGdmlE->NewAttr(childN, nullptr, "zOrder", TString::Format("%i", it));
1640 fGdmlE->NewAttr(childN, nullptr, "zPosition", TString::Format(fltPrecision.Data(), geoShape->GetZ(it)));
1641 fGdmlE->NewAttr(childN, nullptr, "xOffset", TString::Format(fltPrecision.Data(), geoShape->GetXOffset(it)));
1642 fGdmlE->NewAttr(childN, nullptr, "yOffset", TString::Format(fltPrecision.Data(), geoShape->GetYOffset(it)));
1643 fGdmlE->NewAttr(childN, nullptr, "scalingFactor", TString::Format(fltPrecision.Data(), geoShape->GetScale(it)));
1645 }
1646 return mainN;
1647}
1648
1649////////////////////////////////////////////////////////////////////////////////
1650/// Creates "ellipsoid" node for GDML
1651/// this is a special case, because ellipsoid is not defined in ROOT
1652/// so when intersection of scaled sphere and TGeoBBox is found,
1653/// it is considered as an ellipsoid
1654
1656{
1657 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "ellipsoid", nullptr);
1659 TGeoScaledShape *leftS = (TGeoScaledShape *)geoShape->GetBoolNode()->GetLeftShape(); // ScaledShape
1660 TGeoBBox *rightS = (TGeoBBox *)geoShape->GetBoolNode()->GetRightShape(); // BBox
1661
1662 fGdmlE->NewAttr(mainN, nullptr, "name", elName.Data());
1663 Double_t sx = leftS->GetScale()->GetScale()[0];
1664 Double_t sy = leftS->GetScale()->GetScale()[1];
1665 Double_t radius = ((TGeoSphere *)leftS->GetShape())->GetRmax();
1666
1667 Double_t ax, by, cz;
1668 cz = radius;
1669 ax = sx * radius;
1670 by = sy * radius;
1671
1672 Double_t dz = rightS->GetDZ();
1673 Double_t zorig = rightS->GetOrigin()[2];
1674 Double_t zcut2 = dz + zorig;
1675 Double_t zcut1 = 2 * zorig - zcut2;
1676
1677 fGdmlE->NewAttr(mainN, nullptr, "ax", TString::Format(fltPrecision.Data(), ax));
1678 fGdmlE->NewAttr(mainN, nullptr, "by", TString::Format(fltPrecision.Data(), by));
1679 fGdmlE->NewAttr(mainN, nullptr, "cz", TString::Format(fltPrecision.Data(), cz));
1680 fGdmlE->NewAttr(mainN, nullptr, "zcut1", TString::Format(fltPrecision.Data(), zcut1));
1681 fGdmlE->NewAttr(mainN, nullptr, "zcut2", TString::Format(fltPrecision.Data(), zcut2));
1682 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1683
1684 return mainN;
1685}
1686
1687////////////////////////////////////////////////////////////////////////////////
1688/// Creates "elcone" (elliptical cone) node for GDML
1689/// this is a special case, because elliptical cone is not defined in ROOT
1690/// so when scaled cone is found, it is considered as a elliptical cone
1691
1693{
1694 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "elcone", nullptr);
1696 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1697 Double_t zcut = ((TGeoCone *)geoShape->GetShape())->GetDz();
1698 Double_t rx1 = ((TGeoCone *)geoShape->GetShape())->GetRmax1();
1699 Double_t rx2 = ((TGeoCone *)geoShape->GetShape())->GetRmax2();
1700 Double_t zmax = zcut * ((rx1 + rx2) / (rx1 - rx2));
1701 Double_t z = zcut + zmax;
1702
1703 Double_t sy = geoShape->GetScale()->GetScale()[1];
1704 Double_t ry1 = sy * rx1;
1705
1706 std::string format(TString::Format("%s/%s", fltPrecision.Data(), fltPrecision.Data()).Data());
1707 fGdmlE->NewAttr(mainN, nullptr, "dx", TString::Format(format.c_str(), rx1, z));
1708 fGdmlE->NewAttr(mainN, nullptr, "dy", TString::Format(format.c_str(), ry1, z));
1709 fGdmlE->NewAttr(mainN, nullptr, "zmax", TString::Format(fltPrecision.Data(), zmax));
1710 fGdmlE->NewAttr(mainN, nullptr, "zcut", TString::Format(fltPrecision.Data(), zcut));
1711 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1712
1713 return mainN;
1714}
1715
1716////////////////////////////////////////////////////////////////////////////////
1717/// Creates "tessellated" (tessellated shape) node for GDML
1718
1720{
1721 // add all vertices to the define section
1723 for (int i = 0; i < geoShape->GetNvertices(); ++i) {
1724 auto vertex = geoShape->GetVertex(i);
1725 TString posName = TString::Format("%s_%d", genname.Data(), i);
1726 Xyz nodPos;
1727 nodPos.x = vertex[0];
1728 nodPos.y = vertex[1];
1729 nodPos.z = vertex[2];
1730 auto childN = CreatePositionN(posName.Data(), nodPos, "position", fDefault_lunit);
1731 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
1732 }
1733 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tessellated", nullptr);
1734 fGdmlE->NewAttr(mainN, nullptr, "name", genname.Data());
1735 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1736
1738 for (Int_t it = 0; it < geoShape->GetNfacets(); it++) {
1739 // add section child node
1740 auto facet = geoShape->GetFacet(it);
1741 bool triangular = facet.GetNvert() == 3;
1742 TString ntype = (triangular) ? "triangular" : "quadrangular";
1743 childN = fGdmlE->NewChild(nullptr, nullptr, ntype.Data(), nullptr);
1744 fGdmlE->NewAttr(childN, nullptr, "vertex1", TString::Format("%s_%d", genname.Data(), facet[0]));
1745 fGdmlE->NewAttr(childN, nullptr, "vertex2", TString::Format("%s_%d", genname.Data(), facet[1]));
1746 fGdmlE->NewAttr(childN, nullptr, "vertex3", TString::Format("%s_%d", genname.Data(), facet[2]));
1747 if (!triangular)
1748 fGdmlE->NewAttr(childN, nullptr, "vertex4", TString::Format("%s_%d", genname.Data(), facet[3]));
1749 fGdmlE->NewAttr(childN, nullptr, "type", "ABSOLUTE");
1751 }
1752 return mainN;
1753}
1754
1755////////////////////////////////////////////////////////////////////////////////
1756/// Creates a scaled node for GDML
1757
1759{
1763 auto const scale = geoShape->GetScale()->GetScale();
1764 TGeoShape *unscaled = geoShape->GetShape();
1765 // Create the unscaled shape record
1767 // retrieve node name by their pointer to make reference
1769
1770 // the unplaced node appended to main structure of nodes (if they are not already there)
1771 if (unscaledN != nullptr) {
1773 fSolCnt++;
1774 } else {
1775 if (uname.Contains("missing_") || uname == "") {
1776 Info("CreateScaledN", "ERROR! Unscaled node is NULL - Scaled shape will be skipped");
1777 return nullptr;
1778 }
1779 }
1780
1781 // create union node and its child nodes (or intersection or subtraction)
1782 /* <scaledSolid name="...">
1783 * <solidref ref="unscaled solid name" />
1784 * <scale name="scale name" x="sx" y="sy" z="sz">
1785 * </scaledSolid>
1786 */
1787 mainN = fGdmlE->NewChild(nullptr, nullptr, "scaledSolid", nullptr);
1788 fGdmlE->NewAttr(mainN, nullptr, "name", nodeName);
1789
1790 //<solidref> (solid)
1791 childN = fGdmlE->NewChild(nullptr, nullptr, "solidref", nullptr);
1792 fGdmlE->NewAttr(childN, nullptr, "ref", uname);
1794
1795 //<scale ame="scale name" x="sx" y="sy" z="sz">
1796 childN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
1797 fGdmlE->NewAttr(childN, nullptr, "name", (nodeName + "scl").Data());
1798 fGdmlE->NewAttr(childN, nullptr, "x", TString::Format(fltPrecision.Data(), scale[0]));
1799 fGdmlE->NewAttr(childN, nullptr, "y", TString::Format(fltPrecision.Data(), scale[1]));
1800 fGdmlE->NewAttr(childN, nullptr, "z", TString::Format(fltPrecision.Data(), scale[2]));
1802
1803 return mainN;
1804}
1805
1806////////////////////////////////////////////////////////////////////////////////
1807/// Creates common part of union intersection and subtraction nodes
1808
1810{
1814 TGeoBoolNode::EGeoBoolType boolType = geoShape->GetBoolNode()->GetBooleanOperator();
1815 switch (boolType) {
1816 case TGeoBoolNode::kGeoUnion: lboolType = "union"; break;
1817 case TGeoBoolNode::kGeoSubtraction: lboolType = "subtraction"; break;
1818 case TGeoBoolNode::kGeoIntersection: lboolType = "intersection"; break;
1819 }
1820
1821 TGDMLWrite::Xyz lrot = GetXYZangles(geoShape->GetBoolNode()->GetLeftMatrix()->Inverse().GetRotationMatrix());
1822 const Double_t *ltr = geoShape->GetBoolNode()->GetLeftMatrix()->GetTranslation();
1823 TGDMLWrite::Xyz rrot = GetXYZangles(geoShape->GetBoolNode()->GetRightMatrix()->Inverse().GetRotationMatrix());
1824 const Double_t *rtr = geoShape->GetBoolNode()->GetRightMatrix()->GetTranslation();
1825
1826 // specific case!
1827 // Ellipsoid tag preparing
1828 // if left == TGeoScaledShape AND right == TGeoBBox
1829 // AND if TGeoScaledShape->GetShape == TGeoSphere
1830 TGeoShape *leftS = geoShape->GetBoolNode()->GetLeftShape();
1831 TGeoShape *rightS = geoShape->GetBoolNode()->GetRightShape();
1832 if (strcmp(leftS->ClassName(), "TGeoScaledShape") == 0 && strcmp(rightS->ClassName(), "TGeoBBox") == 0) {
1833 if (strcmp(((TGeoScaledShape *)leftS)->GetShape()->ClassName(), "TGeoSphere") == 0) {
1834 if (lboolType == "intersection") {
1836 return mainN;
1837 }
1838 }
1839 }
1840
1842 // translation
1843 translL.x = ltr[0];
1844 translL.y = ltr[1];
1845 translL.z = ltr[2];
1846 translR.x = rtr[0];
1847 translR.y = rtr[1];
1848 translR.z = rtr[2];
1849
1850 // left and right nodes are created here also their names are created
1851 ndL = ChooseObject(geoShape->GetBoolNode()->GetLeftShape());
1852 ndR = ChooseObject(geoShape->GetBoolNode()->GetRightShape());
1853
1854 // retrieve left and right node names by their pointer to make reference
1855 TString lname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetLeftShape())];
1856 TString rname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetRightShape())];
1857
1858 // left and right nodes appended to main structure of nodes (if they are not already there)
1859 if (ndL != nullptr) {
1861 fSolCnt++;
1862 } else {
1863 if (lname.Contains("missing_") || lname == "") {
1864 Info("CreateCommonBoolN", "ERROR! Left node is NULL - Boolean Shape will be skipped");
1865 return nullptr;
1866 }
1867 }
1868 if (ndR != nullptr) {
1870 fSolCnt++;
1871 } else {
1872 if (rname.Contains("missing_") || rname == "") {
1873 Info("CreateCommonBoolN", "ERROR! Right node is NULL - Boolean Shape will be skipped");
1874 return nullptr;
1875 }
1876 }
1877
1878 // create union node and its child nodes (or intersection or subtraction)
1879 /* <union name="...">
1880 * <first ref="left name" />
1881 * <second ref="right name" />
1882 * <firstposition .../>
1883 * <firstrotation .../>
1884 * <position .../>
1885 * <rotation .../>
1886 * </union>
1887 */
1888 mainN = fGdmlE->NewChild(nullptr, nullptr, lboolType.Data(), nullptr);
1889 fGdmlE->NewAttr(mainN, nullptr, "name", nodeName);
1890
1891 //<first> (left)
1892 childN = fGdmlE->NewChild(nullptr, nullptr, "first", nullptr);
1893 fGdmlE->NewAttr(childN, nullptr, "ref", lname);
1895
1896 //<second> (right)
1897 childN = fGdmlE->NewChild(nullptr, nullptr, "second", nullptr);
1898 fGdmlE->NewAttr(childN, nullptr, "ref", rname);
1900
1901 //<firstposition> (left)
1902 if ((translL.x != 0.0) || (translL.y != 0.0) || (translL.z != 0.0)) {
1903 childN = CreatePositionN((nodeName + lname + "pos").Data(), translL, "firstposition", fDefault_lunit);
1905 }
1906 //<firstrotation> (left)
1907 if ((lrot.x != 0.0) || (lrot.y != 0.0) || (lrot.z != 0.0)) {
1908 childN = CreateRotationN((nodeName + lname + "rot").Data(), lrot, "firstrotation");
1910 }
1911 //<position> (right)
1912 if ((translR.x != 0.0) || (translR.y != 0.0) || (translR.z != 0.0)) {
1913 childN = CreatePositionN((nodeName + rname + "pos").Data(), translR, "position", fDefault_lunit);
1915 }
1916 //<rotation> (right)
1917 if ((rrot.x != 0.0) || (rrot.y != 0.0) || (rrot.z != 0.0)) {
1918 childN = CreateRotationN((nodeName + rname + "rot").Data(), rrot, "rotation");
1920 }
1921
1922 return mainN;
1923}
1924
1925////////////////////////////////////////////////////////////////////////////////
1926/// Creates "opticalsurface" node for GDML
1927
1929{
1930 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "opticalsurface", nullptr);
1932 std::string name = make_NCName(geoSurf->GetName());
1933 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1934 fGdmlE->NewAttr(mainN, nullptr, "model", TGeoOpticalSurface::ModelToString(geoSurf->GetModel()));
1935 fGdmlE->NewAttr(mainN, nullptr, "finish", TGeoOpticalSurface::FinishToString(geoSurf->GetFinish()));
1936 fGdmlE->NewAttr(mainN, nullptr, "type", TGeoOpticalSurface::TypeToString(geoSurf->GetType()));
1937 fGdmlE->NewAttr(mainN, nullptr, "value", TString::Format(fltPrecision.Data(), geoSurf->GetValue()));
1938
1939 // Write properties
1940 TList const &properties = geoSurf->GetProperties();
1941 if (properties.GetSize()) {
1942 TIter next(&properties);
1944 while ((property = (TNamed *)next()))
1946 }
1947 // Write CONST properties
1948 TList const &const_properties = geoSurf->GetConstProperties();
1949 if (const_properties.GetSize()) {
1950 TIter next(&const_properties);
1952 while ((property = (TNamed *)next()))
1954 }
1955 return mainN;
1956}
1957
1958////////////////////////////////////////////////////////////////////////////////
1959/// Creates "skinsurface" node for GDML
1960
1962{
1963 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "skinsurface", nullptr);
1964 std::string name = make_NCName(geoSurf->GetName());
1965 std::string prop = make_NCName(geoSurf->GetTitle());
1966
1967 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1968 fGdmlE->NewAttr(mainN, nullptr, "surfaceproperty", prop.c_str());
1969
1970 // Cretate the logical volume reference node
1971 XMLNodePointer_t childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
1972 const TString &volname = fNameList->fLst[TString::Format("%p", geoSurf->GetVolume())];
1973 fGdmlE->NewAttr(childN, nullptr, "ref", volname.Data());
1975 return mainN;
1976}
1977
1978////////////////////////////////////////////////////////////////////////////////
1979/// Creates "bordersurface" node for GDML
1980
1982{
1983 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "bordersurface", nullptr);
1984 std::string name = make_NCName(geoSurf->GetName());
1985 std::string prop = make_NCName(geoSurf->GetTitle());
1986
1987 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1988 fGdmlE->NewAttr(mainN, nullptr, "surfaceproperty", prop.c_str());
1989
1990 // Cretate the logical volume reference nodes
1991 XMLNodePointer_t childN = fGdmlE->NewChild(nullptr, nullptr, "physvolref", nullptr);
1993 fGdmlE->NewAttr(childN, nullptr, "ref", physvolname);
1995
1996 childN = fGdmlE->NewChild(nullptr, nullptr, "physvolref", nullptr);
1997 physvolname = fNameList->fLst[TString::Format("%p", geoSurf->GetNode2())];
1998 fGdmlE->NewAttr(childN, nullptr, "ref", physvolname);
2000 return mainN;
2001}
2002
2003////////////////////////////////////////////////////////////////////////////////
2004/// Creates "position" kind of node for GDML
2005
2006XMLNodePointer_t TGDMLWrite::CreatePositionN(const char *name, Xyz position, const char *type, const char *unit)
2007{
2008 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, type, nullptr);
2010 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2011 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), position.x));
2012 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), position.y));
2013 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), position.z));
2014 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2015 return mainN;
2016}
2017
2018////////////////////////////////////////////////////////////////////////////////
2019/// Creates "rotation" kind of node for GDML
2020
2021XMLNodePointer_t TGDMLWrite::CreateRotationN(const char *name, Xyz rotation, const char *type, const char *unit)
2022{
2023 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, type, nullptr);
2025 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2026 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), rotation.x));
2027 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), rotation.y));
2028 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), rotation.z));
2029 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2030 return mainN;
2031}
2032
2033////////////////////////////////////////////////////////////////////////////////
2034/// Creates "matrix" kind of node for GDML
2035
2037{
2038 std::stringstream vals;
2039 size_t cols = matrix->GetCols();
2040 size_t rows = matrix->GetRows();
2041 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "matrix", nullptr);
2042 fGdmlE->NewAttr(mainN, nullptr, "name", matrix->GetName());
2043 fGdmlE->NewAttr(mainN, nullptr, "coldim", TString::Format("%zu", cols));
2044 for (size_t i = 0; i < rows; ++i) {
2045 for (size_t j = 0; j < cols; ++j) {
2046 vals << matrix->Get(i, j);
2047 if (j < cols - 1)
2048 vals << ' ';
2049 }
2050 if (i < rows - 1)
2051 vals << '\n';
2052 }
2053 fGdmlE->NewAttr(mainN, nullptr, "values", vals.str().c_str());
2054 return mainN;
2055}
2056
2057////////////////////////////////////////////////////////////////////////////////
2058/// Creates "constant" kind of node for GDML
2059
2061{
2062 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "matrix", nullptr);
2064 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2065 fGdmlE->NewAttr(mainN, nullptr, "coldim", "1");
2066 fGdmlE->NewAttr(mainN, nullptr, "values", TString::Format(fltPrecision.Data(), value));
2067 return mainN;
2068}
2069
2070////////////////////////////////////////////////////////////////////////////////
2071/// Creates "setup" node for GDML
2072
2074{
2075 XMLNodePointer_t setupN = fGdmlE->NewChild(nullptr, nullptr, "setup", nullptr);
2076 fGdmlE->NewAttr(setupN, nullptr, "name", name);
2077 fGdmlE->NewAttr(setupN, nullptr, "version", version);
2078 XMLNodePointer_t fworldN = fGdmlE->NewChild(setupN, nullptr, "world", nullptr);
2079 fGdmlE->NewAttr(fworldN, nullptr, "ref", topVolName);
2080 return setupN;
2081}
2082
2083////////////////////////////////////////////////////////////////////////////////
2084/// Creates "volume" node for GDML
2085
2086XMLNodePointer_t TGDMLWrite::StartVolumeN(const char *name, const char *solid, const char *material)
2087{
2089 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "volume", nullptr);
2090 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2091
2092 childN = fGdmlE->NewChild(nullptr, nullptr, "materialref", nullptr);
2093 fGdmlE->NewAttr(childN, nullptr, "ref", material);
2095
2096 childN = fGdmlE->NewChild(nullptr, nullptr, "solidref", nullptr);
2097 fGdmlE->NewAttr(childN, nullptr, "ref", solid);
2099
2100 return mainN;
2101}
2102
2103////////////////////////////////////////////////////////////////////////////////
2104/// Creates "assembly" node for GDML
2105
2107{
2108 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "assembly", nullptr);
2109 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2110
2111 return mainN;
2112}
2113
2114////////////////////////////////////////////////////////////////////////////////
2115/// Creates "physvol" node for GDML
2116
2118 const char *rotref, XMLNodePointer_t scaleN)
2119{
2120 fPhysVolCnt++;
2122 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "physvol", nullptr);
2123 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2124 fGdmlE->NewAttr(mainN, nullptr, "copynumber", TString::Format("%d", copyno));
2125
2126 childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
2127 fGdmlE->NewAttr(childN, nullptr, "ref", volref);
2129
2130 childN = fGdmlE->NewChild(nullptr, nullptr, "positionref", nullptr);
2131 fGdmlE->NewAttr(childN, nullptr, "ref", posref);
2133
2134 // if is not empty string add this node
2135 if (strcmp(rotref, "") != 0) {
2136 childN = fGdmlE->NewChild(nullptr, nullptr, "rotationref", nullptr);
2137 fGdmlE->NewAttr(childN, nullptr, "ref", rotref);
2139 }
2140 if (scaleN != nullptr) {
2142 }
2143
2144 return mainN;
2145}
2146
2147////////////////////////////////////////////////////////////////////////////////
2148/// Creates "divisionvol" node for GDML
2149
2151 const char *unit, const char *volref)
2152{
2153 XMLNodePointer_t childN = nullptr;
2154 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "divisionvol", nullptr);
2155 fGdmlE->NewAttr(mainN, nullptr, "axis", axis);
2156 fGdmlE->NewAttr(mainN, nullptr, "number", TString::Format("%i", number));
2158 if (fgG4Compatibility == kTRUE) {
2159 // if eg. full length is 20 and width * number = 20,0001 problem in geant4
2160 // unit is either in cm or degrees nothing else
2161 width = (floor(width * 1E4)) * 1E-4;
2162 if ((offset >= 0.) && (strcmp(axis, "kPhi") == 0)) {
2165 // put to range from 0 to 360 add decimals and then put to range 0 -> -360
2166 offset = (offsetI % 360) + decimals - 360;
2167 }
2168 }
2169 fGdmlE->NewAttr(mainN, nullptr, "width", TString::Format(fltPrecision.Data(), width));
2170
2171 fGdmlE->NewAttr(mainN, nullptr, "offset", TString::Format(fltPrecision.Data(), offset));
2172 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2173 if (strcmp(volref, "") != 0) {
2174 childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
2175 fGdmlE->NewAttr(childN, nullptr, "ref", volref);
2176 }
2178
2179 return mainN;
2180}
2181
2182////////////////////////////////////////////////////////////////////////////////
2183/// Chooses the object and method that should be used for processing object
2184
2186{
2187 const char *clsname = geoShape->ClassName();
2189
2190 if (CanProcess((TObject *)geoShape) == kFALSE) {
2191 return nullptr;
2192 }
2193
2194 // process different shapes
2195 if (strcmp(clsname, "TGeoBBox") == 0) {
2197 } else if (strcmp(clsname, "TGeoParaboloid") == 0) {
2199 } else if (strcmp(clsname, "TGeoSphere") == 0) {
2201 } else if (strcmp(clsname, "TGeoArb8") == 0) {
2203 } else if (strcmp(clsname, "TGeoConeSeg") == 0) {
2205 } else if (strcmp(clsname, "TGeoCone") == 0) {
2207 } else if (strcmp(clsname, "TGeoPara") == 0) {
2209 } else if (strcmp(clsname, "TGeoTrap") == 0) {
2211 } else if (strcmp(clsname, "TGeoGtra") == 0) {
2213 } else if (strcmp(clsname, "TGeoTrd1") == 0) {
2215 } else if (strcmp(clsname, "TGeoTrd2") == 0) {
2217 } else if (strcmp(clsname, "TGeoTubeSeg") == 0) {
2219 } else if (strcmp(clsname, "TGeoCtub") == 0) {
2221 } else if (strcmp(clsname, "TGeoTube") == 0) {
2223 } else if (strcmp(clsname, "TGeoPcon") == 0) {
2225 } else if (strcmp(clsname, "TGeoTorus") == 0) {
2227 } else if (strcmp(clsname, "TGeoPgon") == 0) {
2229 } else if (strcmp(clsname, "TGeoEltu") == 0) {
2231 } else if (strcmp(clsname, "TGeoHype") == 0) {
2233 } else if (strcmp(clsname, "TGeoXtru") == 0) {
2235 } else if (strcmp(clsname, "TGeoTessellated") == 0) {
2237 } else if (strcmp(clsname, "TGeoScaledShape") == 0) {
2239 } else if (strcmp(clsname, "TGeoCompositeShape") == 0) {
2241 } else if (strcmp(clsname, "TGeoUnion") == 0) {
2243 } else if (strcmp(clsname, "TGeoIntersection") == 0) {
2245 } else if (strcmp(clsname, "TGeoSubtraction") == 0) {
2247 } else {
2248 Info("ChooseObject", "ERROR! %s Solid CANNOT be processed, solid is NOT supported", clsname);
2249 solidN = nullptr;
2250 }
2251 if (solidN == nullptr) {
2252 if (fNameList->fLst[TString::Format("%p", geoShape)] == "") {
2253 TString missingName = geoShape->GetName();
2254 GenName("missing_" + missingName, TString::Format("%p", geoShape));
2255 } else {
2257 "missing_" + fNameList->fLst[TString::Format("%p", geoShape)];
2258 }
2259 }
2260
2261 return solidN;
2262}
2263
2264////////////////////////////////////////////////////////////////////////////////
2265/// Retrieves X Y Z angles from rotation matrix
2266
2268{
2270 Double_t a, b, c;
2271 Double_t rad = 180.0 / TMath::ACos(-1.0);
2272 const Double_t *r = rotationMatrix;
2273 Double_t cosb = TMath::Sqrt(r[0] * r[0] + r[1] * r[1]);
2274 if (cosb > 0.00001) {
2275 a = TMath::ATan2(r[5], r[8]) * rad;
2276 b = TMath::ATan2(-r[2], cosb) * rad;
2277 c = TMath::ATan2(r[1], r[0]) * rad;
2278 } else {
2279 a = TMath::ATan2(-r[7], r[4]) * rad;
2280 b = TMath::ATan2(-r[2], cosb) * rad;
2281 c = 0;
2282 }
2283 lxyz.x = a;
2284 lxyz.y = b;
2285 lxyz.z = c;
2286 return lxyz;
2287}
2288
2289////////////////////////////////////////////////////////////////////////////////
2290/// Method creating cutTube as an intersection of tube and two boxes
2291/// - not used anymore because cutube is supported in Geant4 9.5
2292
2294{
2295 Double_t rmin = geoShape->GetRmin();
2296 Double_t rmax = geoShape->GetRmax();
2297 Double_t z = geoShape->GetDz();
2298 Double_t startphi = geoShape->GetPhi1();
2299 Double_t deltaphi = geoShape->GetPhi2();
2300 Double_t x1 = geoShape->GetNlow()[0];
2301 Double_t y1 = geoShape->GetNlow()[1];
2302 Double_t z1 = geoShape->GetNlow()[2];
2303 Double_t x2 = geoShape->GetNhigh()[0];
2304 Double_t y2 = geoShape->GetNhigh()[1];
2305 Double_t z2 = geoShape->GetNhigh()[2];
2306 TString xname = geoShape->GetName();
2307
2308 Double_t h0 = 2. * ((TGeoBBox *)geoShape)->GetDZ();
2309 Double_t h1 = 2 * z;
2310 Double_t h2 = 2 * z;
2311 Double_t boxdx = 1E8 * (2 * rmax) + (2 * z);
2312
2313 TGeoTubeSeg *T = new TGeoTubeSeg((xname + "T").Data(), rmin, rmax, h0, startphi, deltaphi);
2314 TGeoBBox *B1 = new TGeoBBox((xname + "B1").Data(), boxdx, boxdx, h1);
2315 TGeoBBox *B2 = new TGeoBBox((xname + "B2").Data(), boxdx, boxdx, h2);
2316
2317 // first box position parameters
2319 Double_t theta1 = 360 - TMath::ATan2(sqrt(x1 * x1 + y1 * y1), z1) * TMath::RadToDeg();
2320
2322 Double_t theta11 = TMath::ATan2(z1, sqrt(x1 * x1 + y1 * y1)) * TMath::RadToDeg();
2323
2327
2328 // second box position parameters
2330 Double_t theta2 = 360 - TMath::ATan2(sqrt(x2 * x2 + y2 * y2), z2) * TMath::RadToDeg();
2331
2333 Double_t theta21 = TMath::ATan2(z2, sqrt(x2 * x2 + y2 * y2)) * TMath::RadToDeg();
2334
2337 Double_t zpos2 = h2 * TMath::Sin((theta21)*TMath::DegToRad()) * (-1);
2338
2339 // positioning
2340 TGeoTranslation *t0 = new TGeoTranslation(0, 0, 0);
2341 TGeoTranslation *t1 = new TGeoTranslation(0 + xpos1, 0 + ypos1, 0 + (zpos1 - z));
2342 TGeoTranslation *t2 = new TGeoTranslation(0 + xpos2, 0 + ypos2, 0 + (zpos2 + z));
2343 TGeoRotation *r0 = new TGeoRotation((xname + "r0").Data());
2344 TGeoRotation *r1 = new TGeoRotation((xname + "r1").Data());
2345 TGeoRotation *r2 = new TGeoRotation((xname + "r2").Data());
2346
2347 r1->SetAngles(phi1, theta1, 0);
2348 r2->SetAngles(phi2, theta2, 0);
2349
2350 TGeoMatrix *m0 = new TGeoCombiTrans(*t0, *r0);
2351 TGeoMatrix *m1 = new TGeoCombiTrans(*t1, *r1);
2352 TGeoMatrix *m2 = new TGeoCombiTrans(*t2, *r2);
2353
2354 TGeoCompositeShape *CS1 = new TGeoCompositeShape((xname + "CS1").Data(), new TGeoIntersection(T, B1, m0, m1));
2355 TGeoCompositeShape *cs = new TGeoCompositeShape((xname + "CS").Data(), new TGeoIntersection(CS1, B2, m0, m2));
2356 delete t0;
2357 delete t1;
2358 delete t2;
2359 delete r0;
2360 delete r1;
2361 delete r2;
2362 return cs;
2363}
2364
2365////////////////////////////////////////////////////////////////////////////////
2366/// Checks whether name2check is in (NameList) list
2367
2369{
2370 Bool_t isIN = list[name2check];
2371 return isIN;
2372}
2373
2374////////////////////////////////////////////////////////////////////////////////
2375/// NCNAME basic restrictions
2376/// Replace "$" character with empty character etc.
2377
2379{
2380 TString newname = oldname.ReplaceAll("$", "");
2381 newname = newname.ReplaceAll(" ", "_");
2382 // :, @, $, %, &, /, +, ,, ;, whitespace characters or different parenthesis
2383 newname = newname.ReplaceAll(":", "");
2384 newname = newname.ReplaceAll("@", "");
2385 newname = newname.ReplaceAll("%", "");
2386 newname = newname.ReplaceAll("&", "");
2387 newname = newname.ReplaceAll("/", "");
2388 newname = newname.ReplaceAll("+", "");
2389 newname = newname.ReplaceAll(";", "");
2390 newname = newname.ReplaceAll("{", "");
2391 newname = newname.ReplaceAll("}", "");
2392 newname = newname.ReplaceAll("(", "");
2393 newname = newname.ReplaceAll(")", "");
2394 newname = newname.ReplaceAll("[", "");
2395 newname = newname.ReplaceAll("]", "");
2396 newname = newname.ReplaceAll("_refl", "");
2397 // workaround if first letter is digit than replace it to "O" (ou character)
2398 TString fstLet = newname(0, 1);
2399 if (fstLet.IsDigit()) {
2400 newname = "O" + newname(1, newname.Length());
2401 }
2402 return newname;
2403}
2404
2405////////////////////////////////////////////////////////////////////////////////
2406/// Important function which is responsible for naming volumes, solids and materials
2407
2409{
2411 if (newname != oldname) {
2412 if (fgkMaxNameErr > fActNameErr) {
2413 Info("GenName", "WARNING! Name of the object was changed because it failed to comply with NCNAME xml datatype "
2414 "restrictions.");
2415 } else if ((fgkMaxNameErr == fActNameErr)) {
2416 Info("GenName", "WARNING! Probably more names are going to be changed to comply with NCNAME xml datatype "
2417 "restriction, but it will not be displayed on the screen.");
2418 }
2419 fActNameErr++;
2420 }
2422 Int_t iter = 0;
2423 switch (fgNamingSpeed) {
2424 case kfastButUglySufix: newname = newname + "0x" + objPointer; break;
2425 case kelegantButSlow:
2426 // 0 means not in the list
2427 iter = fNameList->fLstIter[newname];
2428 if (iter == 0) {
2429 nameIter = "";
2430 } else {
2431 nameIter = TString::Format("0x%i", iter);
2432 }
2435 break;
2437 // no change
2438 break;
2439 }
2440 // store the name (mapped to pointer)
2442 return newname;
2443}
2444
2445////////////////////////////////////////////////////////////////////////////////
2446/// Method which tests whether solids can be processed
2447
2449{
2451 isProcessed = pointer->TestBit(fgkProcBit);
2452 pointer->SetBit(fgkProcBit, kTRUE);
2453 return !(isProcessed);
2454}
2455
2456////////////////////////////////////////////////////////////////////////////////
2457/// Method that retrieves axis and unit along which object is divided
2458
2460{
2462 unit = fDefault_lunit;
2463 switch (divAxis) {
2464 case 1:
2465 if (strcmp(pattName, "TGeoPatternX") == 0) {
2466 return "kXAxis";
2467 } else if (strcmp(pattName, "TGeoPatternCylR") == 0) {
2468 return "kRho";
2469 }
2470 break;
2471 case 2:
2472 if (strcmp(pattName, "TGeoPatternY") == 0) {
2473 return "kYAxis";
2474 } else if (strcmp(pattName, "TGeoPatternCylPhi") == 0) {
2475 unit = "deg";
2476 return "kPhi";
2477 }
2478 break;
2479 case 3:
2480 if (strcmp(pattName, "TGeoPatternZ") == 0) {
2481 return "kZAxis";
2482 }
2483 break;
2484 default: return "kUndefined"; break;
2485 }
2486 return "kUndefined";
2487}
2488
2489////////////////////////////////////////////////////////////////////////////////
2490/// Check for null parameter to skip the NULL objects
2491
2493{
2494 if (parValue == 0.) {
2495 Info("IsNullParam", "ERROR! %s is NULL due to %s = %.12g, Volume based on this shape will be skipped",
2496 objName.Data(), parName.Data(), parValue);
2497 return kTRUE;
2498 }
2499 return kFALSE;
2500}
2501
2502////////////////////////////////////////////////////////////////////////////////
2503/// Unsetting bits that were changed in gGeoManager during export so that export
2504/// can be run more times with the same instance of gGeoManager.
2505
2507{
2508 TIter next(geoMng->GetListOfVolumes());
2509 TGeoVolume *vol;
2510 while ((vol = (TGeoVolume *)next())) {
2511 ((TObject *)vol->GetShape())->SetBit(fgkProcBit, kFALSE);
2513 }
2514}
2515
2516////////////////////////////////////////////////////////////////////////////////
2517//
2518// Backwards compatibility for old DD4hep version (to be removed in the future)
2519//
2520////////////////////////////////////////////////////////////////////////////////
2521
2522////////////////////////////////////////////////////////////////////////////////
2523// Backwards compatibility (to be removed in the future): Wrapper to only selectively write one branch
2525{
2526 TList materials, volumes, nodes;
2527 MaterialExtractor extract;
2528 if (!volume) {
2529 Info("WriteGDMLfile", "Invalid Volume reference to extract GDML information!");
2530 return;
2531 }
2532 extract(volume);
2533 for (TGeoMaterial *m : extract.materials)
2534 materials.Add(m);
2535 fTopVolumeName = volume->GetName();
2536 fTopVolume = volume;
2537 fSurfaceList.clear();
2538 fVolumeList.clear();
2539 fNodeList.clear();
2540 WriteGDMLfile(geomanager, volume, &materials, filename, option);
2541 materials.Clear("nodelete");
2542 volumes.Clear("nodelete");
2543 nodes.Clear("nodelete");
2544}
2545
2546////////////////////////////////////////////////////////////////////////////////
2547/// Wrapper of all exporting methods
2548/// Creates blank GDML file and fills it with gGeoManager structure converted
2549/// to GDML structure of xml nodes
2550
2553{
2554 // option processing
2555 option.ToLower();
2556 if (option.Contains("g")) {
2558 Info("WriteGDMLfile", "Geant4 compatibility mode set");
2559 } else {
2561 }
2562 if (option.Contains("f")) {
2564 Info("WriteGDMLfile", "Fast naming convention with pointer suffix set");
2565 } else if (option.Contains("n")) {
2567 Info("WriteGDMLfile", "Naming without prefix set - be careful uniqness of name is not ensured");
2568 } else {
2570 Info("WriteGDMLfile", "Potentially slow with incremental suffix naming convention set");
2571 }
2572
2573 // local variables
2574 Int_t outputLayout = 1;
2575 const char *krootNodeName = "gdml";
2576 const char *knsRefGeneral = "http://www.w3.org/2001/XMLSchema-instance";
2577 const char *knsNameGeneral = "xsi";
2578 const char *knsRefGdml = "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd";
2579 const char *knsNameGdml = "xsi:noNamespaceSchemaLocation";
2580
2581 // First create engine
2582 fGdmlE = new TXMLEngine;
2584
2585 // create blank GDML file
2586 fGdmlFile = fGdmlE->NewDoc();
2587
2588 // create root node and add it to blank GDML file
2589 XMLNodePointer_t rootNode = fGdmlE->NewChild(nullptr, nullptr, krootNodeName, nullptr);
2591
2592 // add namespaces to root node
2595
2596 // initialize general lists and <define>, <solids>, <structure> nodes
2597 fIsotopeList = new StructLst;
2598 fElementList = new StructLst;
2599
2600 fNameList = new NameLst;
2601
2602 fDefineNode = fGdmlE->NewChild(nullptr, nullptr, "define", nullptr);
2603 fSolidsNode = fGdmlE->NewChild(nullptr, nullptr, "solids", nullptr);
2604 fStructureNode = fGdmlE->NewChild(nullptr, nullptr, "structure", nullptr);
2605 //========================
2606
2607 // initialize list of accepted patterns for divisions (in ExtractVolumes)
2608 fAccPatt = new StructLst;
2609 fAccPatt->fLst["TGeoPatternX"] = kTRUE;
2610 fAccPatt->fLst["TGeoPatternY"] = kTRUE;
2611 fAccPatt->fLst["TGeoPatternZ"] = kTRUE;
2612 fAccPatt->fLst["TGeoPatternCylR"] = kTRUE;
2613 fAccPatt->fLst["TGeoPatternCylPhi"] = kTRUE;
2614 //========================
2615
2616 // initialize list of rejected shapes for divisions (in ExtractVolumes)
2617 fRejShape = new StructLst;
2618 // this shapes are rejected because, it is not possible to divide trd2
2619 // in Y axis and while only trd2 object is imported from GDML
2620 // it causes a problem when TGeoTrd1 is divided in Y axis
2621 fRejShape->fLst["TGeoTrd1"] = kTRUE;
2622 fRejShape->fLst["TGeoTrd2"] = kTRUE;
2623 //=========================
2624
2625 // Initialize global counters
2626 fActNameErr = 0;
2627 fVolCnt = 0;
2628 fPhysVolCnt = 0;
2629 fSolCnt = 0;
2630
2631 // calling main extraction functions (with measuring time)
2632 time_t startT, endT;
2633 startT = time(nullptr);
2634 ExtractMatrices(geomanager->GetListOfGDMLMatrices());
2637
2638 Info("WriteGDMLfile", "Extracting volumes");
2639 ExtractVolumes(volume);
2640 Info("WriteGDMLfile", "%i solids added", fSolCnt);
2641 Info("WriteGDMLfile", "%i volumes added", fVolCnt);
2642 Info("WriteGDMLfile", "%i physvolumes added", fPhysVolCnt);
2643 ExtractSkinSurfaces(geomanager->GetListOfSkinSurfaces());
2644 ExtractBorderSurfaces(geomanager->GetListOfBorderSurfaces());
2645 ExtractOpticalSurfaces(geomanager->GetListOfOpticalSurfaces());
2646 endT = time(nullptr);
2647 //<gdml>
2648 fGdmlE->AddChild(rootNode, fDefineNode); // <define>...</define>
2649 fGdmlE->AddChild(rootNode, fMaterialsNode); // <materials>...</materials>
2650 fGdmlE->AddChild(rootNode, fSolidsNode); // <solids>...</solids>
2651 fGdmlE->AddChild(rootNode, fStructureNode); // <structure>...</structure>
2652 fGdmlE->AddChild(rootNode, CreateSetupN(fTopVolumeName.Data())); // <setup>...</setup>
2653 //</gdml>
2655 TString tdiffS = (tdiffI == 0 ? TString("< 1 s") : TString::Format("%.0lf s", tdiffI));
2656 Info("WriteGDMLfile", "Exporting time: %s", tdiffS.Data());
2657 //=========================
2658
2659 // Saving document
2661 Info("WriteGDMLfile", "File %s saved", filename);
2662 // cleaning
2664 // unset processing bits:
2666 delete fGdmlE;
2667}
2668
2669////////////////////////////////////////////////////////////////////////////////
2670/// Method extracting geometry structure recursively
2671
2673{
2676 TGeoPatternFinder *pattFinder = nullptr;
2679
2680 // create the name for volume/assembly
2681 if (volume == fTopVolume) {
2682 // not needed a special function for generating name
2683 volname = volume->GetName();
2685 // register name to the pointer
2686 fNameList->fLst[TString::Format("%p", volume)] = volname;
2687 } else {
2688 volname = GenName(volume->GetName(), TString::Format("%p", volume));
2689 }
2690
2691 // start to create main volume/assembly node
2692 if (volume->IsAssembly()) {
2694 } else {
2695 // get reference material and add solid to <solids> + get name
2696 matname = fNameList->fLst[TString::Format("%p", volume->GetMaterial())];
2697 solname = ExtractSolid(volume->GetShape());
2698 // If solid is not supported or corrupted
2699 if (solname == "-1") {
2700 Info("ExtractVolumes", "ERROR! %s volume was not added, because solid is either not supported or corrupted",
2701 volname.Data());
2702 // set volume as missing volume
2703 fNameList->fLst[TString::Format("%p", volume)] = "missing_" + volname;
2704 return;
2705 }
2707
2708 // divisionvol can't be in assembly
2709 pattFinder = volume->GetFinder();
2710 // if found pattern
2711 if (pattFinder) {
2712 pattClsName = TString::Format("%s", pattFinder->ClassName());
2713 TString shapeCls = TString::Format("%s", volume->GetShape()->ClassName());
2714 // if pattern in accepted pattern list and not in shape rejected list
2715 if ((fAccPatt->fLst[pattClsName] == kTRUE) && (fRejShape->fLst[shapeCls] != kTRUE)) {
2716 isPattern = kTRUE;
2717 }
2718 }
2719 }
2720 // get all nodes in volume
2721 TObjArray *nodeLst = volume->GetNodes();
2722 TIter next(nodeLst);
2725 Int_t nCnt = 0;
2726 // loop through all nodes
2727 while ((geoNode = (TGeoNode *)next())) {
2728 // get volume of current node and if not processed then process it
2729 TGeoVolume *subvol = geoNode->GetVolume();
2730 if (subvol->TestAttBit(fgkProcBitVol) == kFALSE) {
2731 subvol->SetAttBit(fgkProcBitVol);
2733 }
2734
2735 // volume of this node has to exist because it was processed recursively
2736 TString nodevolname = fNameList->fLst[TString::Format("%p", geoNode->GetVolume())];
2737 if (nodevolname.Contains("missing_")) {
2738 continue;
2739 }
2740 if (nCnt == 0) { // save name of the first node for divisionvol
2742 }
2743
2744 if (isPattern == kFALSE) {
2745 // create name for node
2746 TString nodename, posname, rotname;
2747 nodename = GenName(geoNode->GetName(), TString::Format("%p", geoNode));
2748 nodename = nodename + "in" + volname;
2749
2750 // create name for position and clear rotation
2751 posname = nodename + "pos";
2752 rotname = "";
2753
2754 // position
2755 const Double_t *pos = geoNode->GetMatrix()->GetTranslation();
2756 Xyz nodPos;
2757 nodPos.x = pos[0];
2758 nodPos.y = pos[1];
2759 nodPos.z = pos[2];
2760 childN = CreatePositionN(posname.Data(), nodPos, "position", fDefault_lunit);
2761 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
2762 // Deal with reflection
2763 XMLNodePointer_t scaleN = nullptr;
2764 Double_t lx, ly, lz;
2765 Double_t xangle = 0;
2766 Double_t zangle = 0;
2767 lx = geoNode->GetMatrix()->GetRotationMatrix()[0];
2768 ly = geoNode->GetMatrix()->GetRotationMatrix()[4];
2769 lz = geoNode->GetMatrix()->GetRotationMatrix()[8];
2770 if (geoNode->GetMatrix()->IsReflection() && TMath::Abs(lx) == 1 && TMath::Abs(ly) == 1 &&
2771 TMath::Abs(lz) == 1) {
2772 scaleN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
2773 fGdmlE->NewAttr(scaleN, nullptr, "name", (nodename + "scl").Data());
2774 fGdmlE->NewAttr(scaleN, nullptr, "x", TString::Format(fltPrecision.Data(), lx));
2775 fGdmlE->NewAttr(scaleN, nullptr, "y", TString::Format(fltPrecision.Data(), ly));
2776 fGdmlE->NewAttr(scaleN, nullptr, "z", TString::Format(fltPrecision.Data(), lz));
2777 // experimentally found out, that rotation should be updated like this
2778 if (lx == -1) {
2779 zangle = 180;
2780 }
2781 if (lz == -1) {
2782 xangle = 180;
2783 }
2784 }
2785
2786 // rotation
2787 TGDMLWrite::Xyz lxyz = GetXYZangles(geoNode->GetMatrix()->GetRotationMatrix());
2788 lxyz.x -= xangle;
2789 lxyz.z -= zangle;
2790 if ((lxyz.x != 0.0) || (lxyz.y != 0.0) || (lxyz.z != 0.0)) {
2791 rotname = nodename + "rot";
2792 childN = CreateRotationN(rotname.Data(), lxyz);
2793 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
2794 }
2795
2796 // create physvol for main volume/assembly node
2798 childN = CreatePhysVolN(physvolname, geoNode->GetNumber(), nodevolname.Data(), posname.Data(), rotname.Data(),
2799 scaleN);
2801 }
2802 nCnt++;
2803 }
2804 // create only one divisionvol node
2805 if (isPattern && pattFinder) {
2806 // retrieve attributes of division
2807 Int_t ndiv, divaxis;
2808 Double_t offset, width, xlo, xhi;
2809 TString axis, unit;
2810
2811 ndiv = pattFinder->GetNdiv();
2812 width = pattFinder->GetStep();
2813
2814 divaxis = pattFinder->GetDivAxis();
2815 volume->GetShape()->GetAxisRange(divaxis, xlo, xhi);
2816
2817 // compute relative start (not positional)
2818 offset = pattFinder->GetStart() - xlo;
2819 axis = GetPattAxis(divaxis, pattClsName, unit);
2820
2821 // create division node
2822 childN = CreateDivisionN(offset, width, ndiv, axis.Data(), unit.Data(), nodeVolNameBak.Data());
2824 }
2825
2826 fVolCnt++;
2827 // add volume/assembly node into the <structure> node
2829}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
#define ClassImp(name)
Definition Rtypes.h:374
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t atom
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h prop
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
Option_t Option_t width
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t property
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
TRObject operator()(const T1 &t1) const
void * XMLNodePointer_t
Definition TXMLEngine.h:17
const_iterator begin() const
const_iterator end() const
This class is used in the process of reading and writing the GDML "matrix" tag.
Definition TGDMLMatrix.h:33
This class contains implementation of converting ROOT's gGeoManager geometry to GDML file.
Definition TGDMLWrite.h:56
XMLNodePointer_t fSolidsNode
Definition TGDMLWrite.h:132
StructLst * fElementList
Definition TGDMLWrite.h:110
void UnsetTemporaryBits(TGeoManager *geoMng)
Unsetting bits that were changed in gGeoManager during export so that export can be run more times wi...
std::map< TString, Bool_t > NameList
Definition TGDMLWrite.h:96
TString fTopVolumeName
Definition TGDMLWrite.h:126
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.
static const UInt_t fgkProcBitVol
Definition TGDMLWrite.h:141
XMLNodePointer_t CreatePolyconeN(TGeoPcon *geoShape)
Creates "polycone" node for GDML.
XMLNodePointer_t CreateFractionN(Double_t percentage, const char *refName)
Creates "fraction" node for GDML.
void ExtractMatrices(TObjArray *matrices)
Method exporting GDML matrices.
XMLDocPointer_t fGdmlFile
Definition TGDMLWrite.h:124
SurfaceList fSurfaceList
Definition TGDMLWrite.h:113
void SetG4Compatibility(Bool_t G4Compatible)
Definition TGDMLWrite.h:84
XMLNodePointer_t CreateParaboloidN(TGeoParaboloid *geoShape)
Creates "paraboloid" node for GDML.
TGeoCompositeShape * CreateFakeCtub(TGeoCtub *geoShape)
Method creating cutTube as an intersection of tube and two boxes.
TGDMLWrite()
Default constructor.
Int_t fIgnoreDummyMaterial
Definition TGDMLWrite.h:122
XMLNodePointer_t CreateBoxN(TGeoBBox *geoShape)
Creates "box" node for GDML.
std::map< TString, Float_t > NameListF
Definition TGDMLWrite.h:99
XMLNodePointer_t CreateMaterialN(TGeoMaterial *material, TString mname)
Creates "material" node for GDML.
VolList fVolumeList
Definition TGDMLWrite.h:114
XMLNodePointer_t CreateSphereN(TGeoSphere *geoShape)
Creates "sphere" node for GDML.
XMLNodePointer_t CreateTwistedTrapN(TGeoGtra *geoShape)
Creates "twistedtrap" node for GDML.
TXMLEngine * fGdmlE
Definition TGDMLWrite.h:128
XMLNodePointer_t CreateZplaneN(Double_t z, Double_t rmin, Double_t rmax)
Creates "zplane" node for GDML.
static const UInt_t fgkMaxNameErr
Definition TGDMLWrite.h:142
Bool_t IsNullParam(Double_t parValue, TString parName, TString objName)
Check for null parameter to skip the NULL objects.
static TGDMLWrite * fgGDMLWrite
Definition TGDMLWrite.h:120
TString ExtractSolid(TGeoShape *volShape)
Method creating solid to xml file and returning its name.
XMLNodePointer_t CreateHypeN(TGeoHype *geoShape)
Creates "hype" node for GDML.
XMLNodePointer_t CreateElConeN(TGeoScaledShape *geoShape)
Creates "elcone" (elliptical cone) node for GDML this is a special case, because elliptical cone is n...
XMLNodePointer_t CreateConstantN(const char *name, Double_t value)
Creates "constant" kind of node for GDML.
XMLNodePointer_t CreateMatrixN(TGDMLMatrix const *matrix)
Creates "matrix" kind of node for GDML.
XMLNodePointer_t CreateOpticalSurfaceN(TGeoOpticalSurface *geoSurf)
Creates "opticalsurface" node for GDML.
StructLst * fRejShape
Definition TGDMLWrite.h:112
UInt_t fSolCnt
Definition TGDMLWrite.h:137
Bool_t fgG4Compatibility
Definition TGDMLWrite.h:123
XMLNodePointer_t CreateMixtureN(TGeoMixture *mixture, XMLNodePointer_t materials, TString mname)
Creates "material" node for GDML with references to other sub elements.
std::map< TString, Int_t > NameListI
Definition TGDMLWrite.h:98
void ExtractBorderSurfaces(TObjArray *surfaces)
Method exporting border surfaces.
static const UInt_t fgkProcBit
floating point precision when writing
Definition TGDMLWrite.h:140
Bool_t CanProcess(TObject *pointer)
Method which tests whether solids can be processed.
UInt_t fActNameErr
Definition TGDMLWrite.h:136
XMLNodePointer_t CreatePolyhedraN(TGeoPgon *geoShape)
Creates "polyhedra" node for GDML.
void WriteGDMLfile(TGeoManager *geomanager, const char *filename="test.gdml", TString option="")
XMLNodePointer_t CreateRotationN(const char *name, Xyz rotation, const char *type="rotation", const char *unit="deg")
Creates "rotation" kind of node for GDML.
NodeList fNodeList
Definition TGDMLWrite.h:115
@ kwithoutSufixNotUniq
Definition TGDMLWrite.h:80
@ kelegantButSlow
Definition TGDMLWrite.h:80
@ kfastButUglySufix
Definition TGDMLWrite.h:80
XMLNodePointer_t CreatePropertyN(TNamed const &property)
Creates "property" node for GDML.
XMLNodePointer_t CreateElementN(TGeoElement *element, XMLNodePointer_t materials, const char *name)
Creates "element" node for GDML element node and attribute.
XMLNodePointer_t CreateConeN(TGeoConeSeg *geoShape)
Creates "cone" node for GDML from TGeoConeSeg object.
XMLNodePointer_t CreateSkinSurfaceN(TGeoSkinSurface *geoSurf)
Creates "skinsurface" node for GDML.
XMLNodePointer_t CreateTessellatedN(TGeoTessellated *geoShape)
Creates "tessellated" (tessellated shape) node for GDML.
XMLNodePointer_t CreateCutTubeN(TGeoCtub *geoShape)
Creates "cutTube" node for GDML.
XMLNodePointer_t fStructureNode
Definition TGDMLWrite.h:133
void ExtractSkinSurfaces(TObjArray *surfaces)
Method exporting skin surfaces.
XMLNodePointer_t fDefineNode
Definition TGDMLWrite.h:130
XMLNodePointer_t CreateTorusN(TGeoTorus *geoShape)
Creates "torus" node for GDML.
Int_t fgNamingSpeed
Definition TGDMLWrite.h:121
TGeoVolume * fTopVolume
Definition TGDMLWrite.h:127
XMLNodePointer_t CreatePositionN(const char *name, Xyz position, const char *type, const char *unit)
Creates "position" kind of node for GDML.
Int_t fVolCnt
Definition TGDMLWrite.h:134
XMLNodePointer_t ExtractMaterials(TList *materialsLst)
Method exporting materials.
XMLNodePointer_t fMaterialsNode
Definition TGDMLWrite.h:131
void ExtractOpticalSurfaces(TObjArray *surfaces)
Method exporting optical surfaces.
XMLNodePointer_t CreateDN(Double_t density, const char *unit="g/cm3")
Creates "D" density node for GDML.
void ExtractConstants(TGeoManager *geom)
Method exporting GDML matrices.
XMLNodePointer_t CreateArb8N(TGeoArb8 *geoShape)
Creates "arb8" node for GDML.
void SetFltPrecision(UInt_t prec)
Definition TGDMLWrite.h:233
XMLNodePointer_t StartVolumeN(const char *name, const char *solid, const char *material)
Creates "volume" node for GDML.
void SetNamingSpeed(ENamingType naming)
Set convention of naming solids and volumes.
TString GetPattAxis(Int_t divAxis, const char *pattName, TString &unit)
Method that retrieves axis and unit along which object is divided.
TString GenName(TString oldname)
NCNAME basic restrictions Replace "$" character with empty character etc.
XMLNodePointer_t CreateXtrusionN(TGeoXtru *geoShape)
Creates "xtru" node for GDML.
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.
Bool_t IsInList(NameList list, TString name2check)
Checks whether name2check is in (NameList) list.
XMLNodePointer_t StartAssemblyN(const char *name)
Creates "assembly" node for GDML.
TString fDefault_lunit
Definition TGDMLWrite.h:125
void SetIgnoreDummyMaterial(bool value)
Ignore dummy material instance, which causes trouble reading GDML in Geant4.
Int_t fPhysVolCnt
Definition TGDMLWrite.h:135
~TGDMLWrite() override
Destructor.
XMLNodePointer_t CreateParaN(TGeoPara *geoShape)
Creates "para" node for GDML.
XMLNodePointer_t CreateSetupN(const char *topVolName, const char *name="default", const char *version="1.0")
Creates "setup" node for GDML.
XMLNodePointer_t CreateTrapN(TGeoTrap *geoShape)
Creates "trap" node for GDML.
Xyz GetXYZangles(const Double_t *rotationMatrix)
Retrieves X Y Z angles from rotation matrix.
UInt_t fFltPrecision
Definition TGDMLWrite.h:138
XMLNodePointer_t CreateTubeN(TGeoTubeSeg *geoShape)
Creates "tube" node for GDML from object TGeoTubeSeg.
XMLNodePointer_t CreateBorderSurfaceN(TGeoBorderSurface *geoSurf)
Creates "bordersurface" node for GDML.
XMLNodePointer_t CreateAtomN(Double_t atom, const char *unit="g/mole")
Creates "atom" node for GDML.
StructLst * fIsotopeList
Definition TGDMLWrite.h:109
StructLst * fAccPatt
Definition TGDMLWrite.h:111
XMLNodePointer_t CreateEltubeN(TGeoEltu *geoShape)
Creates "eltube" node for GDML.
void ExtractVolumes(TGeoNode *topNode)
Method extracting geometry structure recursively.
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...
XMLNodePointer_t CreateIsotopN(TGeoIsotope *isotope, const char *name)
Creates "isotope" node for GDML.
XMLNodePointer_t CreateCommonBoolN(TGeoCompositeShape *geoShape)
Creates common part of union intersection and subtraction nodes.
XMLNodePointer_t CreateTrdN(TGeoTrd1 *geoShape)
Creates "trd" node for GDML from object TGeoTrd1.
XMLNodePointer_t CreateScaledN(TGeoScaledShape *geoShape)
Creates a scaled node for GDML.
NameLst * fNameList
Definition TGDMLWrite.h:117
XMLNodePointer_t ChooseObject(TGeoShape *geoShape)
Chooses the object and method that should be used for processing object.
An arbitrary trapezoid with less than 8 vertices standing on two parallel planes perpendicular to Z a...
Definition TGeoArb8.h:17
void SetAttBit(UInt_t f)
Definition TGeoAtt.h:61
Box class.
Definition TGeoBBox.h:17
Class describing rotation + translation.
Definition TGeoMatrix.h:317
Composite shapes are Boolean combinations of two or more shape components.
A cone segment is a cone having a range in phi.
Definition TGeoCone.h:99
The cones are defined by 5 parameters:
Definition TGeoCone.h:17
The cut tubes constructor has the form:
Definition TGeoTube.h:173
Base class for chemical elements.
Definition TGeoElement.h:31
An elliptical tube is defined by the two semi-axes A and B.
Definition TGeoEltu.h:17
A twisted trapezoid.
Definition TGeoArb8.h:149
A hyperboloid is represented as a solid limited by two planes perpendicular to the Z axis (top and bo...
Definition TGeoHype.h:17
Boolean node representing an intersection between two components.
an isotope defined by the atomic number, number of nucleons and atomic weight (g/mole)
Definition TGeoElement.h:92
The manager class for any TGeo geometry.
Definition TGeoManager.h:44
static EDefaultUnits GetDefaultUnits()
static UInt_t GetExportPrecision()
Base class describing materials.
TList const & GetConstProperties() const
TList const & GetProperties() const
virtual Double_t GetA() const
virtual Double_t GetDensity() const
virtual Double_t GetZ() const
Geometrical transformation package.
Definition TGeoMatrix.h:38
Media are used to store properties related to tracking and which are useful only when using geometry ...
Definition TGeoMedium.h:23
Mixtures of elements.
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:39
TGeoVolume * GetVolume() const
Definition TGeoNode.h:99
This is a wrapper class to G4OpticalSurface.
static const char * ModelToString(ESurfaceModel model)
static const char * TypeToString(ESurfaceType type)
static const char * FinishToString(ESurfaceFinish finish)
Parallelepiped class.
Definition TGeoPara.h:17
A paraboloid is defined by the revolution surface generated by a parabola and is bounded by two plane...
base finder class for patterns. A pattern is specifying a division type
A polycone is represented by a sequence of tubes/cones, glued together at defined Z planes.
Definition TGeoPcon.h:17
Polygons are defined in the same way as polycones, the difference being just that the segments betwee...
Definition TGeoPgon.h:20
Class describing rotations.
Definition TGeoMatrix.h:168
A shape scaled by a TGeoScale transformation.
Base abstract class for all shapes.
Definition TGeoShape.h:25
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
TGeoSphere are not just balls having internal and external radii, but sectors of a sphere having defi...
Definition TGeoSphere.h:17
Tessellated solid class.
The torus is defined by its axial radius, its inner and outer radius.
Definition TGeoTorus.h:17
Class describing translations.
Definition TGeoMatrix.h:116
A general trapezoid.
Definition TGeoArb8.h:96
A trapezoid with only X varying with Z.
Definition TGeoTrd1.h:17
A trapezoid with only X varying with Z.
Definition TGeoTrd2.h:17
A tube segment is a tube having a range in phi.
Definition TGeoTube.h:94
Cylindrical tube class.
Definition TGeoTube.h:17
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
TGeoMaterial * GetMaterial() const
Definition TGeoVolume.h:174
TObjArray * GetNodes()
Definition TGeoVolume.h:169
TGeoPatternFinder * GetFinder() const
Definition TGeoVolume.h:177
static TGeoMedium * DummyMedium()
TGeoShape * GetShape() const
Definition TGeoVolume.h:190
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
A TGeoXtru shape is represented by the extrusion of an arbitrary polygon with fixed outline between s...
Definition TGeoXtru.h:22
A doubly linked list.
Definition TList.h:38
void Clear(Option_t *option="") override
Remove all objects from the list.
Definition TList.cxx:400
void Add(TObject *obj) override
Definition TList.h:81
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
An array of TObjects.
Definition TObjArray.h:31
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:205
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1057
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:864
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1099
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
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:2378
XMLNodePointer_t NewChild(XMLNodePointer_t parent, XMLNsPointer_t ns, const char *name, const char *content=nullptr)
create new child element for parent node
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
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...
void FreeDoc(XMLDocPointer_t xmldoc)
frees allocated document data and deletes document itself
void AddChild(XMLNodePointer_t parent, XMLNodePointer_t child)
add child element to xmlnode
XMLNsPointer_t NewNS(XMLNodePointer_t xmlnode, const char *reference, const char *name=nullptr)
create namespace attribute for xmlnode.
XMLDocPointer_t NewDoc(const char *version="1.0")
creates new xml document with provided version
void SetSkipComments(Bool_t on=kTRUE)
Definition TXMLEngine.h:48
void DocSetRootElement(XMLDocPointer_t xmldoc, XMLNodePointer_t xmlnode)
set main (root) node for document
std::ostream & Info()
Definition hadd.cxx:171
TH1F * h1
Definition legend1.C:5
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:636
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:650
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:598
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:592
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:72
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
TMarker m
Definition textangle.C:8
auto * t1
Definition textangle.C:20