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
167
169
170namespace {
171
172// Helper to replace string patterns
173std::string str_replace(const std::string &str, const std::string &pattern, const std::string &replacement)
174{
175 std::string res = str;
176 for (size_t id = res.find(pattern); id != std::string::npos; id = res.find(pattern))
177 res.replace(id, pattern.length(), replacement);
178 return res;
179}
180
181// Create a NCN compliant name (no '/' and no '#')
182std::string make_NCName(const std::string &in)
183{
184 std::string res = str_replace(in, "/", "_");
185 res = str_replace(res, "#", "_");
186 return res;
187}
188
189// Materials extractor from a volume tree
190struct MaterialExtractor {
191 std::set<TGeoMaterial *> materials;
192 void operator()(const TGeoVolume *v)
193 {
194 materials.insert(v->GetMaterial());
195 for (Int_t i = 0; i < v->GetNdaughters(); ++i)
196 (*this)(v->GetNode(i)->GetVolume());
197 }
198};
199} // namespace
200
201////////////////////////////////////////////////////////////////////////////////
202/// Default constructor.
203
205 : TObject(),
206 fIsotopeList(nullptr),
207 fElementList(nullptr),
208 fAccPatt(nullptr),
209 fRejShape(nullptr),
210 fNameList(nullptr),
211 fgNamingSpeed(0),
212 fgG4Compatibility(false),
213 fGdmlFile(nullptr),
214 fTopVolumeName(0),
215 fGdmlE(nullptr),
216 fDefineNode(nullptr),
217 fMaterialsNode(nullptr),
218 fSolidsNode(nullptr),
219 fStructureNode(nullptr),
220 fVolCnt(0),
221 fPhysVolCnt(0),
222 fActNameErr(0),
223 fSolCnt(0),
224 fFltPrecision(17) // %.17g
225{
226 if (fgGDMLWrite)
227 delete fgGDMLWrite;
228 fgGDMLWrite = this;
229}
230
231////////////////////////////////////////////////////////////////////////////////
232/// Destructor.
233
235{
236 delete fIsotopeList;
237 delete fElementList;
238 delete fAccPatt;
239 delete fRejShape;
240 delete fNameList;
241
242 fgGDMLWrite = nullptr;
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Set convention of naming solids and volumes
247
252
253////////////////////////////////////////////////////////////////////////////////
254/// Ignore dummy material instance, which causes trouble reading GDML in Geant4
255
260
261////////////////////////////////////////////////////////////////////////////////
262// wrapper of all main methods for extraction
264{
265 TList *materials = geomanager->GetListOfMaterials();
266 TGeoNode *node = geomanager->GetTopNode();
267 if (!node) {
268 Info("WriteGDMLfile", "Top volume does not exist!");
269 return;
270 }
271 fTopVolume = node->GetVolume();
273 WriteGDMLfile(geomanager, node, materials, filename, option);
274}
275
276////////////////////////////////////////////////////////////////////////////////
277// Wrapper to only selectively write one branch of the volume hierarchy to file
279{
280 TGeoVolume *volume = node->GetVolume();
281 TList materials, volumes, nodes;
282 MaterialExtractor extract;
283 if (!volume) {
284 Info("WriteGDMLfile", "Invalid Volume reference to extract GDML information!");
285 return;
286 }
287 extract(volume);
288 for (TGeoMaterial *m : extract.materials)
289 materials.Add(m);
290 fTopVolumeName = volume->GetName();
291 fTopVolume = volume;
292 fSurfaceList.clear();
293 fVolumeList.clear();
294 fNodeList.clear();
295 WriteGDMLfile(geomanager, node, &materials, filename, option);
296 materials.Clear("nodelete");
297 volumes.Clear("nodelete");
298 nodes.Clear("nodelete");
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// Wrapper of all exporting methods
303/// Creates blank GDML file and fills it with gGeoManager structure converted
304/// to GDML structure of xml nodes
305
308{
309 // option processing
310 option.ToLower();
311 if (option.Contains("g")) {
313 Info("WriteGDMLfile", "Geant4 compatibility mode set");
314 } else {
316 }
317 if (option.Contains("f")) {
319 Info("WriteGDMLfile", "Fast naming convention with pointer suffix set");
320 } else if (option.Contains("n")) {
322 Info("WriteGDMLfile", "Naming without prefix set - be careful uniqness of name is not ensured");
323 } else {
325 Info("WriteGDMLfile", "Potentially slow with incremental suffix naming convention set");
326 }
327
329
331 switch (def_units) {
332 case TGeoManager::kG4Units: fDefault_lunit = "mm"; break;
333 case TGeoManager::kRootUnits: fDefault_lunit = "cm"; break;
334 default: // G4 units
335 fDefault_lunit = "mm";
336 break;
337 }
338
339 // local variables
341 const char *krootNodeName = "gdml";
342 const char *knsRefGeneral = "http://www.w3.org/2001/XMLSchema-instance";
343 const char *knsNameGeneral = "xsi";
344 const char *knsRefGdml = "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd";
345 const char *knsNameGdml = "xsi:noNamespaceSchemaLocation";
346
347 // First create engine
348 fGdmlE = new TXMLEngine;
350
351 // create blank GDML file
353
354 // create root node and add it to blank GDML file
355 XMLNodePointer_t rootNode = fGdmlE->NewChild(nullptr, nullptr, krootNodeName, nullptr);
357
358 // add namespaces to root node
361
362 // initialize general lists and <define>, <solids>, <structure> nodes
365
366 fNameList = new NameLst;
367
368 fDefineNode = fGdmlE->NewChild(nullptr, nullptr, "define", nullptr);
369 fSolidsNode = fGdmlE->NewChild(nullptr, nullptr, "solids", nullptr);
370 fStructureNode = fGdmlE->NewChild(nullptr, nullptr, "structure", nullptr);
371 //========================
372
373 // initialize list of accepted patterns for divisions (in ExtractVolumes)
374 fAccPatt = new StructLst;
375 fAccPatt->fLst["TGeoPatternX"] = kTRUE;
376 fAccPatt->fLst["TGeoPatternY"] = kTRUE;
377 fAccPatt->fLst["TGeoPatternZ"] = kTRUE;
378 fAccPatt->fLst["TGeoPatternCylR"] = kTRUE;
379 fAccPatt->fLst["TGeoPatternCylPhi"] = kTRUE;
380 //========================
381
382 // initialize list of rejected shapes for divisions (in ExtractVolumes)
383 fRejShape = new StructLst;
384 // this shapes are rejected because, it is not possible to divide trd2
385 // in Y axis and while only trd2 object is imported from GDML
386 // it causes a problem when TGeoTrd1 is divided in Y axis
387 fRejShape->fLst["TGeoTrd1"] = kTRUE;
388 fRejShape->fLst["TGeoTrd2"] = kTRUE;
389 //=========================
390
391 // Initialize global counters
392 fActNameErr = 0;
393 fVolCnt = 0;
394 fPhysVolCnt = 0;
395 fSolCnt = 0;
396
397 // calling main extraction functions (with measuring time)
398 time_t startT, endT;
399 startT = time(nullptr);
400 ExtractMatrices(geomanager->GetListOfGDMLMatrices());
403
404 Info("WriteGDMLfile", "Extracting volumes");
405 ExtractVolumes(node);
406 Info("WriteGDMLfile", "%i solids added", fSolCnt);
407 Info("WriteGDMLfile", "%i volumes added", fVolCnt);
408 Info("WriteGDMLfile", "%i physvolumes added", fPhysVolCnt);
409 ExtractSkinSurfaces(geomanager->GetListOfSkinSurfaces());
410 ExtractBorderSurfaces(geomanager->GetListOfBorderSurfaces());
411 ExtractOpticalSurfaces(geomanager->GetListOfOpticalSurfaces());
412 endT = time(nullptr);
413 //<gdml>
414 fGdmlE->AddChild(rootNode, fDefineNode); // <define>...</define>
415 fGdmlE->AddChild(rootNode, fMaterialsNode); // <materials>...</materials>
416 fGdmlE->AddChild(rootNode, fSolidsNode); // <solids>...</solids>
417 fGdmlE->AddChild(rootNode, fStructureNode); // <structure>...</structure>
418 fGdmlE->AddChild(rootNode, CreateSetupN(fTopVolumeName.Data())); // <setup>...</setup>
419 //</gdml>
421 TString tdiffS = (tdiffI == 0 ? TString("< 1 s") : TString::Format("%.0lf s", tdiffI));
422 Info("WriteGDMLfile", "Exporting time: %s", tdiffS.Data());
423 //=========================
424
425 // Saving document
427 Info("WriteGDMLfile", "File %s saved", filename);
428 // cleaning
430 // unset processing bits:
432 delete fGdmlE;
433}
434
435////////////////////////////////////////////////////////////////////////////////
436/// Method exporting GDML matrices
437
439{
440 if (!matrixList->GetEntriesFast())
441 return;
443 TIter next(matrixList);
445 while ((matrix = (TGDMLMatrix *)next())) {
448 }
449}
450
451////////////////////////////////////////////////////////////////////////////////
452/// Method exporting GDML matrices
453
455{
456 if (!geom->GetNproperties())
457 return;
461 for (Int_t i = 0; i < geom->GetNproperties(); ++i) {
462 value = geom->GetProperty(i, property);
465 }
466}
467
468////////////////////////////////////////////////////////////////////////////////
469/// Method exporting optical surfaces
470
472{
473 if (!surfaces->GetEntriesFast())
474 return;
476 TIter next(surfaces);
478 while ((surf = (TGeoOpticalSurface *)next())) {
479 if (fSurfaceList.find(surf) == fSurfaceList.end())
480 continue;
483 // Info("ExtractSkinSurfaces", "Extracted optical surface: %s",surf->GetName());
484 }
485}
486
487////////////////////////////////////////////////////////////////////////////////
488/// Method exporting skin surfaces
489
491{
492 if (!surfaces->GetEntriesFast())
493 return;
495 TIter next(surfaces);
497 while ((surf = (TGeoSkinSurface *)next())) {
498 if (fVolumeList.find(surf->GetVolume()) == fVolumeList.end())
499 continue;
502 fSurfaceList.insert(surf->GetSurface());
503 // Info("ExtractSkinSurfaces", "Extracted skin surface: %s",surf->GetName());
504 }
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// Method exporting border surfaces
509
511{
512 if (!surfaces->GetEntriesFast())
513 return;
515 TIter next(surfaces);
517 while ((surf = (TGeoBorderSurface *)next())) {
518 auto ia = fNodeList.find(surf->GetNode1());
519 auto ib = fNodeList.find(surf->GetNode2());
520 if (ia == fNodeList.end() && ib == fNodeList.end()) {
521 continue;
522 } else if (ia == fNodeList.end() && ib != fNodeList.end()) {
523 Warning("ExtractBorderSurfaces",
524 "Inconsistent border surface extraction %s: Node %s"
525 " is not part of GDML!",
526 surf->GetName(), surf->GetNode1()->GetName());
527 continue;
528 } else if (ia != fNodeList.end() && ib == fNodeList.end()) {
529 Warning("ExtractBorderSurfaces",
530 "Inconsistent border surface extraction %s: Node %s"
531 " is not part of GDML!",
532 surf->GetName(), surf->GetNode2()->GetName());
533 continue;
534 }
537 fSurfaceList.insert(surf->GetSurface());
538 // Info("ExtractBorderSurfaces", "Extracted border surface: %s",surf->GetName());
539 }
540}
541
542////////////////////////////////////////////////////////////////////////////////
543/// Method exporting materials
544
546{
547 Info("ExtractMaterials", "Extracting materials");
548 // crate main <materials> node
549 XMLNodePointer_t materialsN = fGdmlE->NewChild(nullptr, nullptr, "materials", nullptr);
550 Int_t matcnt = 0;
551
552 // go through materials - iterator and object declaration
553 TIter next(materialsLst);
556 TGeoMaterial *dummy_mat = dummy_med ? dummy_med->GetMaterial() : nullptr;
557 std::string dummy_nam = dummy_mat ? dummy_mat->GetName() : "dummy";
558
559 while ((lmaterial = (TGeoMaterial *)next())) {
560 // check for dummy material: if requested, ignore it
561 std::string mname = lmaterial->GetName();
563 Info("ExtractMaterials", "Skip dummy material: %s", dummy_nam.c_str());
564 continue;
565 }
566 // generate uniq name
568
569 if (lmaterial->IsMixture()) {
573 } else {
576 }
577 matcnt++;
578 }
579 Info("ExtractMaterials", "%i materials added", matcnt);
580 return materialsN;
581}
582
583////////////////////////////////////////////////////////////////////////////////
584/// Method creating solid to xml file and returning its name
585
587{
589 TString solname = "";
590 solidN = ChooseObject(volShape); // volume->GetShape()
592 if (solidN != nullptr)
593 fSolCnt++;
595 if (solname.Contains("missing_")) {
596 solname = "-1";
597 }
598 return solname;
599}
600
601////////////////////////////////////////////////////////////////////////////////
602/// Method extracting geometry structure recursively
603
605{
607 TGeoVolume *volume = node->GetVolume();
609 TGeoPatternFinder *pattFinder = nullptr;
612
613 fNodeList.insert(node);
614 fVolumeList.insert(volume);
615 // create the name for volume/assembly
616 if (volume == fTopVolume) {
617 // not needed a special function for generating name
618 volname = volume->GetName();
620 // register name to the pointer
621 fNameList->fLst[TString::Format("%p", volume)] = volname;
622 } else {
623 volname = GenName(volume->GetName(), TString::Format("%p", volume));
624 }
625
626 // start to create main volume/assembly node
627 if (volume->IsAssembly()) {
629 } else {
630 // get reference material and add solid to <solids> + get name
631 matname = fNameList->fLst[TString::Format("%p", volume->GetMaterial())];
632 solname = ExtractSolid(volume->GetShape());
633 // If solid is not supported or corrupted
634 if (solname == "-1") {
635 Info("ExtractVolumes", "ERROR! %s volume was not added, because solid is either not supported or corrupted",
636 volname.Data());
637 // set volume as missing volume
638 fNameList->fLst[TString::Format("%p", volume)] = "missing_" + volname;
639 return;
640 }
642
643 // divisionvol can't be in assembly
644 pattFinder = volume->GetFinder();
645 // if found pattern
646 if (pattFinder) {
647 pattClsName = TString::Format("%s", pattFinder->ClassName());
648 TString shapeCls = TString::Format("%s", volume->GetShape()->ClassName());
649 // if pattern in accepted pattern list and not in shape rejected list
650 if ((fAccPatt->fLst[pattClsName] == kTRUE) && (fRejShape->fLst[shapeCls] != kTRUE)) {
652 }
653 }
654 }
655 // get all nodes in volume
656 TObjArray *nodeLst = volume->GetNodes();
657 TIter next(nodeLst);
660 Int_t nCnt = 0;
661 // loop through all nodes
662 while ((geoNode = (TGeoNode *)next())) {
663 // get volume of current node and if not processed then process it
664 TGeoVolume *subvol = geoNode->GetVolume();
665 fNodeList.insert(geoNode);
666 if (subvol->TestAttBit(fgkProcBitVol) == kFALSE) {
667 subvol->SetAttBit(fgkProcBitVol);
669 }
670
671 // volume of this node has to exist because it was processed recursively
673 if (nodevolname.Contains("missing_")) {
674 continue;
675 }
676 if (nCnt == 0) { // save name of the first node for divisionvol
678 }
679
680 if (isPattern == kFALSE) {
681 // create name for node
682 TString nodename, posname, rotname;
683 nodename = GenName(geoNode->GetName(), TString::Format("%p", geoNode));
684 nodename = nodename + "in" + volname;
685
686 // create name for position and clear rotation
687 posname = nodename + "pos";
688 rotname = "";
689
690 // position
691 const Double_t *pos = geoNode->GetMatrix()->GetTranslation();
692 Xyz nodPos;
693 nodPos.x = pos[0];
694 nodPos.y = pos[1];
695 nodPos.z = pos[2];
696 childN = CreatePositionN(posname.Data(), nodPos, "position", fDefault_lunit);
697 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
698 // Deal with reflection
699 XMLNodePointer_t scaleN = nullptr;
700 Double_t lx, ly, lz;
701 Double_t xangle = 0;
702 Double_t zangle = 0;
703 lx = geoNode->GetMatrix()->GetRotationMatrix()[0];
704 ly = geoNode->GetMatrix()->GetRotationMatrix()[4];
705 lz = geoNode->GetMatrix()->GetRotationMatrix()[8];
706 if (geoNode->GetMatrix()->IsReflection() && TMath::Abs(lx) == 1 && TMath::Abs(ly) == 1 &&
707 TMath::Abs(lz) == 1) {
708 scaleN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
709 fGdmlE->NewAttr(scaleN, nullptr, "name", (nodename + "scl").Data());
710 fGdmlE->NewAttr(scaleN, nullptr, "x", TString::Format(fltPrecision.Data(), lx));
711 fGdmlE->NewAttr(scaleN, nullptr, "y", TString::Format(fltPrecision.Data(), ly));
712 fGdmlE->NewAttr(scaleN, nullptr, "z", TString::Format(fltPrecision.Data(), lz));
713 // experimentally found out, that rotation should be updated like this
714 if (lx == -1) {
715 zangle = 180;
716 }
717 if (lz == -1) {
718 xangle = 180;
719 }
720 }
721
722 // rotation
723 TGDMLWrite::Xyz lxyz = GetXYZangles(geoNode->GetMatrix()->GetRotationMatrix());
724 lxyz.x -= xangle;
725 lxyz.z -= zangle;
726 if ((lxyz.x != 0.0) || (lxyz.y != 0.0) || (lxyz.z != 0.0)) {
727 rotname = nodename + "rot";
728 childN = CreateRotationN(rotname.Data(), lxyz);
729 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
730 }
731
732 // create physvol for main volume/assembly node
734 childN = CreatePhysVolN(physvolname, geoNode->GetNumber(), nodevolname.Data(), posname.Data(), rotname.Data(),
735 scaleN);
737 }
738 nCnt++;
739 }
740 // create only one divisionvol node
741 if (isPattern && pattFinder) {
742 // retrieve attributes of division
743 Int_t ndiv, divaxis;
744 Double_t offset, width, xlo, xhi;
745 TString axis, unit;
746
747 ndiv = pattFinder->GetNdiv();
748 width = pattFinder->GetStep();
749
750 divaxis = pattFinder->GetDivAxis();
751 volume->GetShape()->GetAxisRange(divaxis, xlo, xhi);
752
753 // compute relative start (not positional)
754 offset = pattFinder->GetStart() - xlo;
755 axis = GetPattAxis(divaxis, pattClsName, unit);
756
757 // create division node
758 childN = CreateDivisionN(offset, width, ndiv, axis.Data(), unit.Data(), nodeVolNameBak.Data());
760 }
761
762 fVolCnt++;
763 // add volume/assembly node into the <structure> node
765}
766
767////////////////////////////////////////////////////////////////////////////////
768/// Creates "atom" node for GDML
769
771{
773 XMLNodePointer_t atomN = fGdmlE->NewChild(nullptr, nullptr, "atom", nullptr);
774 fGdmlE->NewAttr(atomN, nullptr, "unit", unit);
775 fGdmlE->NewAttr(atomN, nullptr, "value", TString::Format(fltPrecision.Data(), atom));
776 return atomN;
777}
778
779////////////////////////////////////////////////////////////////////////////////
780/// Creates "D" density node for GDML
781
783{
785 XMLNodePointer_t densN = fGdmlE->NewChild(nullptr, nullptr, "D", nullptr);
786 fGdmlE->NewAttr(densN, nullptr, "unit", unit);
787 fGdmlE->NewAttr(densN, nullptr, "value", TString::Format(fltPrecision.Data(), density));
788 return densN;
789}
790
791////////////////////////////////////////////////////////////////////////////////
792/// Creates "fraction" node for GDML
793
795{
797 XMLNodePointer_t fractN = fGdmlE->NewChild(nullptr, nullptr, "fraction", nullptr);
799 fGdmlE->NewAttr(fractN, nullptr, "ref", refName);
800 return fractN;
801}
802
803////////////////////////////////////////////////////////////////////////////////
804/// Creates "property" node for GDML
805
807{
808 XMLNodePointer_t propertyN = fGdmlE->NewChild(nullptr, nullptr, "property", nullptr);
809 fGdmlE->NewAttr(propertyN, nullptr, "name", property.GetName());
810 fGdmlE->NewAttr(propertyN, nullptr, "ref", property.GetTitle());
811 return propertyN;
812}
813
814////////////////////////////////////////////////////////////////////////////////
815/// Creates "isotope" node for GDML
816
818{
819 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "isotope", nullptr);
820 fGdmlE->NewAttr(mainN, nullptr, "name", name);
821 fGdmlE->NewAttr(mainN, nullptr, "N", TString::Format("%i", isotope->GetN()));
822 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format("%i", isotope->GetZ()));
824 return mainN;
825}
826
827////////////////////////////////////////////////////////////////////////////////
828/// Creates "element" node for GDML
829/// element node and attribute
830
832{
833 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "element", nullptr);
834 fGdmlE->NewAttr(mainN, nullptr, "name", name);
835 // local associative arrays for saving isotopes and their weight
836 // inside element
839
840 if (element->HasIsotopes()) {
841 Int_t nOfIso = element->GetNisotopes();
842 // go through isotopes
843 for (Int_t idx = 0; idx < nOfIso; idx++) {
844 TGeoIsotope *myIsotope = element->GetIsotope(idx);
845 if (!myIsotope) {
846 Fatal("CreateElementN", "Missing isotopes for element %s", element->GetName());
847 return mainN;
848 }
849
850 // Get name of the Isotope (
851 TString lname = myIsotope->GetName();
852 //_iso suffix is added to avoid problems with same names
853 // for material, element and isotopes
854 lname = TString::Format("%s_iso", lname.Data());
855
856 // cumulates abundance, in case 2 isotopes with same names
857 // within one element
858 wPercentage[lname] += element->GetRelativeAbundance(idx);
859 wCounter[lname]++;
860
861 // check whether isotope name is not in list of isotopes
863 continue;
864 }
865 // add isotope to list of isotopes and to main <materials> node
868 fGdmlE->AddChild(materials, isoNode);
869 }
870 // loop through asoc array of isotopes
871 for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
872 if (itr->second > 1) {
873 Info("CreateMixtureN", "WARNING! 2 equal isotopes in one element. Check: %s isotope of %s element",
874 itr->first.Data(), name);
875 }
876 // add fraction child to element with reference to isotope
877 fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
878 }
879 } else {
880 fGdmlE->NewAttr(mainN, nullptr, "formula", element->GetName());
881 Int_t valZ = element->Z();
882 // Z can't be <1 in Geant4 and Z is optional parameter
883 if (valZ >= 1) {
884 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format("%i", valZ));
885 }
886 Int_t valN = element->N();
887 fGdmlE->NewAttr(mainN, nullptr, "N", TString::Format("%i", valN));
889 }
890 return mainN;
891}
892
893////////////////////////////////////////////////////////////////////////////////
894/// Creates "material" node for GDML with references to other sub elements
895
897{
898 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "material", nullptr);
899 fGdmlE->NewAttr(mainN, nullptr, "name", mname);
900
901 // Write properties
902 TList const &properties = mixture->GetProperties();
903 if (properties.GetSize()) {
904 TIter next(&properties);
906 while ((property = (TNamed *)next()))
908 }
909 // Write CONST properties
910 TList const &const_properties = mixture->GetConstProperties();
911 if (const_properties.GetSize()) {
912 TIter next(&const_properties);
914 while ((property = (TNamed *)next()))
916 }
917
918 fGdmlE->AddChild(mainN, CreateDN(mixture->GetDensity()));
919 // local associative arrays for saving elements and their weight
920 // inside mixture
923
924 Int_t nOfElm = mixture->GetNelements();
925 // go through elements
926 for (Int_t idx = 0; idx < nOfElm; idx++) {
927 TGeoElement *myElement = mixture->GetElement(idx);
928
929 // Get name of the element
930 // NOTE: that for element - GetTitle() returns the "name" tag
931 // and GetName() returns "formula" tag (see createElementN)
932 TString lname = myElement->GetTitle();
933 //_elm suffix is added to avoid problems with same names
934 // for material and element
935 lname = TString::Format("%s_elm", lname.Data());
936
937 // cumulates percentage, in case 2 elements with same names within one mixture
938 wPercentage[lname] += mixture->GetWmixt()[idx];
939 wCounter[lname]++;
940
941 // check whether element name is not in list of elements already created
943 continue;
944 }
945
946 // add element to list of elements and to main <materials> node
949 fGdmlE->AddChild(materials, elmNode);
950 }
951 // loop through asoc array
952 for (NameListI::iterator itr = wCounter.begin(); itr != wCounter.end(); ++itr) {
953 if (itr->second > 1) {
954 Info("CreateMixtureN", "WARNING! 2 equal elements in one material. Check: %s element of %s material",
955 itr->first.Data(), mname.Data());
956 }
957 // add fraction child to material with reference to element
958 fGdmlE->AddChild(mainN, CreateFractionN(wPercentage[itr->first], itr->first.Data()));
959 }
960
961 return mainN;
962}
963
964////////////////////////////////////////////////////////////////////////////////
965/// Creates "material" node for GDML
966
968{
969 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "material", nullptr);
970 fGdmlE->NewAttr(mainN, nullptr, "name", mname);
971 Double_t valZ = material->GetZ();
972 // Z can't be zero in Geant4 so this is workaround for vacuum
975 tmpname.ToLower();
976 if (valZ < 1) {
977 if (tmpname == "vacuum") {
978 valZ = 1;
979 } else {
980 if (fgG4Compatibility == kTRUE) {
981 Info("CreateMaterialN",
982 "WARNING! value of Z in %s material can't be < 1 in Geant4, that is why it was changed to 1, please "
983 "check it manually! ",
984 mname.Data());
985 valZ = 1;
986 } else {
987 Info("CreateMaterialN", "WARNING! value of Z in %s material can't be < 1 in Geant4", mname.Data());
988 }
989 }
990 }
991 fGdmlE->NewAttr(mainN, nullptr, "Z", TString::Format(fltPrecision.Data(), valZ)); // material->GetZ()));
992
993 // Create properties if any: Properties according to the GDML schema MUST come first!
994 TList const &properties = material->GetProperties();
995 if (properties.GetSize()) {
996 TIter next(&properties);
998 while ((property = (TNamed *)next()))
1000 }
1001 // Write CONST properties
1002 TList const &const_properties = material->GetConstProperties();
1003 if (const_properties.GetSize()) {
1004 TIter next(&const_properties);
1006 while ((property = (TNamed *)next()))
1008 }
1009
1010 // Now add the other children
1011 fGdmlE->AddChild(mainN, CreateDN(material->GetDensity()));
1012 fGdmlE->AddChild(mainN, CreateAtomN(material->GetA()));
1013 return mainN;
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// Creates "box" node for GDML
1018
1020{
1021 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "box", nullptr);
1023 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1024 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1025 if (IsNullParam(geoShape->GetDX(), "DX", lname) || IsNullParam(geoShape->GetDY(), "DY", lname) ||
1026 IsNullParam(geoShape->GetDZ(), "DZ", lname)) {
1027 return nullptr;
1028 }
1029 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDX()));
1030 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDY()));
1031 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDZ()));
1032
1033 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1034 return mainN;
1035}
1036
1037////////////////////////////////////////////////////////////////////////////////
1038/// Creates "paraboloid" node for GDML
1039
1041{
1042 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "paraboloid", nullptr);
1044 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1045 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1046 if (IsNullParam(geoShape->GetRhi(), "Rhi", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1047 return nullptr;
1048 }
1049 fGdmlE->NewAttr(mainN, nullptr, "rlo", TString::Format(fltPrecision.Data(), geoShape->GetRlo()));
1050 fGdmlE->NewAttr(mainN, nullptr, "rhi", TString::Format(fltPrecision.Data(), geoShape->GetRhi()));
1051 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1052
1053 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1054 return mainN;
1055}
1056
1057////////////////////////////////////////////////////////////////////////////////
1058/// Creates "sphere" node for GDML
1059
1061{
1062 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "sphere", nullptr);
1064 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1065 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1066 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
1067 return nullptr;
1068 }
1069
1070 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1071 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1072 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1073 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1074 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1075 fGdmlE->NewAttr(mainN, nullptr, "starttheta", TString::Format(fltPrecision.Data(), geoShape->GetTheta1()));
1076 fGdmlE->NewAttr(mainN, nullptr, "deltatheta",
1077 TString::Format(fltPrecision.Data(), geoShape->GetTheta2() - geoShape->GetTheta1()));
1078
1079 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1080 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1081 return mainN;
1082}
1083
1084////////////////////////////////////////////////////////////////////////////////
1085/// Creates "arb8" node for GDML
1086
1088{
1089 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "arb8", nullptr);
1091 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1092 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1093 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1094 return nullptr;
1095 }
1096
1097 fGdmlE->NewAttr(mainN, nullptr, "v1x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[0]));
1098 fGdmlE->NewAttr(mainN, nullptr, "v1y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[1]));
1099 fGdmlE->NewAttr(mainN, nullptr, "v2x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[2]));
1100 fGdmlE->NewAttr(mainN, nullptr, "v2y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[3]));
1101 fGdmlE->NewAttr(mainN, nullptr, "v3x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[4]));
1102 fGdmlE->NewAttr(mainN, nullptr, "v3y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[5]));
1103 fGdmlE->NewAttr(mainN, nullptr, "v4x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[6]));
1104 fGdmlE->NewAttr(mainN, nullptr, "v4y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[7]));
1105 fGdmlE->NewAttr(mainN, nullptr, "v5x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[8]));
1106 fGdmlE->NewAttr(mainN, nullptr, "v5y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[9]));
1107 fGdmlE->NewAttr(mainN, nullptr, "v6x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[10]));
1108 fGdmlE->NewAttr(mainN, nullptr, "v6y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[11]));
1109 fGdmlE->NewAttr(mainN, nullptr, "v7x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[12]));
1110 fGdmlE->NewAttr(mainN, nullptr, "v7y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[13]));
1111 fGdmlE->NewAttr(mainN, nullptr, "v8x", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[14]));
1112 fGdmlE->NewAttr(mainN, nullptr, "v8y", TString::Format(fltPrecision.Data(), geoShape->GetVertices()[15]));
1113 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1114
1115 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1116 return mainN;
1117}
1118
1119////////////////////////////////////////////////////////////////////////////////
1120/// Creates "cone" node for GDML from TGeoConeSeg object
1121
1123{
1124 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "cone", nullptr);
1126 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1127 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1128 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1129 return nullptr;
1130 }
1131
1132 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1133 fGdmlE->NewAttr(mainN, nullptr, "rmin1", TString::Format(fltPrecision.Data(), geoShape->GetRmin1()));
1134 fGdmlE->NewAttr(mainN, nullptr, "rmin2", TString::Format(fltPrecision.Data(), geoShape->GetRmin2()));
1135 fGdmlE->NewAttr(mainN, nullptr, "rmax1", TString::Format(fltPrecision.Data(), geoShape->GetRmax1()));
1136 fGdmlE->NewAttr(mainN, nullptr, "rmax2", TString::Format(fltPrecision.Data(), geoShape->GetRmax2()));
1137 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1138 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1139 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1140
1141 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1142 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1143 return mainN;
1144}
1145
1146////////////////////////////////////////////////////////////////////////////////
1147/// Creates "cone" node for GDML from TGeoCone object
1148
1150{
1151 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "cone", nullptr);
1153 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1154 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1155 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1156 return nullptr;
1157 }
1158
1159 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1160 fGdmlE->NewAttr(mainN, nullptr, "rmin1", TString::Format(fltPrecision.Data(), geoShape->GetRmin1()));
1161 fGdmlE->NewAttr(mainN, nullptr, "rmin2", TString::Format(fltPrecision.Data(), geoShape->GetRmin2()));
1162 fGdmlE->NewAttr(mainN, nullptr, "rmax1", TString::Format(fltPrecision.Data(), geoShape->GetRmax1()));
1163 fGdmlE->NewAttr(mainN, nullptr, "rmax2", TString::Format(fltPrecision.Data(), geoShape->GetRmax2()));
1164 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format("%i", 0));
1165 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format("%i", 360));
1166
1167 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1168 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1169 return mainN;
1170}
1171
1172////////////////////////////////////////////////////////////////////////////////
1173/// Creates "para" node for GDML
1174
1176{
1177 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "para", nullptr);
1179 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1180
1181 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), 2 * geoShape->GetX()));
1182 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), 2 * geoShape->GetY()));
1183 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetZ()));
1184 fGdmlE->NewAttr(mainN, nullptr, "alpha", TString::Format(fltPrecision.Data(), geoShape->GetAlpha()));
1185 fGdmlE->NewAttr(mainN, nullptr, "theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1186 fGdmlE->NewAttr(mainN, nullptr, "phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1187
1188 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1189 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1190 return mainN;
1191}
1192
1193////////////////////////////////////////////////////////////////////////////////
1194/// Creates "trap" node for GDML
1195
1197{
1200
1201 // if one base equals 0 create Arb8 instead of trap
1202 if ((geoShape->GetBl1() == 0 || geoShape->GetTl1() == 0 || geoShape->GetH1() == 0) ||
1203 (geoShape->GetBl2() == 0 || geoShape->GetTl2() == 0 || geoShape->GetH2() == 0)) {
1205 return mainN;
1206 }
1207
1208 // if is twisted then create arb8
1209 if (geoShape->IsTwisted()) {
1211 return mainN;
1212 }
1213
1214 mainN = fGdmlE->NewChild(nullptr, nullptr, "trap", nullptr);
1215 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1216 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1217 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1218 return nullptr;
1219 }
1220
1221 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1222 fGdmlE->NewAttr(mainN, nullptr, "theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1223 fGdmlE->NewAttr(mainN, nullptr, "phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1224 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl1()));
1225 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl1()));
1226 fGdmlE->NewAttr(mainN, nullptr, "x3", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl2()));
1227 fGdmlE->NewAttr(mainN, nullptr, "x4", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl2()));
1228 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH1()));
1229 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH2()));
1230
1231 fGdmlE->NewAttr(mainN, nullptr, "alpha1", TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()));
1232 fGdmlE->NewAttr(mainN, nullptr, "alpha2", TString::Format(fltPrecision.Data(), geoShape->GetAlpha2()));
1233
1234 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1235 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1236 return mainN;
1237}
1238
1239////////////////////////////////////////////////////////////////////////////////
1240/// Creates "twistedtrap" node for GDML
1241
1243{
1246
1247 // if one base equals 0 create Arb8 instead of twisted trap
1248 if ((geoShape->GetBl1() == 0 && geoShape->GetTl1() == 0 && geoShape->GetH1() == 0) ||
1249 (geoShape->GetBl2() == 0 && geoShape->GetTl2() == 0 && geoShape->GetH2() == 0)) {
1251 return mainN;
1252 }
1253
1254 // if is twisted then create arb8
1255 if (geoShape->IsTwisted()) {
1257 return mainN;
1258 }
1259
1260 // if parameter twistAngle (PhiTwist) equals zero create trap node
1261 if (geoShape->GetTwistAngle() == 0) {
1263 return mainN;
1264 }
1265
1266 mainN = fGdmlE->NewChild(nullptr, nullptr, "twistedtrap", nullptr);
1267 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1268 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1269 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1270 return nullptr;
1271 }
1272
1273 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1274 fGdmlE->NewAttr(mainN, nullptr, "Theta", TString::Format(fltPrecision.Data(), geoShape->GetTheta()));
1275 fGdmlE->NewAttr(mainN, nullptr, "Phi", TString::Format(fltPrecision.Data(), geoShape->GetPhi()));
1276 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl1()));
1277 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl1()));
1278 fGdmlE->NewAttr(mainN, nullptr, "x3", TString::Format(fltPrecision.Data(), 2 * geoShape->GetBl2()));
1279 fGdmlE->NewAttr(mainN, nullptr, "x4", TString::Format(fltPrecision.Data(), 2 * geoShape->GetTl2()));
1280 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH1()));
1281 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetH2()));
1282
1283 fGdmlE->NewAttr(mainN, nullptr, "Alph", TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()));
1284
1285 // check if alpha1 equals to alpha2 (converting to string - to avoid problems with floats)
1286 if (TString::Format(fltPrecision.Data(), geoShape->GetAlpha1()) !=
1287 TString::Format(fltPrecision.Data(), geoShape->GetAlpha2())) {
1288 Info("CreateTwistedTrapN",
1289 "ERROR! Object %s is not exported correctly because parameter Alpha2 is not declared in GDML schema",
1290 lname.Data());
1291 }
1292 // fGdmlE->NewAttr(mainN,0, "alpha2", TString::Format(fltPrecision.Data(), geoShape->GetAlpha2()));
1293 fGdmlE->NewAttr(mainN, nullptr, "PhiTwist", TString::Format(fltPrecision.Data(), geoShape->GetTwistAngle()));
1294
1295 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1296 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1297 return mainN;
1298}
1299
1300////////////////////////////////////////////////////////////////////////////////
1301/// Creates "trd" node for GDML from object TGeoTrd1
1302
1304{
1305 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "trd", nullptr);
1307 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1308 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1309 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1310 return nullptr;
1311 }
1312
1313 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx1()));
1314 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx2()));
1315 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy()));
1316 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy()));
1317 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1318
1319 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1320 return mainN;
1321}
1322
1323////////////////////////////////////////////////////////////////////////////////
1324/// Creates "trd" node for GDML from object TGeoTrd2
1325
1327{
1328 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "trd", nullptr);
1330 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1331 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1332 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1333 return nullptr;
1334 }
1335
1336 fGdmlE->NewAttr(mainN, nullptr, "x1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx1()));
1337 fGdmlE->NewAttr(mainN, nullptr, "x2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDx2()));
1338 fGdmlE->NewAttr(mainN, nullptr, "y1", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy1()));
1339 fGdmlE->NewAttr(mainN, nullptr, "y2", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDy2()));
1340 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1341
1342 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1343 return mainN;
1344}
1345
1346////////////////////////////////////////////////////////////////////////////////
1347/// Creates "tube" node for GDML from object TGeoTubeSeg
1348
1350{
1351 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tube", nullptr);
1353 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1354 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1355 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1356 return nullptr;
1357 }
1358
1359 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1360 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1361 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1362 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1363 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1364 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1365
1366 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1367 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1368 return mainN;
1369}
1370
1371////////////////////////////////////////////////////////////////////////////////
1372/// Creates "cutTube" node for GDML
1373
1375{
1378
1379 mainN = fGdmlE->NewChild(nullptr, nullptr, "cutTube", nullptr);
1380 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1381 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1382 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1383 return nullptr;
1384 }
1385 // This is not needed, because cutTube is already supported by Geant4 9.5
1386 if (fgG4Compatibility == kTRUE && kFALSE) {
1389
1390 // register name for cuttube shape (so it will be found during volume export)
1393 Info("CreateCutTubeN", "WARNING! %s - CutTube was replaced by intersection of TGeoTubSeg and two TGeoBBoxes",
1394 lname.Data());
1395 return mainN;
1396 }
1397 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1398 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1399 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1400 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1401 fGdmlE->NewAttr(mainN, nullptr, "deltaphi",
1402 TString::Format(fltPrecision.Data(), geoShape->GetPhi2() - geoShape->GetPhi1()));
1403 fGdmlE->NewAttr(mainN, nullptr, "lowX", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[0]));
1404 fGdmlE->NewAttr(mainN, nullptr, "lowY", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[1]));
1405 fGdmlE->NewAttr(mainN, nullptr, "lowZ", TString::Format(fltPrecision.Data(), geoShape->GetNlow()[2]));
1406 fGdmlE->NewAttr(mainN, nullptr, "highX", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[0]));
1407 fGdmlE->NewAttr(mainN, nullptr, "highY", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[1]));
1408 fGdmlE->NewAttr(mainN, nullptr, "highZ", TString::Format(fltPrecision.Data(), geoShape->GetNhigh()[2]));
1409
1410 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1411 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1412
1413 return mainN;
1414}
1415
1416////////////////////////////////////////////////////////////////////////////////
1417/// Creates "tube" node for GDML from object TGeoTube
1418
1420{
1421 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tube", nullptr);
1423 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1424 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1425 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname) || IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1426 return nullptr;
1427 }
1428
1429 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1430 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1431 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1432 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format("%i", 0));
1433 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format("%i", 360));
1434
1435 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1436 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1437 return mainN;
1438}
1439
1440////////////////////////////////////////////////////////////////////////////////
1441/// Creates "zplane" node for GDML
1442
1444{
1445 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "zplane", nullptr);
1447
1448 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), z));
1449 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), rmin));
1450 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), rmax));
1451
1452 return mainN;
1453}
1454
1455////////////////////////////////////////////////////////////////////////////////
1456/// Creates "polycone" node for GDML
1457
1459{
1460 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "polycone", nullptr);
1462 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1463 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1464
1465 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1466 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1467
1468 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1469 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1470 Int_t nZPlns = geoShape->GetNz();
1471 for (Int_t it = 0; it < nZPlns; it++) {
1472 // add zplane child node
1473 fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1474 // compare actual plane and next plane
1475 if ((it < nZPlns - 1) && (geoShape->GetZ(it) == geoShape->GetZ(it + 1))) {
1476 // rmin of actual is greater then rmax of next one
1477 // | |rmax next
1478 // __ ...| |... __ < rmin actual
1479 // | | | |
1480 if (geoShape->GetRmin(it) > geoShape->GetRmax(it + 1)) {
1481 // adding plane from rmax next to rmin actual at the same z position
1482 if (fgG4Compatibility == kTRUE) {
1484 CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it + 1), geoShape->GetRmin(it)));
1485 Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4",
1486 lname.Data());
1487 } else {
1488 Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4",
1489 lname.Data());
1490 }
1491 }
1492 // rmin of next is greater then rmax of actual
1493 // | | | |
1494 // | |...___...| | rmin next
1495 // | | > rmax act
1496 if (geoShape->GetRmin(it + 1) > geoShape->GetRmax(it)) {
1497 // adding plane from rmax act to rmin next at the same z position
1498 if (fgG4Compatibility == kTRUE) {
1500 CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmax(it), geoShape->GetRmin(it + 1)));
1501 Info("CreatePolyconeN", "WARNING! One plane was added to %s solid to be compatible with Geant4",
1502 lname.Data());
1503 } else {
1504 Info("CreatePolyconeN", "WARNING! Solid %s definition seems not contiguous may cause problems in Geant4",
1505 lname.Data());
1506 }
1507 }
1508 }
1509 }
1510 return mainN;
1511}
1512
1513////////////////////////////////////////////////////////////////////////////////
1514/// Creates "torus" node for GDML
1515
1517{
1518 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "torus", nullptr);
1520 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1521 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1522 if (IsNullParam(geoShape->GetRmax(), "Rmax", lname)) {
1523 return nullptr;
1524 }
1525
1526 fGdmlE->NewAttr(mainN, nullptr, "rtor", TString::Format(fltPrecision.Data(), geoShape->GetR()));
1527 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1528 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1529 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1530 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1531
1532 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1533 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1534
1535 return mainN;
1536}
1537
1538////////////////////////////////////////////////////////////////////////////////
1539/// Creates "polyhedra" node for GDML
1540
1542{
1543 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "polyhedra", nullptr);
1545 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1546
1547 fGdmlE->NewAttr(mainN, nullptr, "startphi", TString::Format(fltPrecision.Data(), geoShape->GetPhi1()));
1548 fGdmlE->NewAttr(mainN, nullptr, "deltaphi", TString::Format(fltPrecision.Data(), geoShape->GetDphi()));
1549 fGdmlE->NewAttr(mainN, nullptr, "numsides", TString::Format("%i", geoShape->GetNedges()));
1550
1551 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1552 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1553 for (Int_t it = 0; it < geoShape->GetNz(); it++) {
1554 // add zplane child node
1555 fGdmlE->AddChild(mainN, CreateZplaneN(geoShape->GetZ(it), geoShape->GetRmin(it), geoShape->GetRmax(it)));
1556 }
1557 return mainN;
1558}
1559
1560////////////////////////////////////////////////////////////////////////////////
1561/// Creates "eltube" node for GDML
1562
1564{
1565 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "eltube", nullptr);
1567 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1568 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1569 if (IsNullParam(geoShape->GetA(), "A", lname) || IsNullParam(geoShape->GetB(), "B", lname) ||
1570 IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1571 return nullptr;
1572 }
1573
1574 fGdmlE->NewAttr(mainN, nullptr, "dx", TString::Format(fltPrecision.Data(), geoShape->GetA()));
1575 fGdmlE->NewAttr(mainN, nullptr, "dy", TString::Format(fltPrecision.Data(), geoShape->GetB()));
1576 fGdmlE->NewAttr(mainN, nullptr, "dz", TString::Format(fltPrecision.Data(), geoShape->GetDz()));
1577
1578 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1579
1580 return mainN;
1581}
1582
1583////////////////////////////////////////////////////////////////////////////////
1584/// Creates "hype" node for GDML
1585
1587{
1588 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "hype", nullptr);
1590 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1591 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1592 if (IsNullParam(geoShape->GetDz(), "Dz", lname)) {
1593 return nullptr;
1594 }
1595
1596 fGdmlE->NewAttr(mainN, nullptr, "rmin", TString::Format(fltPrecision.Data(), geoShape->GetRmin()));
1597 fGdmlE->NewAttr(mainN, nullptr, "rmax", TString::Format(fltPrecision.Data(), geoShape->GetRmax()));
1598 fGdmlE->NewAttr(mainN, nullptr, "inst", TString::Format(fltPrecision.Data(), geoShape->GetStIn()));
1599 fGdmlE->NewAttr(mainN, nullptr, "outst", TString::Format(fltPrecision.Data(), geoShape->GetStOut()));
1600 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), 2 * geoShape->GetDz()));
1601
1602 fGdmlE->NewAttr(mainN, nullptr, "aunit", "deg");
1603 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1604
1605 return mainN;
1606}
1607
1608////////////////////////////////////////////////////////////////////////////////
1609/// Creates "xtru" node for GDML
1610
1612{
1613 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "xtru", nullptr);
1615 TString lname = GenName(geoShape->GetName(), TString::Format("%p", geoShape));
1616 fGdmlE->NewAttr(mainN, nullptr, "name", lname);
1617
1618 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1620 Int_t vertNum = geoShape->GetNvert();
1621 Int_t secNum = geoShape->GetNz();
1622 if (vertNum < 3 || secNum < 2) {
1623 Info("CreateXtrusionN", "ERROR! TGeoXtru %s has only %i vertices and %i sections. It was not exported",
1624 lname.Data(), vertNum, secNum);
1625 mainN = nullptr;
1626 return mainN;
1627 }
1628 for (Int_t it = 0; it < vertNum; it++) {
1629 // add twoDimVertex child node
1630 childN = fGdmlE->NewChild(nullptr, nullptr, "twoDimVertex", nullptr);
1631 fGdmlE->NewAttr(childN, nullptr, "x", TString::Format(fltPrecision.Data(), geoShape->GetX(it)));
1632 fGdmlE->NewAttr(childN, nullptr, "y", TString::Format(fltPrecision.Data(), geoShape->GetY(it)));
1634 }
1635 for (Int_t it = 0; it < secNum; it++) {
1636 // add section child node
1637 childN = fGdmlE->NewChild(nullptr, nullptr, "section", nullptr);
1638 fGdmlE->NewAttr(childN, nullptr, "zOrder", TString::Format("%i", it));
1639 fGdmlE->NewAttr(childN, nullptr, "zPosition", TString::Format(fltPrecision.Data(), geoShape->GetZ(it)));
1640 fGdmlE->NewAttr(childN, nullptr, "xOffset", TString::Format(fltPrecision.Data(), geoShape->GetXOffset(it)));
1641 fGdmlE->NewAttr(childN, nullptr, "yOffset", TString::Format(fltPrecision.Data(), geoShape->GetYOffset(it)));
1642 fGdmlE->NewAttr(childN, nullptr, "scalingFactor", TString::Format(fltPrecision.Data(), geoShape->GetScale(it)));
1644 }
1645 return mainN;
1646}
1647
1648////////////////////////////////////////////////////////////////////////////////
1649/// Creates "ellipsoid" node for GDML
1650/// this is a special case, because ellipsoid is not defined in ROOT
1651/// so when intersection of scaled sphere and TGeoBBox is found,
1652/// it is considered as an ellipsoid
1653
1655{
1656 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "ellipsoid", nullptr);
1658 TGeoScaledShape *leftS = (TGeoScaledShape *)geoShape->GetBoolNode()->GetLeftShape(); // ScaledShape
1659 TGeoBBox *rightS = (TGeoBBox *)geoShape->GetBoolNode()->GetRightShape(); // BBox
1660
1661 fGdmlE->NewAttr(mainN, nullptr, "name", elName.Data());
1662 Double_t sx = leftS->GetScale()->GetScale()[0];
1663 Double_t sy = leftS->GetScale()->GetScale()[1];
1664 Double_t radius = ((TGeoSphere *)leftS->GetShape())->GetRmax();
1665
1666 Double_t ax, by, cz;
1667 cz = radius;
1668 ax = sx * radius;
1669 by = sy * radius;
1670
1671 Double_t dz = rightS->GetDZ();
1672 Double_t zorig = rightS->GetOrigin()[2];
1673 Double_t zcut2 = dz + zorig;
1674 Double_t zcut1 = 2 * zorig - zcut2;
1675
1676 fGdmlE->NewAttr(mainN, nullptr, "ax", TString::Format(fltPrecision.Data(), ax));
1677 fGdmlE->NewAttr(mainN, nullptr, "by", TString::Format(fltPrecision.Data(), by));
1678 fGdmlE->NewAttr(mainN, nullptr, "cz", TString::Format(fltPrecision.Data(), cz));
1679 fGdmlE->NewAttr(mainN, nullptr, "zcut1", TString::Format(fltPrecision.Data(), zcut1));
1680 fGdmlE->NewAttr(mainN, nullptr, "zcut2", TString::Format(fltPrecision.Data(), zcut2));
1681 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1682
1683 return mainN;
1684}
1685
1686////////////////////////////////////////////////////////////////////////////////
1687/// Creates "elcone" (elliptical cone) node for GDML
1688/// this is a special case, because elliptical cone is not defined in ROOT
1689/// so when scaled cone is found, it is considered as a elliptical cone
1690
1692{
1693 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "elcone", nullptr);
1695 fGdmlE->NewAttr(mainN, nullptr, "name", GenName(geoShape->GetName(), TString::Format("%p", geoShape)));
1696 Double_t zcut = ((TGeoCone *)geoShape->GetShape())->GetDz();
1697 Double_t rx1 = ((TGeoCone *)geoShape->GetShape())->GetRmax1();
1698 Double_t rx2 = ((TGeoCone *)geoShape->GetShape())->GetRmax2();
1699 Double_t zmax = zcut * ((rx1 + rx2) / (rx1 - rx2));
1700 Double_t z = zcut + zmax;
1701
1702 Double_t sy = geoShape->GetScale()->GetScale()[1];
1703 Double_t ry1 = sy * rx1;
1704
1705 std::string format(TString::Format("%s/%s", fltPrecision.Data(), fltPrecision.Data()).Data());
1706 fGdmlE->NewAttr(mainN, nullptr, "dx", TString::Format(format.c_str(), rx1, z));
1707 fGdmlE->NewAttr(mainN, nullptr, "dy", TString::Format(format.c_str(), ry1, z));
1708 fGdmlE->NewAttr(mainN, nullptr, "zmax", TString::Format(fltPrecision.Data(), zmax));
1709 fGdmlE->NewAttr(mainN, nullptr, "zcut", TString::Format(fltPrecision.Data(), zcut));
1710 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1711
1712 return mainN;
1713}
1714
1715////////////////////////////////////////////////////////////////////////////////
1716/// Creates "tessellated" (tessellated shape) node for GDML
1717
1719{
1720 // add all vertices to the define section
1722 for (int i = 0; i < geoShape->GetNvertices(); ++i) {
1723 auto vertex = geoShape->GetVertex(i);
1724 TString posName = TString::Format("%s_%d", genname.Data(), i);
1725 Xyz nodPos;
1726 nodPos.x = vertex[0];
1727 nodPos.y = vertex[1];
1728 nodPos.z = vertex[2];
1729 auto childN = CreatePositionN(posName.Data(), nodPos, "position", fDefault_lunit);
1730 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
1731 }
1732 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "tessellated", nullptr);
1733 fGdmlE->NewAttr(mainN, nullptr, "name", genname.Data());
1734 fGdmlE->NewAttr(mainN, nullptr, "lunit", fDefault_lunit);
1735
1737 for (Int_t it = 0; it < geoShape->GetNfacets(); it++) {
1738 // add section child node
1739 auto facet = geoShape->GetFacet(it);
1740 bool triangular = facet.GetNvert() == 3;
1741 TString ntype = (triangular) ? "triangular" : "quadrangular";
1742 childN = fGdmlE->NewChild(nullptr, nullptr, ntype.Data(), nullptr);
1743 fGdmlE->NewAttr(childN, nullptr, "vertex1", TString::Format("%s_%d", genname.Data(), facet[0]));
1744 fGdmlE->NewAttr(childN, nullptr, "vertex2", TString::Format("%s_%d", genname.Data(), facet[1]));
1745 fGdmlE->NewAttr(childN, nullptr, "vertex3", TString::Format("%s_%d", genname.Data(), facet[2]));
1746 if (!triangular)
1747 fGdmlE->NewAttr(childN, nullptr, "vertex4", TString::Format("%s_%d", genname.Data(), facet[3]));
1748 fGdmlE->NewAttr(childN, nullptr, "type", "ABSOLUTE");
1750 }
1751 return mainN;
1752}
1753
1754////////////////////////////////////////////////////////////////////////////////
1755/// Creates a scaled node for GDML
1756
1758{
1762 auto const scale = geoShape->GetScale()->GetScale();
1763 TGeoShape *unscaled = geoShape->GetShape();
1764 // Create the unscaled shape record
1766 // retrieve node name by their pointer to make reference
1768
1769 // the unplaced node appended to main structure of nodes (if they are not already there)
1770 if (unscaledN != nullptr) {
1772 fSolCnt++;
1773 } else {
1774 if (uname.Contains("missing_") || uname == "") {
1775 Info("CreateScaledN", "ERROR! Unscaled node is NULL - Scaled shape will be skipped");
1776 return nullptr;
1777 }
1778 }
1779
1780 // create union node and its child nodes (or intersection or subtraction)
1781 /* <scaledSolid name="...">
1782 * <solidref ref="unscaled solid name" />
1783 * <scale name="scale name" x="sx" y="sy" z="sz">
1784 * </scaledSolid>
1785 */
1786 mainN = fGdmlE->NewChild(nullptr, nullptr, "scaledSolid", nullptr);
1787 fGdmlE->NewAttr(mainN, nullptr, "name", nodeName);
1788
1789 //<solidref> (solid)
1790 childN = fGdmlE->NewChild(nullptr, nullptr, "solidref", nullptr);
1791 fGdmlE->NewAttr(childN, nullptr, "ref", uname);
1793
1794 //<scale ame="scale name" x="sx" y="sy" z="sz">
1795 childN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
1796 fGdmlE->NewAttr(childN, nullptr, "name", (nodeName + "scl").Data());
1797 fGdmlE->NewAttr(childN, nullptr, "x", TString::Format(fltPrecision.Data(), scale[0]));
1798 fGdmlE->NewAttr(childN, nullptr, "y", TString::Format(fltPrecision.Data(), scale[1]));
1799 fGdmlE->NewAttr(childN, nullptr, "z", TString::Format(fltPrecision.Data(), scale[2]));
1801
1802 return mainN;
1803}
1804
1805////////////////////////////////////////////////////////////////////////////////
1806/// Creates common part of union intersection and subtraction nodes
1807
1809{
1813 TGeoBoolNode::EGeoBoolType boolType = geoShape->GetBoolNode()->GetBooleanOperator();
1814 switch (boolType) {
1815 case TGeoBoolNode::kGeoUnion: lboolType = "union"; break;
1816 case TGeoBoolNode::kGeoSubtraction: lboolType = "subtraction"; break;
1817 case TGeoBoolNode::kGeoIntersection: lboolType = "intersection"; break;
1818 }
1819
1820 TGDMLWrite::Xyz lrot = GetXYZangles(geoShape->GetBoolNode()->GetLeftMatrix()->Inverse().GetRotationMatrix());
1821 const Double_t *ltr = geoShape->GetBoolNode()->GetLeftMatrix()->GetTranslation();
1822 TGDMLWrite::Xyz rrot = GetXYZangles(geoShape->GetBoolNode()->GetRightMatrix()->Inverse().GetRotationMatrix());
1823 const Double_t *rtr = geoShape->GetBoolNode()->GetRightMatrix()->GetTranslation();
1824
1825 // specific case!
1826 // Ellipsoid tag preparing
1827 // if left == TGeoScaledShape AND right == TGeoBBox
1828 // AND if TGeoScaledShape->GetShape == TGeoSphere
1829 TGeoShape *leftS = geoShape->GetBoolNode()->GetLeftShape();
1830 TGeoShape *rightS = geoShape->GetBoolNode()->GetRightShape();
1831 if (strcmp(leftS->ClassName(), "TGeoScaledShape") == 0 && strcmp(rightS->ClassName(), "TGeoBBox") == 0) {
1832 if (strcmp(((TGeoScaledShape *)leftS)->GetShape()->ClassName(), "TGeoSphere") == 0) {
1833 if (lboolType == "intersection") {
1835 return mainN;
1836 }
1837 }
1838 }
1839
1841 // translation
1842 translL.x = ltr[0];
1843 translL.y = ltr[1];
1844 translL.z = ltr[2];
1845 translR.x = rtr[0];
1846 translR.y = rtr[1];
1847 translR.z = rtr[2];
1848
1849 // left and right nodes are created here also their names are created
1850 ndL = ChooseObject(geoShape->GetBoolNode()->GetLeftShape());
1851 ndR = ChooseObject(geoShape->GetBoolNode()->GetRightShape());
1852
1853 // retrieve left and right node names by their pointer to make reference
1854 TString lname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetLeftShape())];
1855 TString rname = fNameList->fLst[TString::Format("%p", geoShape->GetBoolNode()->GetRightShape())];
1856
1857 // left and right nodes appended to main structure of nodes (if they are not already there)
1858 if (ndL != nullptr) {
1860 fSolCnt++;
1861 } else {
1862 if (lname.Contains("missing_") || lname == "") {
1863 Info("CreateCommonBoolN", "ERROR! Left node is NULL - Boolean Shape will be skipped");
1864 return nullptr;
1865 }
1866 }
1867 if (ndR != nullptr) {
1869 fSolCnt++;
1870 } else {
1871 if (rname.Contains("missing_") || rname == "") {
1872 Info("CreateCommonBoolN", "ERROR! Right node is NULL - Boolean Shape will be skipped");
1873 return nullptr;
1874 }
1875 }
1876
1877 // create union node and its child nodes (or intersection or subtraction)
1878 /* <union name="...">
1879 * <first ref="left name" />
1880 * <second ref="right name" />
1881 * <firstposition .../>
1882 * <firstrotation .../>
1883 * <position .../>
1884 * <rotation .../>
1885 * </union>
1886 */
1887 mainN = fGdmlE->NewChild(nullptr, nullptr, lboolType.Data(), nullptr);
1888 fGdmlE->NewAttr(mainN, nullptr, "name", nodeName);
1889
1890 //<first> (left)
1891 childN = fGdmlE->NewChild(nullptr, nullptr, "first", nullptr);
1892 fGdmlE->NewAttr(childN, nullptr, "ref", lname);
1894
1895 //<second> (right)
1896 childN = fGdmlE->NewChild(nullptr, nullptr, "second", nullptr);
1897 fGdmlE->NewAttr(childN, nullptr, "ref", rname);
1899
1900 //<firstposition> (left)
1901 if ((translL.x != 0.0) || (translL.y != 0.0) || (translL.z != 0.0)) {
1902 childN = CreatePositionN((nodeName + lname + "pos").Data(), translL, "firstposition", fDefault_lunit);
1904 }
1905 //<firstrotation> (left)
1906 if ((lrot.x != 0.0) || (lrot.y != 0.0) || (lrot.z != 0.0)) {
1907 childN = CreateRotationN((nodeName + lname + "rot").Data(), lrot, "firstrotation");
1909 }
1910 //<position> (right)
1911 if ((translR.x != 0.0) || (translR.y != 0.0) || (translR.z != 0.0)) {
1912 childN = CreatePositionN((nodeName + rname + "pos").Data(), translR, "position", fDefault_lunit);
1914 }
1915 //<rotation> (right)
1916 if ((rrot.x != 0.0) || (rrot.y != 0.0) || (rrot.z != 0.0)) {
1917 childN = CreateRotationN((nodeName + rname + "rot").Data(), rrot, "rotation");
1919 }
1920
1921 return mainN;
1922}
1923
1924////////////////////////////////////////////////////////////////////////////////
1925/// Creates "opticalsurface" node for GDML
1926
1928{
1929 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "opticalsurface", nullptr);
1931 std::string name = make_NCName(geoSurf->GetName());
1932 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1933 fGdmlE->NewAttr(mainN, nullptr, "model", TGeoOpticalSurface::ModelToString(geoSurf->GetModel()));
1934 fGdmlE->NewAttr(mainN, nullptr, "finish", TGeoOpticalSurface::FinishToString(geoSurf->GetFinish()));
1935 fGdmlE->NewAttr(mainN, nullptr, "type", TGeoOpticalSurface::TypeToString(geoSurf->GetType()));
1936 fGdmlE->NewAttr(mainN, nullptr, "value", TString::Format(fltPrecision.Data(), geoSurf->GetValue()));
1937
1938 // Write properties
1939 TList const &properties = geoSurf->GetProperties();
1940 if (properties.GetSize()) {
1941 TIter next(&properties);
1943 while ((property = (TNamed *)next()))
1945 }
1946 // Write CONST properties
1947 TList const &const_properties = geoSurf->GetConstProperties();
1948 if (const_properties.GetSize()) {
1949 TIter next(&const_properties);
1951 while ((property = (TNamed *)next()))
1953 }
1954 return mainN;
1955}
1956
1957////////////////////////////////////////////////////////////////////////////////
1958/// Creates "skinsurface" node for GDML
1959
1961{
1962 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "skinsurface", nullptr);
1963 std::string name = make_NCName(geoSurf->GetName());
1964 std::string prop = make_NCName(geoSurf->GetTitle());
1965
1966 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1967 fGdmlE->NewAttr(mainN, nullptr, "surfaceproperty", prop.c_str());
1968
1969 // Cretate the logical volume reference node
1970 XMLNodePointer_t childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
1971 const TString &volname = fNameList->fLst[TString::Format("%p", geoSurf->GetVolume())];
1972 fGdmlE->NewAttr(childN, nullptr, "ref", volname.Data());
1974 return mainN;
1975}
1976
1977////////////////////////////////////////////////////////////////////////////////
1978/// Creates "bordersurface" node for GDML
1979
1981{
1982 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "bordersurface", nullptr);
1983 std::string name = make_NCName(geoSurf->GetName());
1984 std::string prop = make_NCName(geoSurf->GetTitle());
1985
1986 fGdmlE->NewAttr(mainN, nullptr, "name", name.c_str());
1987 fGdmlE->NewAttr(mainN, nullptr, "surfaceproperty", prop.c_str());
1988
1989 // Cretate the logical volume reference nodes
1990 XMLNodePointer_t childN = fGdmlE->NewChild(nullptr, nullptr, "physvolref", nullptr);
1992 fGdmlE->NewAttr(childN, nullptr, "ref", physvolname);
1994
1995 childN = fGdmlE->NewChild(nullptr, nullptr, "physvolref", nullptr);
1996 physvolname = fNameList->fLst[TString::Format("%p", geoSurf->GetNode2())];
1997 fGdmlE->NewAttr(childN, nullptr, "ref", physvolname);
1999 return mainN;
2000}
2001
2002////////////////////////////////////////////////////////////////////////////////
2003/// Creates "position" kind of node for GDML
2004
2005XMLNodePointer_t TGDMLWrite::CreatePositionN(const char *name, Xyz position, const char *type, const char *unit)
2006{
2007 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, type, nullptr);
2009 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2010 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), position.x));
2011 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), position.y));
2012 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), position.z));
2013 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2014 return mainN;
2015}
2016
2017////////////////////////////////////////////////////////////////////////////////
2018/// Creates "rotation" kind of node for GDML
2019
2020XMLNodePointer_t TGDMLWrite::CreateRotationN(const char *name, Xyz rotation, const char *type, const char *unit)
2021{
2022 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, type, nullptr);
2024 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2025 fGdmlE->NewAttr(mainN, nullptr, "x", TString::Format(fltPrecision.Data(), rotation.x));
2026 fGdmlE->NewAttr(mainN, nullptr, "y", TString::Format(fltPrecision.Data(), rotation.y));
2027 fGdmlE->NewAttr(mainN, nullptr, "z", TString::Format(fltPrecision.Data(), rotation.z));
2028 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2029 return mainN;
2030}
2031
2032////////////////////////////////////////////////////////////////////////////////
2033/// Creates "matrix" kind of node for GDML
2034
2036{
2037 std::stringstream vals;
2038 size_t cols = matrix->GetCols();
2039 size_t rows = matrix->GetRows();
2040 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "matrix", nullptr);
2041 fGdmlE->NewAttr(mainN, nullptr, "name", matrix->GetName());
2042 fGdmlE->NewAttr(mainN, nullptr, "coldim", TString::Format("%zu", cols));
2043 for (size_t i = 0; i < rows; ++i) {
2044 for (size_t j = 0; j < cols; ++j) {
2045 vals << matrix->Get(i, j);
2046 if (j < cols - 1)
2047 vals << ' ';
2048 }
2049 if (i < rows - 1)
2050 vals << '\n';
2051 }
2052 fGdmlE->NewAttr(mainN, nullptr, "values", vals.str().c_str());
2053 return mainN;
2054}
2055
2056////////////////////////////////////////////////////////////////////////////////
2057/// Creates "constant" kind of node for GDML
2058
2060{
2061 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "matrix", nullptr);
2063 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2064 fGdmlE->NewAttr(mainN, nullptr, "coldim", "1");
2065 fGdmlE->NewAttr(mainN, nullptr, "values", TString::Format(fltPrecision.Data(), value));
2066 return mainN;
2067}
2068
2069////////////////////////////////////////////////////////////////////////////////
2070/// Creates "setup" node for GDML
2071
2073{
2074 XMLNodePointer_t setupN = fGdmlE->NewChild(nullptr, nullptr, "setup", nullptr);
2075 fGdmlE->NewAttr(setupN, nullptr, "name", name);
2076 fGdmlE->NewAttr(setupN, nullptr, "version", version);
2077 XMLNodePointer_t fworldN = fGdmlE->NewChild(setupN, nullptr, "world", nullptr);
2078 fGdmlE->NewAttr(fworldN, nullptr, "ref", topVolName);
2079 return setupN;
2080}
2081
2082////////////////////////////////////////////////////////////////////////////////
2083/// Creates "volume" node for GDML
2084
2085XMLNodePointer_t TGDMLWrite::StartVolumeN(const char *name, const char *solid, const char *material)
2086{
2088 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "volume", nullptr);
2089 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2090
2091 childN = fGdmlE->NewChild(nullptr, nullptr, "materialref", nullptr);
2092 fGdmlE->NewAttr(childN, nullptr, "ref", material);
2094
2095 childN = fGdmlE->NewChild(nullptr, nullptr, "solidref", nullptr);
2096 fGdmlE->NewAttr(childN, nullptr, "ref", solid);
2098
2099 return mainN;
2100}
2101
2102////////////////////////////////////////////////////////////////////////////////
2103/// Creates "assembly" node for GDML
2104
2106{
2107 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "assembly", nullptr);
2108 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2109
2110 return mainN;
2111}
2112
2113////////////////////////////////////////////////////////////////////////////////
2114/// Creates "physvol" node for GDML
2115
2117 const char *rotref, XMLNodePointer_t scaleN)
2118{
2119 fPhysVolCnt++;
2121 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "physvol", nullptr);
2122 fGdmlE->NewAttr(mainN, nullptr, "name", name);
2123 fGdmlE->NewAttr(mainN, nullptr, "copynumber", TString::Format("%d", copyno));
2124
2125 childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
2126 fGdmlE->NewAttr(childN, nullptr, "ref", volref);
2128
2129 childN = fGdmlE->NewChild(nullptr, nullptr, "positionref", nullptr);
2130 fGdmlE->NewAttr(childN, nullptr, "ref", posref);
2132
2133 // if is not empty string add this node
2134 if (strcmp(rotref, "") != 0) {
2135 childN = fGdmlE->NewChild(nullptr, nullptr, "rotationref", nullptr);
2136 fGdmlE->NewAttr(childN, nullptr, "ref", rotref);
2138 }
2139 if (scaleN != nullptr) {
2141 }
2142
2143 return mainN;
2144}
2145
2146////////////////////////////////////////////////////////////////////////////////
2147/// Creates "divisionvol" node for GDML
2148
2150 const char *unit, const char *volref)
2151{
2152 XMLNodePointer_t childN = nullptr;
2153 XMLNodePointer_t mainN = fGdmlE->NewChild(nullptr, nullptr, "divisionvol", nullptr);
2154 fGdmlE->NewAttr(mainN, nullptr, "axis", axis);
2155 fGdmlE->NewAttr(mainN, nullptr, "number", TString::Format("%i", number));
2157 if (fgG4Compatibility == kTRUE) {
2158 // if eg. full length is 20 and width * number = 20,0001 problem in geant4
2159 // unit is either in cm or degrees nothing else
2160 width = (floor(width * 1E4)) * 1E-4;
2161 if ((offset >= 0.) && (strcmp(axis, "kPhi") == 0)) {
2164 // put to range from 0 to 360 add decimals and then put to range 0 -> -360
2165 offset = (offsetI % 360) + decimals - 360;
2166 }
2167 }
2168 fGdmlE->NewAttr(mainN, nullptr, "width", TString::Format(fltPrecision.Data(), width));
2169
2170 fGdmlE->NewAttr(mainN, nullptr, "offset", TString::Format(fltPrecision.Data(), offset));
2171 fGdmlE->NewAttr(mainN, nullptr, "unit", unit);
2172 if (strcmp(volref, "") != 0) {
2173 childN = fGdmlE->NewChild(nullptr, nullptr, "volumeref", nullptr);
2174 fGdmlE->NewAttr(childN, nullptr, "ref", volref);
2175 }
2177
2178 return mainN;
2179}
2180
2181////////////////////////////////////////////////////////////////////////////////
2182/// Chooses the object and method that should be used for processing object
2183
2185{
2186 const char *clsname = geoShape->ClassName();
2188
2189 if (CanProcess((TObject *)geoShape) == kFALSE) {
2190 return nullptr;
2191 }
2192
2193 // process different shapes
2194 if (strcmp(clsname, "TGeoBBox") == 0) {
2196 } else if (strcmp(clsname, "TGeoParaboloid") == 0) {
2198 } else if (strcmp(clsname, "TGeoSphere") == 0) {
2200 } else if (strcmp(clsname, "TGeoArb8") == 0) {
2202 } else if (strcmp(clsname, "TGeoConeSeg") == 0) {
2204 } else if (strcmp(clsname, "TGeoCone") == 0) {
2206 } else if (strcmp(clsname, "TGeoPara") == 0) {
2208 } else if (strcmp(clsname, "TGeoTrap") == 0) {
2210 } else if (strcmp(clsname, "TGeoGtra") == 0) {
2212 } else if (strcmp(clsname, "TGeoTrd1") == 0) {
2214 } else if (strcmp(clsname, "TGeoTrd2") == 0) {
2216 } else if (strcmp(clsname, "TGeoTubeSeg") == 0) {
2218 } else if (strcmp(clsname, "TGeoCtub") == 0) {
2220 } else if (strcmp(clsname, "TGeoTube") == 0) {
2222 } else if (strcmp(clsname, "TGeoPcon") == 0) {
2224 } else if (strcmp(clsname, "TGeoTorus") == 0) {
2226 } else if (strcmp(clsname, "TGeoPgon") == 0) {
2228 } else if (strcmp(clsname, "TGeoEltu") == 0) {
2230 } else if (strcmp(clsname, "TGeoHype") == 0) {
2232 } else if (strcmp(clsname, "TGeoXtru") == 0) {
2234 } else if (strcmp(clsname, "TGeoTessellated") == 0) {
2236 } else if (strcmp(clsname, "TGeoScaledShape") == 0) {
2238 } else if (strcmp(clsname, "TGeoCompositeShape") == 0) {
2240 } else if (strcmp(clsname, "TGeoUnion") == 0) {
2242 } else if (strcmp(clsname, "TGeoIntersection") == 0) {
2244 } else if (strcmp(clsname, "TGeoSubtraction") == 0) {
2246 } else {
2247 Info("ChooseObject", "ERROR! %s Solid CANNOT be processed, solid is NOT supported", clsname);
2248 solidN = nullptr;
2249 }
2250 if (solidN == nullptr) {
2251 if (fNameList->fLst[TString::Format("%p", geoShape)] == "") {
2252 TString missingName = geoShape->GetName();
2253 GenName("missing_" + missingName, TString::Format("%p", geoShape));
2254 } else {
2256 "missing_" + fNameList->fLst[TString::Format("%p", geoShape)];
2257 }
2258 }
2259
2260 return solidN;
2261}
2262
2263////////////////////////////////////////////////////////////////////////////////
2264/// Retrieves X Y Z angles from rotation matrix
2265
2267{
2269 Double_t a, b, c;
2270 Double_t rad = 180.0 / TMath::ACos(-1.0);
2271 const Double_t *r = rotationMatrix;
2272 Double_t cosb = TMath::Sqrt(r[0] * r[0] + r[1] * r[1]);
2273 if (cosb > 0.00001) {
2274 a = TMath::ATan2(r[5], r[8]) * rad;
2275 b = TMath::ATan2(-r[2], cosb) * rad;
2276 c = TMath::ATan2(r[1], r[0]) * rad;
2277 } else {
2278 a = TMath::ATan2(-r[7], r[4]) * rad;
2279 b = TMath::ATan2(-r[2], cosb) * rad;
2280 c = 0;
2281 }
2282 lxyz.x = a;
2283 lxyz.y = b;
2284 lxyz.z = c;
2285 return lxyz;
2286}
2287
2288////////////////////////////////////////////////////////////////////////////////
2289/// Method creating cutTube as an intersection of tube and two boxes
2290/// - not used anymore because cutube is supported in Geant4 9.5
2291
2293{
2294 Double_t rmin = geoShape->GetRmin();
2295 Double_t rmax = geoShape->GetRmax();
2296 Double_t z = geoShape->GetDz();
2297 Double_t startphi = geoShape->GetPhi1();
2298 Double_t deltaphi = geoShape->GetPhi2();
2299 Double_t x1 = geoShape->GetNlow()[0];
2300 Double_t y1 = geoShape->GetNlow()[1];
2301 Double_t z1 = geoShape->GetNlow()[2];
2302 Double_t x2 = geoShape->GetNhigh()[0];
2303 Double_t y2 = geoShape->GetNhigh()[1];
2304 Double_t z2 = geoShape->GetNhigh()[2];
2305 TString xname = geoShape->GetName();
2306
2307 Double_t h0 = 2. * ((TGeoBBox *)geoShape)->GetDZ();
2308 Double_t h1 = 2 * z;
2309 Double_t h2 = 2 * z;
2310 Double_t boxdx = 1E8 * (2 * rmax) + (2 * z);
2311
2312 TGeoTubeSeg *T = new TGeoTubeSeg((xname + "T").Data(), rmin, rmax, h0, startphi, deltaphi);
2313 TGeoBBox *B1 = new TGeoBBox((xname + "B1").Data(), boxdx, boxdx, h1);
2314 TGeoBBox *B2 = new TGeoBBox((xname + "B2").Data(), boxdx, boxdx, h2);
2315
2316 // first box position parameters
2318 Double_t theta1 = 360 - TMath::ATan2(sqrt(x1 * x1 + y1 * y1), z1) * TMath::RadToDeg();
2319
2321 Double_t theta11 = TMath::ATan2(z1, sqrt(x1 * x1 + y1 * y1)) * TMath::RadToDeg();
2322
2326
2327 // second box position parameters
2329 Double_t theta2 = 360 - TMath::ATan2(sqrt(x2 * x2 + y2 * y2), z2) * TMath::RadToDeg();
2330
2332 Double_t theta21 = TMath::ATan2(z2, sqrt(x2 * x2 + y2 * y2)) * TMath::RadToDeg();
2333
2336 Double_t zpos2 = h2 * TMath::Sin((theta21)*TMath::DegToRad()) * (-1);
2337
2338 // positioning
2339 TGeoTranslation *t0 = new TGeoTranslation(0, 0, 0);
2340 TGeoTranslation *t1 = new TGeoTranslation(0 + xpos1, 0 + ypos1, 0 + (zpos1 - z));
2341 TGeoTranslation *t2 = new TGeoTranslation(0 + xpos2, 0 + ypos2, 0 + (zpos2 + z));
2342 TGeoRotation *r0 = new TGeoRotation((xname + "r0").Data());
2343 TGeoRotation *r1 = new TGeoRotation((xname + "r1").Data());
2344 TGeoRotation *r2 = new TGeoRotation((xname + "r2").Data());
2345
2346 r1->SetAngles(phi1, theta1, 0);
2347 r2->SetAngles(phi2, theta2, 0);
2348
2349 TGeoMatrix *m0 = new TGeoCombiTrans(*t0, *r0);
2350 TGeoMatrix *m1 = new TGeoCombiTrans(*t1, *r1);
2351 TGeoMatrix *m2 = new TGeoCombiTrans(*t2, *r2);
2352
2353 TGeoCompositeShape *CS1 = new TGeoCompositeShape((xname + "CS1").Data(), new TGeoIntersection(T, B1, m0, m1));
2354 TGeoCompositeShape *cs = new TGeoCompositeShape((xname + "CS").Data(), new TGeoIntersection(CS1, B2, m0, m2));
2355 delete t0;
2356 delete t1;
2357 delete t2;
2358 delete r0;
2359 delete r1;
2360 delete r2;
2361 return cs;
2362}
2363
2364////////////////////////////////////////////////////////////////////////////////
2365/// Checks whether name2check is in (NameList) list
2366
2368{
2369 Bool_t isIN = list[name2check];
2370 return isIN;
2371}
2372
2373////////////////////////////////////////////////////////////////////////////////
2374/// NCNAME basic restrictions
2375/// Replace "$" character with empty character etc.
2376
2378{
2379 TString newname = oldname.ReplaceAll("$", "");
2380 newname = newname.ReplaceAll(" ", "_");
2381 // :, @, $, %, &, /, +, ,, ;, whitespace characters or different parenthesis
2382 newname = newname.ReplaceAll(":", "");
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("_refl", "");
2396 // workaround if first letter is digit than replace it to "O" (ou character)
2397 TString fstLet = newname(0, 1);
2398 if (fstLet.IsDigit()) {
2399 newname = "O" + newname(1, newname.Length());
2400 }
2401 return newname;
2402}
2403
2404////////////////////////////////////////////////////////////////////////////////
2405/// Important function which is responsible for naming volumes, solids and materials
2406
2408{
2410 if (newname != oldname) {
2411 if (fgkMaxNameErr > fActNameErr) {
2412 Info("GenName", "WARNING! Name of the object was changed because it failed to comply with NCNAME xml datatype "
2413 "restrictions.");
2414 } else if ((fgkMaxNameErr == fActNameErr)) {
2415 Info("GenName", "WARNING! Probably more names are going to be changed to comply with NCNAME xml datatype "
2416 "restriction, but it will not be displayed on the screen.");
2417 }
2418 fActNameErr++;
2419 }
2421 Int_t iter = 0;
2422 switch (fgNamingSpeed) {
2423 case kfastButUglySufix: newname = newname + "0x" + objPointer; break;
2424 case kelegantButSlow:
2425 // 0 means not in the list
2426 iter = fNameList->fLstIter[newname];
2427 if (iter == 0) {
2428 nameIter = "";
2429 } else {
2430 nameIter = TString::Format("0x%i", iter);
2431 }
2434 break;
2436 // no change
2437 break;
2438 }
2439 // store the name (mapped to pointer)
2441 return newname;
2442}
2443
2444////////////////////////////////////////////////////////////////////////////////
2445/// Method which tests whether solids can be processed
2446
2448{
2450 isProcessed = pointer->TestBit(fgkProcBit);
2451 pointer->SetBit(fgkProcBit, kTRUE);
2452 return !(isProcessed);
2453}
2454
2455////////////////////////////////////////////////////////////////////////////////
2456/// Method that retrieves axis and unit along which object is divided
2457
2459{
2461 unit = fDefault_lunit;
2462 switch (divAxis) {
2463 case 1:
2464 if (strcmp(pattName, "TGeoPatternX") == 0) {
2465 return "kXAxis";
2466 } else if (strcmp(pattName, "TGeoPatternCylR") == 0) {
2467 return "kRho";
2468 }
2469 break;
2470 case 2:
2471 if (strcmp(pattName, "TGeoPatternY") == 0) {
2472 return "kYAxis";
2473 } else if (strcmp(pattName, "TGeoPatternCylPhi") == 0) {
2474 unit = "deg";
2475 return "kPhi";
2476 }
2477 break;
2478 case 3:
2479 if (strcmp(pattName, "TGeoPatternZ") == 0) {
2480 return "kZAxis";
2481 }
2482 break;
2483 default: return "kUndefined"; break;
2484 }
2485 return "kUndefined";
2486}
2487
2488////////////////////////////////////////////////////////////////////////////////
2489/// Check for null parameter to skip the NULL objects
2490
2492{
2493 if (parValue == 0.) {
2494 Info("IsNullParam", "ERROR! %s is NULL due to %s = %.12g, Volume based on this shape will be skipped",
2495 objName.Data(), parName.Data(), parValue);
2496 return kTRUE;
2497 }
2498 return kFALSE;
2499}
2500
2501////////////////////////////////////////////////////////////////////////////////
2502/// Unsetting bits that were changed in gGeoManager during export so that export
2503/// can be run more times with the same instance of gGeoManager.
2504
2506{
2507 TIter next(geoMng->GetListOfVolumes());
2508 TGeoVolume *vol;
2509 while ((vol = (TGeoVolume *)next())) {
2510 ((TObject *)vol->GetShape())->SetBit(fgkProcBit, kFALSE);
2512 }
2513}
2514
2515////////////////////////////////////////////////////////////////////////////////
2516//
2517// Backwards compatibility for old DD4hep version (to be removed in the future)
2518//
2519////////////////////////////////////////////////////////////////////////////////
2520
2521////////////////////////////////////////////////////////////////////////////////
2522// Backwards compatibility (to be removed in the future): Wrapper to only selectively write one branch
2524{
2525 TList materials, volumes, nodes;
2526 MaterialExtractor extract;
2527 if (!volume) {
2528 Info("WriteGDMLfile", "Invalid Volume reference to extract GDML information!");
2529 return;
2530 }
2531 extract(volume);
2532 for (TGeoMaterial *m : extract.materials)
2533 materials.Add(m);
2534 fTopVolumeName = volume->GetName();
2535 fTopVolume = volume;
2536 fSurfaceList.clear();
2537 fVolumeList.clear();
2538 fNodeList.clear();
2539 WriteGDMLfile(geomanager, volume, &materials, filename, option);
2540 materials.Clear("nodelete");
2541 volumes.Clear("nodelete");
2542 nodes.Clear("nodelete");
2543}
2544
2545////////////////////////////////////////////////////////////////////////////////
2546/// Wrapper of all exporting methods
2547/// Creates blank GDML file and fills it with gGeoManager structure converted
2548/// to GDML structure of xml nodes
2549
2552{
2553 // option processing
2554 option.ToLower();
2555 if (option.Contains("g")) {
2557 Info("WriteGDMLfile", "Geant4 compatibility mode set");
2558 } else {
2560 }
2561 if (option.Contains("f")) {
2563 Info("WriteGDMLfile", "Fast naming convention with pointer suffix set");
2564 } else if (option.Contains("n")) {
2566 Info("WriteGDMLfile", "Naming without prefix set - be careful uniqness of name is not ensured");
2567 } else {
2569 Info("WriteGDMLfile", "Potentially slow with incremental suffix naming convention set");
2570 }
2571
2572 // local variables
2573 Int_t outputLayout = 1;
2574 const char *krootNodeName = "gdml";
2575 const char *knsRefGeneral = "http://www.w3.org/2001/XMLSchema-instance";
2576 const char *knsNameGeneral = "xsi";
2577 const char *knsRefGdml = "http://service-spi.web.cern.ch/service-spi/app/releases/GDML/schema/gdml.xsd";
2578 const char *knsNameGdml = "xsi:noNamespaceSchemaLocation";
2579
2580 // First create engine
2581 fGdmlE = new TXMLEngine;
2583
2584 // create blank GDML file
2585 fGdmlFile = fGdmlE->NewDoc();
2586
2587 // create root node and add it to blank GDML file
2588 XMLNodePointer_t rootNode = fGdmlE->NewChild(nullptr, nullptr, krootNodeName, nullptr);
2590
2591 // add namespaces to root node
2594
2595 // initialize general lists and <define>, <solids>, <structure> nodes
2596 fIsotopeList = new StructLst;
2597 fElementList = new StructLst;
2598
2599 fNameList = new NameLst;
2600
2601 fDefineNode = fGdmlE->NewChild(nullptr, nullptr, "define", nullptr);
2602 fSolidsNode = fGdmlE->NewChild(nullptr, nullptr, "solids", nullptr);
2603 fStructureNode = fGdmlE->NewChild(nullptr, nullptr, "structure", nullptr);
2604 //========================
2605
2606 // initialize list of accepted patterns for divisions (in ExtractVolumes)
2607 fAccPatt = new StructLst;
2608 fAccPatt->fLst["TGeoPatternX"] = kTRUE;
2609 fAccPatt->fLst["TGeoPatternY"] = kTRUE;
2610 fAccPatt->fLst["TGeoPatternZ"] = kTRUE;
2611 fAccPatt->fLst["TGeoPatternCylR"] = kTRUE;
2612 fAccPatt->fLst["TGeoPatternCylPhi"] = kTRUE;
2613 //========================
2614
2615 // initialize list of rejected shapes for divisions (in ExtractVolumes)
2616 fRejShape = new StructLst;
2617 // this shapes are rejected because, it is not possible to divide trd2
2618 // in Y axis and while only trd2 object is imported from GDML
2619 // it causes a problem when TGeoTrd1 is divided in Y axis
2620 fRejShape->fLst["TGeoTrd1"] = kTRUE;
2621 fRejShape->fLst["TGeoTrd2"] = kTRUE;
2622 //=========================
2623
2624 // Initialize global counters
2625 fActNameErr = 0;
2626 fVolCnt = 0;
2627 fPhysVolCnt = 0;
2628 fSolCnt = 0;
2629
2630 // calling main extraction functions (with measuring time)
2631 time_t startT, endT;
2632 startT = time(nullptr);
2633 ExtractMatrices(geomanager->GetListOfGDMLMatrices());
2636
2637 Info("WriteGDMLfile", "Extracting volumes");
2638 ExtractVolumes(volume);
2639 Info("WriteGDMLfile", "%i solids added", fSolCnt);
2640 Info("WriteGDMLfile", "%i volumes added", fVolCnt);
2641 Info("WriteGDMLfile", "%i physvolumes added", fPhysVolCnt);
2642 ExtractSkinSurfaces(geomanager->GetListOfSkinSurfaces());
2643 ExtractBorderSurfaces(geomanager->GetListOfBorderSurfaces());
2644 ExtractOpticalSurfaces(geomanager->GetListOfOpticalSurfaces());
2645 endT = time(nullptr);
2646 //<gdml>
2647 fGdmlE->AddChild(rootNode, fDefineNode); // <define>...</define>
2648 fGdmlE->AddChild(rootNode, fMaterialsNode); // <materials>...</materials>
2649 fGdmlE->AddChild(rootNode, fSolidsNode); // <solids>...</solids>
2650 fGdmlE->AddChild(rootNode, fStructureNode); // <structure>...</structure>
2651 fGdmlE->AddChild(rootNode, CreateSetupN(fTopVolumeName.Data())); // <setup>...</setup>
2652 //</gdml>
2654 TString tdiffS = (tdiffI == 0 ? TString("< 1 s") : TString::Format("%.0lf s", tdiffI));
2655 Info("WriteGDMLfile", "Exporting time: %s", tdiffS.Data());
2656 //=========================
2657
2658 // Saving document
2660 Info("WriteGDMLfile", "File %s saved", filename);
2661 // cleaning
2663 // unset processing bits:
2665 delete fGdmlE;
2666}
2667
2668////////////////////////////////////////////////////////////////////////////////
2669/// Method extracting geometry structure recursively
2670
2672{
2675 TGeoPatternFinder *pattFinder = nullptr;
2678
2679 // create the name for volume/assembly
2680 if (volume == fTopVolume) {
2681 // not needed a special function for generating name
2682 volname = volume->GetName();
2684 // register name to the pointer
2685 fNameList->fLst[TString::Format("%p", volume)] = volname;
2686 } else {
2687 volname = GenName(volume->GetName(), TString::Format("%p", volume));
2688 }
2689
2690 // start to create main volume/assembly node
2691 if (volume->IsAssembly()) {
2693 } else {
2694 // get reference material and add solid to <solids> + get name
2695 matname = fNameList->fLst[TString::Format("%p", volume->GetMaterial())];
2696 solname = ExtractSolid(volume->GetShape());
2697 // If solid is not supported or corrupted
2698 if (solname == "-1") {
2699 Info("ExtractVolumes", "ERROR! %s volume was not added, because solid is either not supported or corrupted",
2700 volname.Data());
2701 // set volume as missing volume
2702 fNameList->fLst[TString::Format("%p", volume)] = "missing_" + volname;
2703 return;
2704 }
2706
2707 // divisionvol can't be in assembly
2708 pattFinder = volume->GetFinder();
2709 // if found pattern
2710 if (pattFinder) {
2711 pattClsName = TString::Format("%s", pattFinder->ClassName());
2712 TString shapeCls = TString::Format("%s", volume->GetShape()->ClassName());
2713 // if pattern in accepted pattern list and not in shape rejected list
2714 if ((fAccPatt->fLst[pattClsName] == kTRUE) && (fRejShape->fLst[shapeCls] != kTRUE)) {
2715 isPattern = kTRUE;
2716 }
2717 }
2718 }
2719 // get all nodes in volume
2720 TObjArray *nodeLst = volume->GetNodes();
2721 TIter next(nodeLst);
2724 Int_t nCnt = 0;
2725 // loop through all nodes
2726 while ((geoNode = (TGeoNode *)next())) {
2727 // get volume of current node and if not processed then process it
2728 TGeoVolume *subvol = geoNode->GetVolume();
2729 if (subvol->TestAttBit(fgkProcBitVol) == kFALSE) {
2730 subvol->SetAttBit(fgkProcBitVol);
2732 }
2733
2734 // volume of this node has to exist because it was processed recursively
2735 TString nodevolname = fNameList->fLst[TString::Format("%p", geoNode->GetVolume())];
2736 if (nodevolname.Contains("missing_")) {
2737 continue;
2738 }
2739 if (nCnt == 0) { // save name of the first node for divisionvol
2741 }
2742
2743 if (isPattern == kFALSE) {
2744 // create name for node
2745 TString nodename, posname, rotname;
2746 nodename = GenName(geoNode->GetName(), TString::Format("%p", geoNode));
2747 nodename = nodename + "in" + volname;
2748
2749 // create name for position and clear rotation
2750 posname = nodename + "pos";
2751 rotname = "";
2752
2753 // position
2754 const Double_t *pos = geoNode->GetMatrix()->GetTranslation();
2755 Xyz nodPos;
2756 nodPos.x = pos[0];
2757 nodPos.y = pos[1];
2758 nodPos.z = pos[2];
2759 childN = CreatePositionN(posname.Data(), nodPos, "position", fDefault_lunit);
2760 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
2761 // Deal with reflection
2762 XMLNodePointer_t scaleN = nullptr;
2763 Double_t lx, ly, lz;
2764 Double_t xangle = 0;
2765 Double_t zangle = 0;
2766 lx = geoNode->GetMatrix()->GetRotationMatrix()[0];
2767 ly = geoNode->GetMatrix()->GetRotationMatrix()[4];
2768 lz = geoNode->GetMatrix()->GetRotationMatrix()[8];
2769 if (geoNode->GetMatrix()->IsReflection() && TMath::Abs(lx) == 1 && TMath::Abs(ly) == 1 &&
2770 TMath::Abs(lz) == 1) {
2771 scaleN = fGdmlE->NewChild(nullptr, nullptr, "scale", nullptr);
2772 fGdmlE->NewAttr(scaleN, nullptr, "name", (nodename + "scl").Data());
2773 fGdmlE->NewAttr(scaleN, nullptr, "x", TString::Format(fltPrecision.Data(), lx));
2774 fGdmlE->NewAttr(scaleN, nullptr, "y", TString::Format(fltPrecision.Data(), ly));
2775 fGdmlE->NewAttr(scaleN, nullptr, "z", TString::Format(fltPrecision.Data(), lz));
2776 // experimentally found out, that rotation should be updated like this
2777 if (lx == -1) {
2778 zangle = 180;
2779 }
2780 if (lz == -1) {
2781 xangle = 180;
2782 }
2783 }
2784
2785 // rotation
2786 TGDMLWrite::Xyz lxyz = GetXYZangles(geoNode->GetMatrix()->GetRotationMatrix());
2787 lxyz.x -= xangle;
2788 lxyz.z -= zangle;
2789 if ((lxyz.x != 0.0) || (lxyz.y != 0.0) || (lxyz.z != 0.0)) {
2790 rotname = nodename + "rot";
2791 childN = CreateRotationN(rotname.Data(), lxyz);
2792 fGdmlE->AddChild(fDefineNode, childN); // adding node to <define> node
2793 }
2794
2795 // create physvol for main volume/assembly node
2797 childN = CreatePhysVolN(physvolname, geoNode->GetNumber(), nodevolname.Data(), posname.Data(), rotname.Data(),
2798 scaleN);
2800 }
2801 nCnt++;
2802 }
2803 // create only one divisionvol node
2804 if (isPattern && pattFinder) {
2805 // retrieve attributes of division
2806 Int_t ndiv, divaxis;
2807 Double_t offset, width, xlo, xhi;
2808 TString axis, unit;
2809
2810 ndiv = pattFinder->GetNdiv();
2811 width = pattFinder->GetStep();
2812
2813 divaxis = pattFinder->GetDivAxis();
2814 volume->GetShape()->GetAxisRange(divaxis, xlo, xhi);
2815
2816 // compute relative start (not positional)
2817 offset = pattFinder->GetStart() - xlo;
2818 axis = GetPattAxis(divaxis, pattClsName, unit);
2819
2820 // create division node
2821 childN = CreateDivisionN(offset, width, ndiv, axis.Data(), unit.Data(), nodeVolNameBak.Data());
2823 }
2824
2825 fVolCnt++;
2826 // add volume/assembly node into the <structure> node
2828}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
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:19
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:151
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:45
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:98
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:399
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:202
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
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1045
Basic string class.
Definition TString.h:138
const char * Data() const
Definition TString.h:384
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:2384
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
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:643
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:657
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:82
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:75
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
TMarker m
Definition textangle.C:8
auto * t1
Definition textangle.C:20