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