ROOT  6.06/09
Reference Guide
TGDMLParse.cxx
Go to the documentation of this file.
1 /* @(#)root/gdml:$Id$ */
2 // Author: Ben Lloyd 09/11/06
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 
13 /*************************************************************************
14 
15 ____________________________________________________________________
16 
17 TGDMLParse Class
18 
19 --------------------------------------------------------------------
20 
21  This class contains the implementation of the GDML parser associated to
22  all the supported GDML elements. User should never need to explicitly
23  instaciate this class. It is internally used by the TGeoManager.
24 
25  Each element process has a 'Binding' to ROOT. The 'binding' is specific
26  mapping of GDML elements (materials, solids, etc) to specific objects which
27  should be instanciated by the converted. In the present case (ROOT) the
28  binding is implemented at the near the end of each process function. Most
29  bindings follow similar format, dependent on what is being added to the
30  geometry.
31 
32  This file also contains the implementation of the TGDMLRefl class. This is
33  just a small helper class used internally by the 'reflection' method (for
34  reflected solids).
35 
36  The presently supported list of TGeo classes is the following:
37 
38  Materials:
39  TGeoElement
40  TGeoMaterial
41  TGeoMixture
42 
43  Solids:
44  TGeoBBox
45  TGeoArb8
46  TGeoTubeSeg
47  TGeoConeSeg
48  TGeoCtub
49  TGeoPcon
50  TGeoTrap
51  TGeoGtra
52  TGeoTrd2
53  TGeoSphere
54  TGeoPara
55  TGeoTorus
56  TGeoHype
57  TGeoPgon
58  TGeoXtru
59  TGeoEltu
60  TGeoParaboloid
61  TGeoCompositeShape (subtraction, union, intersection)
62 
63  Approximated Solids:
64  Ellipsoid (approximated to a TGeoBBox)
65  Elliptical cone (approximated to a TGeoCone)
66 
67  Geometry:
68  TGeoVolume
69  TGeoVolumeAssembly
70  divisions
71  reflection
72 
73 When most solids or volumes are added to the geometry they
74 
75 
76  Whenever a new element is added to GDML schema, this class needs to be extended.
77  The appropriate method (process) needs to be implemented, as well as the new
78  element process then needs to be linked thru the function TGDMLParse
79 
80  For any question or remarks concerning this code, please send an email to
81  ben.lloyd@cern.ch
82 
83 ****************************************************************************/
84 
85 #include "TGeoManager.h"
86 #include "TGeoMatrix.h"
87 #include "TXMLEngine.h"
88 #include "TGeoVolume.h"
89 #include "TGeoBBox.h"
90 #include "TGeoParaboloid.h"
91 #include "TGeoArb8.h"
92 #include "TGeoTube.h"
93 #include "TGeoCone.h"
94 #include "TGeoTrd2.h"
95 #include "TGeoPcon.h"
96 #include "TGeoPgon.h"
97 #include "TGeoSphere.h"
98 #include "TGeoTorus.h"
99 #include "TGeoPara.h"
100 #include "TGeoHype.h"
101 #include "TGeoEltu.h"
102 #include "TGeoXtru.h"
103 #include "TGeoScaledShape.h"
104 #include "TGeoVolume.h"
105 #include "TROOT.h"
106 #include "TMath.h"
107 #include "TMap.h"
108 #include "TObjString.h"
109 #include "TGeoExtension.h"
110 #include "TGeoMaterial.h"
111 #include "TGeoBoolNode.h"
112 #include "TGeoMedium.h"
113 #include "TGeoElement.h"
114 #include "TGeoShape.h"
115 #include "TGeoCompositeShape.h"
116 #include "TGDMLParse.h"
117 #include <stdlib.h>
118 #include <string>
119 
121 
122 ////////////////////////////////////////////////////////////////////////////////
123 ///creates the new instance of the XMLEngine called 'gdml', using the filename >>
124 ///then parses the file and creates the DOM tree. Then passes the DOM to the
125 ///next function to translate it.
126 
127 TGeoVolume* TGDMLParse::GDMLReadFile(const char* filename)
128 {
129  // First create engine
130  TXMLEngine* gdml = new TXMLEngine;
131  gdml->SetSkipComments(kTRUE);
132 
133  // Now try to parse xml file
134  XMLDocPointer_t gdmldoc = gdml->ParseFile(filename);
135  if (gdmldoc == 0) {
136  delete gdml;
137  return 0;
138  } else {
139 
140  // take access to main node
141  XMLNodePointer_t mainnode = gdml->DocGetRootElement(gdmldoc);
142 
143  fFileEngine[fFILENO] = gdml;
144  fStartFile = filename;
145  fCurrentFile = filename;
146 
147  // display recursively all nodes and subnodes
148  ParseGDML(gdml, mainnode);
149 
150  // Release memory before exit
151  gdml->FreeDoc(gdmldoc);
152  delete gdml;
153 
154  }
155  return fWorld;
156 
157 }
158 
159 ////////////////////////////////////////////////////////////////////////////////
160 ///this function recursively moves thru the DOM tree of the GDML file. It checks for
161 ///key words along the way and if a key word is found it calls the corresponding
162 ///function to interpret the node.
163 
165 {
166  XMLAttrPointer_t attr = gdml->GetFirstAttr(node);
167  const char* name = gdml->GetNodeName(node);
168  XMLNodePointer_t parentn = gdml->GetParent(node);
169  const char* parent = gdml->GetNodeName(parentn);
170  XMLNodePointer_t childtmp = 0;
171 
172  const char* posistr = "position";
173  const char* setustr = "setup";
174  const char* consstr = "constant";
175  const char* varistr = "variable";
176  const char* rotastr = "rotation";
177  const char* scalstr = "scale";
178  const char* elemstr = "element";
179  const char* istpstr = "isotope";
180  const char* matestr = "material";
181  const char* volustr = "volume";
182  const char* assestr = "assembly";
183  const char* twtrstr = "twistedtrap"; //name changed according to schema
184  const char* cutTstr = "cutTube";
185  const char* bboxstr = "box";
186  const char* xtrustr = "xtru";
187  const char* arb8str = "arb8";
188  const char* tubestr = "tube";
189  const char* conestr = "cone";
190  const char* polystr = "polycone";
191  const char* hypestr = "hype";
192  const char* trapstr = "trap";
193  const char* trdstr = "trd";
194  const char* sphestr = "sphere";
195  const char* orbstr = "orb";
196  const char* parastr = "para";
197  const char* torustr = "torus";
198  const char* hedrstr = "polyhedra";
199  const char* eltustr = "eltube";
200  const char* subtstr = "subtraction";
201  const char* uniostr = "union";
202  const char* parbstr = "paraboloid";
203  const char* intestr = "intersection";
204  const char* reflstr = "reflectedSolid";
205  const char* ellistr = "ellipsoid";
206  const char* elcnstr = "elcone";
207  const char* usrstr = "userinfo";
208  Bool_t hasIsotopes;
209  Bool_t hasIsotopesExtended;
210 
211  if ((strcmp(name, posistr)) == 0) {
212  node = PosProcess(gdml, node, attr);
213  } else if ((strcmp(name, rotastr)) == 0) {
214  node = RotProcess(gdml, node, attr);
215  } else if ((strcmp(name, scalstr)) == 0) {
216  node = SclProcess(gdml, node, attr);
217  } else if ((strcmp(name, setustr)) == 0) {
218  node = TopProcess(gdml, node);
219  } else if ((strcmp(name, consstr)) == 0) {
220  node = ConProcess(gdml, node, attr);
221  } else if ((strcmp(name, varistr)) == 0) {
222  node = ConProcess(gdml, node, attr);
223  }
224  //*************eleprocess********************************
225 
226  else if (((strcmp(name, "atom")) == 0) && ((strcmp(parent, elemstr)) == 0)) {
227  hasIsotopes = kFALSE;
228  hasIsotopesExtended = kFALSE;
229  node = EleProcess(gdml, node, parentn, hasIsotopes, hasIsotopesExtended);
230  }
231  else if ((strcmp(name, elemstr) == 0) && !gdml->HasAttr(node, "Z")) {
232  hasIsotopes = kTRUE;
233  hasIsotopesExtended = kFALSE;
234  node = EleProcess(gdml, node, parentn, hasIsotopes, hasIsotopesExtended);
235  }
236 
237  else if ((strcmp(name, elemstr) == 0) && gdml->HasAttr(node, "Z")) {
238  childtmp = gdml->GetChild(node);
239  if ((strcmp(gdml->GetNodeName(childtmp), "fraction") == 0) ){
240  hasIsotopes = kFALSE;
241  hasIsotopesExtended = kTRUE;
242  node = EleProcess(gdml, node, parentn, hasIsotopes, hasIsotopesExtended);}
243  }
244 
245  //********isoprocess******************************
246 
247  else if (((strcmp(name, "atom")) == 0) && ((strcmp(parent, istpstr)) == 0)) {
248  node = IsoProcess(gdml, node, parentn);
249  }
250 
251  //********matprocess***********************************
252  else if ((strcmp(name, matestr)) == 0 && gdml->HasAttr(node, "Z")) {
253  childtmp = gdml->GetChild(node);
254 // if ((strcmp(gdml->GetNodeName(childtmp), "fraction") == 0) || (strcmp(gdml->GetNodeName(childtmp), "D") == 0)){
255  // Bool_t frac = kFALSE;
256  Bool_t atom = kFALSE;
257  while(childtmp) {
258  // frac = strcmp(gdml->GetNodeName(childtmp),"fraction")==0;
259  atom = strcmp(gdml->GetNodeName(childtmp),"atom")==0;
260  gdml->ShiftToNext(childtmp);
261  }
262  int z = (atom) ? 1 : 0;
263  node = MatProcess(gdml, node, attr, z);
264  }
265  else if ((strcmp(name, matestr)) == 0 && !gdml->HasAttr(node, "Z")) {
266  int z = 0;
267  node = MatProcess(gdml, node, attr, z);
268  }
269 
270  //*********************************************
271  else if ((strcmp(name, volustr)) == 0) {
272  node = VolProcess(gdml, node);
273  } else if ((strcmp(name, bboxstr)) == 0) {
274  node = Box(gdml, node, attr);
275  } else if ((strcmp(name, ellistr)) == 0) {
276  node = Ellipsoid(gdml, node, attr);
277  } else if ((strcmp(name, elcnstr)) == 0) {
278  node = ElCone(gdml, node, attr);
279  } else if ((strcmp(name, cutTstr)) == 0) {
280  node = CutTube(gdml, node, attr);
281  } else if ((strcmp(name, arb8str)) == 0) {
282  node = Arb8(gdml, node, attr);
283  } else if ((strcmp(name, tubestr)) == 0) {
284  node = Tube(gdml, node, attr);
285  } else if ((strcmp(name, conestr)) == 0) {
286  node = Cone(gdml, node, attr);
287  } else if ((strcmp(name, polystr)) == 0) {
288  node = Polycone(gdml, node, attr);
289  } else if ((strcmp(name, trapstr)) == 0) {
290  node = Trap(gdml, node, attr);
291  } else if ((strcmp(name, trdstr)) == 0) {
292  node = Trd(gdml, node, attr);
293  } else if ((strcmp(name, sphestr)) == 0) {
294  node = Sphere(gdml, node, attr);
295  } else if ((strcmp(name, xtrustr)) == 0) {
296  node = Xtru(gdml, node, attr);
297  } else if ((strcmp(name, twtrstr)) == 0) {
298  node = TwistTrap(gdml, node, attr);
299  } else if ((strcmp(name, hypestr)) == 0) {
300  node = Hype(gdml, node, attr);
301  } else if ((strcmp(name, orbstr)) == 0) {
302  node = Orb(gdml, node, attr);
303  } else if ((strcmp(name, parastr)) == 0) {
304  node = Para(gdml, node, attr);
305  } else if ((strcmp(name, torustr)) == 0) {
306  node = Torus(gdml, node, attr);
307  } else if ((strcmp(name, eltustr)) == 0) {
308  node = ElTube(gdml, node, attr);
309  } else if ((strcmp(name, hedrstr)) == 0) {
310  node = Polyhedra(gdml, node, attr);
311  } else if ((strcmp(name, parbstr)) == 0) {
312  node = Paraboloid(gdml, node, attr);
313  } else if ((strcmp(name, subtstr)) == 0) {
314  node = BooSolid(gdml, node, attr, 1);
315  } else if ((strcmp(name, intestr)) == 0) {
316  node = BooSolid(gdml, node, attr, 2);
317  } else if ((strcmp(name, uniostr)) == 0) {
318  node = BooSolid(gdml, node, attr, 3);
319  } else if ((strcmp(name, reflstr)) == 0) {
320  node = Reflection(gdml, node, attr);
321  } else if ((strcmp(name, assestr)) == 0) {
322  node = AssProcess(gdml, node);
323  } else if ((strcmp(name, usrstr)) == 0) {
324  node = UsrProcess(gdml, node);
325  //CHECK FOR TAGS NOT SUPPORTED
326  } else if (((strcmp(name, "gdml")) != 0) && ((strcmp(name, "define")) != 0) &&
327  ((strcmp(name, "element")) != 0) && ((strcmp(name, "materials")) != 0) &&
328  ((strcmp(name, "solids")) != 0) && ((strcmp(name, "structure")) != 0) &&
329  ((strcmp(name, "zplane")) != 0) && ((strcmp(name, "first")) != 0) &&
330  ((strcmp(name, "second")) != 0) && ((strcmp(name, "twoDimVertex")) != 0) &&
331  ((strcmp(name, "firstposition")) != 0) && ((strcmp(name, "firstpositionref")) != 0) &&
332  ((strcmp(name, "firstrotation")) != 0) && ((strcmp(name, "firstrotationref")) != 0) &&
333  ((strcmp(name, "section")) != 0) && ((strcmp(name, "world")) != 0) &&
334  ((strcmp(name, "isotope")) != 0)) {
335  std::cout << "Error: Unsupported GDML Tag Used :" << name << ". Please Check Geometry/Schema." << std::endl;
336  }
337 
338  // Check for Child node - if present call this funct. recursively until no more
339 
340  XMLNodePointer_t child = gdml->GetChild(node);
341  while (child != 0) {
342  ParseGDML(gdml, child);
343  child = gdml->GetNext(child);
344  }
345 
346  return fWorldName;
347 
348 }
349 
350 ////////////////////////////////////////////////////////////////////////////////
351 
352 double TGDMLParse::Evaluate(const char* evalline)
353 {
354  //takes a string containing a mathematical expression and returns the value of the expression
355 
356  return TFormula("TFormula", evalline).Eval(0);
357 }
358 
359 ////////////////////////////////////////////////////////////////////////////////
360 ///When using the 'divide' process in the geometry this function
361 ///sets the variable 'axis' depending on what is specified.
362 
363 Int_t TGDMLParse::SetAxis(const char* axisString)
364 {
365  Int_t axis = 0;
366 
367  if ((strcmp(axisString, "kXAxis")) == 0) {
368  axis = 1;
369  } else if ((strcmp(axisString, "kYAxis")) == 0) {
370  axis = 2;
371  } else if ((strcmp(axisString, "kZAxis")) == 0) {
372  axis = 3;
373  } else if ((strcmp(axisString, "kRho")) == 0) {
374  axis = 1;
375  } else if ((strcmp(axisString, "kPhi")) == 0) {
376  axis = 2;
377  }
378 
379  return axis;
380 }
381 
382 ////////////////////////////////////////////////////////////////////////////////
383 ///this function looks thru a string for the chars '0x' next to
384 ///each other, when it finds this, it calls another function to strip
385 ///the hex address. It does this recursively until the end of the
386 ///string is reached, returning a string without any hex addresses.
387 
388 const char* TGDMLParse::NameShort(const char* name)
389 {
390  static TString stripped;
391  stripped = name;
392  Int_t index = stripped.Index("0x");
393  if (index >= 0) stripped = stripped(0, index);
394  return stripped.Data();
395 }
396 
397 ////////////////////////////////////////////////////////////////////////////////
398 ///In the define section of the GDML file, constants can be declared.
399 ///when the constant keyword is found, this function is called, and the
400 ///name and value of the constant is stored in the "fformvec" vector as
401 ///a TFormula class, representing a constant function
402 
404 {
405  TString name = "";
406  TString value = "";
407  TString tempattr;
408 
409  while (attr != 0) {
410  tempattr = gdml->GetAttrName(attr);
411  tempattr.ToLower();
412 
413  if (tempattr == "name") {
414  name = gdml->GetAttrValue(attr);
415  }
416  if (tempattr == "value") {
417  value = gdml->GetAttrValue(attr);
418  }
419  attr = gdml->GetNextAttr(attr);
420  }
421 
422  //if ((strcmp(fCurrentFile, fStartFile)) != 0) {
423  // name = TString::Format("%s_%s", name.Data(), fCurrentFile);
424  //}
425 
426  fformvec.push_back(new TFormula(name, value));
427 
428  return node;
429 }
430 ////////////////////////////////////////////////////////////////////////////////
431 ///Throughout the GDML file, a unit can de specified. Whether it be
432 ///angular or linear, values can be used as well as abbreviations such as
433 /// 'mm' or 'deg'. This function is passed the specified unit and if it is
434 ///found, replaces it with the appropriate value.
435 
436 TString TGDMLParse::GetScale(const char* unit)
437 {
438  TString retunit = "";
439 
440  if (strcmp(unit, "mm") == 0) {
441  retunit = "0.1";
442  } else if (strcmp(unit, "milimeter") == 0) {
443  retunit = "0.1";
444  } else if (strcmp(unit, "cm") == 0) {
445  retunit = "1.0";
446  } else if (strcmp(unit, "centimeter") == 0) {
447  retunit = "1.0";
448  } else if (strcmp(unit, "m") == 0) {
449  retunit = "100.0";
450  } else if (strcmp(unit, "meter") == 0) {
451  retunit = "100.0";
452  } else if (strcmp(unit, "km") == 0) {
453  retunit = "100000.0";
454  } else if (strcmp(unit, "kilometer") == 0) {
455  retunit = "100000.0";
456  } else if (strcmp(unit, "rad") == 0) {
457  retunit = TString::Format("%.12f", TMath::RadToDeg());
458  } else if (strcmp(unit, "radian") == 0) {
459  retunit = TString::Format("%.12f", TMath::RadToDeg());
460  } else if (strcmp(unit, "deg") == 0) {
461  retunit = "1.0";
462  } else if (strcmp(unit, "degree") == 0) {
463  retunit = "1.0";
464  } else if (strcmp(unit, "pi") == 0) {
465  retunit = "pi";
466  } else if (strcmp(unit, "avogadro") == 0) {
467  retunit = TString::Format("%.12g", TMath::Na());
468  } else {
469  Fatal("GetScale", "Unit <%s> not known", unit);
470  retunit = "0";
471  }
472  return retunit;
473 
474 }
475 
476 ////////////////////////////////////////////////////////////////////////////////
477 ///Throughout the GDML file, a unit can de specified. Whether it be
478 ///angular or linear, values can be used as well as abbreviations such as
479 /// 'mm' or 'deg'. This function is passed the specified unit and if it is
480 ///found, replaces it with the appropriate value.
481 
483 {
484  Double_t retunit = 0.;
485  TString unit(sunit);
486  unit.ToLower();
487 
488  if ((unit == "mm") || (unit == "milimeter")) {
489  retunit = 0.1;
490  } else if ((unit == "cm") || (unit == "centimeter")) {
491  retunit = 1.0;
492  } else if ((unit == "m") || (unit == "meter")) {
493  retunit = 100.0;
494  } else if ((unit == "km") || (unit == "kilometer")) {
495  retunit = 100000.0;
496  } else if ((unit == "rad") || (unit == "radian")) {
497  retunit = TMath::RadToDeg();
498  } else if ((unit == "deg") || (unit == "degree")) {
499  retunit = 1.0;
500  } else if ((unit == "ev") || (unit == "electronvolt")) {
501  retunit = 0.000000001;
502  } else if ((unit == "kev") || (unit == "kiloelectronvolt")) {
503  retunit = 0.000001;
504  } else if ((unit == "mev") || (unit == "megaelectronvolt")) {
505  retunit = 0.001;
506  } else if ((unit == "gev") || (unit == "gigaelectronvolt")) {
507  retunit = 1;
508  } else if (unit == "pi") {
509  retunit = TMath::Pi();
510  } else if (unit == "avogadro") {
511  retunit = TMath::Na();
512  } else {
513  Fatal("GetScaleVal", "Unit <%s> not known", sunit);
514  retunit = 0;
515  }
516  return retunit;
517 }
518 
519 ////////////////////////////////////////////////////////////////////////////////
520 /// Convert number in string format to double value.
521 
522 Double_t TGDMLParse::Value(const char* svalue) const
523 {
524  static TString s;
525  s = svalue;
526  return s.Atof();
527 }
528 
529 ////////////////////////////////////////////////////////////////////////////////
530 ///In the define section of the GDML file, positions can be declared.
531 ///when the position keyword is found, this function is called, and the
532 ///name and values of the position are converted into type TGeoPosition
533 ///and stored in fposmap map using the name as its key. This function
534 ///can also be called when declaring solids.
535 
537 {
538  TString lunit = "mm";
539  TString xpos = "0";
540  TString ypos = "0";
541  TString zpos = "0";
542  TString name = "0";
543  TString tempattr;
544 
545  while (attr != 0) {
546 
547  tempattr = gdml->GetAttrName(attr);
548  tempattr.ToLower();
549 
550  if (tempattr == "name") {
551  name = gdml->GetAttrValue(attr);
552  } else if (tempattr == "x") {
553  xpos = gdml->GetAttrValue(attr);
554  } else if (tempattr == "y") {
555  ypos = gdml->GetAttrValue(attr);
556  } else if (tempattr == "z") {
557  zpos = gdml->GetAttrValue(attr);
558  } else if (tempattr == "unit") {
559  lunit = gdml->GetAttrValue(attr);
560  }
561 
562  attr = gdml->GetNextAttr(attr);
563  }
564 
565  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
566  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
567  }
568 
569  Double_t retunit = GetScaleVal(lunit);
570  Double_t xline = Value(xpos)*retunit;
571  Double_t yline = Value(ypos)*retunit;
572  Double_t zline = Value(zpos)*retunit;
573 
574  TGeoTranslation* pos = new TGeoTranslation(xline, yline, zline);
575 
576  fposmap[name.Data()] = pos;
577 
578  return node;
579 
580 }
581 
582 ////////////////////////////////////////////////////////////////////////////////
583 ///In the define section of the GDML file, rotations can be declared.
584 ///when the rotation keyword is found, this function is called, and the
585 ///name and values of the rotation are converted into type TGeoRotation
586 ///and stored in frotmap map using the name as its key. This function
587 ///can also be called when declaring solids.
588 
590 {
591  TString aunit = "rad";
592  TString xpos = "0";
593  TString ypos = "0";
594  TString zpos = "0";
595  TString name = "";
596  TString tempattr;
597 
598  while (attr != 0) {
599 
600  tempattr = gdml->GetAttrName(attr);
601  tempattr.ToLower();
602 
603  if (tempattr == "name") {
604  name = gdml->GetAttrValue(attr);
605  } else if (tempattr == "x") {
606  xpos = gdml->GetAttrValue(attr);
607  } else if (tempattr == "y") {
608  ypos = gdml->GetAttrValue(attr);
609  } else if (tempattr == "z") {
610  zpos = gdml->GetAttrValue(attr);
611  } else if (tempattr == "unit") {
612  aunit = gdml->GetAttrValue(attr);
613  }
614 
615  attr = gdml->GetNextAttr(attr);
616  }
617 
618  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
619  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
620  }
621 
622  Double_t retunit = GetScaleVal(aunit);
623 
624  Double_t xline = Value(xpos)*retunit;
625  Double_t yline = Value(ypos)*retunit;
626  Double_t zline = Value(zpos)*retunit;
627 
628  TGeoRotation* rot = new TGeoRotation();
629 
630  rot->RotateZ(-zline);
631  rot->RotateY(-yline);
632  rot->RotateX(-xline);
633 
634  frotmap[name.Data()] = rot;
635 
636  return node;
637 
638 }
639 
640 ////////////////////////////////////////////////////////////////////////////////
641 ///In the define section of the GDML file, rotations can be declared.
642 ///when the scale keyword is found, this function is called, and the
643 ///name and values of the scale are converted into type TGeoScale
644 ///and stored in fsclmap map using the name as its key. This function
645 ///can also be called when declaring solids.
646 
648 {
649  TString xpos = "0";
650  TString ypos = "0";
651  TString zpos = "0";
652  TString name = "";
653  TString tempattr;
654 
655  while (attr != 0) {
656 
657  tempattr = gdml->GetAttrName(attr);
658  tempattr.ToLower();
659 
660  if (tempattr == "name") {
661  name = gdml->GetAttrValue(attr);
662  } else if (tempattr == "x") {
663  xpos = gdml->GetAttrValue(attr);
664  } else if (tempattr == "y") {
665  ypos = gdml->GetAttrValue(attr);
666  } else if (tempattr == "z") {
667  zpos = gdml->GetAttrValue(attr);
668  }
669 
670  attr = gdml->GetNextAttr(attr);
671  }
672 
673  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
674  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
675  }
676 
677  TGeoScale* scl = new TGeoScale(Value(xpos), Value(ypos), Value(zpos));
678 
679  fsclmap[name.Data()] = scl;
680 
681  return node;
682 }
683 
684 ////////////////////////////////////////////////////////////////////////////////
685 ///In the material section of the GDML file, an isotope may be declared.
686 ///when the isotope keyword is found, this function is called, and the
687 ///required parameters are taken and stored, these are then bound and
688 ///converted to type TGeoIsotope and stored in fisomap map using the name
689 ///as its key.
690 
692 {
693  TString z = "0";
694  TString name = "";
695  TString n = "0";
696  TString atom = "0";
697  TString tempattr;
698 
699  //obtain attributes for the element
700 
701  XMLAttrPointer_t attr = gdml->GetFirstAttr(parentn);
702 
703  while (attr != 0) {
704 
705  tempattr = gdml->GetAttrName(attr);
706  tempattr.ToLower();
707 
708  if (tempattr == "name") {
709  name = gdml->GetAttrValue(attr);
710  } else if (tempattr == "z") {
711  z = gdml->GetAttrValue(attr);
712  } else if (tempattr == "n") {
713  n = gdml->GetAttrValue(attr);
714  }
715 
716  attr = gdml->GetNextAttr(attr);
717  }
718 
719  //get the atom value for the element
720 
721  attr = gdml->GetFirstAttr(node);
722 
723  while (attr != 0) {
724 
725  tempattr = gdml->GetAttrName(attr);
726 
727  if (tempattr == "value") {
728  atom = gdml->GetAttrValue(attr);
729  }
730 
731  attr = gdml->GetNextAttr(attr);
732  }
733 
734  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
735  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
736  }
737 
738  Int_t z2 = (Int_t)Value(z);
739  Int_t n2 = (Int_t)Value(n);
740  Double_t atom2 = Value(atom);
741 
742  TGeoIsotope* iso = new TGeoIsotope(NameShort(name), z2 , n2, atom2);
743  fisomap[name.Data()] = iso;
744 
745  return node;
746 
747 }
748 
749 //___________________________________________________________
750 //XMLNodePointer_t TGDMLParse::EleProcess(TXMLEngine* gdml, XMLNodePointer_t node, XMLNodePointer_t parentn, Bool_t hasIsotopes)
751 XMLNodePointer_t TGDMLParse::EleProcess(TXMLEngine* gdml, XMLNodePointer_t node, XMLNodePointer_t parentn, Bool_t hasIsotopes, Bool_t hasIsotopesExtended)
752 
753 {
754 
755  //when the element keyword is found, this function is called, and the
756  //name and values of the element are converted into type TGeoElement and
757  //stored in felemap map using the name as its key.
758 
759  TString z = "0";
760  TString name = "";
761  TString formula = "";
762  TString atom = "0";
763  TString tempattr;
764  Int_t ncompo = 0;
765  typedef FracMap::iterator fractions;
766  FracMap fracmap;
767 
768  XMLNodePointer_t child = 0;
769 
770  //obtain attributes for the element
771 
772  XMLAttrPointer_t attr = gdml->GetFirstAttr(node);
773 
774  if (hasIsotopes) {
775 
776  // Get the name of the element
777  while (attr != 0) {
778  tempattr = gdml->GetAttrName(attr);
779  if (tempattr == "name") {
780  name = gdml->GetAttrValue(attr);
781 
782  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
783  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
784  }
785  break;
786  }
787  attr = gdml->GetNextAttr(attr);
788  }
789  // Get component isotopes. Loop all children.
790  child = gdml->GetChild(node);
791  while (child != 0) {
792 
793  // Check for fraction node name
794  if ((strcmp(gdml->GetNodeName(child), "fraction")) == 0) {
795  Double_t n = 0;
796  TString ref = "";
797  ncompo = ncompo + 1;
798  attr = gdml->GetFirstAttr(child);
799  while (attr != 0) {
800  tempattr = gdml->GetAttrName(attr);
801  tempattr.ToLower();
802  if (tempattr == "n") {
803  n = Value(gdml->GetAttrValue(attr));
804  } else if (tempattr == "ref") {
805  ref = gdml->GetAttrValue(attr);
806  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
807  ref = TString::Format("%s_%s", ref.Data(), fCurrentFile);
808  }
809  }
810  attr = gdml->GetNextAttr(attr);
811  } // loop on child attributes
812  fracmap[ref.Data()] = n;
813  }
814  child = gdml->GetNext(child);
815  } // loop on children
816  // Create TGeoElement - note: Object(name, title) corresponds to Element(formula, name)
817  TGeoElement *ele = new TGeoElement(NameShort(name), NameShort(name), ncompo);
818  for (fractions f = fracmap.begin(); f != fracmap.end(); f++) {
819  if (fisomap.find(f->first) != fisomap.end()) {
820  ele->AddIsotope((TGeoIsotope*)fisomap[f->first], f->second);
821  }
822  }
823  felemap[name.Data()] = ele;
824  return child;
825  } // hasisotopes end loop
826 
827  //*************************
828 
829 
830  if (hasIsotopesExtended) {
831 
832  while (attr != 0) {
833  tempattr = gdml->GetAttrName(attr);
834 
835  if (tempattr == "name") {
836  name = gdml->GetAttrValue(attr);
837 
838  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
839  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
840  }
841  break;
842  }
843  attr = gdml->GetNextAttr(attr);
844  }
845  // Get component isotopes. Loop all children.
846  child = gdml->GetChild(node);
847  while (child != 0) {
848 
849  // Check for fraction node name
850  if ((strcmp(gdml->GetNodeName(child), "fraction")) == 0) {
851  Double_t n = 0;
852  TString ref = "";
853  ncompo = ncompo + 1;
854  attr = gdml->GetFirstAttr(child);
855  while (attr != 0) {
856  tempattr = gdml->GetAttrName(attr);
857  tempattr.ToLower();
858  if (tempattr == "n") {
859  n = Value(gdml->GetAttrValue(attr));
860  } else if (tempattr == "ref") {
861  ref = gdml->GetAttrValue(attr);
862  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
863  ref = TString::Format("%s_%s", ref.Data(), fCurrentFile);
864  }
865  }
866  attr = gdml->GetNextAttr(attr);
867  } // loop on child attributes
868  fracmap[ref.Data()] = n;
869  }
870  child = gdml->GetNext(child);
871  } // loop on children
872  // Create TGeoElement - note: Object(name, title) corresponds to Element(formula, name)
873  TGeoElement *ele = new TGeoElement(NameShort(name), NameShort(name), ncompo);
874  for (fractions f = fracmap.begin(); f != fracmap.end(); f++) {
875  if (fisomap.find(f->first) != fisomap.end()) {
876  ele->AddIsotope((TGeoIsotope*)fisomap[f->first], f->second);
877  }
878  }
879  felemap[name.Data()] = ele;
880  return child;
881  } // hasisotopesExtended end loop
882 
883 
884 
885  //***************************
886 
887  attr = gdml->GetFirstAttr(parentn);
888  while (attr != 0) {
889 
890  tempattr = gdml->GetAttrName(attr);
891  tempattr.ToLower();
892 
893  if (tempattr == "name") {
894  name = gdml->GetAttrValue(attr);
895 
896  } else if (tempattr == "z") {
897  z = gdml->GetAttrValue(attr);
898  } else if (tempattr == "formula") {
899  formula = gdml->GetAttrValue(attr);
900  }
901 
902  attr = gdml->GetNextAttr(attr);
903  }
904 
905  //get the atom value for the element
906 
907  attr = gdml->GetFirstAttr(node);
908 
909  while (attr != 0) {
910 
911  tempattr = gdml->GetAttrName(attr);
912  tempattr.ToLower();
913 
914  if (tempattr == "value") {
915  atom = gdml->GetAttrValue(attr);
916  }
917 
918  attr = gdml->GetNextAttr(attr);
919  }
920 
921  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
922  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
923  }
924 
925  Int_t z2 = (Int_t)Value(z);
926  Double_t atom2 = Value(atom);
927 
928  TGeoElement* ele = new TGeoElement(formula, NameShort(name), z2 , atom2);
929 
930  felemap[name.Data()] = ele;
931 
932  return node;
933 
934 }
935 
936 ////////////////////////////////////////////////////////////////////////////////
937 ///In the materials section of the GDML file, materials can be declared.
938 ///when the material keyword is found, this function is called, and the
939 ///name and values of the material are converted into type TGeoMaterial
940 ///and stored in fmatmap map using the name as its key. Mixtures can also
941 /// be declared, and they are converted to TGeoMixture and stored in
942 ///fmixmap. These mixtures and materials are then all converted into one
943 ///common type - TGeoMedium. The map fmedmap is then built up of all the
944 ///mixtures and materials.
945 
947 {
948  //!Map to hold fractions while being processed
949  typedef FracMap::iterator fractions;
950 // typedef FracMap::iterator i;
951  FracMap fracmap;
952 
953  static int medid = 0;
954  XMLNodePointer_t child = gdml->GetChild(node);
955  TString tempattr = "";
956  Int_t ncompo = 0, mixflag = 2;
957  Double_t density = 0;
958  TString name = "";
959  TGeoMixture* mix = 0;
960  TGeoMaterial* mat = 0;
961  TString tempconst = "";
962  TString matname;
963  Bool_t composite = kFALSE;
964 
965  if (z == 1) {
966  Double_t a = 0;
967  Double_t d = 0;
968 
969  while (child != 0) {
970  attr = gdml->GetFirstAttr(child);
971 
972  if ((strcmp(gdml->GetNodeName(child), "atom")) == 0) {
973  while (attr != 0) {
974  tempattr = gdml->GetAttrName(attr);
975  tempattr.ToLower();
976 
977  if (tempattr == "value") {
978  a = Value(gdml->GetAttrValue(attr));
979  }
980  attr = gdml->GetNextAttr(attr);
981  }
982  }
983 
984  if ((strcmp(gdml->GetNodeName(child), "D")) == 0) {
985  while (attr != 0) {
986  tempattr = gdml->GetAttrName(attr);
987  tempattr.ToLower();
988 
989  if (tempattr == "value") {
990  d = Value(gdml->GetAttrValue(attr));
991  }
992  attr = gdml->GetNextAttr(attr);
993  }
994  }
995  child = gdml->GetNext(child);
996  }
997  //still in the is Z else...but not in the while..
998 
999  name = gdml->GetAttr(node, "name");
1000 
1001  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1002  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1003  }
1004 
1005  //CHECK FOR CONSTANTS
1006  tempconst = gdml->GetAttr(node, "Z");
1007 
1008  Double_t valZ = Value(tempconst);
1009 
1010  TString tmpname = name;
1011  //deal with special case - Z of vacuum is always 0
1012  tmpname.ToLower();
1013  if (tmpname == "vacuum") {
1014  valZ = 0;
1015  }
1016  mat = new TGeoMaterial(NameShort(name), a, valZ, d);
1017  mixflag = 0;
1018  //Note: Object(name, title) corresponds to Element(formula, name)
1019  TGeoElement* mat_ele = new TGeoElement(NameShort(name), NameShort(name), atoi(tempconst), a);
1020  felemap[name.Data()] = mat_ele;
1021 
1022  }
1023 
1024  else if (z == 0) {
1025  while (child != 0) {
1026  attr = gdml->GetFirstAttr(child);
1027 
1028  if ((strcmp(gdml->GetNodeName(child), "fraction")) == 0) {
1029  Double_t n = 0;
1030  TString ref = "";
1031  ncompo = ncompo + 1;
1032 
1033  while (attr != 0) {
1034  tempattr = gdml->GetAttrName(attr);
1035  tempattr.ToLower();
1036 
1037  if (tempattr == "n") {
1038  n = Value(gdml->GetAttrValue(attr));
1039  } else if (tempattr == "ref") {
1040  ref = gdml->GetAttrValue(attr);
1041  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1042  ref = TString::Format("%s_%s", ref.Data(), fCurrentFile);
1043  }
1044 
1045  }
1046 
1047  attr = gdml->GetNextAttr(attr);
1048  }
1049  fracmap[ref.Data()] = n;
1050 
1051  }
1052 
1053  else if ((strcmp(gdml->GetNodeName(child), "composite")) == 0) {
1054  composite = kTRUE;
1055  Double_t n = 0;
1056  TString ref = "";
1057  ncompo = ncompo + 1;
1058 
1059  while (attr != 0) {
1060  tempattr = gdml->GetAttrName(attr);
1061  tempattr.ToLower();
1062 
1063  if (tempattr == "n") {
1064  n = Value(gdml->GetAttrValue(attr));
1065  } else if (tempattr == "ref") {
1066  ref = gdml->GetAttrValue(attr);
1067  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1068  ref = TString::Format("%s_%s", ref.Data(), fCurrentFile);
1069  }
1070  }
1071 
1072  attr = gdml->GetNextAttr(attr);
1073  }
1074 
1075  fracmap[ref.Data()] = n;
1076 
1077  }
1078  else if ((strcmp(gdml->GetNodeName(child), "D")) == 0) {
1079  while (attr != 0) {
1080  tempattr = gdml->GetAttrName(attr);
1081  tempattr.ToLower();
1082 
1083  if (tempattr == "value") {
1084  density = Value(gdml->GetAttrValue(attr));
1085  }
1086 
1087  attr = gdml->GetNextAttr(attr);
1088  }
1089  }
1090 
1091  child = gdml->GetNext(child);
1092  }
1093  //still in the not Z else...but not in the while..
1094 
1095  name = gdml->GetAttr(node, "name");
1096  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1097  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1098  }
1099  //mix = new TGeoMixture(NameShort(name), 0 /*ncompo*/, density);
1100  mix = new TGeoMixture(NameShort(name), ncompo, density);
1101  mixflag = 1;
1102  Int_t natoms;
1103  Double_t weight;
1104 
1105  for (fractions f = fracmap.begin(); f != fracmap.end(); f++) {
1106  matname = f->first;
1107  matname = NameShort(matname);
1108 
1110 
1111  if (mattmp || (felemap.find(f->first) != felemap.end())) {
1112  if (composite) {
1113  natoms = (Int_t)f->second;
1114 
1115  mix->AddElement(felemap[f->first], natoms);
1116 
1117  }
1118 
1119  else {
1120  weight = f->second;
1121  if (mattmp){
1122 
1123  mix->AddElement(mattmp, weight);
1124  }
1125  else{
1126 
1127  mix->AddElement(felemap[f->first], weight);
1128  }
1129  }
1130  }
1131  }
1132 
1133  }//end of not Z else
1134 
1135  medid = medid + 1;
1136 
1137  TGeoMedium* med = 0;
1138 
1139  if (mixflag == 1) {
1140  fmixmap[name.Data()] = mix;
1141  med = new TGeoMedium(NameShort(name), medid, mix);
1142  } else if (mixflag == 0) {
1143  fmatmap[name.Data()] = mat;
1144  med = new TGeoMedium(NameShort(name), medid, mat);
1145  }
1146 
1147  fmedmap[name.Data()] = med;
1148 
1149  return child;
1150 
1151 }
1152 
1153 ////////////////////////////////////////////////////////////////////////////////
1154 ///In the structure section of the GDML file, volumes can be declared.
1155 ///when the volume keyword is found, this function is called, and the
1156 ///name and values of the volume are converted into type TGeoVolume and
1157 ///stored in fvolmap map using the name as its key. Volumes reference to
1158 ///a solid declared higher up in the solids section of the GDML file.
1159 ///Some volumes reference to other physical volumes to contain inside
1160 ///that volume, declaring positions and rotations within that volume.
1161 ///When each 'physvol' is declared, a matrix for its rotation and
1162 ///translation is built and the 'physvol node' is added to the original
1163 ///volume using TGeoVolume->AddNode.
1164 ///volume division is also declared within the volume node, and once the
1165 ///values for the division have been collected, using TGeoVolume->divide,
1166 ///the division can be applied.
1167 
1169 {
1171  XMLNodePointer_t subchild;
1172  XMLNodePointer_t subsubchild;
1173 
1174  XMLNodePointer_t child = gdml->GetChild(node);
1175  TString name;
1176  TString solidname = "";
1177  TString tempattr = "";
1178  TGeoShape* solid = 0;
1179  TGeoMedium* medium = 0;
1180  TGeoVolume* vol = 0;
1181  TGeoVolume* lv = 0;
1182  TGeoShape* reflex = 0;
1183  const Double_t* parentrot = 0;
1184  int yesrefl = 0;
1185  TString reftemp = "";
1186  TMap *auxmap = 0;
1187 
1188  while (child != 0) {
1189  if ((strcmp(gdml->GetNodeName(child), "solidref")) == 0) {
1190 
1191  reftemp = gdml->GetAttr(child, "ref");
1192  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1193  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1194  }
1195  if (fsolmap.find(reftemp.Data()) != fsolmap.end()) {
1196  solid = fsolmap[reftemp.Data()];
1197  } else if (freflectmap.find(reftemp.Data()) != freflectmap.end()) {
1198  solidname = reftemp;
1199  reflex = fsolmap[freflectmap[reftemp.Data()]];
1200  } else {
1201  printf("Solid: %s, Not Yet Defined!\n", reftemp.Data());
1202  }
1203  }
1204 
1205  if ((strcmp(gdml->GetNodeName(child), "materialref")) == 0) {
1206  reftemp = gdml->GetAttr(child, "ref");
1207  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1208  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1209  }
1210  if (fmedmap.find(reftemp.Data()) != fmedmap.end()) {
1211  medium = fmedmap[reftemp.Data()];
1212  } else {
1213  printf("Medium: %s, Not Yet Defined!\n", gdml->GetAttr(child, "ref"));
1214  }
1215  }
1216 
1217  child = gdml->GetNext(child);
1218  }
1219 
1220  name = gdml->GetAttr(node, "name");
1221 
1222  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1223  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1224  }
1225 
1226  if (reflex == 0) {
1227  vol = new TGeoVolume(NameShort(name), solid, medium);
1228  } else {
1229  vol = new TGeoVolume(NameShort(name), reflex, medium);
1230  freflvolmap[name.Data()] = solidname;
1231  TGDMLRefl* parentrefl = freflsolidmap[solidname.Data()];
1232  parentrot = parentrefl->GetMatrix()->GetRotationMatrix();
1233  yesrefl = 1;
1234  }
1235 
1236  fvolmap[name.Data()] = vol;
1237 
1238  //PHYSVOL - run through child nodes of VOLUME again..
1239 
1240  child = gdml->GetChild(node);
1241 
1242  while (child != 0) {
1243  if ((strcmp(gdml->GetNodeName(child), "physvol")) == 0) {
1244 
1245  TString volref = "";
1246 
1247  TGeoTranslation* pos = 0;
1248  TGeoRotation* rot = 0;
1249  TGeoScale* scl = 0;
1250  TString pnodename = gdml->GetAttr(child, "name");
1251  TString scopynum = gdml->GetAttr(child, "copynumber");
1252  Int_t copynum = (scopynum.IsNull()) ? 0 : (Int_t)Value(scopynum);
1253 
1254  subchild = gdml->GetChild(child);
1255 
1256  while (subchild != 0) {
1257  tempattr = gdml->GetNodeName(subchild);
1258  tempattr.ToLower();
1259 
1260  if (tempattr == "volumeref") {
1261  reftemp = gdml->GetAttr(subchild, "ref");
1262  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1263  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1264  }
1265  lv = fvolmap[reftemp.Data()];
1266  volref = reftemp;
1267  }
1268  else if (tempattr == "file") {
1269  const char* filevol;
1270  const char* prevfile = fCurrentFile;
1271 
1272  fCurrentFile = gdml->GetAttr(subchild, "name");
1273  filevol = gdml->GetAttr(subchild, "volname");
1274 
1275  TXMLEngine* gdml2 = new TXMLEngine;
1276  gdml2->SetSkipComments(kTRUE);
1277  XMLDocPointer_t filedoc1 = gdml2->ParseFile(fCurrentFile);
1278  if (filedoc1 == 0) {
1279  Fatal("VolProcess", "Bad filename given %s", fCurrentFile);
1280  }
1281  // take access to main node
1282  XMLNodePointer_t mainnode2 = gdml2->DocGetRootElement(filedoc1);
1283  //increase depth counter + add DOM pointer
1284  fFILENO = fFILENO + 1;
1285  fFileEngine[fFILENO] = gdml2;
1286 
1287  if (ffilemap.find(fCurrentFile) != ffilemap.end()) {
1288  volref = ffilemap[fCurrentFile];
1289  } else {
1290  volref = ParseGDML(gdml2, mainnode2);
1291  ffilemap[fCurrentFile] = volref;
1292  }
1293 
1294  if (filevol) {
1295  volref = filevol;
1296  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1297  volref = TString::Format("%s_%s", volref.Data(), fCurrentFile);
1298  }
1299  }
1300 
1301  fFILENO = fFILENO - 1;
1302  gdml = fFileEngine[fFILENO];
1303  fCurrentFile = prevfile;
1304 
1305  lv = fvolmap[volref.Data()];
1306  //File tree complete - Release memory before exit
1307 
1308  gdml->FreeDoc(filedoc1);
1309  delete gdml2;
1310  }
1311  else if (tempattr == "position") {
1312  attr = gdml->GetFirstAttr(subchild);
1313  PosProcess(gdml, subchild, attr);
1314  reftemp = gdml->GetAttr(subchild, "name");
1315  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1316  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1317  }
1318  pos = fposmap[reftemp.Data()];
1319  } else if (tempattr == "positionref") {
1320  reftemp = gdml->GetAttr(subchild, "ref");
1321  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1322  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1323  }
1324  if (fposmap.find(reftemp.Data()) != fposmap.end()) pos = fposmap[reftemp.Data()];
1325  else std::cout << "ERROR! Physvol's position " << reftemp << " not found!" << std::endl;
1326  } else if (tempattr == "rotation") {
1327  attr = gdml->GetFirstAttr(subchild);
1328  RotProcess(gdml, subchild, attr);
1329  reftemp = gdml->GetAttr(subchild, "name");
1330  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1331  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1332  }
1333  rot = frotmap[reftemp.Data()];
1334  } else if (tempattr == "rotationref") {
1335  reftemp = gdml->GetAttr(subchild, "ref");
1336  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1337  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1338  }
1339  if (frotmap.find(reftemp.Data()) != frotmap.end()) rot = frotmap[reftemp.Data()];
1340  else std::cout << "ERROR! Physvol's rotation " << reftemp << " not found!" << std::endl;
1341  } else if (tempattr == "scale") {
1342  attr = gdml->GetFirstAttr(subchild);
1343  SclProcess(gdml, subchild, attr);
1344  reftemp = gdml->GetAttr(subchild, "name");
1345  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1346  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1347  }
1348  scl = fsclmap[reftemp.Data()];
1349  } else if (tempattr == "scaleref") {
1350  reftemp = gdml->GetAttr(subchild, "ref");
1351  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1352  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1353  }
1354  if (fsclmap.find(reftemp.Data()) != fsclmap.end()) scl = fsclmap[reftemp.Data()];
1355  else std::cout << "ERROR! Physvol's scale " << reftemp << " not found!" << std::endl;
1356  }
1357 
1358  subchild = gdml->GetNext(subchild);
1359  }
1360 
1361  //ADD PHYSVOL TO GEOMETRY
1362  fVolID = fVolID + 1;
1363 
1364  TGeoHMatrix *transform = new TGeoHMatrix();
1365 
1366  if (pos != 0) transform->SetTranslation(pos->GetTranslation());
1367  if (rot != 0) transform->SetRotation(rot->GetRotationMatrix());
1368 
1369  if (scl != 0) { // Scaling must be added to the rotation matrix!
1370 
1371  Double_t scale3x3[9];
1372  memset(scale3x3, 0, 9 * sizeof(Double_t));
1373  const Double_t *diagonal = scl->GetScale();
1374 
1375  scale3x3[0] = diagonal[0];
1376  scale3x3[4] = diagonal[1];
1377  scale3x3[8] = diagonal[2];
1378 
1379  TGeoRotation scaleMatrix;
1380  scaleMatrix.SetMatrix(scale3x3);
1381  transform->Multiply(&scaleMatrix);
1382  }
1383 
1384 // BEGIN: reflectedSolid. Remove lines between if reflectedSolid will be removed from GDML!!!
1385 
1386  if (freflvolmap.find(volref.Data()) != freflvolmap.end()) {
1387  // if the volume is a reflected volume the matrix needs to be CHANGED
1388  TGDMLRefl* temprefl = freflsolidmap[freflvolmap[volref.Data()]];
1389  transform->Multiply(temprefl->GetMatrix());
1390  }
1391 
1392  if (yesrefl == 1) {
1393  // reflection is done per solid so that we cancel it if exists in mother volume!!!
1394  TGeoRotation prot;
1395  prot.SetMatrix(parentrot);
1396  transform->MultiplyLeft(&prot);
1397  }
1398 
1399 // END: reflectedSolid
1400 
1401  vol->AddNode(lv, copynum, transform);
1402  if (!pnodename.IsNull())
1403  ((TNamed*)vol->GetNodes()->Last())->SetName(pnodename);
1404  } else if ((strcmp(gdml->GetNodeName(child), "divisionvol")) == 0) {
1405 
1406  TString divVolref = "";
1407  Int_t axis = 0;
1408  TString number = "";
1409  TString width = "";
1410  TString offset = "";
1411  TString lunit = "mm";
1412 
1413  attr = gdml->GetFirstAttr(child);
1414 
1415  while (attr != 0) {
1416 
1417  tempattr = gdml->GetAttrName(attr);
1418  tempattr.ToLower();
1419 
1420  if (tempattr == "axis") {
1421  axis = SetAxis(gdml->GetAttrValue(attr));
1422  } else if (tempattr == "number") {
1423  number = gdml->GetAttrValue(attr);
1424  } else if (tempattr == "width") {
1425  width = gdml->GetAttrValue(attr);
1426  } else if (tempattr == "offset") {
1427  offset = gdml->GetAttrValue(attr);
1428  } else if (tempattr == "unit") {
1429  lunit = gdml->GetAttrValue(attr);
1430  }
1431 
1432  attr = gdml->GetNextAttr(attr);
1433 
1434  }
1435 
1436  subchild = gdml->GetChild(child);
1437 
1438  while (subchild != 0) {
1439  tempattr = gdml->GetNodeName(subchild);
1440  tempattr.ToLower();
1441 
1442  if (tempattr == "volumeref") {
1443  reftemp = gdml->GetAttr(subchild, "ref");
1444  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1445  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1446  }
1447  divVolref = reftemp;
1448  }
1449 
1450  subchild = gdml->GetNext(subchild);
1451  }
1452 
1453 
1454  Double_t numberline = Value(number);
1455  Double_t retunit = GetScaleVal(lunit);
1456  Double_t step = Value(width) * retunit;
1457  Double_t offsetline = Value(offset) * retunit;
1458 
1459  fVolID = fVolID + 1;
1460  Double_t xlo, xhi;
1461  vol->GetShape()->GetAxisRange(axis, xlo, xhi);
1462 
1463  Int_t ndiv = (Int_t)numberline;
1464  Double_t start = xlo + offsetline;
1465 
1466  Int_t numed = 0;
1467  TGeoVolume *old = fvolmap[NameShort(reftemp)];
1468  if (old) {
1469  // We need to recreate the content of the divided volume
1470  old = fvolmap[NameShort(reftemp)];
1471  // medium id
1472  numed = old->GetMedium()->GetId();
1473  }
1474  TGeoVolume *divvol = vol->Divide(NameShort(reftemp), axis, ndiv, start, step, numed);
1475  if (!divvol) {
1476  Fatal("VolProcess", "Cannot divide volume %s", vol->GetName());
1477  return child;
1478  }
1479  if (old && old->GetNdaughters()) {
1480  divvol->ReplayCreation(old);
1481  }
1482  fvolmap[NameShort(reftemp)] = divvol;
1483 
1484  }//end of Division else if
1485 
1486 
1487  else if ((strcmp(gdml->GetNodeName(child), "replicavol")) == 0) {
1488 
1489  TString divVolref = "";
1490  Int_t axis = 0;
1491  TString number = "";
1492  TString width = "";
1493  TString offset = "";
1494  TString wunit = "mm";
1495  TString ounit = "mm";
1496  Double_t wvalue = 0;
1497  Double_t ovalue = 0;
1498 
1499 
1500  attr = gdml->GetFirstAttr(child);
1501 
1502  while (attr != 0) {
1503 
1504  tempattr = gdml->GetAttrName(attr);
1505  tempattr.ToLower();
1506 
1507  if (tempattr == "number") {
1508  number = gdml->GetAttrValue(attr);
1509  }
1510  attr = gdml->GetNextAttr(attr);
1511  }
1512 
1513  subchild = gdml->GetChild(child);
1514 
1515  while (subchild != 0) {
1516  tempattr = gdml->GetNodeName(subchild);
1517  tempattr.ToLower();
1518 
1519  if (tempattr == "volumeref") {
1520  reftemp = gdml->GetAttr(subchild, "ref");
1521  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1522  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1523  }
1524  divVolref = reftemp;
1525  }
1526 
1527  if (tempattr == "replicate_along_axis") {
1528  subsubchild = gdml->GetChild(subchild);
1529 
1530  while (subsubchild != 0) {
1531  if ((strcmp(gdml->GetNodeName(subsubchild), "width")) == 0) {
1532  attr = gdml->GetFirstAttr(subsubchild);
1533  while (attr != 0) {
1534  tempattr = gdml->GetAttrName(attr);
1535  tempattr.ToLower();
1536  if (tempattr == "value") {
1537  wvalue = Value(gdml->GetAttrValue(attr));
1538  }
1539  else if (tempattr == "unit"){
1540  wunit = gdml->GetAttrValue(attr);
1541  }
1542 
1543  attr = gdml->GetNextAttr(attr);
1544  }
1545  }
1546  else if ((strcmp(gdml->GetNodeName(subsubchild), "offset")) == 0) {
1547  attr = gdml->GetFirstAttr(subsubchild);
1548  while (attr != 0) {
1549  tempattr = gdml->GetAttrName(attr);
1550  tempattr.ToLower();
1551  if (tempattr == "value") {
1552  ovalue = Value(gdml->GetAttrValue(attr));
1553  }
1554  else if (tempattr == "unit"){
1555  ounit = gdml->GetAttrValue(attr);
1556  }
1557  attr = gdml->GetNextAttr(attr);
1558  }
1559  }
1560  else if ((strcmp(gdml->GetNodeName(subsubchild), "direction")) == 0) {
1561  attr = gdml->GetFirstAttr(subsubchild);
1562  while (attr != 0) {
1563  tempattr = gdml->GetAttrName(attr);
1564  tempattr.ToLower();
1565  if (tempattr == "x") {
1566  axis = 1;
1567  }
1568  else if (tempattr == "y"){
1569  axis = 2;
1570  }
1571  else if (tempattr == "z"){
1572  axis = 3;
1573  }
1574  else if (tempattr == "rho"){
1575  axis = 1;
1576  }
1577  else if (tempattr == "phi"){
1578  axis = 2;
1579  }
1580 
1581  attr = gdml->GetNextAttr(attr);
1582  }
1583  }
1584 
1585  subsubchild = gdml->GetNext(subsubchild);
1586  }
1587 
1588  }
1589 
1590  subchild = gdml->GetNext(subchild);
1591  }
1592 
1593 
1594  Double_t retwunit = GetScaleVal(wunit);
1595  Double_t retounit = GetScaleVal(ounit);
1596 
1597  Double_t numberline = Value(number);
1598  Double_t widthline = wvalue*retwunit;
1599  Double_t offsetline = ovalue*retounit;
1600 
1601  fVolID = fVolID + 1;
1602  Double_t xlo, xhi;
1603  vol->GetShape()->GetAxisRange(axis, xlo, xhi);
1604 
1605  Int_t ndiv = (Int_t)numberline;
1606  Double_t start = xlo + offsetline;
1607 
1608  Double_t step = widthline;
1609  Int_t numed = 0;
1610  TGeoVolume *old = fvolmap[NameShort(reftemp)];
1611  if (old) {
1612  // We need to recreate the content of the divided volume
1613  old = fvolmap[NameShort(reftemp)];
1614  // medium id
1615  numed = old->GetMedium()->GetId();
1616  }
1617  TGeoVolume *divvol = vol->Divide(NameShort(reftemp), axis, ndiv, start, step, numed);
1618  if (!divvol) {
1619  Fatal("VolProcess", "Cannot divide volume %s", vol->GetName());
1620  return child;
1621  }
1622  if (old && old->GetNdaughters()) {
1623  divvol->ReplayCreation(old);
1624  }
1625  fvolmap[NameShort(reftemp)] = divvol;
1626 
1627  } //End of replicavol
1628  else if (strcmp(gdml->GetNodeName(child), "auxiliary") == 0) {
1629  TString auxType, auxUnit, auxValue;
1630  if(!auxmap) {
1631  printf("Auxiliary values for volume %s\n",vol->GetName());
1632  auxmap = new TMap();
1633  vol->SetUserExtension(new TGeoRCExtension(auxmap));
1634  }
1635  attr = gdml->GetFirstAttr(child);
1636  while(attr) {
1637  if (!strcmp(gdml->GetAttrName(attr),"auxtype")) auxType = gdml->GetAttrValue(attr);
1638  else if (!strcmp(gdml->GetAttrName(attr),"auxvalue")) auxValue = gdml->GetAttrValue(attr);
1639  else if (!strcmp(gdml->GetAttrName(attr),"auxunit")) auxUnit = gdml->GetAttrValue(attr);
1640  attr = gdml->GetNextAttr(attr);
1641  }
1642  if (!auxUnit.IsNull()) auxValue = TString::Format("%s*%s", auxValue.Data(), auxUnit.Data());
1643  auxmap->Add(new TObjString(auxType),new TObjString(auxValue));
1644  printf(" %s: %s\n", auxType.Data(), auxValue.Data());
1645  }
1646 
1647  child = gdml->GetNext(child);
1648  }
1649 
1650  return child;
1651 
1652 }
1653 
1654 ////////////////////////////////////////////////////////////////////////////////
1655 ///In the solid section of the GDML file, boolean solids can be
1656 ///declared. when the subtraction, intersection or union keyword
1657 ///is found, this function is called, and the values (rotation and
1658 ///translation) of the solid are converted into type TGeoCompositeShape
1659 ///and stored in fsolmap map using the name as its key.
1660 ///
1661 ///1 = SUBTRACTION
1662 ///2 = INTERSECTION
1663 ///3 = UNION
1664 
1666 {
1667  TString reftemp = "";
1668  TString tempattr = "";
1669  XMLNodePointer_t child = gdml->GetChild(node);
1670 
1671  TGeoShape* first = 0;
1672  TGeoShape* second = 0;
1673 
1674  TGeoTranslation* firstPos = new TGeoTranslation(0, 0, 0);
1675  TGeoTranslation* secondPos = new TGeoTranslation(0, 0, 0);
1676 
1677  TGeoRotation* firstRot = new TGeoRotation();
1678  TGeoRotation* secondRot = new TGeoRotation();
1679 
1680  firstRot->RotateZ(0);
1681  firstRot->RotateY(0);
1682  firstRot->RotateX(0);
1683 
1684  secondRot->RotateZ(0);
1685  secondRot->RotateY(0);
1686  secondRot->RotateX(0);
1687 
1688  TString name = gdml->GetAttr(node, "name");
1689 
1690  if ((strcmp(fCurrentFile, fStartFile)) != 0)
1691  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1692 
1693  while (child != 0) {
1694  tempattr = gdml->GetNodeName(child);
1695  tempattr.ToLower();
1696 
1697  if (tempattr == "first") {
1698  reftemp = gdml->GetAttr(child, "ref");
1699  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1700  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1701  }
1702  if (fsolmap.find(reftemp.Data()) != fsolmap.end()) {
1703  first = fsolmap[reftemp.Data()];
1704  }
1705  } else if (tempattr == "second") {
1706  reftemp = gdml->GetAttr(child, "ref");
1707  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1708  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1709  }
1710  if (fsolmap.find(reftemp.Data()) != fsolmap.end()) {
1711  second = fsolmap[reftemp.Data()];
1712  }
1713  } else if (tempattr == "position") {
1714  attr = gdml->GetFirstAttr(child);
1715  PosProcess(gdml, child, attr);
1716  reftemp = gdml->GetAttr(child, "name");
1717  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1718  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1719  }
1720  secondPos = fposmap[reftemp.Data()];
1721  } else if (tempattr == "positionref") {
1722  reftemp = gdml->GetAttr(child, "ref");
1723  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1724  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1725  }
1726  if (fposmap.find(reftemp.Data()) != fposmap.end()) {
1727  secondPos = fposmap[reftemp.Data()];
1728  }
1729  } else if (tempattr == "rotation") {
1730  attr = gdml->GetFirstAttr(child);
1731  RotProcess(gdml, child, attr);
1732  reftemp = gdml->GetAttr(child, "name");
1733  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1734  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1735  }
1736  secondRot = frotmap[reftemp.Data()];
1737  } else if (tempattr == "rotationref") {
1738  reftemp = gdml->GetAttr(child, "ref");
1739  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1740  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1741  }
1742  if (frotmap.find(reftemp.Data()) != frotmap.end()) {
1743  secondRot = frotmap[reftemp.Data()];
1744  }
1745  } else if (tempattr == "firstposition") {
1746  attr = gdml->GetFirstAttr(child);
1747  PosProcess(gdml, child, attr);
1748  reftemp = gdml->GetAttr(child, "name");
1749  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1750  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1751  }
1752  firstPos = fposmap[reftemp.Data()];
1753  } else if (tempattr == "firstpositionref") {
1754  reftemp = gdml->GetAttr(child, "ref");
1755  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1756  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1757  }
1758  if (fposmap.find(reftemp.Data()) != fposmap.end()) {
1759  firstPos = fposmap[reftemp.Data()];
1760  }
1761  } else if (tempattr == "firstrotation") {
1762  attr = gdml->GetFirstAttr(child);
1763  RotProcess(gdml, child, attr);
1764  reftemp = gdml->GetAttr(child, "name");
1765  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1766  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1767  }
1768  firstRot = frotmap[reftemp.Data()];
1769  } else if (tempattr == "firstrotationref") {
1770  reftemp = gdml->GetAttr(child, "ref");
1771  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1772  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1773  }
1774  if (frotmap.find(reftemp.Data()) != frotmap.end()) {
1775  firstRot = frotmap[reftemp.Data()];
1776  }
1777  }
1778  child = gdml->GetNext(child);
1779  }
1780 
1781  TGeoMatrix* firstMatrix = new TGeoCombiTrans(*firstPos, firstRot->Inverse());
1782  TGeoMatrix* secondMatrix = new TGeoCombiTrans(*secondPos, secondRot->Inverse());
1783 
1784  TGeoCompositeShape* boolean = 0;
1785  if (!first || !second) {
1786  Fatal("BooSolid", "Incomplete solid %s, missing shape components", name.Data());
1787  return child;
1788  }
1789  switch (num) {
1790  case 1:
1791  boolean = new TGeoCompositeShape(NameShort(name), new TGeoSubtraction(first, second, firstMatrix, secondMatrix));
1792  break; // SUBTRACTION
1793  case 2:
1794  boolean = new TGeoCompositeShape(NameShort(name), new TGeoIntersection(first, second, firstMatrix, secondMatrix));
1795  break; // INTERSECTION
1796  case 3:
1797  boolean = new TGeoCompositeShape(NameShort(name), new TGeoUnion(first, second, firstMatrix, secondMatrix));
1798  break; // UNION
1799  default:
1800  break;
1801  }
1802 
1803  fsolmap[name.Data()] = boolean;
1804 
1805  return child;
1806 }
1807 
1808 ////////////////////////////////////////////////////////////////////////////////
1809 ///User data to be processed
1811 {
1812  Warning("ParseGDML", "<userinfo> not supported yet. Skipping.");
1813  XMLNodePointer_t child = gdml->GetChild(node);
1814  while (child != 0) {
1815  child = gdml->GetNext(child);
1816  }
1817  return child;
1818 }
1819 
1820 ////////////////////////////////////////////////////////////////////////////////
1821 ///In the structure section of the GDML file, assembly volumes can be
1822 ///declared. when the assembly keyword is found, this function is called,
1823 ///and the name is converted into type TGeoVolumeAssembly and
1824 ///stored in fvolmap map using the name as its key. Some assembly volumes
1825 ///reference to other physical volumes to contain inside that assembly,
1826 ///declaring positions and rotations within that volume. When each 'physvol'
1827 ///is declared, a matrix for its rotation and translation is built and the
1828 ///'physvol node' is added to the original assembly using TGeoVolume->AddNode.
1829 
1831 {
1832  TString name = gdml->GetAttr(node, "name");
1833  TString reftemp = "";
1834 
1835  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1836  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1837  }
1838 
1840  XMLNodePointer_t subchild;
1841  XMLNodePointer_t child = gdml->GetChild(node);
1842  TString tempattr = "";
1843  TGeoVolume* lv = 0;
1844  TGeoTranslation* pos = 0;
1845  TGeoRotation* rot = 0;
1846  TGeoCombiTrans* matr;
1847 
1848  TGeoVolumeAssembly* assem = new TGeoVolumeAssembly(NameShort(name));
1849 
1850 
1851  //PHYSVOL - run through child nodes of VOLUME again..
1852 
1853 // child = gdml->GetChild(node);
1854 
1855  while (child != 0) {
1856  if ((strcmp(gdml->GetNodeName(child), "physvol")) == 0) {
1857  TString pnodename = gdml->GetAttr(child, "name");
1858  TString scopynum = gdml->GetAttr(child, "copynumber");
1859  Int_t copynum = (scopynum.IsNull()) ? 0 : (Int_t)Value(scopynum);
1860 
1861  subchild = gdml->GetChild(child);
1862  pos = new TGeoTranslation(0, 0, 0);
1863  rot = new TGeoRotation();
1864 
1865  while (subchild != 0) {
1866  tempattr = gdml->GetNodeName(subchild);
1867  tempattr.ToLower();
1868 
1869  if (tempattr == "volumeref") {
1870  reftemp = gdml->GetAttr(subchild, "ref");
1871  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1872  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1873  }
1874  lv = fvolmap[reftemp.Data()];
1875  } else if (tempattr == "positionref") {
1876  reftemp = gdml->GetAttr(subchild, "ref");
1877  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1878  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1879  }
1880  if (fposmap.find(reftemp.Data()) != fposmap.end()) {
1881  pos = fposmap[reftemp.Data()];
1882  }
1883  } else if (tempattr == "position") {
1884  attr = gdml->GetFirstAttr(subchild);
1885  PosProcess(gdml, subchild, attr);
1886  reftemp = gdml->GetAttr(subchild, "name");
1887  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1888  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1889  }
1890  pos = fposmap[reftemp.Data()];
1891  } else if (tempattr == "rotationref") {
1892  reftemp = gdml->GetAttr(subchild, "ref");
1893  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1894  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1895  }
1896  if (frotmap.find(reftemp.Data()) != frotmap.end()) {
1897  rot = frotmap[reftemp.Data()];
1898  }
1899  } else if (tempattr == "rotation") {
1900  attr = gdml->GetFirstAttr(subchild);
1901  RotProcess(gdml, subchild, attr);
1902  reftemp = gdml->GetAttr(subchild, "name");
1903  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1904  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1905  }
1906  rot = frotmap[reftemp.Data()];
1907  }
1908 
1909  subchild = gdml->GetNext(subchild);
1910  }
1911 
1912  //ADD PHYSVOL TO GEOMETRY
1913  fVolID = fVolID + 1;
1914  matr = new TGeoCombiTrans(*pos, *rot);
1915  assem->AddNode(lv, copynum, matr);
1916  if (!pnodename.IsNull())
1917  ((TNamed*)assem->GetNodes()->Last())->SetName(pnodename);
1918 
1919  }
1920  child = gdml->GetNext(child);
1921  }
1922 
1923  fvolmap[name.Data()] = assem;
1924  return child;
1925 }
1926 
1927 ////////////////////////////////////////////////////////////////////////////////
1928 ///In the setup section of the GDML file, the top volume need to be
1929 ///declared. when the setup keyword is found, this function is called,
1930 ///and the top volume ref is taken and 'world' is set
1931 
1933 {
1934  const char* name = gdml->GetAttr(node, "name");
1935  gGeoManager->SetName(name);
1936  XMLNodePointer_t child = gdml->GetChild(node);
1937  TString reftemp = "";
1938 
1939  while (child != 0) {
1940 
1941  if ((strcmp(gdml->GetNodeName(child), "world") == 0)) {
1942  //const char* reftemp;
1943  //TString reftemp = "";
1944  reftemp = gdml->GetAttr(child, "ref");
1945 
1946 
1947  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1948  reftemp = TString::Format("%s_%s", reftemp.Data(), fCurrentFile);
1949 
1950  }
1951  fWorld = fvolmap[reftemp.Data()];
1952  fWorldName = reftemp.Data();
1953 
1954  }
1955  child = gdml->GetNext(child);
1956  }
1957  return node;
1958 }
1959 
1960 ////////////////////////////////////////////////////////////////////////////////
1961 ///In the solids section of the GDML file, a box may be declared.
1962 ///when the box keyword is found, this function is called, and the
1963 ///dimensions required are taken and stored, these are then bound and
1964 ///converted to type TGeoBBox and stored in fsolmap map using the name
1965 ///as its key.
1966 
1968 {
1969  TString lunit = "mm";
1970  TString xpos = "0";
1971  TString ypos = "0";
1972  TString zpos = "0";
1973  TString name = "";
1974  TString tempattr;
1975 
1976  while (attr != 0) {
1977 
1978  tempattr = gdml->GetAttrName(attr);
1979  tempattr.ToLower();
1980 
1981  if (tempattr == "name") {
1982  name = gdml->GetAttrValue(attr);
1983  } else if (tempattr == "x") {
1984  xpos = gdml->GetAttrValue(attr);
1985  } else if (tempattr == "y") {
1986  ypos = gdml->GetAttrValue(attr);
1987  } else if (tempattr == "z") {
1988  zpos = gdml->GetAttrValue(attr);
1989  } else if (tempattr == "lunit") {
1990  lunit = gdml->GetAttrValue(attr);
1991  }
1992 
1993  attr = gdml->GetNextAttr(attr);
1994  }
1995 
1996  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
1997  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
1998  }
1999 
2000  Double_t retunit = GetScaleVal(lunit);
2001 
2002  Double_t xline = 0.5*Value(xpos)*retunit;
2003  Double_t yline = 0.5*Value(ypos)*retunit;
2004  Double_t zline = 0.5*Value(zpos)*retunit;
2005 
2006 
2007  TGeoBBox* box = new TGeoBBox(NameShort(name), xline, yline, zline);
2008 
2009  fsolmap[name.Data()] = box;
2010 
2011  return node;
2012 
2013 }
2014 
2015 ////////////////////////////////////////////////////////////////////////////////
2016 ///In the solids section of the GDML file, an ellipsoid may be declared.
2017 ///Unfortunately, the ellipsoid is not supported under ROOT so,
2018 ///when the ellipsoid keyword is found, this function is called
2019 ///to convert it to a simple box with similar dimensions, and the
2020 ///dimensions required are taken and stored, these are then bound and
2021 ///converted to type TGeoBBox and stored in fsolmap map using the name
2022 ///as its key.
2023 
2025 {
2026  TString lunit = "mm";
2027  TString ax = "0";
2028  TString by = "0";
2029  TString cz = "0";
2030  //initialization to empty string
2031  TString zcut1 = "";
2032  TString zcut2 = "";
2033  TString name = "";
2034  TString tempattr;
2035 
2036  while (attr != 0) {
2037 
2038  tempattr = gdml->GetAttrName(attr);
2039  tempattr.ToLower();
2040 
2041  if (tempattr == "name") {
2042  name = gdml->GetAttrValue(attr);
2043  } else if (tempattr == "ax") {
2044  ax = gdml->GetAttrValue(attr);
2045  } else if (tempattr == "by") {
2046  by = gdml->GetAttrValue(attr);
2047  } else if (tempattr == "cz") {
2048  cz = gdml->GetAttrValue(attr);
2049  } else if (tempattr == "zcut1") {
2050  zcut1 = gdml->GetAttrValue(attr);
2051  } else if (tempattr == "zcut2") {
2052  zcut2 = gdml->GetAttrValue(attr);
2053  } else if (tempattr == "lunit") {
2054  lunit = gdml->GetAttrValue(attr);
2055  }
2056 
2057  attr = gdml->GetNextAttr(attr);
2058  }
2059 
2060  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2061  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2062  }
2063 
2064  Double_t retunit = GetScaleVal(lunit);
2065 
2066  Double_t dx = Value(ax)*retunit;
2067  Double_t dy = Value(by)*retunit;
2068  Double_t radius = Value(cz)*retunit;
2069  Double_t sx = dx / radius;
2070  Double_t sy = dy / radius;
2071  Double_t sz = 1.;
2072  Double_t z1, z2;
2073  //Initialization of cutting
2074  if (zcut1 == "") {
2075  z1 = -radius;
2076  } else {
2077  z1 = Value(zcut1)*retunit;
2078  }
2079  if (zcut2 == "") {
2080  z2 = radius;
2081  } else {
2082  z2 = Value(zcut2)*retunit;
2083  }
2084 
2085  TGeoSphere *sph = new TGeoSphere(0, radius);
2086  TGeoScale *scl = new TGeoScale("", sx, sy, sz);
2087  TGeoScaledShape *shape = new TGeoScaledShape(NameShort(name), sph, scl);
2088 
2089  Double_t origin[3] = {0., 0., 0.};
2090  origin[2] = 0.5 * (z1 + z2);
2091  Double_t dz = 0.5 * (z2 - z1);
2092  TGeoBBox *pCutBox = new TGeoBBox("cutBox", dx, dy, dz, origin);
2093  TGeoBoolNode *pBoolNode = new TGeoIntersection(shape, pCutBox, 0, 0);
2094  TGeoCompositeShape *cs = new TGeoCompositeShape(NameShort(name), pBoolNode);
2095  fsolmap[name.Data()] = cs;
2096 
2097  return node;
2098 
2099 }
2100 
2101 ////////////////////////////////////////////////////////////////////////////////
2102 ///In the solids section of the GDML file, an elliptical cone may be declared.
2103 ///Unfortunately, the elliptical cone is not supported under ROOT so,
2104 ///when the elcone keyword is found, this function is called
2105 ///to convert it to a simple box with similar dimensions, and the
2106 ///dimensions required are taken and stored, these are then bound and
2107 ///converted to type TGeoBBox and stored in fsolmap map using the name
2108 ///as its key.
2109 
2111 {
2112  TString lunit = "mm";
2113  TString dx = "0";
2114  TString dy = "0";
2115  TString zmax = "0";
2116  TString zcut = "0";
2117  TString name = "";
2118  TString tempattr;
2119 
2120  while (attr != 0) {
2121 
2122  tempattr = gdml->GetAttrName(attr);
2123  tempattr.ToLower();
2124 
2125  if (tempattr == "name") {
2126  name = gdml->GetAttrValue(attr);
2127  } else if (tempattr == "dx") {
2128  dx = gdml->GetAttrValue(attr);
2129  } else if (tempattr == "dy") {
2130  dy = gdml->GetAttrValue(attr);
2131  } else if (tempattr == "zmax") {
2132  zmax = gdml->GetAttrValue(attr);
2133  } else if (tempattr == "zcut") {
2134  zcut = gdml->GetAttrValue(attr);
2135  } else if (tempattr == "lunit") {
2136  lunit = gdml->GetAttrValue(attr);
2137  }
2138 
2139  attr = gdml->GetNextAttr(attr);
2140  }
2141 
2142  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2143  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2144  }
2145 
2146  //semiaxises of elliptical cone (elcone) are different then ellipsoid
2147 
2148  Double_t retunit = GetScaleVal(lunit);
2149 
2150  //dxline and dyline are without units because they are as a ration
2151  Double_t dxratio = Value(dx);
2152  Double_t dyratio = Value(dy);
2153  Double_t z = Value(zmax)*retunit;
2154  Double_t z1 = Value(zcut)*retunit;
2155 
2156  if (z1 <= 0) {
2157  Info("ElCone", "ERROR! Parameter zcut = %.12g is not set properly, elcone will not be imported.", z1);
2158  return node;
2159  }
2160  if (z1 > z){
2161  z1 = z;
2162  }
2163  Double_t rx1 = (z + z1) * dxratio;
2164  Double_t ry1 = (z + z1) * dyratio;
2165  Double_t rx2 = (z - z1) * dxratio;
2166  Double_t sx = 1.;
2167  Double_t sy = ry1 / rx1;
2168  Double_t sz = 1.;
2169 
2170  TGeoCone *con = new TGeoCone(z1, 0, rx1, 0, rx2);
2171  TGeoScale *scl = new TGeoScale("", sx, sy, sz);
2172  TGeoScaledShape *shape = new TGeoScaledShape(NameShort(name), con, scl);
2173 
2174  fsolmap[name.Data()] = shape;
2175 
2176  return node;
2177 
2178 }
2179 
2180 ////////////////////////////////////////////////////////////////////////////////
2181 ///In the solids section of the GDML file, a Paraboloid may be declared.
2182 ///when the paraboloid keyword is found, this function is called, and the
2183 ///dimensions required are taken and stored, these are then bound and
2184 ///converted to type TGeoParaboloid and stored in fsolmap map using the name
2185 ///as its key.
2186 
2188 {
2189  TString lunit = "mm";
2190  TString rlopos = "0";
2191  TString rhipos = "0";
2192  TString dzpos = "0";
2193  TString name = "";
2194  TString tempattr;
2195 
2196  while (attr != 0) {
2197 
2198  tempattr = gdml->GetAttrName(attr);
2199  tempattr.ToLower();
2200 
2201  if (tempattr == "name") {
2202  name = gdml->GetAttrValue(attr);
2203  } else if (tempattr == "rlo") {
2204  rlopos = gdml->GetAttrValue(attr);
2205  } else if (tempattr == "rhi") {
2206  rhipos = gdml->GetAttrValue(attr);
2207  } else if (tempattr == "dz") {
2208  dzpos = gdml->GetAttrValue(attr);
2209  } else if (tempattr == "lunit") {
2210  lunit = gdml->GetAttrValue(attr);
2211  }
2212 
2213  attr = gdml->GetNextAttr(attr);
2214  }
2215 
2216  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2217  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2218  }
2219 
2220  Double_t retunit = GetScaleVal(lunit);
2221 
2222  Double_t rlo = Value(rlopos)*retunit;
2223  Double_t rhi = Value(rhipos)*retunit;
2224  Double_t dz = Value(dzpos)*retunit;
2225 
2226  TGeoParaboloid* paraboloid = new TGeoParaboloid(NameShort(name), rlo, rhi, dz);
2227 
2228  fsolmap[name.Data()] = paraboloid;
2229 
2230  return node;
2231 
2232 }
2233 
2234 ////////////////////////////////////////////////////////////////////////////////
2235 ///In the solids section of the GDML file, an Arb8 may be declared.
2236 ///when the arb8 keyword is found, this function is called, and the
2237 ///dimensions required are taken and stored, these are then bound and
2238 ///converted to type TGeoArb8 and stored in fsolmap map using the name
2239 ///as its key.
2240 
2242 {
2243  TString lunit = "mm";
2244  TString v1xpos = "0";
2245  TString v1ypos = "0";
2246  TString v2xpos = "0";
2247  TString v2ypos = "0";
2248  TString v3xpos = "0";
2249  TString v3ypos = "0";
2250  TString v4xpos = "0";
2251  TString v4ypos = "0";
2252  TString v5xpos = "0";
2253  TString v5ypos = "0";
2254  TString v6xpos = "0";
2255  TString v6ypos = "0";
2256  TString v7xpos = "0";
2257  TString v7ypos = "0";
2258  TString v8xpos = "0";
2259  TString v8ypos = "0";
2260  TString dzpos = "0";
2261  TString name = "";
2262  TString tempattr;
2263 
2264  while (attr != 0) {
2265 
2266  tempattr = gdml->GetAttrName(attr);
2267  tempattr.ToLower();
2268 
2269  if (tempattr == "name") {
2270  name = gdml->GetAttrValue(attr);
2271  } else if (tempattr == "v1x") {
2272  v1xpos = gdml->GetAttrValue(attr);
2273  } else if (tempattr == "v1y") {
2274  v1ypos = gdml->GetAttrValue(attr);
2275  } else if (tempattr == "v2x") {
2276  v2xpos = gdml->GetAttrValue(attr);
2277  } else if (tempattr == "v2y") {
2278  v2ypos = gdml->GetAttrValue(attr);
2279  } else if (tempattr == "v3x") {
2280  v3xpos = gdml->GetAttrValue(attr);
2281  } else if (tempattr == "v3y") {
2282  v3ypos = gdml->GetAttrValue(attr);
2283  } else if (tempattr == "v4x") {
2284  v4xpos = gdml->GetAttrValue(attr);
2285  } else if (tempattr == "v4y") {
2286  v4ypos = gdml->GetAttrValue(attr);
2287  } else if (tempattr == "v5x") {
2288  v5xpos = gdml->GetAttrValue(attr);
2289  } else if (tempattr == "v5y") {
2290  v5ypos = gdml->GetAttrValue(attr);
2291  } else if (tempattr == "v6x") {
2292  v6xpos = gdml->GetAttrValue(attr);
2293  } else if (tempattr == "v6y") {
2294  v6ypos = gdml->GetAttrValue(attr);
2295  } else if (tempattr == "v7x") {
2296  v7xpos = gdml->GetAttrValue(attr);
2297  } else if (tempattr == "v7y") {
2298  v7ypos = gdml->GetAttrValue(attr);
2299  } else if (tempattr == "v8x") {
2300  v8xpos = gdml->GetAttrValue(attr);
2301  } else if (tempattr == "v8y") {
2302  v8ypos = gdml->GetAttrValue(attr);
2303  } else if (tempattr == "dz") {
2304  dzpos = gdml->GetAttrValue(attr);
2305  } else if (tempattr == "lunit") {
2306  lunit = gdml->GetAttrValue(attr);
2307  }
2308 
2309  attr = gdml->GetNextAttr(attr);
2310  }
2311 
2312  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2313  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2314  }
2315 
2316  Double_t retunit = GetScaleVal(lunit);
2317 
2318  Double_t v1x = Value(v1xpos)*retunit;
2319  Double_t v1y = Value(v1ypos)*retunit;
2320  Double_t v2x = Value(v2xpos)*retunit;
2321  Double_t v2y = Value(v2ypos)*retunit;
2322  Double_t v3x = Value(v3xpos)*retunit;
2323  Double_t v3y = Value(v3ypos)*retunit;
2324  Double_t v4x = Value(v4xpos)*retunit;
2325  Double_t v4y = Value(v4ypos)*retunit;
2326  Double_t v5x = Value(v5xpos)*retunit;
2327  Double_t v5y = Value(v5ypos)*retunit;
2328  Double_t v6x = Value(v6xpos)*retunit;
2329  Double_t v6y = Value(v6ypos)*retunit;
2330  Double_t v7x = Value(v7xpos)*retunit;
2331  Double_t v7y = Value(v7ypos)*retunit;
2332  Double_t v8x = Value(v8xpos)*retunit;
2333  Double_t v8y = Value(v8ypos)*retunit;
2334  Double_t dz = Value(dzpos)*retunit;
2335 
2336 
2337  TGeoArb8* arb8 = new TGeoArb8(NameShort(name), dz);
2338 
2339  arb8->SetVertex(0, v1x, v1y);
2340  arb8->SetVertex(1, v2x, v2y);
2341  arb8->SetVertex(2, v3x, v3y);
2342  arb8->SetVertex(3, v4x, v4y);
2343  arb8->SetVertex(4, v5x, v5y);
2344  arb8->SetVertex(5, v6x, v6y);
2345  arb8->SetVertex(6, v7x, v7y);
2346  arb8->SetVertex(7, v8x, v8y);
2347 
2348  fsolmap[name.Data()] = arb8;
2349 
2350  return node;
2351 
2352 }
2353 
2354 ////////////////////////////////////////////////////////////////////////////////
2355 ///In the solids section of the GDML file, a Tube may be declared.
2356 ///when the tube keyword is found, this function is called, and the
2357 ///dimensions required are taken and stored, these are then bound and
2358 ///converted to type TGeoTubeSeg and stored in fsolmap map using the name
2359 ///as its key.
2360 
2362 {
2363  TString lunit = "mm";
2364  TString aunit = "rad";
2365  TString rmin = "0";
2366  TString rmax = "0";
2367  TString z = "0";
2368  TString startphi = "0";
2369  TString deltaphi = "0";
2370  TString name = "";
2371  TString tempattr;
2372 
2373  while (attr != 0) {
2374 
2375  tempattr = gdml->GetAttrName(attr);
2376  tempattr.ToLower();
2377 
2378  if (tempattr == "name") {
2379  name = gdml->GetAttrValue(attr);
2380  } else if (tempattr == "rmin") {
2381  rmin = gdml->GetAttrValue(attr);
2382  } else if (tempattr == "rmax") {
2383  rmax = gdml->GetAttrValue(attr);
2384  } else if (tempattr == "z") {
2385  z = gdml->GetAttrValue(attr);
2386  } else if (tempattr == "lunit") {
2387  lunit = gdml->GetAttrValue(attr);
2388  } else if (tempattr == "aunit") {
2389  aunit = gdml->GetAttrValue(attr);
2390  } else if (tempattr == "startphi") {
2391  startphi = gdml->GetAttrValue(attr);
2392  } else if (tempattr == "deltaphi") {
2393  deltaphi = gdml->GetAttrValue(attr);
2394  }
2395 
2396  attr = gdml->GetNextAttr(attr);
2397  }
2398 
2399  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2400  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2401  }
2402 
2403  Double_t retlunit = GetScaleVal(lunit);
2404  Double_t retaunit = GetScaleVal(aunit);
2405 
2406  Double_t rminline = Value(rmin)*retlunit;
2407  Double_t rmaxline = Value(rmax)*retlunit;
2408  Double_t zline = Value(z)*retlunit;
2409  Double_t startphideg = Value(startphi)*retaunit;
2410  Double_t deltaphideg = Value(deltaphi)*retaunit;
2411  Double_t endphideg = startphideg + deltaphideg;
2412 
2413  TGeoShape *tube = 0;
2414  if (deltaphideg < 360.)
2415  tube = new TGeoTubeSeg(NameShort(name), rminline,
2416  rmaxline,
2417  zline / 2,
2418  startphideg,
2419  endphideg);
2420  else
2421  tube = new TGeoTube(NameShort(name), rminline,
2422  rmaxline,
2423  zline / 2);
2424  fsolmap[name.Data()] = tube;
2425 
2426  return node;
2427 
2428 }
2429 
2430 ////////////////////////////////////////////////////////////////////////////////
2431 ///In the solids section of the GDML file, a Cut Tube may be declared.
2432 ///when the cutTube keyword is found, this function is called, and the
2433 ///dimensions required are taken and stored, these are then bound and
2434 ///converted to type TGeoCtub and stored in fsolmap map using the name
2435 ///as its key.
2436 
2438 {
2439  TString lunit = "mm";
2440  TString aunit = "rad";
2441  TString rmin = "0";
2442  TString rmax = "0";
2443  TString z = "0";
2444  TString startphi = "0";
2445  TString deltaphi = "0";
2446  TString lowX = "0";
2447  TString lowY = "0";
2448  TString lowZ = "0";
2449  TString highX = "0";
2450  TString highY = "0";
2451  TString highZ = "0";
2452  TString name = "";
2453  TString tempattr;
2454 
2455  while (attr != 0) {
2456 
2457  tempattr = gdml->GetAttrName(attr);
2458  tempattr.ToLower();
2459 
2460  if (tempattr == "name") {
2461  name = gdml->GetAttrValue(attr);
2462  } else if (tempattr == "rmin") {
2463  rmin = gdml->GetAttrValue(attr);
2464  } else if (tempattr == "rmax") {
2465  rmax = gdml->GetAttrValue(attr);
2466  } else if (tempattr == "z") {
2467  z = gdml->GetAttrValue(attr);
2468  } else if (tempattr == "lunit") {
2469  lunit = gdml->GetAttrValue(attr);
2470  } else if (tempattr == "aunit") {
2471  aunit = gdml->GetAttrValue(attr);
2472  } else if (tempattr == "startphi") {
2473  startphi = gdml->GetAttrValue(attr);
2474  } else if (tempattr == "deltaphi") {
2475  deltaphi = gdml->GetAttrValue(attr);
2476  } else if (tempattr == "lowx") {
2477  lowX = gdml->GetAttrValue(attr);
2478  } else if (tempattr == "lowy") {
2479  lowY = gdml->GetAttrValue(attr);
2480  } else if (tempattr == "lowz") {
2481  lowZ = gdml->GetAttrValue(attr);
2482  } else if (tempattr == "highx") {
2483  highX = gdml->GetAttrValue(attr);
2484  } else if (tempattr == "highy") {
2485  highY = gdml->GetAttrValue(attr);
2486  } else if (tempattr == "highz") {
2487  highZ = gdml->GetAttrValue(attr);
2488  }
2489 
2490  attr = gdml->GetNextAttr(attr);
2491  }
2492 
2493  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2494  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2495  }
2496 
2497  Double_t retlunit = GetScaleVal(lunit);
2498  Double_t retaunit = GetScaleVal(aunit);
2499 
2500  Double_t rminline = Value(rmin)*retlunit;
2501  Double_t rmaxline = Value(rmax)*retlunit;
2502  Double_t zline = Value(z)*retlunit;
2503  Double_t startphiline = Value(startphi)*retaunit;
2504  Double_t deltaphiline = Value(deltaphi)*retaunit + startphiline;
2505  Double_t lowXline = Value(lowX)*retlunit;
2506  Double_t lowYline = Value(lowY)*retlunit;
2507  Double_t lowZline = Value(lowZ)*retlunit;
2508  Double_t highXline = Value(highX)*retlunit;
2509  Double_t highYline = Value(highY)*retlunit;
2510  Double_t highZline = Value(highZ)*retlunit;
2511 
2512 
2513  TGeoCtub* cuttube = new TGeoCtub(NameShort(name), rminline,
2514  rmaxline,
2515  zline / 2,
2516  startphiline,
2517  deltaphiline,
2518  lowXline,
2519  lowYline,
2520  lowZline,
2521  highXline,
2522  highYline,
2523  highZline);
2524 
2525 
2526  fsolmap[name.Data()] = cuttube;
2527 
2528  return node;
2529 
2530 }
2531 
2532 ////////////////////////////////////////////////////////////////////////////////
2533 ///In the solids section of the GDML file, a cone may be declared.
2534 ///when the cone keyword is found, this function is called, and the
2535 ///dimensions required are taken and stored, these are then bound and
2536 ///converted to type TGeoConSeg and stored in fsolmap map using the name
2537 ///as its key.
2538 
2540 {
2541  TString lunit = "mm";
2542  TString aunit = "rad";
2543  TString rmin1 = "0";
2544  TString rmax1 = "0";
2545  TString rmin2 = "0";
2546  TString rmax2 = "0";
2547  TString z = "0";
2548  TString startphi = "0";
2549  TString deltaphi = "0";
2550  TString name = "";
2551  TString tempattr;
2552 
2553  while (attr != 0) {
2554 
2555  tempattr = gdml->GetAttrName(attr);
2556  tempattr.ToLower();
2557 
2558  if (tempattr == "name") {
2559  name = gdml->GetAttrValue(attr);
2560  } else if (tempattr == "rmin1") {
2561  rmin1 = gdml->GetAttrValue(attr);
2562  } else if (tempattr == "rmax1") {
2563  rmax1 = gdml->GetAttrValue(attr);
2564  } else if (tempattr == "rmin2") {
2565  rmin2 = gdml->GetAttrValue(attr);
2566  } else if (tempattr == "rmax2") {
2567  rmax2 = gdml->GetAttrValue(attr);
2568  } else if (tempattr == "z") {
2569  z = gdml->GetAttrValue(attr);
2570  } else if (tempattr == "lunit") {
2571  lunit = gdml->GetAttrValue(attr);
2572  } else if (tempattr == "aunit") {
2573  aunit = gdml->GetAttrValue(attr);
2574  } else if (tempattr == "startphi") {
2575  startphi = gdml->GetAttrValue(attr);
2576  } else if (tempattr == "deltaphi") {
2577  deltaphi = gdml->GetAttrValue(attr);
2578  }
2579 
2580  attr = gdml->GetNextAttr(attr);
2581  }
2582 
2583  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2584  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2585  }
2586 
2587  Double_t retlunit = GetScaleVal(lunit);
2588  Double_t retaunit = GetScaleVal(aunit);
2589 
2590  Double_t rmin1line = Value(rmin1)*retlunit;
2591  Double_t rmax1line = Value(rmax1)*retlunit;
2592  Double_t rmin2line = Value(rmin2)*retlunit;
2593  Double_t rmax2line = Value(rmax2)*retlunit;
2594  Double_t zline = Value(z)*retlunit;
2595  Double_t sphi = Value(startphi)*retaunit;
2596  Double_t dphi = Value(deltaphi)*retaunit;
2597  Double_t ephi = sphi + dphi;
2598 
2599  TGeoShape *cone = 0;
2600  if (dphi < 360.)
2601  cone = new TGeoConeSeg(NameShort(name), zline / 2,
2602  rmin1line,
2603  rmax1line,
2604  rmin2line,
2605  rmax2line,
2606  sphi, ephi);
2607  else
2608  cone = new TGeoCone(NameShort(name), zline / 2,
2609  rmin1line,
2610  rmax1line,
2611  rmin2line,
2612  rmax2line);
2613 
2614  fsolmap[name.Data()] = cone;
2615 
2616  return node;
2617 
2618 }
2619 
2620 ////////////////////////////////////////////////////////////////////////////////
2621 ///In the solids section of the GDML file, a Trap may be declared.
2622 ///when the trap keyword is found, this function is called, and the
2623 ///dimensions required are taken and stored, these are then bound and
2624 ///converted to type TGeoTrap and stored in fsolmap map using the name
2625 ///as its key.
2626 
2628 {
2629  TString lunit = "mm";
2630  TString aunit = "rad";
2631  TString x1 = "0";
2632  TString x2 = "0";
2633  TString x3 = "0";
2634  TString x4 = "0";
2635  TString y1 = "0";
2636  TString y2 = "0";
2637  TString z = "0";
2638  TString phi = "0";
2639  TString theta = "0";
2640  TString alpha1 = "0";
2641  TString alpha2 = "0";
2642  TString name = "";
2643  TString tempattr;
2644 
2645  while (attr != 0) {
2646 
2647  tempattr = gdml->GetAttrName(attr);
2648  tempattr.ToLower();
2649 
2650  if (tempattr == "name") {
2651  name = gdml->GetAttrValue(attr);
2652  } else if (tempattr == "x1") {
2653  x1 = gdml->GetAttrValue(attr);
2654  } else if (tempattr == "x2") {
2655  x2 = gdml->GetAttrValue(attr);
2656  } else if (tempattr == "x3") {
2657  x3 = gdml->GetAttrValue(attr);
2658  } else if (tempattr == "x4") {
2659  x4 = gdml->GetAttrValue(attr);
2660  } else if (tempattr == "y1") {
2661  y1 = gdml->GetAttrValue(attr);
2662  } else if (tempattr == "y2") {
2663  y2 = gdml->GetAttrValue(attr);
2664  } else if (tempattr == "z") {
2665  z = gdml->GetAttrValue(attr);
2666  } else if (tempattr == "lunit") {
2667  lunit = gdml->GetAttrValue(attr);
2668  } else if (tempattr == "aunit") {
2669  aunit = gdml->GetAttrValue(attr);
2670  } else if (tempattr == "phi") {
2671  phi = gdml->GetAttrValue(attr);
2672  } else if (tempattr == "theta") {
2673  theta = gdml->GetAttrValue(attr);
2674  } else if (tempattr == "alpha1") {
2675  alpha1 = gdml->GetAttrValue(attr);
2676  } else if (tempattr == "alpha2") {
2677  alpha2 = gdml->GetAttrValue(attr);
2678  }
2679 
2680  attr = gdml->GetNextAttr(attr);
2681  }
2682 
2683  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2684  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2685  }
2686 
2687  Double_t retlunit = GetScaleVal(lunit);
2688  Double_t retaunit = GetScaleVal(aunit);
2689 
2690  Double_t x1line = Value(x1)*retlunit;
2691  Double_t x2line = Value(x2)*retlunit;
2692  Double_t x3line = Value(x3)*retlunit;
2693  Double_t x4line = Value(x4)*retlunit;
2694  Double_t y1line = Value(y1)*retlunit;
2695  Double_t y2line = Value(y2)*retlunit;
2696  Double_t zline = Value(z)*retlunit;
2697  Double_t philine = Value(phi)*retaunit;
2698  Double_t thetaline = Value(theta)*retaunit;
2699  Double_t alpha1line = Value(alpha1)*retaunit;
2700  Double_t alpha2line = Value(alpha2)*retaunit;
2701 
2702  TGeoTrap* trap = new TGeoTrap(NameShort(name), zline / 2,
2703  thetaline,
2704  philine,
2705  y1line / 2,
2706  x1line / 2,
2707  x2line / 2,
2708  alpha1line,
2709  y2line / 2,
2710  x3line / 2,
2711  x4line / 2,
2712  alpha2line);
2713 
2714  fsolmap[name.Data()] = trap;
2715 
2716  return node;
2717 
2718 }
2719 
2720 ////////////////////////////////////////////////////////////////////////////////
2721 ///In the solids section of the GDML file, a Trd may be declared.
2722 ///when the trd keyword is found, this function is called, and the
2723 ///dimensions required are taken and stored, these are then bound and
2724 ///converted to type TGeoTrd2 and stored in fsolmap map using the name
2725 ///as its key.
2726 
2728 {
2729  TString lunit = "mm";
2730  TString x1 = "0";
2731  TString x2 = "0";
2732  TString y1 = "0";
2733  TString y2 = "0";
2734  TString z = "0";
2735  TString name = "";
2736  TString tempattr;
2737 
2738  while (attr != 0) {
2739 
2740  tempattr = gdml->GetAttrName(attr);
2741  tempattr.ToLower();
2742 
2743  if (tempattr == "name") {
2744  name = gdml->GetAttrValue(attr);
2745  } else if (tempattr == "x1") {
2746  x1 = gdml->GetAttrValue(attr);
2747  } else if (tempattr == "x2") {
2748  x2 = gdml->GetAttrValue(attr);
2749  } else if (tempattr == "y1") {
2750  y1 = gdml->GetAttrValue(attr);
2751  } else if (tempattr == "y2") {
2752  y2 = gdml->GetAttrValue(attr);
2753  } else if (tempattr == "z") {
2754  z = gdml->GetAttrValue(attr);
2755  } else if (tempattr == "lunit") {
2756  lunit = gdml->GetAttrValue(attr);
2757  }
2758 
2759  attr = gdml->GetNextAttr(attr);
2760  }
2761 
2762  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2763  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2764  }
2765 
2766  Double_t retlunit = GetScaleVal(lunit);
2767 
2768  Double_t x1line = Value(x1)*retlunit;
2769  Double_t x2line = Value(x2)*retlunit;
2770  Double_t y1line = Value(y1)*retlunit;
2771  Double_t y2line = Value(y2)*retlunit;
2772  Double_t zline = Value(z)*retlunit;
2773 
2774  TGeoTrd2* trd = new TGeoTrd2(NameShort(name),
2775  x1line / 2,
2776  x2line / 2,
2777  y1line / 2,
2778  y2line / 2,
2779  zline / 2);
2780 
2781  fsolmap[name.Data()] = trd;
2782 
2783  return node;
2784 
2785 }
2786 
2787 ////////////////////////////////////////////////////////////////////////////////
2788 ///In the solids section of the GDML file, a Polycone may be declared.
2789 ///when the polycone keyword is found, this function is called, and the
2790 ///dimensions required are taken and stored, these are then bound and
2791 ///converted to type TGeoPCon and stored in fsolmap map using the name
2792 ///as its key. Polycone has Zplanes, planes along the z axis specifying
2793 ///the rmin, rmax dimenstions at that point along z.
2794 
2796 {
2797  TString lunit = "mm";
2798  TString aunit = "rad";
2799  TString rmin = "0";
2800  TString rmax = "0";
2801  TString z = "0";
2802  TString startphi = "0";
2803  TString deltaphi = "0";
2804  TString name = "";
2805  TString tempattr;
2806 
2807  while (attr != 0) {
2808 
2809  tempattr = gdml->GetAttrName(attr);
2810  tempattr.ToLower();
2811 
2812  if (tempattr == "name") {
2813  name = gdml->GetAttrValue(attr);
2814  } else if (tempattr == "lunit") {
2815  lunit = gdml->GetAttrValue(attr);
2816  } else if (tempattr == "aunit") {
2817  aunit = gdml->GetAttrValue(attr);
2818  } else if (tempattr == "startphi") {
2819  startphi = gdml->GetAttrValue(attr);
2820  } else if (tempattr == "deltaphi") {
2821  deltaphi = gdml->GetAttrValue(attr);
2822  }
2823  attr = gdml->GetNextAttr(attr);
2824  }
2825 
2826  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2827  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2828  }
2829 
2830  Double_t retlunit = GetScaleVal(lunit);
2831  Double_t retaunit = GetScaleVal(aunit);
2832 
2833  //START TO LOOK THRU CHILD (ZPLANE) NODES...
2834 
2835  XMLNodePointer_t child = gdml->GetChild(node);
2836  int numplanes = 0;
2837 
2838  while (child != 0) {
2839  numplanes = numplanes + 1;
2840  child = gdml->GetNext(child);
2841  }
2842 
2843  int cols;
2844  int i;
2845  cols = 3;
2846  double ** table = new double*[numplanes];
2847  for (i = 0; i < numplanes; i++) {
2848  table[i] = new double[cols];
2849  }
2850 
2851  child = gdml->GetChild(node);
2852  int planeno = 0;
2853 
2854  while (child != 0) {
2855  if (strcmp(gdml->GetNodeName(child), "zplane") == 0) {
2856  //removed original dec
2857  Double_t rminline = 0;
2858  Double_t rmaxline = 0;
2859  Double_t zline = 0;
2860 
2861  attr = gdml->GetFirstAttr(child);
2862 
2863  while (attr != 0) {
2864  tempattr = gdml->GetAttrName(attr);
2865  tempattr.ToLower();
2866 
2867  if (tempattr == "rmin") {
2868  rmin = gdml->GetAttrValue(attr);
2869  rminline = Value(rmin)*retlunit;
2870  table[planeno][0] = rminline;
2871  } else if (tempattr == "rmax") {
2872  rmax = gdml->GetAttrValue(attr);
2873  rmaxline = Value(rmax)*retlunit;
2874  table[planeno][1] = rmaxline;
2875  } else if (tempattr == "z") {
2876  z = gdml->GetAttrValue(attr);
2877  zline = Value(z)*retlunit;
2878  table[planeno][2] = zline;
2879  }
2880  attr = gdml->GetNextAttr(attr);
2881  }
2882  }
2883  planeno = planeno + 1;
2884  child = gdml->GetNext(child);
2885  }
2886 
2887  Double_t startphiline = Value(startphi)*retaunit;
2888  Double_t deltaphiline = Value(deltaphi)*retaunit;
2889 
2890  TGeoPcon* poly = new TGeoPcon(NameShort(name),
2891  startphiline,
2892  deltaphiline,
2893  numplanes);
2894  Int_t zno = 0;
2895 
2896  for (int j = 0; j < numplanes; j++) {
2897  poly->DefineSection(zno, table[j][2], table[j][0], table[j][1]);
2898  zno = zno + 1;
2899  }
2900 
2901  fsolmap[name.Data()] = poly;
2902  for (i = 0; i < numplanes; i++) {
2903  delete [] table[i];
2904  }
2905  delete [] table;
2906 
2907  return node;
2908 }
2909 
2910 ////////////////////////////////////////////////////////////////////////////////
2911 ///In the solids section of the GDML file, a Polyhedra may be declared.
2912 ///when the polyhedra keyword is found, this function is called, and the
2913 ///dimensions required are taken and stored, these are then bound and
2914 ///converted to type TGeoPgon and stored in fsolmap map using the name
2915 ///as its key. Polycone has Zplanes, planes along the z axis specifying
2916 ///the rmin, rmax dimenstions at that point along z.
2917 
2919 {
2920  TString lunit = "mm";
2921  TString aunit = "rad";
2922  TString rmin = "0";
2923  TString rmax = "0";
2924  TString z = "0";
2925  TString startphi = "0";
2926  TString deltaphi = "0";
2927  TString numsides = "1";
2928  TString name = "";
2929  TString tempattr;
2930 
2931  while (attr != 0) {
2932 
2933  tempattr = gdml->GetAttrName(attr);
2934  tempattr.ToLower();
2935 
2936  if (tempattr == "name") {
2937  name = gdml->GetAttrValue(attr);
2938  } else if (tempattr == "lunit") {
2939  lunit = gdml->GetAttrValue(attr);
2940  } else if (tempattr == "aunit") {
2941  aunit = gdml->GetAttrValue(attr);
2942  } else if (tempattr == "startphi") {
2943  startphi = gdml->GetAttrValue(attr);
2944  } else if (tempattr == "deltaphi") {
2945  deltaphi = gdml->GetAttrValue(attr);
2946  } else if (tempattr == "numsides") {
2947  numsides = gdml->GetAttrValue(attr);
2948  }
2949 
2950  attr = gdml->GetNextAttr(attr);
2951  }
2952 
2953  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
2954  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
2955  }
2956 
2957  Double_t retlunit = GetScaleVal(lunit);
2958  Double_t retaunit = GetScaleVal(aunit);
2959 
2960  //START TO LOOK THRU CHILD (ZPLANE) NODES...
2961 
2962  XMLNodePointer_t child = gdml->GetChild(node);
2963  int numplanes = 0;
2964 
2965  while (child != 0) {
2966  numplanes = numplanes + 1;
2967  child = gdml->GetNext(child);
2968  }
2969 
2970  int cols;
2971  int i;
2972  cols = 3;
2973  double ** table = new double*[numplanes];
2974  for (i = 0; i < numplanes; i++) {
2975  table[i] = new double[cols];
2976  }
2977 
2978  child = gdml->GetChild(node);
2979  int planeno = 0;
2980 
2981  while (child != 0) {
2982  if (strcmp(gdml->GetNodeName(child), "zplane") == 0) {
2983 
2984  Double_t rminline = 0;
2985  Double_t rmaxline = 0;
2986  Double_t zline = 0;
2987  attr = gdml->GetFirstAttr(child);
2988 
2989  while (attr != 0) {
2990  tempattr = gdml->GetAttrName(attr);
2991  tempattr.ToLower();
2992 
2993  if (tempattr == "rmin") {
2994  rmin = gdml->GetAttrValue(attr);
2995  rminline = Value(rmin)*retlunit;
2996  table[planeno][0] = rminline;
2997  } else if (tempattr == "rmax") {
2998  rmax = gdml->GetAttrValue(attr);
2999  rmaxline = Value(rmax)*retlunit;
3000  table[planeno][1] = rmaxline;
3001  } else if (tempattr == "z") {
3002  z = gdml->GetAttrValue(attr);
3003  zline = Value(z)*retlunit;
3004  table[planeno][2] = zline;
3005  }
3006 
3007  attr = gdml->GetNextAttr(attr);
3008  }
3009  }
3010  planeno = planeno + 1;
3011  child = gdml->GetNext(child);
3012  }
3013 
3014  Double_t startphiline = Value(startphi)*retaunit;
3015  Double_t deltaphiline = Value(deltaphi)*retaunit;
3016  Int_t numsidesline = (int)Value(numsides);
3017 
3018  TGeoPgon* polyg = new TGeoPgon(NameShort(name),
3019  startphiline,
3020  deltaphiline,
3021  numsidesline,
3022  numplanes);
3023  Int_t zno = 0;
3024 
3025  for (int j = 0; j < numplanes; j++) {
3026  polyg->DefineSection(zno, table[j][2], table[j][0], table[j][1]);
3027  zno = zno + 1;
3028  }
3029 
3030  fsolmap[name.Data()] = polyg;
3031  for (i = 0; i < numplanes; i++) {
3032  delete [] table[i];
3033  }
3034  delete [] table;
3035 
3036  return node;
3037 
3038 }
3039 
3040 ////////////////////////////////////////////////////////////////////////////////
3041 ///In the solids section of the GDML file, a Sphere may be declared.
3042 ///when the sphere keyword is found, this function is called, and the
3043 ///dimensions required are taken and stored, these are then bound and
3044 ///converted to type TGeoSphere and stored in fsolmap map using the name
3045 ///as its key.
3046 
3048 {
3049  TString lunit = "mm";
3050  TString aunit = "rad";
3051  TString rmin = "0";
3052  TString rmax = "0";
3053  TString startphi = "0";
3054  TString deltaphi = "0";
3055  TString starttheta = "0";
3056  TString deltatheta = "0";
3057  TString name = "";
3058  TString tempattr;
3059 
3060  while (attr != 0) {
3061  tempattr = gdml->GetAttrName(attr);
3062  tempattr.ToLower();
3063 
3064  if (tempattr == "name") {
3065  name = gdml->GetAttrValue(attr);
3066  } else if (tempattr == "rmin") {
3067  rmin = gdml->GetAttrValue(attr);
3068  } else if (tempattr == "rmax") {
3069  rmax = gdml->GetAttrValue(attr);
3070  } else if (tempattr == "lunit") {
3071  lunit = gdml->GetAttrValue(attr);
3072  } else if (tempattr == "aunit") {
3073  aunit = gdml->GetAttrValue(attr);
3074  } else if (tempattr == "startphi") {
3075  startphi = gdml->GetAttrValue(attr);
3076  } else if (tempattr == "deltaphi") {
3077  deltaphi = gdml->GetAttrValue(attr);
3078  } else if (tempattr == "starttheta") {
3079  starttheta = gdml->GetAttrValue(attr);
3080  } else if (tempattr == "deltatheta") {
3081  deltatheta = gdml->GetAttrValue(attr);
3082  }
3083 
3084  attr = gdml->GetNextAttr(attr);
3085  }
3086 
3087  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3088  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3089  }
3090 
3091  Double_t retlunit = GetScaleVal(lunit);
3092  Double_t retaunit = GetScaleVal(aunit);
3093 
3094  Double_t rminline = Value(rmin)*retlunit;
3095  Double_t rmaxline = Value(rmax)*retlunit;
3096  Double_t startphiline = Value(startphi)*retaunit;
3097  Double_t deltaphiline = startphiline+ Value(deltaphi)*retaunit;
3098  Double_t startthetaline = Value(starttheta)*retaunit;
3099  Double_t deltathetaline = startthetaline + Value(deltatheta)*retaunit;
3100 
3101  TGeoSphere* sphere = new TGeoSphere(NameShort(name),
3102  rminline,
3103  rmaxline,
3104  startthetaline,
3105  deltathetaline,
3106  startphiline,
3107  deltaphiline);
3108 
3109  fsolmap[name.Data()] = sphere;
3110 
3111  return node;
3112 
3113 }
3114 
3115 ////////////////////////////////////////////////////////////////////////////////
3116 ///In the solids section of the GDML file, a Torus may be declared.
3117 ///when the torus keyword is found, this function is called, and the
3118 ///dimensions required are taken and stored, these are then bound and
3119 ///converted to type TGeoTorus and stored in fsolmap map using the name
3120 ///as its key.
3121 
3123 {
3124  TString lunit = "mm";
3125  TString aunit = "rad";
3126  TString rmin = "0";
3127  TString rmax = "0";
3128  TString rtor = "0";
3129  TString startphi = "0";
3130  TString deltaphi = "0";
3131  TString name = "";
3132  TString tempattr;
3133 
3134  while (attr != 0) {
3135 
3136  tempattr = gdml->GetAttrName(attr);
3137  tempattr.ToLower();
3138 
3139  if (tempattr == "name") {
3140  name = gdml->GetAttrValue(attr);
3141  } else if (tempattr == "rmin") {
3142  rmin = gdml->GetAttrValue(attr);
3143  } else if (tempattr == "rmax") {
3144  rmax = gdml->GetAttrValue(attr);
3145  } else if (tempattr == "rtor") {
3146  rtor = gdml->GetAttrValue(attr);
3147  } else if (tempattr == "lunit") {
3148  lunit = gdml->GetAttrValue(attr);
3149  } else if (tempattr == "aunit") {
3150  aunit = gdml->GetAttrValue(attr);
3151  } else if (tempattr == "startphi") {
3152  startphi = gdml->GetAttrValue(attr);
3153  } else if (tempattr == "deltaphi") {
3154  deltaphi = gdml->GetAttrValue(attr);
3155  }
3156 
3157  attr = gdml->GetNextAttr(attr);
3158  }
3159 
3160  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3161  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3162  }
3163 
3164  Double_t retlunit = GetScaleVal(lunit);
3165  Double_t retaunit = GetScaleVal(aunit);
3166 
3167  Double_t rminline = Value(rmin)*retlunit;
3168  Double_t rmaxline = Value(rmax)*retlunit;
3169  Double_t rtorline = Value(rtor)*retlunit;
3170  Double_t startphiline = Value(startphi)*retaunit;
3171  Double_t deltaphiline = Value(deltaphi)*retaunit;
3172 
3173 
3174  TGeoTorus* torus = new TGeoTorus(NameShort(name), rtorline,
3175  rminline,
3176  rmaxline,
3177  startphiline,
3178  deltaphiline);
3179 
3180  fsolmap[name.Data()] = torus;
3181 
3182  return node;
3183 
3184 }
3185 
3186 ////////////////////////////////////////////////////////////////////////////////
3187 ///In the solids section of the GDML file, a Hype may be declared.
3188 ///when the hype keyword is found, this function is called, and the
3189 ///dimensions required are taken and stored, these are then bound and
3190 ///converted to type TGeoHype and stored in fsolmap map using the name
3191 ///as its key.
3192 
3194 {
3195  TString lunit = "mm";
3196  TString aunit = "rad";
3197  TString rmin = "0";
3198  TString rmax = "0";
3199  TString z = "0";
3200  TString inst = "0";
3201  TString outst = "0";
3202  TString name = "";
3203  TString tempattr;
3204 
3205  while (attr != 0) {
3206  tempattr = gdml->GetAttrName(attr);
3207  tempattr.ToLower();
3208 
3209  if (tempattr == "name") {
3210  name = gdml->GetAttrValue(attr);
3211  } else if (tempattr == "rmin") {
3212  rmin = gdml->GetAttrValue(attr);
3213  } else if (tempattr == "rmax") {
3214  rmax = gdml->GetAttrValue(attr);
3215  } else if (tempattr == "z") {
3216  z = gdml->GetAttrValue(attr);
3217  } else if (tempattr == "lunit") {
3218  lunit = gdml->GetAttrValue(attr);
3219  } else if (tempattr == "aunit") {
3220  aunit = gdml->GetAttrValue(attr);
3221  } else if (tempattr == "inst") {
3222  inst = gdml->GetAttrValue(attr);
3223  } else if (tempattr == "outst") {
3224  outst = gdml->GetAttrValue(attr);
3225  }
3226 
3227  attr = gdml->GetNextAttr(attr);
3228  }
3229 
3230  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3231  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3232  }
3233 
3234  Double_t retlunit = GetScaleVal(lunit);
3235  Double_t retaunit = GetScaleVal(aunit);
3236 
3237  Double_t rminline = Value(rmin)*retlunit;
3238  Double_t rmaxline = Value(rmax)*retlunit;
3239  Double_t zline = Value(z)*retlunit;
3240  Double_t instline = Value(inst)*retaunit;
3241  Double_t outstline = Value(outst)*retaunit;
3242 
3243 
3244  TGeoHype* hype = new TGeoHype(NameShort(name),
3245  rminline,
3246  instline,
3247  rmaxline,
3248  outstline,
3249  zline / 2);
3250 
3251  fsolmap[name.Data()] = hype;
3252 
3253  return node;
3254 
3255 }
3256 
3257 ////////////////////////////////////////////////////////////////////////////////
3258 ///In the solids section of the GDML file, a Para may be declared.
3259 ///when the para keyword is found, this function is called, and the
3260 ///dimensions required are taken and stored, these are then bound and
3261 ///converted to type TGeoPara and stored in fsolmap map using the name
3262 ///as its key.
3263 
3265 {
3266  TString lunit = "mm";
3267  TString aunit = "rad";
3268  TString x = "0";
3269  TString y = "0";
3270  TString z = "0";
3271  TString phi = "0";
3272  TString theta = "0";
3273  TString alpha = "0";
3274  TString name = "";
3275  TString tempattr;
3276 
3277  while (attr != 0) {
3278 
3279  tempattr = gdml->GetAttrName(attr);
3280  tempattr.ToLower();
3281 
3282  if (tempattr == "name") {
3283  name = gdml->GetAttrValue(attr);
3284  } else if (tempattr == "x") {
3285  x = gdml->GetAttrValue(attr);
3286  } else if (tempattr == "y") {
3287  y = gdml->GetAttrValue(attr);
3288  } else if (tempattr == "z") {
3289  z = gdml->GetAttrValue(attr);
3290  } else if (tempattr == "lunit") {
3291  lunit = gdml->GetAttrValue(attr);
3292  } else if (tempattr == "aunit") {
3293  aunit = gdml->GetAttrValue(attr);
3294  } else if (tempattr == "phi") {
3295  phi = gdml->GetAttrValue(attr);
3296  } else if (tempattr == "theta") {
3297  theta = gdml->GetAttrValue(attr);
3298  } else if (tempattr == "alpha") {
3299  alpha = gdml->GetAttrValue(attr);
3300  }
3301 
3302  attr = gdml->GetNextAttr(attr);
3303  }
3304 
3305  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3306  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3307  }
3308 
3309  Double_t retlunit = GetScaleVal(lunit);
3310  Double_t retaunit = GetScaleVal(aunit);
3311 
3312  Double_t xline = Value(x)*retlunit;
3313  Double_t yline = Value(y)*retlunit;
3314  Double_t zline = Value(z)*retlunit;
3315  Double_t philine = Value(phi)*retaunit;
3316  Double_t alphaline = Value(alpha)*retaunit;
3317  Double_t thetaline = Value(theta)*retaunit;
3318 
3319 
3320  TGeoPara* para = new TGeoPara(NameShort(name),
3321  xline / 2,
3322  yline / 2,
3323  zline / 2,
3324  alphaline,
3325  thetaline,
3326  philine);
3327 
3328  fsolmap[name.Data()] = para;
3329 
3330  return node;
3331 
3332 }
3333 
3334 ////////////////////////////////////////////////////////////////////////////////
3335 ///In the solids section of the GDML file, a TwistTrap may be declared.
3336 ///when the twistedtrap keyword is found, this function is called, and the
3337 ///dimensions required are taken and stored, these are then bound and
3338 ///converted to type TGeoGTra and stored in fsolmap map using the name
3339 ///as its key.
3340 
3342 {
3343  TString lunit = "mm";
3344  TString aunit = "rad";
3345  TString x1 = "0";
3346  TString x2 = "0";
3347  TString x3 = "0";
3348  TString x4 = "0";
3349  TString y1 = "0";
3350  TString y2 = "0";
3351  TString z = "0";
3352  TString phi = "0";
3353  TString theta = "0";
3354  TString alpha1 = "0";
3355  TString alpha2 = "0";
3356  TString twist = "0";
3357  TString name = "";
3358  TString tempattr;
3359 
3360  while (attr != 0) {
3361 
3362  tempattr = gdml->GetAttrName(attr);
3363  tempattr.ToLower();
3364 
3365  if (tempattr == "name") {
3366  name = gdml->GetAttrValue(attr);
3367  } else if (tempattr == "x1") {
3368  x1 = gdml->GetAttrValue(attr);
3369  } else if (tempattr == "x2") {
3370  x2 = gdml->GetAttrValue(attr);
3371  } else if (tempattr == "x3") {
3372  x3 = gdml->GetAttrValue(attr);
3373  } else if (tempattr == "x4") {
3374  x4 = gdml->GetAttrValue(attr);
3375  } else if (tempattr == "y1") {
3376  y1 = gdml->GetAttrValue(attr);
3377  } else if (tempattr == "y2") {
3378  y2 = gdml->GetAttrValue(attr);
3379  } else if (tempattr == "z") {
3380  z = gdml->GetAttrValue(attr);
3381  } else if (tempattr == "lunit") {
3382  lunit = gdml->GetAttrValue(attr);
3383  } else if (tempattr == "aunit") {
3384  aunit = gdml->GetAttrValue(attr);
3385  } else if (tempattr == "phi") {
3386  phi = gdml->GetAttrValue(attr);
3387  } else if (tempattr == "theta") {
3388  theta = gdml->GetAttrValue(attr);
3389  } else if (tempattr == "alph") { //gdml schema knows only alph attribute
3390  alpha1 = gdml->GetAttrValue(attr);
3391  alpha2 = alpha1;
3392  //} else if (tempattr == "alpha2") {
3393  // alpha2 = gdml->GetAttrValue(attr);
3394  } else if (tempattr == "phitwist") {
3395  twist = gdml->GetAttrValue(attr);
3396  }
3397 
3398  attr = gdml->GetNextAttr(attr);
3399  }
3400 
3401  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3402  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3403  }
3404 
3405  Double_t retlunit = GetScaleVal(lunit);
3406  Double_t retaunit = GetScaleVal(aunit);
3407 
3408  Double_t x1line = Value(x1)*retlunit;
3409  Double_t x2line = Value(x2)*retlunit;
3410  Double_t x3line = Value(x3)*retlunit;
3411  Double_t x4line = Value(x4)*retlunit;
3412  Double_t y1line = Value(y1)*retlunit;
3413  Double_t y2line = Value(y2)*retlunit;
3414  Double_t zline = Value(z)*retlunit;
3415  Double_t philine = Value(phi)*retaunit;
3416  Double_t thetaline = Value(theta)*retaunit;
3417  Double_t alpha1line = Value(alpha1)*retaunit;
3418  Double_t alpha2line = Value(alpha2)*retaunit;
3419  Double_t twistline = Value(twist)*retaunit;
3420 
3421 
3422  TGeoGtra* twtrap = new TGeoGtra(NameShort(name), zline / 2,
3423  thetaline,
3424  philine,
3425  twistline,
3426  y1line / 2,
3427  x1line / 2,
3428  x2line / 2,
3429  alpha1line,
3430  y2line / 2,
3431  x3line / 2,
3432  x4line / 2,
3433  alpha2line);
3434 
3435  fsolmap[name.Data()] = twtrap;
3436 
3437  return node;
3438 
3439 }
3440 
3441 
3442 ////////////////////////////////////////////////////////////////////////////////
3443 ///In the solids section of the GDML file, a ElTube may be declared.
3444 ///when the eltube keyword is found, this function is called, and the
3445 ///dimensions required are taken and stored, these are then bound and
3446 ///converted to type TGeoEltu and stored in fsolmap map using the name
3447 ///as its key.
3448 
3450 {
3451  TString lunit = "mm";
3452  TString xpos = "0";
3453  TString ypos = "0";
3454  TString zpos = "0";
3455  TString name = "";
3456  TString tempattr;
3457 
3458  while (attr != 0) {
3459 
3460  tempattr = gdml->GetAttrName(attr);
3461  tempattr.ToLower();
3462 
3463  if (tempattr == "name") {
3464  name = gdml->GetAttrValue(attr);
3465  } else if (tempattr == "dx") {
3466  xpos = gdml->GetAttrValue(attr);
3467  } else if (tempattr == "dy") {
3468  ypos = gdml->GetAttrValue(attr);
3469  } else if (tempattr == "dz") {
3470  zpos = gdml->GetAttrValue(attr);
3471  } else if (tempattr == "lunit") {
3472  lunit = gdml->GetAttrValue(attr);
3473  }
3474 
3475  attr = gdml->GetNextAttr(attr);
3476  }
3477 
3478  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3479  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3480  }
3481 
3482  Double_t retunit = GetScaleVal(lunit);
3483 
3484  Double_t xline = Value(xpos)*retunit;
3485  Double_t yline = Value(ypos)*retunit;
3486  Double_t zline = Value(zpos)*retunit;
3487 
3488  TGeoEltu* eltu = new TGeoEltu(NameShort(name), xline,
3489  yline,
3490  zline);
3491 
3492  fsolmap[name.Data()] = eltu;
3493 
3494  return node;
3495 
3496 }
3497 ////////////////////////////////////////////////////////////////////////////////
3498 ///In the solids section of the GDML file, an Orb may be declared.
3499 ///when the orb keyword is found, this function is called, and the
3500 ///dimensions required are taken and stored, these are then bound and
3501 ///converted to type TGeoSphere and stored in fsolmap map using the name
3502 ///as its key.
3503 
3505 {
3506  TString lunit = "mm";
3507  TString r = "0";
3508  TString name = "";
3509  TString tempattr;
3510 
3511  while (attr != 0) {
3512 
3513  tempattr = gdml->GetAttrName(attr);
3514  tempattr.ToLower();
3515 
3516  if (tempattr == "name") {
3517  name = gdml->GetAttrValue(attr);
3518  } else if (tempattr == "r") {
3519  r = gdml->GetAttrValue(attr);
3520  } else if (tempattr == "lunit") {
3521  lunit = gdml->GetAttrValue(attr);
3522  }
3523 
3524  attr = gdml->GetNextAttr(attr);
3525  }
3526 
3527  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3528  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3529  }
3530 
3531  Double_t retunit = GetScaleVal(lunit);
3532 
3533  Double_t rline = Value(r)*retunit;
3534 
3535  TGeoSphere* orb = new TGeoSphere(NameShort(name), 0, rline, 0, 180, 0, 360);
3536 
3537  fsolmap[name.Data()] = orb;
3538 
3539  return node;
3540 
3541 }
3542 
3543 
3544 ////////////////////////////////////////////////////////////////////////////////
3545 ///In the solids section of the GDML file, an Xtru may be declared.
3546 ///when the xtru keyword is found, this function is called, and the
3547 ///dimensions required are taken and stored, these are then bound and
3548 ///converted to type TGeoXtru and stored in fsolmap map using the name
3549 ///as its key. The xtru has child nodes of either 'twoDimVertex'or
3550 ///'section'. These two nodes define the real structure of the shape.
3551 ///The twoDimVertex's define the x,y sizes of a vertice. The section links
3552 ///the vertice to a position within the xtru.
3553 
3555 {
3556  TString lunit = "mm";
3557 // TString aunit = "rad";
3558  TString x = "0";
3559  TString y = "0";
3560  TString zorder = "0";
3561  TString zpos = "0";
3562  TString xoff = "0";
3563  TString yoff = "0";
3564  TString scale = "0";
3565  TString name = "";
3566  TString tempattr;
3567 
3568  while (attr != 0) {
3569 
3570  tempattr = gdml->GetAttrName(attr);
3571  tempattr.ToLower();
3572 
3573  if (tempattr == "name") {
3574  name = gdml->GetAttrValue(attr);
3575  } else if (tempattr == "lunit") {
3576  lunit = gdml->GetAttrValue(attr);
3577  }
3578 
3579  attr = gdml->GetNextAttr(attr);
3580  }
3581 
3582  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3583  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3584  }
3585 
3586  Double_t retlunit = GetScaleVal(lunit);
3587 
3588  //START TO LOOK THRU CHILD NODES...
3589 
3590  XMLNodePointer_t child = gdml->GetChild(node);
3591  int nosects = 0;
3592  int noverts = 0;
3593 
3594  while (child != 0) {
3595  tempattr = gdml->GetNodeName(child);
3596 
3597  if (tempattr == "twoDimVertex") {
3598  noverts = noverts + 1;
3599  } else if (tempattr == "section") {
3600  nosects = nosects + 1;
3601  }
3602 
3603  child = gdml->GetNext(child);
3604  }
3605 
3606  //Build the dynamic arrays..
3607  int cols;
3608  int i;
3609  double *vertx = new double[noverts];
3610  double *verty = new double[noverts];
3611  cols = 5;
3612  double ** section = new double*[nosects];
3613  for (i = 0; i < nosects; i++) {
3614  section[i] = new double[cols];
3615  }
3616 
3617  child = gdml->GetChild(node);
3618  int sect = 0;
3619  int vert = 0;
3620 
3621  while (child != 0) {
3622  if (strcmp(gdml->GetNodeName(child), "twoDimVertex") == 0) {
3623  Double_t xline = 0;
3624  Double_t yline = 0;
3625 
3626  attr = gdml->GetFirstAttr(child);
3627 
3628  while (attr != 0) {
3629  tempattr = gdml->GetAttrName(attr);
3630 
3631  if (tempattr == "x") {
3632  x = gdml->GetAttrValue(attr);
3633  xline = Value(x)*retlunit;
3634  vertx[vert] = xline;
3635  } else if (tempattr == "y") {
3636  y = gdml->GetAttrValue(attr);
3637  yline = Value(y)*retlunit;
3638  verty[vert] = yline;
3639  }
3640 
3641  attr = gdml->GetNextAttr(attr);
3642  }
3643 
3644  vert = vert + 1;
3645  }
3646 
3647  else if (strcmp(gdml->GetNodeName(child), "section") == 0) {
3648 
3649  Double_t zposline = 0;
3650  Double_t xoffline = 0;
3651  Double_t yoffline = 0;
3652 
3653  attr = gdml->GetFirstAttr(child);
3654 
3655  while (attr != 0) {
3656  tempattr = gdml->GetAttrName(attr);
3657 
3658  if (tempattr == "zOrder") {
3659  zorder = gdml->GetAttrValue(attr);
3660  section[sect][0] = Value(zorder);
3661  } else if (tempattr == "zPosition") {
3662  zpos = gdml->GetAttrValue(attr);
3663  zposline = Value(zpos)*retlunit;
3664  section[sect][1] = zposline;
3665  } else if (tempattr == "xOffset") {
3666  xoff = gdml->GetAttrValue(attr);
3667  xoffline = Value(xoff)*retlunit;
3668  section[sect][2] = xoffline;
3669  } else if (tempattr == "yOffset") {
3670  yoff = gdml->GetAttrValue(attr);
3671  yoffline = Value(yoff)*retlunit;
3672  section[sect][3] = yoffline;
3673  } else if (tempattr == "scalingFactor") {
3674  scale = gdml->GetAttrValue(attr);
3675  section[sect][4] = Value(scale);
3676  }
3677 
3678  attr = gdml->GetNextAttr(attr);
3679  }
3680 
3681  sect = sect + 1;
3682  }
3683  child = gdml->GetNext(child);
3684  }
3685 
3686  TGeoXtru* xtru = new TGeoXtru(nosects);
3687  xtru->SetName(NameShort(name));
3688  xtru->DefinePolygon(vert, vertx, verty);
3689 
3690  for (int j = 0; j < sect; j++) {
3691  xtru->DefineSection((int)section[j][0], section[j][1], section[j][2], section[j][3], section[j][4]);
3692  }
3693 
3694  fsolmap[name.Data()] = xtru;
3695  delete [] vertx;
3696  delete [] verty;
3697  for (i = 0; i < nosects; i++) {
3698  delete [] section[i];
3699  }
3700  delete [] section;
3701  return node;
3702 }
3703 
3704 ////////////////////////////////////////////////////////////////////////////////
3705 ///In the solids section of the GDML file, a Reflected Solid may be
3706 ///declared when the ReflectedSolid keyword is found, this function
3707 ///is called. The rotation, position and scale for the reflection are
3708 ///applied to a matrix that is then stored in the class object
3709 ///TGDMLRefl. This is then stored in the map freflsolidmap, with
3710 ///the reflection name as a reference. also the name of the solid to
3711 ///be reflected is stored in a map called freflectmap with the reflection
3712 ///name as a reference.
3713 
3715 {
3716  std::cout << "WARNING! The reflectedSolid is obsolete! Use scale transformation instead!" << std::endl;
3717 
3718  TString sx = "0";
3719  TString sy = "0";
3720  TString sz = "0";
3721  TString rx = "0";
3722  TString ry = "0";
3723  TString rz = "0";
3724  TString dx = "0";
3725  TString dy = "0";
3726  TString dz = "0";
3727  TString name = "0";
3728  TString solid = "0";
3729  TString tempattr;
3730 
3731  while (attr != 0) {
3732 
3733  tempattr = gdml->GetAttrName(attr);
3734  tempattr.ToLower();
3735 
3736  if (tempattr == "name") {
3737  name = gdml->GetAttrValue(attr);
3738  } else if (tempattr == "sx") {
3739  sx = gdml->GetAttrValue(attr);
3740  } else if (tempattr == "sy") {
3741  sy = gdml->GetAttrValue(attr);
3742  } else if (tempattr == "sz") {
3743  sz = gdml->GetAttrValue(attr);
3744  } else if (tempattr == "rx") {
3745  rx = gdml->GetAttrValue(attr);
3746  } else if (tempattr == "ry") {
3747  ry = gdml->GetAttrValue(attr);
3748  } else if (tempattr == "rz") {
3749  rz = gdml->GetAttrValue(attr);
3750  } else if (tempattr == "dx") {
3751  dx = gdml->GetAttrValue(attr);
3752  } else if (tempattr == "dy") {
3753  dy = gdml->GetAttrValue(attr);
3754  } else if (tempattr == "dz") {
3755  dz = gdml->GetAttrValue(attr);
3756  } else if (tempattr == "solid") {
3757  solid = gdml->GetAttrValue(attr);
3758  }
3759  attr = gdml->GetNextAttr(attr);
3760  }
3761 
3762  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3763  name = TString::Format("%s_%s", name.Data(), fCurrentFile);
3764  }
3765  if ((strcmp(fCurrentFile, fStartFile)) != 0) {
3766  solid = TString::Format("%s_%s", solid.Data(), fCurrentFile);
3767  }
3768 
3769  TGeoRotation* rot = new TGeoRotation();
3770  rot->RotateZ(-(Value(rz)));
3771  rot->RotateY(-(Value(ry)));
3772  rot->RotateX(-(Value(rx)));
3773 
3774  if (atoi(sx) == -1) {
3775  rot->ReflectX(kTRUE);
3776  }
3777  if (atoi(sy) == -1) {
3778  rot->ReflectY(kTRUE);
3779  }
3780  if (atoi(sz) == -1) {
3781  rot->ReflectZ(kTRUE);
3782  }
3783 
3784  TGeoCombiTrans* relf_matx = new TGeoCombiTrans(Value(dx), Value(dy), Value(dz), rot);
3785 
3786  TGDMLRefl* reflsol = new TGDMLRefl(NameShort(name), solid, relf_matx);
3787  freflsolidmap[name.Data()] = reflsol;
3788  freflectmap[name.Data()] = solid;
3789 
3790  return node;
3791 }
3792 
3793 
3794 
3795 //===================================================================
3796 
3798 
3799 /******************************************************************
3800 ____________________________________________________________
3801 
3802 TGDMLRefl Class
3803 
3804 ------------------------------------------------------------
3805 
3806 This class is a helper class for TGDMLParse. It assists in the
3807 reflection process. This process takes a previously defined solid
3808 and can reflect the matrix of it. This class stores the name of the
3809 reflected solid, along with the name of the solid that is being
3810 reflected, and finally the reflected solid's matrix. This is then
3811 recalled when the volume is used in the structure part of the gdml
3812 file.
3813 
3814 ******************************************************************/
3815 
3816 ////////////////////////////////////////////////////////////////////////////////
3817 ///this constructor method stores the values brought in as params.
3818 
3819 TGDMLRefl::TGDMLRefl(const char* name, const char* solid, TGeoMatrix* matrix)
3820 {
3821  fNameS = name;
3822  fSolid = solid;
3823  fMatrix = matrix;
3824 }
3825 
3826 ////////////////////////////////////////////////////////////////////////////////
3827 ///this accessor method returns the matrix.
3828 
3830 {
3831  return fMatrix;
3832 }
XMLNodePointer_t UsrProcess(TXMLEngine *gdml, XMLNodePointer_t node)
User data to be processed.
XMLNodePointer_t ConProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the define section of the GDML file, constants can be declared.
Definition: TGDMLParse.cxx:403
XMLNodePointer_t TwistTrap(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a TwistTrap may be declared.
TList * GetListOfMaterials() const
Definition: TGeoManager.h:463
RotMap frotmap
Map containing position names and the TGeoTranslation for it.
Definition: TGDMLParse.h:211
virtual TGeoMatrix & Inverse() const
Return a temporary inverse of this.
Definition: TGeoMatrix.cxx:882
XMLNodePointer_t Trd(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Trd may be declared.
MatMap fmatmap
Map containing element names and the TGeoElement for it.
Definition: TGDMLParse.h:215
TGeoVolume * fWorld
Definition: TGDMLParse.h:106
Double_t Eval(Double_t x) const
Definition: TFormula.cxx:2538
Collectable string class.
Definition: TObjString.h:32
void MultiplyLeft(const TGeoMatrix *left)
multiply to the left with an other transformation if right is identity matrix, just return ...
XMLNodePointer_t Paraboloid(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Paraboloid may be declared.
XMLNodePointer_t Xtru(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, an Xtru may be declared.
virtual const Double_t * GetRotationMatrix() const =0
virtual void SetName(const char *name)
Change (i.e.
Definition: TNamed.cxx:128
XMLNodePointer_t Reflection(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Reflected Solid may be declared when the ReflectedSolid key...
EleMap felemap
Map containing isotope names and the TGeoIsotope for it.
Definition: TGDMLParse.h:214
Double_t Atof() const
Return floating-point value contained in string.
Definition: TString.cxx:2030
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition: TObject.cxx:892
virtual void RotateX(Double_t angle)
Rotate about X axis of the master frame with angle expressed in degrees.
Definition: TGeoMatrix.cxx:993
XMLNodePointer_t Cone(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a cone may be declared.
void SetTranslation(const Double_t *vect)
Definition: TGeoMatrix.h:450
static const char * filename()
void Add(TObject *obj)
This function may not be used (but we need to provide it since it is a pure virtual in TCollection)...
Definition: TMap.cxx:52
virtual void DefineSection(Int_t snum, Double_t z, Double_t x0=0., Double_t y0=0., Double_t scale=1.)
defines z position of a section plane, rmin and rmax at this z.
Definition: TGeoXtru.cxx:753
TString GetScale(const char *unit)
Throughout the GDML file, a unit can de specified.
Definition: TGDMLParse.cxx:436
XMLNodePointer_t GetNext(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
return next to xmlnode node if realnode==kTRUE, any special nodes in between will be skipped ...
Double_t RadToDeg()
Definition: TMath.h:49
Basic string class.
Definition: TString.h:137
void Multiply(const TGeoMatrix *right)
multiply to the right with an other transformation if right is identity matrix, just return ...
XMLNodePointer_t IsoProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLNodePointer_t parentn)
In the material section of the GDML file, an isotope may be declared.
Definition: TGDMLParse.cxx:691
XMLNodePointer_t Polycone(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Polycone may be declared.
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1088
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:496
void SetSkipComments(Bool_t on=kTRUE)
Definition: TXMLEngine.h:50
void ReplayCreation(const TGeoVolume *other)
Recreate the content of the other volume without pointer copying.
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
XMLAttrPointer_t GetNextAttr(XMLAttrPointer_t xmlattr)
return next attribute in the list
Definition: TXMLEngine.cxx:582
void FreeDoc(XMLDocPointer_t xmldoc)
frees allocated document data and deletes document itself
Int_t GetNdaughters() const
Definition: TGeoVolume.h:362
double GetScaleVal(const char *unit)
Throughout the GDML file, a unit can de specified.
Definition: TGDMLParse.cxx:482
TObjArray * GetNodes()
Definition: TGeoVolume.h:183
ReflectionsMap freflectmap
Map containing volume names and the TGeoVolume for it.
Definition: TGDMLParse.h:220
const char * Data() const
Definition: TString.h:349
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition: TObject.cxx:946
const char * fCurrentFile
Definition: TGDMLParse.h:111
XMLNodePointer_t SclProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the define section of the GDML file, rotations can be declared.
Definition: TGDMLParse.cxx:647
TGeoMatrix * fMatrix
solid name being reflected
Definition: TGDMLParse.h:56
static const double x2[5]
XMLNodePointer_t BooSolid(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr, int num)
In the solid section of the GDML file, boolean solids can be declared.
virtual void ReflectY(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to ZX.
Double_t x[n]
Definition: legend1.C:17
double Evaluate(const char *evalline)
Definition: TGDMLParse.cxx:352
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:2334
FormVec fformvec
Map containing files parsed during entire parsing, with their world volume name.
Definition: TGDMLParse.h:224
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
void * XMLDocPointer_t
Definition: TXMLEngine.h:22
TGeoMatrix * GetMatrix()
this accessor method returns the matrix.
TCanvas * fractions()
Definition: fractions.C:1
virtual const Double_t * GetScale() const
Definition: TGeoMatrix.h:275
double Value(const char *svalue) const
Convert number in string format to double value.
Definition: TGDMLParse.cxx:522
void AddIsotope(TGeoIsotope *isotope, Double_t relativeAbundance)
Add an isotope for this element. All isotopes have to be isotopes of the same element.
static const double x4[22]
XMLNodePointer_t Arb8(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, an Arb8 may be declared.
virtual void SetVertex(Int_t vnum, Double_t x, Double_t y)
Set values for a given vertex.
Definition: TGeoArb8.cxx:1209
virtual void RotateY(Double_t angle)
Rotate about Y axis of the master frame with angle expressed in degrees.
FileMap ffilemap
Map containing reflected volume names and the solid ref for it.
Definition: TGDMLParse.h:223
virtual TGeoVolume * Divide(const char *divname, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
Division a la G3.
const char * GetNodeName(XMLNodePointer_t xmlnode)
returns name of xmlnode
Definition: TXMLEngine.cxx:930
XMLNodePointer_t CutTube(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Cut Tube may be declared.
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
Add a TGeoNode to the list of nodes.
Definition: TGeoVolume.cxx:948
XMLNodePointer_t Sphere(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Sphere may be declared.
XMLNodePointer_t Orb(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, an Orb may be declared.
void AddElement(Double_t a, Double_t z, Double_t weight)
add an element to the mixture using fraction by weight Check if the element is already defined ...
XMLNodePointer_t PosProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the define section of the GDML file, positions can be declared.
Definition: TGDMLParse.cxx:536
XMLNodePointer_t Torus(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Torus may be declared.
XMLNodePointer_t Box(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a box may be declared.
XMLNodePointer_t Ellipsoid(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, an ellipsoid may be declared.
void ShiftToNext(XMLNodePointer_t &xmlnode, Bool_t realnode=kTRUE)
shifts specified node to next if realnode==kTRUE, any special nodes in between will be skipped ...
ROOT::R::TRInterface & r
Definition: Object.C:4
std::map< std::string, double > FracMap
Definition: TGDMLParse.h:207
XMLNodePointer_t Para(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Para may be declared.
XMLNodePointer_t Trap(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Trap may be declared.
virtual const Double_t * GetTranslation() const
Definition: TGeoMatrix.h:168
const char * fStartFile
Definition: TGDMLParse.h:110
The F O R M U L A class.
Definition: TFormula.h:89
Int_t GetId() const
Definition: TGeoMedium.h:50
virtual const Double_t * GetRotationMatrix() const
Definition: TGeoMatrix.h:234
const char * GetAttrValue(XMLAttrPointer_t xmlattr)
return value of attribute
Definition: TXMLEngine.cxx:603
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
XMLNodePointer_t RotProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the define section of the GDML file, rotations can be declared.
Definition: TGDMLParse.cxx:589
XMLAttrPointer_t GetFirstAttr(XMLNodePointer_t xmlnode)
return first attribute in the list, namespace (if exists) will be skiped
Definition: TXMLEngine.cxx:568
Bool_t IsNull() const
Definition: TString.h:387
TXMLEngine * fFileEngine[20]
Definition: TGDMLParse.h:109
XMLNodePointer_t Hype(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Hype may be declared.
void * XMLNodePointer_t
Definition: TXMLEngine.h:19
IsoMap fisomap
Map containing scale names and the TGeoScale for it.
Definition: TGDMLParse.h:213
virtual void AddNode(TGeoVolume *vol, Int_t copy_no, TGeoMatrix *mat=0, Option_t *option="")
Add a component to the assembly.
void SetUserExtension(TGeoExtension *ext)
Connect user-defined extension to the volume.
virtual void RotateZ(Double_t angle)
Rotate about Z axis of the master frame with angle expressed in degrees.
XMLNodePointer_t TopProcess(TXMLEngine *gdml, XMLNodePointer_t node)
In the setup section of the GDML file, the top volume need to be declared.
Bool_t HasAttr(XMLNodePointer_t xmlnode, const char *name)
checks if node has attribute of specified name
Definition: TXMLEngine.cxx:446
Double_t Pi()
Definition: TMath.h:44
XMLDocPointer_t ParseFile(const char *filename, Int_t maxbuf=100000)
Parses content of file and tries to produce xml structures.
void * XMLAttrPointer_t
Definition: TXMLEngine.h:21
TGeoShape * GetShape() const
Definition: TGeoVolume.h:204
static const double x1[5]
XMLNodePointer_t Polyhedra(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Polyhedra may be declared.
double f(double x)
virtual void ReflectX(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to YZ.
R__EXTERN TGeoManager * gGeoManager
Definition: TGeoManager.h:556
double Double_t
Definition: RtypesCore.h:55
const char * GetAttr(XMLNodePointer_t xmlnode, const char *name)
returns value of attribute for xmlnode
Definition: TXMLEngine.cxx:460
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
TMap implements an associative array of (key,value) pairs using a THashTable for efficient retrieval ...
Definition: TMap.h:44
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
MixMap fmixmap
Map containing medium names and the TGeoMedium for it.
Definition: TGDMLParse.h:217
Double_t y[n]
Definition: legend1.C:17
Int_t SetAxis(const char *axisString)
When using the 'divide' process in the geometry this function sets the variable 'axis' depending on w...
Definition: TGDMLParse.cxx:363
Double_t Na()
Definition: TMath.h:104
#define name(a, b)
Definition: linkTestLib0.cpp:5
const char * NameShort(const char *name)
this function looks thru a string for the chars '0x' next to each other, when it finds this...
Definition: TGDMLParse.cxx:388
void SetMatrix(const Double_t *rot)
Definition: TGeoMatrix.h:229
MedMap fmedmap
Map containing material names and the TGeoMaterial for it.
Definition: TGDMLParse.h:216
SclMap fsclmap
Map containing rotation names and the TGeoRotation for it.
Definition: TGDMLParse.h:212
SolMap fsolmap
Map containing mixture names and the TGeoMixture for it.
Definition: TGDMLParse.h:218
const char * ParseGDML(TXMLEngine *gdml, XMLNodePointer_t node)
this function recursively moves thru the DOM tree of the GDML file.
Definition: TGDMLParse.cxx:164
ReflVolMap freflvolmap
Map containing reflection names and the TGDMLRefl for it - containing refl matrix.
Definition: TGDMLParse.h:222
TObject * Last() const
Return the object in the last filled slot. Returns 0 if no entries.
Definition: TObjArray.cxx:478
VolMap fvolmap
Map containing solid names and the TGeoShape for it.
Definition: TGDMLParse.h:219
XMLNodePointer_t GetChild(XMLNodePointer_t xmlnode, Bool_t realnode=kTRUE)
returns first child of xml node
Definition: TXMLEngine.cxx:993
XMLNodePointer_t GetParent(XMLNodePointer_t xmlnode)
returns parent of xmlnode
TGeoMedium * GetMedium() const
Definition: TGeoVolume.h:189
void SetRotation(const Double_t *matrix)
Definition: TGeoMatrix.h:451
XMLNodePointer_t DocGetRootElement(XMLDocPointer_t xmldoc)
returns root node of document
XMLNodePointer_t ElTube(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a ElTube may be declared.
TLine * lv
Definition: textalign.C:5
virtual void DefineSection(Int_t snum, Double_t z, Double_t rmin, Double_t rmax)
Defines z position of a section plane, rmin and rmax at this z.
Definition: TGeoPcon.cxx:615
virtual Double_t GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const =0
ClassImp(TGDMLParse) TGeoVolume *TGDMLParse
creates the new instance of the XMLEngine called 'gdml', using the filename >> then parses the file a...
Definition: TGDMLParse.cxx:120
ReflSolidMap freflsolidmap
Map containing reflection names and the Solid name ir references to.
Definition: TGDMLParse.h:221
TGMatrixLayout * fMatrix
XMLNodePointer_t MatProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr, int z)
In the materials section of the GDML file, materials can be declared.
Definition: TGDMLParse.cxx:946
XMLNodePointer_t EleProcess(TXMLEngine *gdml, XMLNodePointer_t node, XMLNodePointer_t parentn, Bool_t hasIsotopes, Bool_t hasIsotopesExtended)
Definition: TGDMLParse.cxx:751
XMLNodePointer_t Tube(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, a Tube may be declared.
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
XMLNodePointer_t AssProcess(TXMLEngine *gdml, XMLNodePointer_t node)
In the structure section of the GDML file, assembly volumes can be declared.
const Bool_t kTRUE
Definition: Rtypes.h:91
XMLNodePointer_t VolProcess(TXMLEngine *gdml, XMLNodePointer_t node)
In the structure section of the GDML file, volumes can be declared.
float value
Definition: math.cpp:443
const char * GetAttrName(XMLAttrPointer_t xmlattr)
return name of the attribute
Definition: TXMLEngine.cxx:592
const Int_t n
Definition: legend1.C:16
TString fWorldName
Definition: TGDMLParse.h:105
gr SetName("gr")
XMLNodePointer_t ElCone(TXMLEngine *gdml, XMLNodePointer_t node, XMLAttrPointer_t attr)
In the solids section of the GDML file, an elliptical cone may be declared.
Bool_t DefinePolygon(Int_t nvert, const Double_t *xv, const Double_t *yv)
Creates the polygon representing the blueprint of any Xtru section.
Definition: TGeoXtru.cxx:720
PosMap fposmap
Definition: TGDMLParse.h:210
static const double x3[11]
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904