Geant4  10.00.p02
G4GDMLReadSolids.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4GDMLReadSolids.cc 81843 2014-06-06 09:11:11Z gcosmo $
27 //
28 // class G4GDMLReadSolids Implementation
29 //
30 // Original author: Zoltan Torzsok, November 2007
31 //
32 // --------------------------------------------------------------------
33 
34 #include "G4GDMLReadSolids.hh"
35 #include "G4Box.hh"
36 #include "G4Cons.hh"
37 #include "G4Ellipsoid.hh"
38 #include "G4EllipticalCone.hh"
39 #include "G4EllipticalTube.hh"
40 #include "G4Hype.hh"
41 #include "G4IntersectionSolid.hh"
42 #include "G4Orb.hh"
43 #include "G4Para.hh"
44 #include "G4Paraboloid.hh"
45 #include "G4Polycone.hh"
46 #include "G4GenericPolycone.hh"
47 #include "G4Polyhedra.hh"
48 #include "G4QuadrangularFacet.hh"
49 #include "G4ReflectedSolid.hh"
50 #include "G4Sphere.hh"
51 #include "G4SolidStore.hh"
52 #include "G4SubtractionSolid.hh"
53 #include "G4GenericTrap.hh"
54 #include "G4TessellatedSolid.hh"
55 #include "G4Tet.hh"
56 #include "G4Torus.hh"
57 #include "G4Transform3D.hh"
58 #include "G4Trap.hh"
59 #include "G4Trd.hh"
60 #include "G4TriangularFacet.hh"
61 #include "G4Tubs.hh"
62 #include "G4CutTubs.hh"
63 #include "G4TwistedBox.hh"
64 #include "G4TwistedTrap.hh"
65 #include "G4TwistedTrd.hh"
66 #include "G4TwistedTubs.hh"
67 #include "G4UnionSolid.hh"
68 #include "G4OpticalSurface.hh"
69 #include "G4SurfaceProperty.hh"
70 
72 {
73 }
74 
76 {
77 }
78 
80 BooleanRead(const xercesc::DOMElement* const booleanElement, const BooleanOp op)
81 {
82  G4String name;
83  G4String first;
84  G4String scnd;
85  G4ThreeVector position(0.0,0.0,0.0);
86  G4ThreeVector rotation(0.0,0.0,0.0);
87  G4ThreeVector firstposition(0.0,0.0,0.0);
88  G4ThreeVector firstrotation(0.0,0.0,0.0);
89 
90  const xercesc::DOMNamedNodeMap* const attributes
91  = booleanElement->getAttributes();
92  XMLSize_t attributeCount = attributes->getLength();
93 
94  for (XMLSize_t attribute_index=0;
95  attribute_index<attributeCount; attribute_index++)
96  {
97  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
98 
99  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
100  { continue; }
101 
102  const xercesc::DOMAttr* const attribute
103  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
104  if (!attribute)
105  {
106  G4Exception("G4GDMLReadSolids::BooleanRead()",
107  "InvalidRead", FatalException, "No attribute found!");
108  return;
109  }
110  const G4String attName = Transcode(attribute->getName());
111  const G4String attValue = Transcode(attribute->getValue());
112 
113  if (attName=="name") { name = GenerateName(attValue); }
114  }
115 
116  for (xercesc::DOMNode* iter = booleanElement->getFirstChild();
117  iter != 0;iter = iter->getNextSibling())
118  {
119  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
120 
121  const xercesc::DOMElement* const child
122  = dynamic_cast<xercesc::DOMElement*>(iter);
123  if (!child)
124  {
125  G4Exception("G4GDMLReadSolids::BooleanRead()",
126  "InvalidRead", FatalException, "No child found!");
127  return;
128  }
129  const G4String tag = Transcode(child->getTagName());
130 
131  if (tag=="first") { first = RefRead(child); } else
132  if (tag=="second") { scnd = RefRead(child); } else
133  if (tag=="position") { VectorRead(child,position); } else
134  if (tag=="rotation") { VectorRead(child,rotation); } else
135  if (tag=="positionref")
136  { position = GetPosition(GenerateName(RefRead(child))); } else
137  if (tag=="rotationref")
138  { rotation = GetRotation(GenerateName(RefRead(child))); } else
139  if (tag=="firstposition") { VectorRead(child,firstposition); } else
140  if (tag=="firstrotation") { VectorRead(child,firstrotation); } else
141  if (tag=="firstpositionref")
142  { firstposition = GetPosition(GenerateName(RefRead(child))); } else
143  if (tag=="firstrotationref")
144  { firstrotation = GetRotation(GenerateName(RefRead(child))); }
145  else
146  {
147  G4String error_msg = "Unknown tag in boolean solid: " + tag;
148  G4Exception("G4GDMLReadSolids::BooleanRead()", "ReadError",
149  FatalException, error_msg);
150  }
151  }
152 
153  G4VSolid* firstSolid = GetSolid(GenerateName(first));
154  G4VSolid* secondSolid = GetSolid(GenerateName(scnd));
155 
156  G4Transform3D transform(GetRotationMatrix(rotation),position);
157 
158  if (( (firstrotation.x()!=0.0) || (firstrotation.y()!=0.0)
159  || (firstrotation.z()!=0.0))
160  || ( (firstposition.x()!=0.0) || (firstposition.y()!=0.0)
161  || (firstposition.z()!=0.0)))
162  {
163  G4Transform3D firsttransform(GetRotationMatrix(firstrotation),
164  firstposition);
165  firstSolid = new G4DisplacedSolid(GenerateName("displaced_"+first),
166  firstSolid, firsttransform);
167  }
168 
169  if (op==UNION)
170  { new G4UnionSolid(name,firstSolid,secondSolid,transform); } else
171  if (op==SUBTRACTION)
172  { new G4SubtractionSolid(name,firstSolid,secondSolid,transform); } else
173  if (op==INTERSECTION)
174  { new G4IntersectionSolid(name,firstSolid,secondSolid,transform); }
175 }
176 
177 void G4GDMLReadSolids::BoxRead(const xercesc::DOMElement* const boxElement)
178 {
179  G4String name;
180  G4double lunit = 1.0;
181  G4double x = 0.0;
182  G4double y = 0.0;
183  G4double z = 0.0;
184 
185  const xercesc::DOMNamedNodeMap* const attributes
186  = boxElement->getAttributes();
187  XMLSize_t attributeCount = attributes->getLength();
188 
189  for (XMLSize_t attribute_index=0;
190  attribute_index<attributeCount; attribute_index++)
191  {
192  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
193 
194  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
195  { continue; }
196 
197  const xercesc::DOMAttr* const attribute
198  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
199  if (!attribute)
200  {
201  G4Exception("G4GDMLReadSolids::BoxRead()",
202  "InvalidRead", FatalException, "No attribute found!");
203  return;
204  }
205  const G4String attName = Transcode(attribute->getName());
206  const G4String attValue = Transcode(attribute->getValue());
207 
208  if (attName=="name") { name = GenerateName(attValue); } else
209  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
210  if (attName=="x") { x = eval.Evaluate(attValue); } else
211  if (attName=="y") { y = eval.Evaluate(attValue); } else
212  if (attName=="z") { z = eval.Evaluate(attValue); }
213  }
214 
215  x *= 0.5*lunit;
216  y *= 0.5*lunit;
217  z *= 0.5*lunit;
218 
219  new G4Box(name,x,y,z);
220 }
221 
222 void G4GDMLReadSolids::ConeRead(const xercesc::DOMElement* const coneElement)
223 {
224  G4String name;
225  G4double lunit = 1.0;
226  G4double aunit = 1.0;
227  G4double rmin1 = 0.0;
228  G4double rmax1 = 0.0;
229  G4double rmin2 = 0.0;
230  G4double rmax2 = 0.0;
231  G4double z = 0.0;
232  G4double startphi = 0.0;
233  G4double deltaphi = 0.0;
234 
235  const xercesc::DOMNamedNodeMap* const attributes
236  = coneElement->getAttributes();
237  XMLSize_t attributeCount = attributes->getLength();
238 
239  for (XMLSize_t attribute_index=0;
240  attribute_index<attributeCount; attribute_index++)
241  {
242  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
243 
244  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
245  { continue; }
246 
247  const xercesc::DOMAttr* const attribute
248  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
249  if (!attribute)
250  {
251  G4Exception("G4GDMLReadSolids::ConeRead()",
252  "InvalidRead", FatalException, "No attribute found!");
253  return;
254  }
255  const G4String attName = Transcode(attribute->getName());
256  const G4String attValue = Transcode(attribute->getValue());
257 
258  if (attName=="name") { name = GenerateName(attValue); } else
259  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
260  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
261  if (attName=="rmin1") { rmin1 = eval.Evaluate(attValue); } else
262  if (attName=="rmax1") { rmax1 = eval.Evaluate(attValue); } else
263  if (attName=="rmin2") { rmin2 = eval.Evaluate(attValue); } else
264  if (attName=="rmax2") { rmax2 = eval.Evaluate(attValue); } else
265  if (attName=="z") { z = eval.Evaluate(attValue); } else
266  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
267  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
268  }
269 
270  rmin1 *= lunit;
271  rmax1 *= lunit;
272  rmin2 *= lunit;
273  rmax2 *= lunit;
274  z *= 0.5*lunit;
275  startphi *= aunit;
276  deltaphi *= aunit;
277 
278  new G4Cons(name,rmin1,rmax1,rmin2,rmax2,z,startphi,deltaphi);
279 }
280 
282 ElconeRead(const xercesc::DOMElement* const elconeElement)
283 {
284  G4String name;
285  G4double lunit = 1.0;
286  G4double dx = 0.0;
287  G4double dy = 0.0;
288  G4double zmax = 0.0;
289  G4double zcut = 0.0;
290 
291  const xercesc::DOMNamedNodeMap* const attributes
292  = elconeElement->getAttributes();
293  XMLSize_t attributeCount = attributes->getLength();
294 
295  for (XMLSize_t attribute_index=0;
296  attribute_index<attributeCount; attribute_index++)
297  {
298  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
299 
300  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
301  { continue; }
302 
303  const xercesc::DOMAttr* const attribute
304  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
305  if (!attribute)
306  {
307  G4Exception("G4GDMLReadSolids::ElconeRead()",
308  "InvalidRead", FatalException, "No attribute found!");
309  return;
310  }
311  const G4String attName = Transcode(attribute->getName());
312  const G4String attValue = Transcode(attribute->getValue());
313 
314  if (attName=="name") { name = GenerateName(attValue); } else
315  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
316  if (attName=="dx") { dx = eval.Evaluate(attValue); } else
317  if (attName=="dy") { dy = eval.Evaluate(attValue); } else
318  if (attName=="zmax") { zmax = eval.Evaluate(attValue); } else
319  if (attName=="zcut") { zcut = eval.Evaluate(attValue); }
320  }
321 
322  zmax *= lunit;
323  zcut *= lunit;
324 
325  new G4EllipticalCone(name,dx,dy,zmax,zcut);
326 }
327 
329 EllipsoidRead(const xercesc::DOMElement* const ellipsoidElement)
330 {
331  G4String name;
332  G4double lunit = 1.0;
333  G4double ax = 0.0;
334  G4double by = 0.0;
335  G4double cz = 0.0;
336  G4double zcut1 = 0.0;
337  G4double zcut2 = 0.0;
338 
339  const xercesc::DOMNamedNodeMap* const attributes
340  = ellipsoidElement->getAttributes();
341  XMLSize_t attributeCount = attributes->getLength();
342 
343  for (XMLSize_t attribute_index=0;
344  attribute_index<attributeCount; attribute_index++)
345  {
346  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
347 
348  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
349  { continue; }
350 
351  const xercesc::DOMAttr* const attribute
352  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
353  if (!attribute)
354  {
355  G4Exception("G4GDMLReadSolids::EllipsoidRead()",
356  "InvalidRead", FatalException, "No attribute found!");
357  return;
358  }
359  const G4String attName = Transcode(attribute->getName());
360  const G4String attValue = Transcode(attribute->getValue());
361 
362  if (attName=="name") { name = GenerateName(attValue); } else
363  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
364  if (attName=="ax") { ax = eval.Evaluate(attValue); } else
365  if (attName=="by") { by = eval.Evaluate(attValue); } else
366  if (attName=="cz") { cz = eval.Evaluate(attValue); } else
367  if (attName=="zcut1") { zcut1 = eval.Evaluate(attValue); } else
368  if (attName=="zcut2") { zcut2 = eval.Evaluate(attValue); }
369  }
370 
371  ax *= lunit;
372  by *= lunit;
373  cz *= lunit;
374  zcut1 *= lunit;
375  zcut2 *= lunit;
376 
377  new G4Ellipsoid(name,ax,by,cz,zcut1,zcut2);
378 }
379 
381 EltubeRead(const xercesc::DOMElement* const eltubeElement)
382 {
383  G4String name;
384  G4double lunit = 1.0;
385  G4double dx = 0.0;
386  G4double dy = 0.0;
387  G4double dz = 0.0;
388 
389  const xercesc::DOMNamedNodeMap* const attributes
390  = eltubeElement->getAttributes();
391  XMLSize_t attributeCount = attributes->getLength();
392 
393  for (XMLSize_t attribute_index=0;
394  attribute_index<attributeCount; attribute_index++)
395  {
396  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
397 
398  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
399  { continue; }
400 
401  const xercesc::DOMAttr* const attribute
402  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
403  if (!attribute)
404  {
405  G4Exception("G4GDMLReadSolids::EltubeRead()",
406  "InvalidRead", FatalException, "No attribute found!");
407  return;
408  }
409  const G4String attName = Transcode(attribute->getName());
410  const G4String attValue = Transcode(attribute->getValue());
411 
412  if (attName=="name") { name = GenerateName(attValue); } else
413  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
414  if (attName=="dx") { dx = eval.Evaluate(attValue); } else
415  if (attName=="dy") { dy = eval.Evaluate(attValue); } else
416  if (attName=="dz") { dz = eval.Evaluate(attValue); }
417  }
418 
419  dx *= lunit;
420  dy *= lunit;
421  dz *= lunit;
422 
423  new G4EllipticalTube(name,dx,dy,dz);
424 }
425 
426 void G4GDMLReadSolids::XtruRead(const xercesc::DOMElement* const xtruElement)
427 {
428  G4String name;
429  G4double lunit = 1.0;
430 
431  const xercesc::DOMNamedNodeMap* const attributes
432  = xtruElement->getAttributes();
433  XMLSize_t attributeCount = attributes->getLength();
434 
435  for (XMLSize_t attribute_index=0;
436  attribute_index<attributeCount; attribute_index++)
437  {
438  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
439 
440  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
441  { continue; }
442 
443  const xercesc::DOMAttr* const attribute
444  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
445  if (!attribute)
446  {
447  G4Exception("G4GDMLReadSolids::XtruRead()",
448  "InvalidRead", FatalException, "No attribute found!");
449  return;
450  }
451  const G4String attName = Transcode(attribute->getName());
452  const G4String attValue = Transcode(attribute->getValue());
453 
454  if (attName=="name") { name = GenerateName(attValue); } else
455  if (attName=="lunit") { lunit = eval.Evaluate(attValue); }
456  }
457 
458  std::vector<G4TwoVector> twoDimVertexList;
459  std::vector<G4ExtrudedSolid::ZSection> sectionList;
460 
461  for (xercesc::DOMNode* iter = xtruElement->getFirstChild();
462  iter != 0; iter = iter->getNextSibling())
463  {
464  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
465 
466  const xercesc::DOMElement* const child
467  = dynamic_cast<xercesc::DOMElement*>(iter);
468  if (!child)
469  {
470  G4Exception("G4GDMLReadSolids::XtruRead()",
471  "InvalidRead", FatalException, "No child found!");
472  return;
473  }
474  const G4String tag = Transcode(child->getTagName());
475 
476  if (tag=="twoDimVertex")
477  { twoDimVertexList.push_back(TwoDimVertexRead(child,lunit)); } else
478  if (tag=="section")
479  { sectionList.push_back(SectionRead(child,lunit)); }
480  }
481 
482  new G4ExtrudedSolid(name,twoDimVertexList,sectionList);
483 }
484 
485 void G4GDMLReadSolids::HypeRead(const xercesc::DOMElement* const hypeElement)
486 {
487  G4String name;
488  G4double lunit = 1.0;
489  G4double aunit = 1.0;
490  G4double rmin = 0.0;
491  G4double rmax = 0.0;
492  G4double inst = 0.0;
493  G4double outst = 0.0;
494  G4double z = 0.0;
495 
496  const xercesc::DOMNamedNodeMap* const attributes
497  = hypeElement->getAttributes();
498  XMLSize_t attributeCount = attributes->getLength();
499 
500  for (XMLSize_t attribute_index=0;
501  attribute_index<attributeCount; attribute_index++)
502  {
503  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
504 
505  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
506  { continue; }
507 
508  const xercesc::DOMAttr* const attribute
509  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
510  if (!attribute)
511  {
512  G4Exception("G4GDMLReadSolids::HypeRead()",
513  "InvalidRead", FatalException, "No attribute found!");
514  return;
515  }
516  const G4String attName = Transcode(attribute->getName());
517  const G4String attValue = Transcode(attribute->getValue());
518 
519  if (attName=="name") { name = GenerateName(attValue); } else
520  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
521  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
522  if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
523  if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
524  if (attName=="inst") { inst = eval.Evaluate(attValue); } else
525  if (attName=="outst") { outst = eval.Evaluate(attValue); } else
526  if (attName=="z") { z = eval.Evaluate(attValue); }
527  }
528 
529  rmin *= lunit;
530  rmax *= lunit;
531  inst *= aunit;
532  outst *= aunit;
533  z *= 0.5*lunit;
534 
535  new G4Hype(name,rmin,rmax,inst,outst,z);
536 }
537 
538 void G4GDMLReadSolids::MultiUnionRead(const xercesc::DOMElement* const)
539 {
540  G4Exception("G4GDMLReadSolids::MultiUnionRead()",
541  "InvalidSetup", JustWarning,
542  "MultiUnion is not supported yet. Sorry!");
543 }
544 
545 void G4GDMLReadSolids::OrbRead(const xercesc::DOMElement* const orbElement)
546 {
547  G4String name;
548  G4double lunit = 1.0;
549  G4double r = 0.0;
550 
551  const xercesc::DOMNamedNodeMap* const attributes
552  = orbElement->getAttributes();
553  XMLSize_t attributeCount = attributes->getLength();
554 
555  for (XMLSize_t attribute_index=0;
556  attribute_index<attributeCount; attribute_index++)
557  {
558  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
559 
560  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
561  { continue; }
562 
563  const xercesc::DOMAttr* const attribute
564  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
565  if (!attribute)
566  {
567  G4Exception("G4GDMLReadSolids::OrbRead()",
568  "InvalidRead", FatalException, "No attribute found!");
569  return;
570  }
571  const G4String attName = Transcode(attribute->getName());
572  const G4String attValue = Transcode(attribute->getValue());
573 
574  if (attName=="name") { name = GenerateName(attValue); } else
575  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
576  if (attName=="r") { r = eval.Evaluate(attValue); }
577  }
578 
579  r *= lunit;
580 
581  new G4Orb(name,r);
582 }
583 
584 void G4GDMLReadSolids::ParaRead(const xercesc::DOMElement* const paraElement)
585 {
586  G4String name;
587  G4double lunit = 1.0;
588  G4double aunit = 1.0;
589  G4double x = 0.0;
590  G4double y = 0.0;
591  G4double z = 0.0;
592  G4double alpha = 0.0;
593  G4double theta = 0.0;
594  G4double phi = 0.0;
595 
596  const xercesc::DOMNamedNodeMap* const attributes
597  = paraElement->getAttributes();
598  XMLSize_t attributeCount = attributes->getLength();
599 
600  for (XMLSize_t attribute_index=0;
601  attribute_index<attributeCount; attribute_index++)
602  {
603  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
604 
605  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
606  { continue; }
607 
608  const xercesc::DOMAttr* const attribute
609  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
610  if (!attribute)
611  {
612  G4Exception("G4GDMLReadSolids::ParaRead()",
613  "InvalidRead", FatalException, "No attribute found!");
614  return;
615  }
616  const G4String attName = Transcode(attribute->getName());
617  const G4String attValue = Transcode(attribute->getValue());
618 
619  if (attName=="name") { name = GenerateName(attValue); } else
620  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
621  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
622  if (attName=="x") { x = eval.Evaluate(attValue); } else
623  if (attName=="y") { y = eval.Evaluate(attValue); } else
624  if (attName=="z") { z = eval.Evaluate(attValue); } else
625  if (attName=="alpha") { alpha = eval.Evaluate(attValue); } else
626  if (attName=="theta") { theta = eval.Evaluate(attValue); } else
627  if (attName=="phi") { phi = eval.Evaluate(attValue); }
628  }
629 
630  x *= 0.5*lunit;
631  y *= 0.5*lunit;
632  z *= 0.5*lunit;
633  alpha *= aunit;
634  theta *= aunit;
635  phi *= aunit;
636 
637  new G4Para(name,x,y,z,alpha,theta,phi);
638 }
639 
641 ParaboloidRead(const xercesc::DOMElement* const paraElement)
642 {
643  G4String name;
644  G4double lunit = 1.0;
645  G4double rlo = 0.0;
646  G4double rhi = 0.0;
647  G4double dz = 0.0;
648 
649  const xercesc::DOMNamedNodeMap* const attributes
650  = paraElement->getAttributes();
651  XMLSize_t attributeCount = attributes->getLength();
652 
653  for (XMLSize_t attribute_index=0;
654  attribute_index<attributeCount; attribute_index++)
655  {
656  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
657 
658  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
659  { continue; }
660 
661  const xercesc::DOMAttr* const attribute
662  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
663  if (!attribute)
664  {
665  G4Exception("G4GDMLReadSolids::ParaboloidRead()",
666  "InvalidRead", FatalException, "No attribute found!");
667  return;
668  }
669  const G4String attName = Transcode(attribute->getName());
670  const G4String attValue = Transcode(attribute->getValue());
671 
672  if (attName=="name") { name = GenerateName(attValue); } else
673  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
674  if (attName=="rlo") { rlo = eval.Evaluate(attValue); } else
675  if (attName=="rhi") { rhi = eval.Evaluate(attValue); } else
676  if (attName=="dz") { dz = eval.Evaluate(attValue); }
677  }
678 
679  rlo *= 1.*lunit;
680  rhi *= 1.*lunit;
681  dz *= 1.*lunit;
682 
683  new G4Paraboloid(name,dz,rlo,rhi);
684 }
685 
687 PolyconeRead(const xercesc::DOMElement* const polyconeElement)
688 {
689  G4String name;
690  G4double lunit = 1.0;
691  G4double aunit = 1.0;
692  G4double startphi = 0.0;
693  G4double deltaphi = 0.0;
694 
695  const xercesc::DOMNamedNodeMap* const attributes
696  = polyconeElement->getAttributes();
697  XMLSize_t attributeCount = attributes->getLength();
698 
699  for (XMLSize_t attribute_index=0;
700  attribute_index<attributeCount; attribute_index++)
701  {
702  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
703 
704  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
705  { continue; }
706 
707  const xercesc::DOMAttr* const attribute
708  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
709  if (!attribute)
710  {
711  G4Exception("G4GDMLReadSolids::PolyconeRead()",
712  "InvalidRead", FatalException, "No attribute found!");
713  return;
714  }
715  const G4String attName = Transcode(attribute->getName());
716  const G4String attValue = Transcode(attribute->getValue());
717 
718  if (attName=="name") { name = GenerateName(attValue); } else
719  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
720  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
721  if (attName=="startphi") { startphi = eval.Evaluate(attValue); }else
722  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
723  }
724 
725  startphi *= aunit;
726  deltaphi *= aunit;
727 
728  std::vector<zplaneType> zplaneList;
729 
730  for (xercesc::DOMNode* iter = polyconeElement->getFirstChild();
731  iter != 0; iter = iter->getNextSibling())
732  {
733  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
734 
735  const xercesc::DOMElement* const child
736  = dynamic_cast<xercesc::DOMElement*>(iter);
737  if (!child)
738  {
739  G4Exception("G4GDMLReadSolids::PolyconeRead()",
740  "InvalidRead", FatalException, "No child found!");
741  return;
742  }
743  const G4String tag = Transcode(child->getTagName());
744 
745  if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
746  }
747 
748  G4int numZPlanes = zplaneList.size();
749 
750  G4double* rmin_array = new G4double[numZPlanes];
751  G4double* rmax_array = new G4double[numZPlanes];
752  G4double* z_array = new G4double[numZPlanes];
753 
754  for (G4int i=0; i<numZPlanes; i++)
755  {
756  rmin_array[i] = zplaneList[i].rmin*lunit;
757  rmax_array[i] = zplaneList[i].rmax*lunit;
758  z_array[i] = zplaneList[i].z*lunit;
759  }
760 
761  new G4Polycone(name,startphi,deltaphi,numZPlanes,
762  z_array,rmin_array,rmax_array);
763 
764  delete [] rmin_array;
765  delete [] rmax_array;
766  delete [] z_array;
767 }
768 
770 GenericPolyconeRead(const xercesc::DOMElement* const polyconeElement)
771 {
772  G4String name;
773  G4double lunit = 1.0;
774  G4double aunit = 1.0;
775  G4double startphi = 0.0;
776  G4double deltaphi = 0.0;
777 
778  const xercesc::DOMNamedNodeMap* const attributes
779  = polyconeElement->getAttributes();
780  XMLSize_t attributeCount = attributes->getLength();
781 
782  for (XMLSize_t attribute_index=0;
783  attribute_index<attributeCount; attribute_index++)
784  {
785  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
786 
787  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
788  { continue; }
789 
790  const xercesc::DOMAttr* const attribute
791  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
792  if (!attribute)
793  {
794  G4Exception("G4GDMLReadSolids::GenericPolyconeRead()",
795  "InvalidRead", FatalException, "No attribute found!");
796  return;
797  }
798  const G4String attName = Transcode(attribute->getName());
799  const G4String attValue = Transcode(attribute->getValue());
800 
801  if (attName=="name") { name = GenerateName(attValue); } else
802  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
803  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
804  if (attName=="startphi") { startphi = eval.Evaluate(attValue); }else
805  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
806  }
807 
808  startphi *= aunit;
809  deltaphi *= aunit;
810 
811  std::vector<rzPointType> rzPointList;
812 
813  for (xercesc::DOMNode* iter = polyconeElement->getFirstChild();
814  iter != 0; iter = iter->getNextSibling())
815  {
816  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
817 
818  const xercesc::DOMElement* const child
819  = dynamic_cast<xercesc::DOMElement*>(iter);
820  if (!child)
821  {
822  G4Exception("G4GDMLReadSolids::GenericPolyconeRead()",
823  "InvalidRead", FatalException, "No child found!");
824  return;
825  }
826  const G4String tag = Transcode(child->getTagName());
827 
828  if (tag=="rzpoint") { rzPointList.push_back(RZPointRead(child)); }
829  }
830 
831  G4int numRZPoints = rzPointList.size();
832 
833  G4double* r_array = new G4double[numRZPoints];
834  G4double* z_array = new G4double[numRZPoints];
835 
836  for (G4int i=0; i<numRZPoints; i++)
837  {
838  r_array[i] = rzPointList[i].r*lunit;
839  z_array[i] = rzPointList[i].z*lunit;
840  }
841  new G4GenericPolycone(name,startphi,deltaphi,numRZPoints,
842  r_array,z_array);
843  delete [] r_array;
844  delete [] z_array;
845 }
846 
848 PolyhedraRead(const xercesc::DOMElement* const polyhedraElement)
849 {
850  G4String name;
851  G4double lunit = 1.0;
852  G4double aunit = 1.0;
853  G4double startphi = 0.0;
854  G4double deltaphi = 0.0;
855  G4int numsides = 0;
856 
857  const xercesc::DOMNamedNodeMap* const attributes
858  = polyhedraElement->getAttributes();
859  XMLSize_t attributeCount = attributes->getLength();
860 
861  for (XMLSize_t attribute_index=0;
862  attribute_index<attributeCount; attribute_index++)
863  {
864  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
865 
866  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
867  { continue; }
868 
869  const xercesc::DOMAttr* const attribute
870  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
871  if (!attribute)
872  {
873  G4Exception("G4GDMLReadSolids::PolyhedraRead()",
874  "InvalidRead", FatalException, "No attribute found!");
875  return;
876  }
877  const G4String attName = Transcode(attribute->getName());
878  const G4String attValue = Transcode(attribute->getValue());
879 
880  if (attName=="name") { name = GenerateName(attValue); } else
881  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
882  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
883  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
884  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
885  if (attName=="numsides") { numsides = eval.EvaluateInteger(attValue); }
886  }
887 
888  startphi *= aunit;
889  deltaphi *= aunit;
890 
891  std::vector<zplaneType> zplaneList;
892 
893  for (xercesc::DOMNode* iter = polyhedraElement->getFirstChild();
894  iter != 0; iter = iter->getNextSibling())
895  {
896  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
897 
898  const xercesc::DOMElement* const child
899  = dynamic_cast<xercesc::DOMElement*>(iter);
900  if (!child)
901  {
902  G4Exception("G4GDMLReadSolids::PolyhedraRead()",
903  "InvalidRead", FatalException, "No child found!");
904  return;
905  }
906  const G4String tag = Transcode(child->getTagName());
907 
908  if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
909  }
910 
911  G4int numZPlanes = zplaneList.size();
912 
913  G4double* rmin_array = new G4double[numZPlanes];
914  G4double* rmax_array = new G4double[numZPlanes];
915  G4double* z_array = new G4double[numZPlanes];
916 
917  for (G4int i=0; i<numZPlanes; i++)
918  {
919  rmin_array[i] = zplaneList[i].rmin*lunit;
920  rmax_array[i] = zplaneList[i].rmax*lunit;
921  z_array[i] = zplaneList[i].z*lunit;
922  }
923 
924  new G4Polyhedra(name,startphi,deltaphi,numsides,numZPlanes,
925  z_array,rmin_array,rmax_array);
926 
927  delete [] rmin_array;
928  delete [] rmax_array;
929  delete [] z_array;
930 }
931 
933 GenericPolyhedraRead(const xercesc::DOMElement* const polyhedraElement)
934 {
935  G4String name;
936  G4double lunit = 1.0;
937  G4double aunit = 1.0;
938  G4double startphi = 0.0;
939  G4double deltaphi = 0.0;
940  G4int numsides = 0;
941 
942  const xercesc::DOMNamedNodeMap* const attributes
943  = polyhedraElement->getAttributes();
944  XMLSize_t attributeCount = attributes->getLength();
945 
946  for (XMLSize_t attribute_index=0;
947  attribute_index<attributeCount; attribute_index++)
948  {
949  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
950 
951  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
952  { continue; }
953 
954  const xercesc::DOMAttr* const attribute
955  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
956  if (!attribute)
957  {
958  G4Exception("G4GDMLReadSolids::GenericPolyhedraRead()",
959  "InvalidRead", FatalException, "No attribute found!");
960  return;
961  }
962  const G4String attName = Transcode(attribute->getName());
963  const G4String attValue = Transcode(attribute->getValue());
964 
965  if (attName=="name") { name = GenerateName(attValue); } else
966  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
967  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
968  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
969  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
970  if (attName=="numsides") { numsides = eval.EvaluateInteger(attValue); }
971  }
972 
973  startphi *= aunit;
974  deltaphi *= aunit;
975 
976  std::vector<rzPointType> rzpointList;
977 
978  for (xercesc::DOMNode* iter = polyhedraElement->getFirstChild();
979  iter != 0; iter = iter->getNextSibling())
980  {
981  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
982 
983  const xercesc::DOMElement* const child
984  = dynamic_cast<xercesc::DOMElement*>(iter);
985  if (!child)
986  {
987  G4Exception("G4GDMLReadSolids::GenericPolyhedraRead()",
988  "InvalidRead", FatalException, "No child found!");
989  return;
990  }
991  const G4String tag = Transcode(child->getTagName());
992 
993  if (tag=="rzpoint") { rzpointList.push_back(RZPointRead(child)); }
994  }
995 
996  G4int numRZPoints = rzpointList.size();
997 
998  G4double* r_array = new G4double[numRZPoints];
999  G4double* z_array = new G4double[numRZPoints];
1000 
1001  for (G4int i=0; i<numRZPoints; i++)
1002  {
1003  r_array[i] = rzpointList[i].r*lunit;
1004  z_array[i] = rzpointList[i].z*lunit;
1005  }
1006 
1007  new G4Polyhedra(name,startphi,deltaphi,numsides,numRZPoints,
1008  r_array,z_array);
1009 
1010  delete [] r_array;
1011  delete [] z_array;
1012 }
1013 
1015 QuadrangularRead(const xercesc::DOMElement* const quadrangularElement)
1016 {
1017  G4ThreeVector vertex1;
1018  G4ThreeVector vertex2;
1019  G4ThreeVector vertex3;
1020  G4ThreeVector vertex4;
1021  G4FacetVertexType type = ABSOLUTE;
1022  G4double lunit = 1.0;
1023 
1024  const xercesc::DOMNamedNodeMap* const attributes
1025  = quadrangularElement->getAttributes();
1026  XMLSize_t attributeCount = attributes->getLength();
1027 
1028  for (XMLSize_t attribute_index=0;
1029  attribute_index<attributeCount; attribute_index++)
1030  {
1031  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1032 
1033  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1034  { continue; }
1035 
1036  const xercesc::DOMAttr* const attribute
1037  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1038  if (!attribute)
1039  {
1040  G4Exception("G4GDMLReadSolids::QuadrangularRead()",
1041  "InvalidRead", FatalException, "No attribute found!");
1042  return 0;
1043  }
1044  const G4String attName = Transcode(attribute->getName());
1045  const G4String attValue = Transcode(attribute->getValue());
1046 
1047  if (attName=="vertex1")
1048  { vertex1 = GetPosition(GenerateName(attValue)); } else
1049  if (attName=="vertex2")
1050  { vertex2 = GetPosition(GenerateName(attValue)); } else
1051  if (attName=="vertex3")
1052  { vertex3 = GetPosition(GenerateName(attValue)); } else
1053  if (attName=="vertex4")
1054  { vertex4 = GetPosition(GenerateName(attValue)); } else
1055  if (attName=="lunit")
1056  { lunit = eval.Evaluate(attValue); } else
1057  if (attName=="type")
1058  { if (attValue=="RELATIVE") { type = RELATIVE; } }
1059  }
1060 
1061  return new G4QuadrangularFacet(vertex1*lunit,vertex2*lunit,
1062  vertex3*lunit,vertex4*lunit,type);
1063 }
1064 
1065 void G4GDMLReadSolids::
1066 ReflectedSolidRead(const xercesc::DOMElement* const reflectedSolidElement)
1067 {
1068  G4String name;
1069  G4double lunit = 1.0;
1070  G4double aunit = 1.0;
1071  G4String solid;
1072  G4ThreeVector scale(1.0,1.0,1.0);
1073  G4ThreeVector rotation;
1075 
1076  const xercesc::DOMNamedNodeMap* const attributes
1077  = reflectedSolidElement->getAttributes();
1078  XMLSize_t attributeCount = attributes->getLength();
1079 
1080  for (XMLSize_t attribute_index=0;
1081  attribute_index<attributeCount; attribute_index++)
1082  {
1083  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1084 
1085  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1086  { continue; }
1087 
1088  const xercesc::DOMAttr* const attribute
1089  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1090  if (!attribute)
1091  {
1092  G4Exception("G4GDMLReadSolids::ReflectedSolidRead()",
1093  "InvalidRead", FatalException, "No attribute found!");
1094  return;
1095  }
1096  const G4String attName = Transcode(attribute->getName());
1097  const G4String attValue = Transcode(attribute->getValue());
1098 
1099  if (attName=="name") { name = GenerateName(attValue); } else
1100  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1101  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1102  if (attName=="solid") { solid = GenerateName(attValue); } else
1103  if (attName=="sx") { scale.setX(eval.Evaluate(attValue)); } else
1104  if (attName=="sy") { scale.setY(eval.Evaluate(attValue)); } else
1105  if (attName=="sz") { scale.setZ(eval.Evaluate(attValue)); } else
1106  if (attName=="rx") { rotation.setX(eval.Evaluate(attValue)); } else
1107  if (attName=="ry") { rotation.setY(eval.Evaluate(attValue)); } else
1108  if (attName=="rz") { rotation.setZ(eval.Evaluate(attValue)); } else
1109  if (attName=="dx") { position.setX(eval.Evaluate(attValue)); } else
1110  if (attName=="dy") { position.setY(eval.Evaluate(attValue)); } else
1111  if (attName=="dz") { position.setZ(eval.Evaluate(attValue)); }
1112  }
1113 
1114  rotation *= aunit;
1115  position *= lunit;
1116 
1117  G4Transform3D transform(GetRotationMatrix(rotation),position);
1118  transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
1119 
1120  new G4ReflectedSolid(name,GetSolid(solid),transform);
1121 }
1122 
1124 SectionRead(const xercesc::DOMElement* const sectionElement,G4double lunit)
1125 {
1126  G4double zPosition = 0.0;
1127  G4TwoVector Offset;
1128  G4double scalingFactor = 1.0;
1129 
1130  const xercesc::DOMNamedNodeMap* const attributes
1131  = sectionElement->getAttributes();
1132  XMLSize_t attributeCount = attributes->getLength();
1133 
1134  for (XMLSize_t attribute_index=0;
1135  attribute_index<attributeCount; attribute_index++)
1136  {
1137  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1138 
1139  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1140  { continue; }
1141 
1142  const xercesc::DOMAttr* const attribute
1143  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1144  if (!attribute)
1145  {
1146  G4Exception("G4GDMLReadSolids::SectionRead()",
1147  "InvalidRead", FatalException, "No attribute found!");
1148  return G4ExtrudedSolid::ZSection(zPosition,Offset,scalingFactor);
1149  }
1150  const G4String attName = Transcode(attribute->getName());
1151  const G4String attValue = Transcode(attribute->getValue());
1152 
1153  if (attName=="zPosition")
1154  { zPosition = eval.Evaluate(attValue)*lunit; } else
1155  if (attName=="xOffset")
1156  { Offset.setX(eval.Evaluate(attValue)*lunit); } else
1157  if (attName=="yOffset")
1158  { Offset.setY(eval.Evaluate(attValue)*lunit); } else
1159  if (attName=="scalingFactor")
1160  { scalingFactor = eval.Evaluate(attValue); }
1161  }
1162 
1163  return G4ExtrudedSolid::ZSection(zPosition,Offset,scalingFactor);
1164 }
1165 
1166 void G4GDMLReadSolids::
1167 SphereRead(const xercesc::DOMElement* const sphereElement)
1168 {
1169  G4String name;
1170  G4double lunit = 1.0;
1171  G4double aunit = 1.0;
1172  G4double rmin = 0.0;
1173  G4double rmax = 0.0;
1174  G4double startphi = 0.0;
1175  G4double deltaphi = 0.0;
1176  G4double starttheta = 0.0;
1177  G4double deltatheta = 0.0;
1178 
1179  const xercesc::DOMNamedNodeMap* const attributes
1180  = sphereElement->getAttributes();
1181  XMLSize_t attributeCount = attributes->getLength();
1182 
1183  for (XMLSize_t attribute_index=0;
1184  attribute_index<attributeCount; attribute_index++)
1185  {
1186  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1187 
1188  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1189  { continue; }
1190 
1191  const xercesc::DOMAttr* const attribute
1192  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1193  if (!attribute)
1194  {
1195  G4Exception("G4GDMLReadSolids::SphereRead()",
1196  "InvalidRead", FatalException, "No attribute found!");
1197  return;
1198  }
1199  const G4String attName = Transcode(attribute->getName());
1200  const G4String attValue = Transcode(attribute->getValue());
1201 
1202  if (attName=="name") { name = GenerateName(attValue); } else
1203  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1204  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1205  if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1206  if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1207  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1208  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
1209  if (attName=="starttheta") { starttheta = eval.Evaluate(attValue); } else
1210  if (attName=="deltatheta") { deltatheta = eval.Evaluate(attValue); }
1211  }
1212 
1213  rmin *= lunit;
1214  rmax *= lunit;
1215  startphi *= aunit;
1216  deltaphi *= aunit;
1217  starttheta *= aunit;
1218  deltatheta *= aunit;
1219 
1220  new G4Sphere(name,rmin,rmax,startphi,deltaphi,starttheta,deltatheta);
1221 }
1222 
1223 void G4GDMLReadSolids::
1224 TessellatedRead(const xercesc::DOMElement* const tessellatedElement)
1225 {
1226  G4String name;
1227 
1228  const xercesc::DOMNamedNodeMap* const attributes
1229  = tessellatedElement->getAttributes();
1230  XMLSize_t attributeCount = attributes->getLength();
1231 
1232  for (XMLSize_t attribute_index=0;
1233  attribute_index<attributeCount; attribute_index++)
1234  {
1235  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1236 
1237  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1238  { continue; }
1239 
1240  const xercesc::DOMAttr* const attribute
1241  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1242  if (!attribute)
1243  {
1244  G4Exception("G4GDMLReadSolids::TessellatedRead()",
1245  "InvalidRead", FatalException, "No attribute found!");
1246  return;
1247  }
1248  const G4String attName = Transcode(attribute->getName());
1249  const G4String attValue = Transcode(attribute->getValue());
1250 
1251  if (attName=="name") { name = GenerateName(attValue); }
1252  }
1253 
1254  G4TessellatedSolid *tessellated = new G4TessellatedSolid(name);
1255 
1256  for (xercesc::DOMNode* iter = tessellatedElement->getFirstChild();
1257  iter != 0; iter = iter->getNextSibling())
1258  {
1259  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
1260 
1261  const xercesc::DOMElement* const child
1262  = dynamic_cast<xercesc::DOMElement*>(iter);
1263  if (!child)
1264  {
1265  G4Exception("G4GDMLReadSolids::TessellatedRead()",
1266  "InvalidRead", FatalException, "No child found!");
1267  return;
1268  }
1269  const G4String tag = Transcode(child->getTagName());
1270 
1271  if (tag=="triangular")
1272  { tessellated->AddFacet(TriangularRead(child)); } else
1273  if (tag=="quadrangular")
1274  { tessellated->AddFacet(QuadrangularRead(child)); }
1275  }
1276 
1277  tessellated->SetSolidClosed(true);
1278 }
1279 
1280 void G4GDMLReadSolids::TetRead(const xercesc::DOMElement* const tetElement)
1281 {
1282  G4String name;
1283  G4ThreeVector vertex1;
1284  G4ThreeVector vertex2;
1285  G4ThreeVector vertex3;
1286  G4ThreeVector vertex4;
1287  G4double lunit = 1.0;
1288 
1289  const xercesc::DOMNamedNodeMap* const attributes
1290  = tetElement->getAttributes();
1291  XMLSize_t attributeCount = attributes->getLength();
1292 
1293  for (XMLSize_t attribute_index=0;
1294  attribute_index<attributeCount;attribute_index++)
1295  {
1296  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1297 
1298  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1299  { continue; }
1300 
1301  const xercesc::DOMAttr* const attribute
1302  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1303  if (!attribute)
1304  {
1305  G4Exception("G4GDMLReadSolids::TetRead()",
1306  "InvalidRead", FatalException, "No attribute found!");
1307  return;
1308  }
1309  const G4String attName = Transcode(attribute->getName());
1310  const G4String attValue = Transcode(attribute->getValue());
1311 
1312  if (attName=="name")
1313  { name = GenerateName(attValue); } else
1314  if (attName=="lunit")
1315  { lunit = eval.Evaluate(attValue); } else
1316  if (attName=="vertex1")
1317  { vertex1 = GetPosition(GenerateName(attValue)); } else
1318  if (attName=="vertex2")
1319  { vertex2 = GetPosition(GenerateName(attValue)); } else
1320  if (attName=="vertex3")
1321  { vertex3 = GetPosition(GenerateName(attValue)); } else
1322  if (attName=="vertex4")
1323  { vertex4 = GetPosition(GenerateName(attValue)); }
1324  }
1325 
1326  new G4Tet(name,vertex1*lunit,vertex2*lunit,vertex3*lunit,vertex4*lunit);
1327 }
1328 
1329 void G4GDMLReadSolids::TorusRead(const xercesc::DOMElement* const torusElement)
1330 {
1331  G4String name;
1332  G4double lunit = 1.0;
1333  G4double aunit = 1.0;
1334  G4double rmin = 0.0;
1335  G4double rmax = 0.0;
1336  G4double rtor = 0.0;
1337  G4double startphi = 0.0;
1338  G4double deltaphi = 0.0;
1339 
1340  const xercesc::DOMNamedNodeMap* const attributes
1341  = torusElement->getAttributes();
1342  XMLSize_t attributeCount = attributes->getLength();
1343 
1344  for (XMLSize_t attribute_index=0;
1345  attribute_index<attributeCount; attribute_index++)
1346  {
1347  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1348 
1349  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1350  { continue; }
1351 
1352  const xercesc::DOMAttr* const attribute
1353  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1354  if (!attribute)
1355  {
1356  G4Exception("G4GDMLReadSolids::TorusRead()",
1357  "InvalidRead", FatalException, "No attribute found!");
1358  return;
1359  }
1360  const G4String attName = Transcode(attribute->getName());
1361  const G4String attValue = Transcode(attribute->getValue());
1362 
1363  if (attName=="name") { name = GenerateName(attValue); } else
1364  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1365  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1366  if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1367  if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1368  if (attName=="rtor") { rtor = eval.Evaluate(attValue); } else
1369  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1370  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1371  }
1372 
1373  rmin *= lunit;
1374  rmax *= lunit;
1375  rtor *= lunit;
1376  startphi *= aunit;
1377  deltaphi *= aunit;
1378 
1379  new G4Torus(name,rmin,rmax,rtor,startphi,deltaphi);
1380 }
1381 
1382 void G4GDMLReadSolids::
1383 GenTrapRead(const xercesc::DOMElement* const gtrapElement)
1384 {
1385  G4String name;
1386  G4double lunit = 1.0;
1387  G4double dz =0.0;
1388  G4double v1x=0.0, v1y=0.0, v2x=0.0, v2y=0.0, v3x=0.0, v3y=0.0,
1389  v4x=0.0, v4y=0.0, v5x=0.0, v5y=0.0, v6x=0.0, v6y=0.0,
1390  v7x=0.0, v7y=0.0, v8x=0.0, v8y=0.0;
1391 
1392  const xercesc::DOMNamedNodeMap* const attributes
1393  = gtrapElement->getAttributes();
1394  XMLSize_t attributeCount = attributes->getLength();
1395 
1396  for (XMLSize_t attribute_index=0;
1397  attribute_index<attributeCount; attribute_index++)
1398  {
1399  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1400 
1401  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1402  { continue; }
1403 
1404  const xercesc::DOMAttr* const attribute
1405  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1406  if (!attribute)
1407  {
1408  G4Exception("G4GDMLReadSolids::GenTrapRead()",
1409  "InvalidRead", FatalException, "No attribute found!");
1410  return;
1411  }
1412  const G4String attName = Transcode(attribute->getName());
1413  const G4String attValue = Transcode(attribute->getValue());
1414 
1415  if (attName=="name") { name = GenerateName(attValue); } else
1416  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1417  if (attName=="dz") { dz = eval.Evaluate(attValue); } else
1418  if (attName=="v1x") { v1x = eval.Evaluate(attValue); } else
1419  if (attName=="v1y") { v1y = eval.Evaluate(attValue); } else
1420  if (attName=="v2x") { v2x = eval.Evaluate(attValue); } else
1421  if (attName=="v2y") { v2y = eval.Evaluate(attValue); } else
1422  if (attName=="v3x") { v3x = eval.Evaluate(attValue); } else
1423  if (attName=="v3y") { v3y = eval.Evaluate(attValue); } else
1424  if (attName=="v4x") { v4x = eval.Evaluate(attValue); } else
1425  if (attName=="v4y") { v4y = eval.Evaluate(attValue); } else
1426  if (attName=="v5x") { v5x = eval.Evaluate(attValue); } else
1427  if (attName=="v5y") { v5y = eval.Evaluate(attValue); } else
1428  if (attName=="v6x") { v6x = eval.Evaluate(attValue); } else
1429  if (attName=="v6y") { v6y = eval.Evaluate(attValue); } else
1430  if (attName=="v7x") { v7x = eval.Evaluate(attValue); } else
1431  if (attName=="v7y") { v7y = eval.Evaluate(attValue); } else
1432  if (attName=="v8x") { v8x = eval.Evaluate(attValue); } else
1433  if (attName=="v8y") { v8y = eval.Evaluate(attValue); }
1434  }
1435 
1436  dz *= lunit;
1437  std::vector<G4TwoVector> vertices;
1438  vertices.push_back(G4TwoVector(v1x*lunit,v1y*lunit));
1439  vertices.push_back(G4TwoVector(v2x*lunit,v2y*lunit));
1440  vertices.push_back(G4TwoVector(v3x*lunit,v3y*lunit));
1441  vertices.push_back(G4TwoVector(v4x*lunit,v4y*lunit));
1442  vertices.push_back(G4TwoVector(v5x*lunit,v5y*lunit));
1443  vertices.push_back(G4TwoVector(v6x*lunit,v6y*lunit));
1444  vertices.push_back(G4TwoVector(v7x*lunit,v7y*lunit));
1445  vertices.push_back(G4TwoVector(v8x*lunit,v8y*lunit));
1446  new G4GenericTrap(name,dz,vertices);
1447 }
1448 
1449 void G4GDMLReadSolids::TrapRead(const xercesc::DOMElement* const trapElement)
1450 {
1451  G4String name;
1452  G4double lunit = 1.0;
1453  G4double aunit = 1.0;
1454  G4double z = 0.0;
1455  G4double theta = 0.0;
1456  G4double phi = 0.0;
1457  G4double y1 = 0.0;
1458  G4double x1 = 0.0;
1459  G4double x2 = 0.0;
1460  G4double alpha1 = 0.0;
1461  G4double y2 = 0.0;
1462  G4double x3 = 0.0;
1463  G4double x4 = 0.0;
1464  G4double alpha2 = 0.0;
1465 
1466  const xercesc::DOMNamedNodeMap* const attributes
1467  = trapElement->getAttributes();
1468  XMLSize_t attributeCount = attributes->getLength();
1469 
1470  for (XMLSize_t attribute_index=0;
1471  attribute_index<attributeCount; attribute_index++)
1472  {
1473  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1474 
1475  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1476  { continue; }
1477 
1478  const xercesc::DOMAttr* const attribute
1479  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1480  if (!attribute)
1481  {
1482  G4Exception("G4GDMLReadSolids::TrapRead()",
1483  "InvalidRead", FatalException, "No attribute found!");
1484  return;
1485  }
1486  const G4String attName = Transcode(attribute->getName());
1487  const G4String attValue = Transcode(attribute->getValue());
1488 
1489  if (attName=="name") { name = GenerateName(attValue); } else
1490  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1491  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1492  if (attName=="z") { z = eval.Evaluate(attValue); } else
1493  if (attName=="theta") { theta = eval.Evaluate(attValue); } else
1494  if (attName=="phi") { phi = eval.Evaluate(attValue); } else
1495  if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1496  if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1497  if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1498  if (attName=="alpha1") { alpha1 = eval.Evaluate(attValue); } else
1499  if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1500  if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1501  if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1502  if (attName=="alpha2") { alpha2 = eval.Evaluate(attValue); }
1503  }
1504 
1505  z *= 0.5*lunit;
1506  theta *= aunit;
1507  phi *= aunit;
1508  y1 *= 0.5*lunit;
1509  x1 *= 0.5*lunit;
1510  x2 *= 0.5*lunit;
1511  alpha1 *= aunit;
1512  y2 *= 0.5*lunit;
1513  x3 *= 0.5*lunit;
1514  x4 *= 0.5*lunit;
1515  alpha2 *= aunit;
1516 
1517  new G4Trap(name,z,theta,phi,y1,x1,x2,alpha1,y2,x3,x4,alpha2);
1518 }
1519 
1520 void G4GDMLReadSolids::TrdRead(const xercesc::DOMElement* const trdElement)
1521 {
1522  G4String name;
1523  G4double lunit = 1.0;
1524  G4double x1 = 0.0;
1525  G4double x2 = 0.0;
1526  G4double y1 = 0.0;
1527  G4double y2 = 0.0;
1528  G4double z = 0.0;
1529 
1530  const xercesc::DOMNamedNodeMap* const attributes = trdElement->getAttributes();
1531  XMLSize_t attributeCount = attributes->getLength();
1532 
1533  for (XMLSize_t attribute_index=0;
1534  attribute_index<attributeCount; attribute_index++)
1535  {
1536  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1537 
1538  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1539  { continue; }
1540 
1541  const xercesc::DOMAttr* const attribute
1542  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1543  if (!attribute)
1544  {
1545  G4Exception("G4GDMLReadSolids::TrdRead()",
1546  "InvalidRead", FatalException, "No attribute found!");
1547  return;
1548  }
1549  const G4String attName = Transcode(attribute->getName());
1550  const G4String attValue = Transcode(attribute->getValue());
1551 
1552  if (attName=="name") { name = GenerateName(attValue); } else
1553  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1554  if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1555  if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1556  if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1557  if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1558  if (attName=="z") { z = eval.Evaluate(attValue); }
1559  }
1560 
1561  x1 *= 0.5*lunit;
1562  x2 *= 0.5*lunit;
1563  y1 *= 0.5*lunit;
1564  y2 *= 0.5*lunit;
1565  z *= 0.5*lunit;
1566 
1567  new G4Trd(name,x1,x2,y1,y2,z);
1568 }
1569 
1571 TriangularRead(const xercesc::DOMElement* const triangularElement)
1572 {
1573  G4ThreeVector vertex1;
1574  G4ThreeVector vertex2;
1575  G4ThreeVector vertex3;
1576  G4FacetVertexType type = ABSOLUTE;
1577  G4double lunit = 1.0;
1578 
1579  const xercesc::DOMNamedNodeMap* const attributes
1580  = triangularElement->getAttributes();
1581  XMLSize_t attributeCount = attributes->getLength();
1582 
1583  for (XMLSize_t attribute_index=0;
1584  attribute_index<attributeCount; attribute_index++)
1585  {
1586  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1587 
1588  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1589  { continue; }
1590 
1591  const xercesc::DOMAttr* const attribute
1592  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1593  if (!attribute)
1594  {
1595  G4Exception("G4GDMLReadSolids::TriangularRead()",
1596  "InvalidRead", FatalException, "No attribute found!");
1597  return 0;
1598  }
1599  const G4String attName = Transcode(attribute->getName());
1600  const G4String attValue = Transcode(attribute->getValue());
1601 
1602  if (attName=="vertex1")
1603  { vertex1 = GetPosition(GenerateName(attValue)); } else
1604  if (attName=="vertex2")
1605  { vertex2 = GetPosition(GenerateName(attValue)); } else
1606  if (attName=="vertex3")
1607  { vertex3 = GetPosition(GenerateName(attValue)); } else
1608  if (attName=="lunit")
1609  { lunit = eval.Evaluate(attValue); } else
1610  if (attName=="type")
1611  { if (attValue=="RELATIVE") { type = RELATIVE; } }
1612  }
1613 
1614  return new G4TriangularFacet(vertex1*lunit,vertex2*lunit,vertex3*lunit,type);
1615 }
1616 
1617 void G4GDMLReadSolids::TubeRead(const xercesc::DOMElement* const tubeElement)
1618 {
1619  G4String name;
1620  G4double lunit = 1.0;
1621  G4double aunit = 1.0;
1622  G4double rmin = 0.0;
1623  G4double rmax = 0.0;
1624  G4double z = 0.0;
1625  G4double startphi = 0.0;
1626  G4double deltaphi = 0.0;
1627 
1628  const xercesc::DOMNamedNodeMap* const attributes
1629  = tubeElement->getAttributes();
1630  XMLSize_t attributeCount = attributes->getLength();
1631 
1632  for (XMLSize_t attribute_index=0;
1633  attribute_index<attributeCount; attribute_index++)
1634  {
1635  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1636 
1637  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1638  { continue; }
1639 
1640  const xercesc::DOMAttr* const attribute
1641  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1642  if (!attribute)
1643  {
1644  G4Exception("G4GDMLReadSolids::TubeRead()",
1645  "InvalidRead", FatalException, "No attribute found!");
1646  return;
1647  }
1648  const G4String attName = Transcode(attribute->getName());
1649  const G4String attValue = Transcode(attribute->getValue());
1650 
1651  if (attName=="name") { name = GenerateName(attValue); } else
1652  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1653  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1654  if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1655  if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1656  if (attName=="z") { z = eval.Evaluate(attValue); } else
1657  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1658  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1659  }
1660 
1661  rmin *= lunit;
1662  rmax *= lunit;
1663  z *= 0.5*lunit;
1664  startphi *= aunit;
1665  deltaphi *= aunit;
1666 
1667  new G4Tubs(name,rmin,rmax,z,startphi,deltaphi);
1668 }
1669 
1670 void G4GDMLReadSolids::CutTubeRead(const xercesc::DOMElement* const cuttubeElement)
1671 {
1672  G4String name;
1673  G4double lunit = 1.0;
1674  G4double aunit = 1.0;
1675  G4double rmin = 0.0;
1676  G4double rmax = 0.0;
1677  G4double z = 0.0;
1678  G4double startphi = 0.0;
1679  G4double deltaphi = 0.0;
1680  G4ThreeVector lowNorm(0);
1681  G4ThreeVector highNorm(0);
1682 
1683  const xercesc::DOMNamedNodeMap* const attributes
1684  = cuttubeElement->getAttributes();
1685  XMLSize_t attributeCount = attributes->getLength();
1686 
1687  for (XMLSize_t attribute_index=0;
1688  attribute_index<attributeCount; attribute_index++)
1689  {
1690  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1691 
1692  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1693  { continue; }
1694 
1695  const xercesc::DOMAttr* const attribute
1696  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1697  if (!attribute)
1698  {
1699  G4Exception("G4GDMLReadSolids::CutTubeRead()",
1700  "InvalidRead", FatalException, "No attribute found!");
1701  return;
1702  }
1703  const G4String attName = Transcode(attribute->getName());
1704  const G4String attValue = Transcode(attribute->getValue());
1705 
1706  if (attName=="name") { name = GenerateName(attValue); } else
1707  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1708  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1709  if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1710  if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1711  if (attName=="z") { z = eval.Evaluate(attValue); } else
1712  if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1713  if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
1714  if (attName=="lowX") { lowNorm.setX (eval.Evaluate(attValue)); } else
1715  if (attName=="lowY") { lowNorm.setY (eval.Evaluate(attValue)); } else
1716  if (attName=="lowZ") { lowNorm.setZ (eval.Evaluate(attValue)); } else
1717  if (attName=="highX") { highNorm.setX (eval.Evaluate(attValue)); } else
1718  if (attName=="highY") { highNorm.setY (eval.Evaluate(attValue)); } else
1719  if (attName=="highZ") { highNorm.setZ (eval.Evaluate(attValue)); }
1720 
1721  }
1722 
1723  rmin *= lunit;
1724  rmax *= lunit;
1725  z *= 0.5*lunit;
1726  startphi *= aunit;
1727  deltaphi *= aunit;
1728 
1729  new G4CutTubs(name,rmin,rmax,z,startphi,deltaphi,lowNorm,highNorm);
1730 }
1731 
1732 void G4GDMLReadSolids::
1733 TwistedboxRead(const xercesc::DOMElement* const twistedboxElement)
1734 {
1735  G4String name;
1736  G4double lunit = 1.0;
1737  G4double aunit = 1.0;
1738  G4double PhiTwist = 0.0;
1739  G4double x = 0.0;
1740  G4double y = 0.0;
1741  G4double z = 0.0;
1742 
1743  const xercesc::DOMNamedNodeMap* const attributes
1744  = twistedboxElement->getAttributes();
1745  XMLSize_t attributeCount = attributes->getLength();
1746 
1747  for (XMLSize_t attribute_index=0;
1748  attribute_index<attributeCount; attribute_index++)
1749  {
1750  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1751 
1752  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1753  { continue; }
1754 
1755  const xercesc::DOMAttr* const attribute
1756  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1757  if (!attribute)
1758  {
1759  G4Exception("G4GDMLReadSolids::TwistedboxRead()",
1760  "InvalidRead", FatalException, "No attribute found!");
1761  return;
1762  }
1763  const G4String attName = Transcode(attribute->getName());
1764  const G4String attValue = Transcode(attribute->getValue());
1765 
1766  if (attName=="name") { name = GenerateName(attValue); } else
1767  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1768  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1769  if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1770  if (attName=="x") { x = eval.Evaluate(attValue); } else
1771  if (attName=="y") { y = eval.Evaluate(attValue); } else
1772  if (attName=="z") { z = eval.Evaluate(attValue); }
1773  }
1774 
1775  PhiTwist *= aunit;
1776  x *= 0.5*lunit;
1777  y *= 0.5*lunit;
1778  z *= 0.5*lunit;
1779 
1780  new G4TwistedBox(name,PhiTwist,x,y,z);
1781 }
1782 
1783 void G4GDMLReadSolids::
1784 TwistedtrapRead(const xercesc::DOMElement* const twistedtrapElement)
1785 {
1786  G4String name;
1787  G4double lunit = 1.0;
1788  G4double aunit = 1.0;
1789  G4double PhiTwist = 0.0;
1790  G4double z = 0.0;
1791  G4double Theta = 0.0;
1792  G4double Phi = 0.0;
1793  G4double y1 = 0.0;
1794  G4double x1 = 0.0;
1795  G4double x2 = 0.0;
1796  G4double y2 = 0.0;
1797  G4double x3 = 0.0;
1798  G4double x4 = 0.0;
1799  G4double Alph = 0.0;
1800 
1801  const xercesc::DOMNamedNodeMap* const attributes
1802  = twistedtrapElement->getAttributes();
1803  XMLSize_t attributeCount = attributes->getLength();
1804 
1805  for (XMLSize_t attribute_index=0;
1806  attribute_index<attributeCount; attribute_index++)
1807  {
1808  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1809 
1810  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1811  { continue; }
1812 
1813  const xercesc::DOMAttr* const attribute
1814  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1815  if (!attribute)
1816  {
1817  G4Exception("G4GDMLReadSolids::TwistedtrapRead()",
1818  "InvalidRead", FatalException, "No attribute found!");
1819  return;
1820  }
1821  const G4String attName = Transcode(attribute->getName());
1822  const G4String attValue = Transcode(attribute->getValue());
1823 
1824  if (attName=="name") { name = GenerateName(attValue); } else
1825  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1826  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1827  if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1828  if (attName=="z") { z = eval.Evaluate(attValue); } else
1829  if (attName=="Theta") { Theta = eval.Evaluate(attValue); } else
1830  if (attName=="Phi") { Phi = eval.Evaluate(attValue); } else
1831  if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1832  if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1833  if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1834  if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1835  if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1836  if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1837  if (attName=="Alph") { Alph = eval.Evaluate(attValue); }
1838  }
1839 
1840 
1841  PhiTwist *= aunit;
1842  z *= 0.5*lunit;
1843  Theta *= aunit;
1844  Phi *= aunit;
1845  Alph *= aunit;
1846  y1 *= 0.5*lunit;
1847  x1 *= 0.5*lunit;
1848  x2 *= 0.5*lunit;
1849  y2 *= 0.5*lunit;
1850  x3 *= 0.5*lunit;
1851  x4 *= 0.5*lunit;
1852 
1853  new G4TwistedTrap(name,PhiTwist,z,Theta,Phi,y1,x1,x2,y2,x3,x4,Alph);
1854 }
1855 
1856 void G4GDMLReadSolids::
1857 TwistedtrdRead(const xercesc::DOMElement* const twistedtrdElement)
1858 {
1859  G4String name;
1860  G4double lunit = 1.0;
1861  G4double aunit = 1.0;
1862  G4double x1 = 0.0;
1863  G4double x2 = 0.0;
1864  G4double y1 = 0.0;
1865  G4double y2 = 0.0;
1866  G4double z = 0.0;
1867  G4double PhiTwist = 0.0;
1868 
1869  const xercesc::DOMNamedNodeMap* const attributes
1870  = twistedtrdElement->getAttributes();
1871  XMLSize_t attributeCount = attributes->getLength();
1872 
1873  for (XMLSize_t attribute_index=0;
1874  attribute_index<attributeCount; attribute_index++)
1875  {
1876  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1877 
1878  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1879  { continue; }
1880 
1881  const xercesc::DOMAttr* const attribute
1882  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1883  if (!attribute)
1884  {
1885  G4Exception("G4GDMLReadSolids::TwistedtrdRead()",
1886  "InvalidRead", FatalException, "No attribute found!");
1887  return;
1888  }
1889  const G4String attName = Transcode(attribute->getName());
1890  const G4String attValue = Transcode(attribute->getValue());
1891 
1892  if (attName=="name") { name = GenerateName(attValue); } else
1893  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1894  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1895  if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1896  if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1897  if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1898  if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1899  if (attName=="z") { z = eval.Evaluate(attValue); } else
1900  if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); }
1901  }
1902 
1903  x1 *= 0.5*lunit;
1904  x2 *= 0.5*lunit;
1905  y1 *= 0.5*lunit;
1906  y2 *= 0.5*lunit;
1907  z *= 0.5*lunit;
1908  PhiTwist *= aunit;
1909 
1910  new G4TwistedTrd(name,x1,x2,y1,y2,z,PhiTwist);
1911 }
1912 
1913 void G4GDMLReadSolids::
1914 TwistedtubsRead(const xercesc::DOMElement* const twistedtubsElement)
1915 {
1916  G4String name;
1917  G4double lunit = 1.0;
1918  G4double aunit = 1.0;
1919  G4double twistedangle = 0.0;
1920  G4double endinnerrad = 0.0;
1921  G4double endouterrad = 0.0;
1922  G4double zlen = 0.0;
1923  G4double phi = 0.0;
1924 
1925  const xercesc::DOMNamedNodeMap* const attributes
1926  = twistedtubsElement->getAttributes();
1927  XMLSize_t attributeCount = attributes->getLength();
1928 
1929  for (XMLSize_t attribute_index=0;
1930  attribute_index<attributeCount; attribute_index++)
1931  {
1932  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1933 
1934  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1935  { continue; }
1936 
1937  const xercesc::DOMAttr* const attribute
1938  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1939  if (!attribute)
1940  {
1941  G4Exception("G4GDMLReadSolids::TwistedtubsRead()",
1942  "InvalidRead", FatalException, "No attribute found!");
1943  return;
1944  }
1945  const G4String attName = Transcode(attribute->getName());
1946  const G4String attValue = Transcode(attribute->getValue());
1947 
1948  if (attName=="name") { name = GenerateName(attValue); } else
1949  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1950  if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1951  if (attName=="twistedangle") { twistedangle=eval.Evaluate(attValue); } else
1952  if (attName=="endinnerrad") { endinnerrad=eval.Evaluate(attValue); } else
1953  if (attName=="endouterrad") { endouterrad=eval.Evaluate(attValue); } else
1954  if (attName=="zlen") { zlen = eval.Evaluate(attValue); } else
1955  if (attName=="phi") { phi = eval.Evaluate(attValue); }
1956  }
1957 
1958  twistedangle *= aunit;
1959  endinnerrad *= lunit;
1960  endouterrad *= lunit;
1961  zlen *= 0.5*lunit;
1962  phi *= aunit;
1963 
1964  new G4TwistedTubs(name,twistedangle,endinnerrad,endouterrad,zlen,phi);
1965 }
1966 
1968 TwoDimVertexRead(const xercesc::DOMElement* const element, G4double lunit)
1969 {
1970  G4TwoVector vec;
1971 
1972  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1973  XMLSize_t attributeCount = attributes->getLength();
1974 
1975  for (XMLSize_t attribute_index=0;
1976  attribute_index<attributeCount; attribute_index++)
1977  {
1978  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1979 
1980  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1981  { continue; }
1982 
1983  const xercesc::DOMAttr* const attribute
1984  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
1985  if (!attribute)
1986  {
1987  G4Exception("G4GDMLReadSolids::TwoDimVertexRead()",
1988  "InvalidRead", FatalException, "No attribute found!");
1989  return vec;
1990  }
1991  const G4String attName = Transcode(attribute->getName());
1992  const G4String attValue = Transcode(attribute->getValue());
1993 
1994  if (attName=="x") { vec.setX(eval.Evaluate(attValue)*lunit); } else
1995  if (attName=="y") { vec.setY(eval.Evaluate(attValue)*lunit); }
1996  }
1997 
1998  return vec;
1999 }
2000 
2002 ZplaneRead(const xercesc::DOMElement* const zplaneElement)
2003 {
2004  zplaneType zplane = {0.,0.,0.};
2005 
2006  const xercesc::DOMNamedNodeMap* const attributes
2007  = zplaneElement->getAttributes();
2008  XMLSize_t attributeCount = attributes->getLength();
2009 
2010  for (XMLSize_t attribute_index=0;
2011  attribute_index<attributeCount; attribute_index++)
2012  {
2013  xercesc::DOMNode* node = attributes->item(attribute_index);
2014 
2015  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
2016 
2017  const xercesc::DOMAttr* const attribute
2018  = dynamic_cast<xercesc::DOMAttr*>(node);
2019  if (!attribute)
2020  {
2021  G4Exception("G4GDMLReadSolids::ZplaneRead()",
2022  "InvalidRead", FatalException, "No attribute found!");
2023  return zplane;
2024  }
2025  const G4String attName = Transcode(attribute->getName());
2026  const G4String attValue = Transcode(attribute->getValue());
2027 
2028  if (attName=="rmin") { zplane.rmin = eval.Evaluate(attValue); } else
2029  if (attName=="rmax") { zplane.rmax = eval.Evaluate(attValue); } else
2030  if (attName=="z") { zplane.z = eval.Evaluate(attValue); }
2031  }
2032 
2033  return zplane;
2034 }
2036 RZPointRead(const xercesc::DOMElement* const zplaneElement)
2037 {
2038  rzPointType rzpoint = {0.,0.};
2039 
2040  const xercesc::DOMNamedNodeMap* const attributes
2041  = zplaneElement->getAttributes();
2042  XMLSize_t attributeCount = attributes->getLength();
2043 
2044  for (XMLSize_t attribute_index=0;
2045  attribute_index<attributeCount; attribute_index++)
2046  {
2047  xercesc::DOMNode* node = attributes->item(attribute_index);
2048 
2049  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
2050 
2051  const xercesc::DOMAttr* const attribute
2052  = dynamic_cast<xercesc::DOMAttr*>(node);
2053  if (!attribute)
2054  {
2055  G4Exception("G4GDMLReadSolids::RZPointRead()",
2056  "InvalidRead", FatalException, "No attribute found!");
2057  return rzpoint;
2058  }
2059  const G4String attName = Transcode(attribute->getName());
2060  const G4String attValue = Transcode(attribute->getValue());
2061 
2062  if (attName=="r") { rzpoint.r = eval.Evaluate(attValue); } else
2063  if (attName=="z") { rzpoint.z = eval.Evaluate(attValue); }
2064  }
2065 
2066  return rzpoint;
2067 
2068 }
2069 
2070 void G4GDMLReadSolids::
2071 OpticalSurfaceRead(const xercesc::DOMElement* const opticalsurfaceElement)
2072 {
2073  G4String name;
2074  G4String smodel;
2075  G4String sfinish;
2076  G4String stype;
2077  G4double value = 0.0;
2078 
2079  const xercesc::DOMNamedNodeMap* const attributes
2080  = opticalsurfaceElement->getAttributes();
2081  XMLSize_t attributeCount = attributes->getLength();
2082 
2083  for (XMLSize_t attribute_index=0;
2084  attribute_index<attributeCount; attribute_index++)
2085  {
2086  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
2087 
2088  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
2089  { continue; }
2090 
2091  const xercesc::DOMAttr* const attribute
2092  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
2093  if (!attribute)
2094  {
2095  G4Exception("G4GDMLReadSolids::OpticalSurfaceRead()",
2096  "InvalidRead", FatalException, "No attribute found!");
2097  return;
2098  }
2099  const G4String attName = Transcode(attribute->getName());
2100  const G4String attValue = Transcode(attribute->getValue());
2101 
2102  if (attName=="name") { name = GenerateName(attValue); } else
2103  if (attName=="model") { smodel = attValue; } else
2104  if (attName=="finish") { sfinish = attValue; } else
2105  if (attName=="type") { stype = attValue; } else
2106  if (attName=="value") { value = eval.Evaluate(attValue); }
2107  }
2108 
2109  G4OpticalSurfaceModel model;
2110  G4OpticalSurfaceFinish finish;
2111  G4SurfaceType type;
2112 
2113  if ((smodel=="glisur") || (smodel=="0")) { model = glisur; } else
2114  if ((smodel=="unified") || (smodel=="1")) { model = unified; }
2115  else { model = LUT; }
2116 
2117  if ((sfinish=="polished") || (sfinish=="0"))
2118  { finish = polished; } else
2119  if ((sfinish=="polishedfrontpainted") || (sfinish=="1"))
2120  { finish = polishedfrontpainted; } else
2121  if ((sfinish=="polishedbackpainted") || (sfinish=="2"))
2122  { finish = polishedbackpainted; } else
2123  if ((sfinish=="ground") || (sfinish=="3"))
2124  { finish = ground; } else
2125  if ((sfinish=="groundfrontpainted") || (sfinish=="4"))
2126  { finish = groundfrontpainted; } else
2127  if ((sfinish=="groundbackpainted") || (sfinish=="5"))
2128  { finish = groundbackpainted; } else
2129  if ((sfinish=="polishedlumirrorair") || (sfinish=="6"))
2130  { finish = polishedlumirrorair; } else
2131  if ((sfinish=="polishedlumirrorglue") || (sfinish=="7"))
2132  { finish = polishedlumirrorglue; } else
2133  if ((sfinish=="polishedair") || (sfinish=="8"))
2134  { finish = polishedair; } else
2135  if ((sfinish=="polishedteflonair") || (sfinish=="9"))
2136  { finish = polishedteflonair; } else
2137  if ((sfinish=="polishedtioair") || (sfinish=="10"))
2138  { finish = polishedtioair; } else
2139  if ((sfinish=="polishedtyvekair") || (sfinish=="11"))
2140  { finish = polishedtyvekair; } else
2141  if ((sfinish=="polishedvm2000air") || (sfinish=="12"))
2142  { finish = polishedvm2000air; } else
2143  if ((sfinish=="polishedvm2000glue") || (sfinish=="13"))
2144  { finish = polishedvm2000glue; } else
2145  if ((sfinish=="etchedlumirrorair") || (sfinish=="14"))
2146  { finish = etchedlumirrorair; } else
2147  if ((sfinish=="etchedlumirrorglue") || (sfinish=="15"))
2148  { finish = etchedlumirrorglue; } else
2149  if ((sfinish=="etchedair") || (sfinish=="16"))
2150  { finish = etchedair; } else
2151  if ((sfinish=="etchedteflonair") || (sfinish=="17"))
2152  { finish = etchedteflonair; } else
2153  if ((sfinish=="etchedtioair") || (sfinish=="18"))
2154  { finish = etchedtioair; } else
2155  if ((sfinish=="etchedtyvekair") || (sfinish=="19"))
2156  { finish = etchedtyvekair; } else
2157  if ((sfinish=="etchedvm2000air") || (sfinish=="20"))
2158  { finish = etchedvm2000air; } else
2159  if ((sfinish=="etchedvm2000glue") || (sfinish=="21"))
2160  { finish = etchedvm2000glue; } else
2161  if ((sfinish=="groundlumirrorair") || (sfinish=="22"))
2162  { finish = groundlumirrorair; } else
2163  if ((sfinish=="groundlumirrorglue") || (sfinish=="23"))
2164  { finish = groundlumirrorglue; } else
2165  if ((sfinish=="groundair") || (sfinish=="24"))
2166  { finish = groundair; } else
2167  if ((sfinish=="groundteflonair") || (sfinish=="25"))
2168  { finish = groundteflonair; } else
2169  if ((sfinish=="groundtioair") || (sfinish=="26"))
2170  { finish = groundtioair; } else
2171  if ((sfinish=="groundtyvekair") || (sfinish=="27"))
2172  { finish = groundtyvekair; } else
2173  if ((sfinish=="groundvm2000air") || (sfinish=="28"))
2174  { finish = groundvm2000air; }
2175  else { finish = groundvm2000glue; }
2176 
2177  if ((stype=="dielectric_metal") || (stype=="0"))
2178  { type = dielectric_metal; } else
2179  if ((stype=="dielectric_dielectric") || (stype=="1"))
2180  { type = dielectric_dielectric; } else
2181  if ((stype=="dielectric_LUT") || (stype=="2"))
2182  { type = dielectric_LUT; } else
2183  if ((stype=="firsov") || (stype=="3"))
2184  { type = firsov; }
2185  else { type = x_ray; }
2186 
2187  new G4OpticalSurface(name,model,finish,type,value);
2188 }
2189 
2190 void G4GDMLReadSolids::SolidsRead(const xercesc::DOMElement* const solidsElement)
2191 {
2192  G4cout << "G4GDML: Reading solids..." << G4endl;
2193 
2194  for (xercesc::DOMNode* iter = solidsElement->getFirstChild();
2195  iter != 0; iter = iter->getNextSibling())
2196  {
2197  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
2198 
2199  const xercesc::DOMElement* const child
2200  = dynamic_cast<xercesc::DOMElement*>(iter);
2201  if (!child)
2202  {
2203  G4Exception("G4GDMLReadSolids::SolidsRead()",
2204  "InvalidRead", FatalException, "No child found!");
2205  return;
2206  }
2207  const G4String tag = Transcode(child->getTagName());
2208  if (tag=="define") { DefineRead(child); } else
2209  if (tag=="box") { BoxRead(child); } else
2210  if (tag=="cone") { ConeRead(child); } else
2211  if (tag=="elcone") { ElconeRead(child); } else
2212  if (tag=="ellipsoid") { EllipsoidRead(child); }else
2213  if (tag=="eltube") { EltubeRead(child); } else
2214  if (tag=="xtru") { XtruRead(child); } else
2215  if (tag=="hype") { HypeRead(child); } else
2216  if (tag=="intersection") { BooleanRead(child,INTERSECTION); } else
2217  if (tag=="multiunion") { MultiUnionRead(child); } else
2218  if (tag=="orb") { OrbRead(child); } else
2219  if (tag=="para") { ParaRead(child); } else
2220  if (tag=="paraboloid") { ParaboloidRead(child); } else
2221  if (tag=="polycone") { PolyconeRead(child); } else
2222  if (tag=="genericPolycone") { GenericPolyconeRead(child); } else
2223  if (tag=="polyhedra") { PolyhedraRead(child); } else
2224  if (tag=="genericPolyhedra") { GenericPolyhedraRead(child); } else
2225  if (tag=="reflectedSolid") { ReflectedSolidRead(child); } else
2226  if (tag=="sphere") { SphereRead(child); } else
2227  if (tag=="subtraction") { BooleanRead(child,SUBTRACTION); } else
2228  if (tag=="tessellated") { TessellatedRead(child); } else
2229  if (tag=="tet") { TetRead(child); } else
2230  if (tag=="torus") { TorusRead(child); } else
2231  if (tag=="arb8") { GenTrapRead(child); } else
2232  if (tag=="trap") { TrapRead(child); } else
2233  if (tag=="trd") { TrdRead(child); } else
2234  if (tag=="tube") { TubeRead(child); } else
2235  if (tag=="cutTube") { CutTubeRead(child); } else
2236  if (tag=="twistedbox") { TwistedboxRead(child); } else
2237  if (tag=="twistedtrap") { TwistedtrapRead(child); } else
2238  if (tag=="twistedtrd") { TwistedtrdRead(child); } else
2239  if (tag=="twistedtubs") { TwistedtubsRead(child); } else
2240  if (tag=="union") { BooleanRead(child,UNION); } else
2241  if (tag=="opticalsurface") { OpticalSurfaceRead(child); } else
2242  if (tag=="loop") { LoopRead(child,&G4GDMLRead::SolidsRead); }
2243  else
2244  {
2245  G4String error_msg = "Unknown tag in solids: " + tag;
2246  G4Exception("G4GDMLReadSolids::SolidsRead()", "ReadError",
2247  FatalException, error_msg);
2248  }
2249  }
2250 }
2251 
2253 {
2254  G4VSolid* solidPtr = G4SolidStore::GetInstance()->GetSolid(ref,false);
2255 
2256  if (!solidPtr)
2257  {
2258  G4String error_msg = "Referenced solid '" + ref + "' was not found!";
2259  G4Exception("G4GDMLReadSolids::GetSolid()", "ReadError",
2260  FatalException, error_msg);
2261  }
2262 
2263  return solidPtr;
2264 }
2265 
2267 GetSurfaceProperty(const G4String& ref) const
2268 {
2269  const G4SurfacePropertyTable* surfaceList
2271  const size_t surfaceCount = surfaceList->size();
2272 
2273  for (size_t i=0; i<surfaceCount; i++)
2274  {
2275  if ((*surfaceList)[i]->GetName() == ref) { return (*surfaceList)[i]; }
2276  }
2277 
2278  G4String error_msg = "Referenced optical surface '" + ref + "' was not found!";
2279  G4Exception("G4GDMLReadSolids::GetSurfaceProperty()", "ReadError",
2280  FatalException, error_msg);
2281 
2282  return 0;
2283 }
zplaneType ZplaneRead(const xercesc::DOMElement *const)
void TubeRead(const xercesc::DOMElement *const)
void TwistedtrapRead(const xercesc::DOMElement *const)
G4int EvaluateInteger(const G4String &)
Definition: G4Para.hh:77
void SetSolidClosed(const G4bool t)
void BoxRead(const xercesc::DOMElement *const)
G4VSolid * GetSolid(const G4String &name, G4bool verbose=true) const
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
void OpticalSurfaceRead(const xercesc::DOMElement *const)
CLHEP::Hep3Vector G4ThreeVector
G4double z
Definition: TRTMaterials.hh:39
virtual void DefineRead(const xercesc::DOMElement *const)
Definition: G4Box.hh:64
void TwistedboxRead(const xercesc::DOMElement *const)
G4String name
Definition: TRTMaterials.hh:40
Definition: G4Tubs.hh:85
Definition: xmlparse.cc:179
void MultiUnionRead(const xercesc::DOMElement *const)
void TorusRead(const xercesc::DOMElement *const)
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void TetRead(const xercesc::DOMElement *const)
G4TwoVector TwoDimVertexRead(const xercesc::DOMElement *const, G4double)
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
void GenericPolyhedraRead(const xercesc::DOMElement *const)
Definition: G4Trd.hh:72
G4SurfaceType
G4TriangularFacet * TriangularRead(const xercesc::DOMElement *const)
int G4int
Definition: G4Types.hh:78
G4VSolid * GetSolid(const G4String &) const
void TwistedtubsRead(const xercesc::DOMElement *const)
void TrapRead(const xercesc::DOMElement *const)
void ConeRead(const xercesc::DOMElement *const)
G4String RefRead(const xercesc::DOMElement *const)
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:179
void SphereRead(const xercesc::DOMElement *const)
void OrbRead(const xercesc::DOMElement *const)
G4ExtrudedSolid::ZSection SectionRead(const xercesc::DOMElement *const, G4double)
G4GLOB_DLL std::ostream G4cout
void BooleanRead(const xercesc::DOMElement *const, const BooleanOp)
G4FacetVertexType
Definition: G4VFacet.hh:56
Definition: G4Tet.hh:65
Definition: G4Hype.hh:67
void ReflectedSolidRead(const xercesc::DOMElement *const)
void EllipsoidRead(const xercesc::DOMElement *const)
Definition: G4Cons.hh:83
G4bool AddFacet(G4VFacet *aFacet)
void GenTrapRead(const xercesc::DOMElement *const)
static G4SolidStore * GetInstance()
void TrdRead(const xercesc::DOMElement *const)
virtual ~G4GDMLReadSolids()
HepGeom::Transform3D G4Transform3D
HepGeom::Scale3D G4Scale3D
Definition: G4Orb.hh:61
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void SolidsRead(const xercesc::DOMElement *const)
void EltubeRead(const xercesc::DOMElement *const)
rzPointType RZPointRead(const xercesc::DOMElement *const)
std::vector< G4SurfaceProperty * > G4SurfacePropertyTable
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
void TessellatedRead(const xercesc::DOMElement *const)
void PolyhedraRead(const xercesc::DOMElement *const)
void XtruRead(const xercesc::DOMElement *const)
void CutTubeRead(const xercesc::DOMElement *const)
G4OpticalSurfaceFinish
int position
Definition: filter.cc:7
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
void ParaboloidRead(const xercesc::DOMElement *const)
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
void GenericPolyconeRead(const xercesc::DOMElement *const)
G4OpticalSurfaceModel
void ElconeRead(const xercesc::DOMElement *const)
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector GetRotation(const G4String &)
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
void ParaRead(const xercesc::DOMElement *const)
double G4double
Definition: G4Types.hh:76
void HypeRead(const xercesc::DOMElement *const)
static const G4double alpha
void PolyconeRead(const xercesc::DOMElement *const)
G4QuadrangularFacet * QuadrangularRead(const xercesc::DOMElement *const)
static const G4SurfacePropertyTable * GetSurfacePropertyTable()
virtual void SolidsRead(const xercesc::DOMElement *const)=0
void TwistedtrdRead(const xercesc::DOMElement *const)
G4double Evaluate(const G4String &)
G4ThreeVector GetPosition(const G4String &)