source: trunk/source/persistency/gdml/src/G4GDMLReadSolids.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 69.8 KB
Line 
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,v 1.32 2010/10/14 16:19:40 gcosmo Exp $
27// GEANT4 tag $Name: gdml-V09-03-09 $
28//
29// class G4GDMLReadSolids Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
35#include "G4GDMLReadSolids.hh"
36
37#include "G4Box.hh"
38#include "G4Cons.hh"
39#include "G4Ellipsoid.hh"
40#include "G4EllipticalCone.hh"
41#include "G4EllipticalTube.hh"
42#include "G4Hype.hh"
43#include "G4IntersectionSolid.hh"
44#include "G4Orb.hh"
45#include "G4Para.hh"
46#include "G4Paraboloid.hh"
47#include "G4Polycone.hh"
48#include "G4Polyhedra.hh"
49#include "G4QuadrangularFacet.hh"
50#include "G4ReflectedSolid.hh"
51#include "G4Sphere.hh"
52#include "G4SolidStore.hh"
53#include "G4SubtractionSolid.hh"
54#include "G4GenericTrap.hh"
55#include "G4TessellatedSolid.hh"
56#include "G4Tet.hh"
57#include "G4Torus.hh"
58#include "G4Transform3D.hh"
59#include "G4Trap.hh"
60#include "G4Trd.hh"
61#include "G4TriangularFacet.hh"
62#include "G4Tubs.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
71G4GDMLReadSolids::G4GDMLReadSolids() : G4GDMLReadMaterials()
72{
73}
74
75G4GDMLReadSolids::~G4GDMLReadSolids()
76{
77}
78
79void G4GDMLReadSolids::
80BooleanRead(const xercesc::DOMElement* const booleanElement, const BooleanOp op)
81{
82   G4String name;
83   G4String first;
84   G4String second;
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") { second = 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(second));
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
177void G4GDMLReadSolids::BoxRead(const xercesc::DOMElement* const boxElement)
178{
179   G4String name;
180   G4double lunit = 1.0;
181   G4double aunit = 1.0;
182   G4double x = 0.0;
183   G4double y = 0.0;
184   G4double z = 0.0;
185
186   const xercesc::DOMNamedNodeMap* const attributes
187         = boxElement->getAttributes();
188   XMLSize_t attributeCount = attributes->getLength();
189
190   for (XMLSize_t attribute_index=0;
191        attribute_index<attributeCount; attribute_index++)
192   {
193      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
194
195      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
196        { continue; }
197
198      const xercesc::DOMAttr* const attribute
199            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
200      if (!attribute)
201      {
202        G4Exception("G4GDMLReadSolids::BoxRead()",
203                    "InvalidRead", FatalException, "No attribute found!");
204        return;
205      }
206      const G4String attName = Transcode(attribute->getName());
207      const G4String attValue = Transcode(attribute->getValue());
208
209      if (attName=="name") { name = GenerateName(attValue); } else
210      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
211      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
212      if (attName=="x") { x = eval.Evaluate(attValue); } else
213      if (attName=="y") { y = eval.Evaluate(attValue); } else
214      if (attName=="z") { z = eval.Evaluate(attValue); }
215   }
216
217   x *= 0.5*lunit;
218   y *= 0.5*lunit;
219   z *= 0.5*lunit;
220
221   new G4Box(name,x,y,z);
222}
223
224void G4GDMLReadSolids::ConeRead(const xercesc::DOMElement* const coneElement)
225{
226   G4String name;
227   G4double lunit = 1.0;
228   G4double aunit = 1.0;
229   G4double rmin1 = 0.0;
230   G4double rmax1 = 0.0;
231   G4double rmin2 = 0.0;
232   G4double rmax2 = 0.0;
233   G4double z = 0.0;
234   G4double startphi = 0.0;
235   G4double deltaphi = 0.0;
236
237   const xercesc::DOMNamedNodeMap* const attributes
238         = coneElement->getAttributes();
239   XMLSize_t attributeCount = attributes->getLength();
240
241   for (XMLSize_t attribute_index=0;
242        attribute_index<attributeCount; attribute_index++)
243   {
244      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
245
246      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
247        { continue; }
248
249      const xercesc::DOMAttr* const attribute
250            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
251      if (!attribute)
252      {
253        G4Exception("G4GDMLReadSolids::ConeRead()",
254                    "InvalidRead", FatalException, "No attribute found!");
255        return;
256      }
257      const G4String attName = Transcode(attribute->getName());
258      const G4String attValue = Transcode(attribute->getValue());
259
260      if (attName=="name") { name = GenerateName(attValue); } else
261      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
262      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
263      if (attName=="rmin1") { rmin1 = eval.Evaluate(attValue); } else
264      if (attName=="rmax1") { rmax1 = eval.Evaluate(attValue); } else
265      if (attName=="rmin2") { rmin2 = eval.Evaluate(attValue); } else
266      if (attName=="rmax2") { rmax2 = eval.Evaluate(attValue); } else
267      if (attName=="z") { z = eval.Evaluate(attValue); } else
268      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
269      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
270   }
271
272   rmin1 *= lunit;
273   rmax1 *= lunit;
274   rmin2 *= lunit;
275   rmax2 *= lunit;
276   z *= 0.5*lunit;
277   startphi *= aunit;
278   deltaphi *= aunit;
279
280   new G4Cons(name,rmin1,rmax1,rmin2,rmax2,z,startphi,deltaphi);
281}
282
283void G4GDMLReadSolids::
284ElconeRead(const xercesc::DOMElement* const elconeElement)
285{
286   G4String name;
287   G4double lunit = 1.0;
288   G4double dx = 0.0;
289   G4double dy = 0.0;
290   G4double zmax = 0.0;
291   G4double zcut = 0.0;
292
293   const xercesc::DOMNamedNodeMap* const attributes
294         = elconeElement->getAttributes();
295   XMLSize_t attributeCount = attributes->getLength();
296
297   for (XMLSize_t attribute_index=0;
298        attribute_index<attributeCount; attribute_index++)
299   {
300      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
301
302      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
303        { continue; }
304
305      const xercesc::DOMAttr* const attribute
306            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
307      if (!attribute)
308      {
309        G4Exception("G4GDMLReadSolids::ElconeRead()",
310                    "InvalidRead", FatalException, "No attribute found!");
311        return;
312      }
313      const G4String attName = Transcode(attribute->getName());
314      const G4String attValue = Transcode(attribute->getValue());
315
316      if (attName=="name") { name = GenerateName(attValue); } else
317      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
318      if (attName=="dx") { dx = eval.Evaluate(attValue); } else
319      if (attName=="dy") { dy = eval.Evaluate(attValue); } else
320      if (attName=="zmax") { zmax = eval.Evaluate(attValue); } else
321      if (attName=="zcut") { zcut = eval.Evaluate(attValue); }
322   }
323
324   dx *= lunit;
325   dy *= lunit;
326   zmax *= lunit;
327   zcut *= lunit;
328
329   new G4EllipticalCone(name,dx,dy,zmax,zcut);
330}
331
332void G4GDMLReadSolids::
333EllipsoidRead(const xercesc::DOMElement* const ellipsoidElement)
334{
335   G4String name;
336   G4double lunit = 1.0;
337   G4double ax = 0.0;
338   G4double by = 0.0;
339   G4double cz = 0.0;
340   G4double zcut1 = 0.0;
341   G4double zcut2 = 0.0; 
342
343   const xercesc::DOMNamedNodeMap* const attributes
344         = ellipsoidElement->getAttributes();
345   XMLSize_t attributeCount = attributes->getLength();
346
347   for (XMLSize_t attribute_index=0;
348        attribute_index<attributeCount; attribute_index++)
349   {
350      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
351
352      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
353        { continue; }
354
355      const xercesc::DOMAttr* const attribute
356            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
357      if (!attribute)
358      {
359        G4Exception("G4GDMLReadSolids::EllipsoidRead()",
360                    "InvalidRead", FatalException, "No attribute found!");
361        return;
362      }
363      const G4String attName = Transcode(attribute->getName());
364      const G4String attValue = Transcode(attribute->getValue());
365
366      if (attName=="name") { name  = GenerateName(attValue); } else
367      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
368      if (attName=="ax") { ax = eval.Evaluate(attValue); } else
369      if (attName=="by") { by = eval.Evaluate(attValue); } else
370      if (attName=="cz") { cz = eval.Evaluate(attValue); } else
371      if (attName=="zcut1") { zcut1 = eval.Evaluate(attValue); } else
372      if (attName=="zcut2") { zcut2 = eval.Evaluate(attValue); }
373   }
374
375   ax *= lunit;
376   by *= lunit;
377   cz *= lunit;
378   zcut1 *= lunit;
379   zcut2 *= lunit; 
380
381   new G4Ellipsoid(name,ax,by,cz,zcut1,zcut2);
382}
383
384void G4GDMLReadSolids::
385EltubeRead(const xercesc::DOMElement* const eltubeElement)
386{
387   G4String name;
388   G4double lunit = 1.0;
389   G4double dx = 0.0;
390   G4double dy = 0.0;
391   G4double dz = 0.0;
392
393   const xercesc::DOMNamedNodeMap* const attributes
394         = eltubeElement->getAttributes();
395   XMLSize_t attributeCount = attributes->getLength();
396
397   for (XMLSize_t attribute_index=0;
398        attribute_index<attributeCount; attribute_index++)
399   {
400      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
401
402      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
403        { continue; }
404
405      const xercesc::DOMAttr* const attribute
406            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
407      if (!attribute)
408      {
409        G4Exception("G4GDMLReadSolids::EltubeRead()",
410                    "InvalidRead", FatalException, "No attribute found!");
411        return;
412      }
413      const G4String attName = Transcode(attribute->getName());
414      const G4String attValue = Transcode(attribute->getValue());
415
416      if (attName=="name") { name = GenerateName(attValue); } else
417      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
418      if (attName=="dx") { dx = eval.Evaluate(attValue); } else
419      if (attName=="dy") { dy = eval.Evaluate(attValue); } else
420      if (attName=="dz") { dz = eval.Evaluate(attValue); }
421   }
422
423   dx *= lunit;
424   dy *= lunit;
425   dz *= lunit;
426
427   new G4EllipticalTube(name,dx,dy,dz);
428}
429
430void G4GDMLReadSolids::XtruRead(const xercesc::DOMElement* const xtruElement)
431{
432   G4String name;
433   G4double lunit = 1.0;
434
435   const xercesc::DOMNamedNodeMap* const attributes
436         = xtruElement->getAttributes();
437   XMLSize_t attributeCount = attributes->getLength();
438
439   for (XMLSize_t attribute_index=0;
440        attribute_index<attributeCount; attribute_index++)
441   {
442      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
443
444      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
445        { continue; } 
446
447      const xercesc::DOMAttr* const attribute
448            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
449      if (!attribute)
450      {
451        G4Exception("G4GDMLReadSolids::XtruRead()",
452                    "InvalidRead", FatalException, "No attribute found!");
453        return;
454      }
455      const G4String attName = Transcode(attribute->getName());
456      const G4String attValue = Transcode(attribute->getValue());
457
458      if (attName=="name") { name = GenerateName(attValue); } else
459      if (attName=="lunit") { lunit = eval.Evaluate(attValue); }
460   }
461
462   std::vector<G4TwoVector> twoDimVertexList;
463   std::vector<G4ExtrudedSolid::ZSection> sectionList;
464
465   for (xercesc::DOMNode* iter = xtruElement->getFirstChild();
466        iter != 0; iter = iter->getNextSibling())
467   {
468      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
469
470      const xercesc::DOMElement* const child
471            = dynamic_cast<xercesc::DOMElement*>(iter);
472      if (!child)
473      {
474        G4Exception("G4GDMLReadSolids::XtruRead()",
475                    "InvalidRead", FatalException, "No child found!");
476        return;
477      }
478      const G4String tag = Transcode(child->getTagName());
479
480      if (tag=="twoDimVertex")
481        { twoDimVertexList.push_back(TwoDimVertexRead(child,lunit)); } else
482      if (tag=="section")
483        { sectionList.push_back(SectionRead(child,lunit)); }
484   }
485
486   new G4ExtrudedSolid(name,twoDimVertexList,sectionList);
487}
488
489void G4GDMLReadSolids::HypeRead(const xercesc::DOMElement* const hypeElement)
490{
491   G4String name;
492   G4double lunit = 1.0;
493   G4double aunit = 1.0;
494   G4double rmin = 0.0;
495   G4double rmax = 0.0;
496   G4double inst = 0.0;
497   G4double outst = 0.0;
498   G4double z = 0.0;
499
500   const xercesc::DOMNamedNodeMap* const attributes
501         = hypeElement->getAttributes();
502   XMLSize_t attributeCount = attributes->getLength();
503
504   for (XMLSize_t attribute_index=0;
505        attribute_index<attributeCount; attribute_index++)
506   {
507      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
508
509      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
510        { continue; }
511
512      const xercesc::DOMAttr* const attribute
513            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
514      if (!attribute)
515      {
516        G4Exception("G4GDMLReadSolids::HypeRead()",
517                    "InvalidRead", FatalException, "No attribute found!");
518        return;
519      }
520      const G4String attName = Transcode(attribute->getName());
521      const G4String attValue = Transcode(attribute->getValue());
522
523      if (attName=="name") { name = GenerateName(attValue); } else
524      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
525      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
526      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
527      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
528      if (attName=="inst") { inst = eval.Evaluate(attValue); } else
529      if (attName=="outst") { outst = eval.Evaluate(attValue); } else
530      if (attName=="z") { z = eval.Evaluate(attValue); }
531   }
532
533   rmin *= lunit;
534   rmax *= lunit;
535   inst *= aunit;
536   outst *= aunit;
537   z *= 0.5*lunit;
538
539   new G4Hype(name,rmin,rmax,inst,outst,z);
540}
541
542void G4GDMLReadSolids::OrbRead(const xercesc::DOMElement* const orbElement)
543{
544   G4String name;
545   G4double lunit = 1.0;
546   G4double r = 0.0;
547
548   const xercesc::DOMNamedNodeMap* const attributes
549         = orbElement->getAttributes();
550   XMLSize_t attributeCount = attributes->getLength();
551
552   for (XMLSize_t attribute_index=0;
553        attribute_index<attributeCount; attribute_index++)
554   {
555      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
556
557      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
558        { continue; }
559
560      const xercesc::DOMAttr* const attribute
561            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
562      if (!attribute)
563      {
564        G4Exception("G4GDMLReadSolids::OrbRead()",
565                    "InvalidRead", FatalException, "No attribute found!");
566        return;
567      }
568      const G4String attName = Transcode(attribute->getName());
569      const G4String attValue = Transcode(attribute->getValue());
570
571      if (attName=="name") { name = GenerateName(attValue); } else
572      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
573      if (attName=="r") { r = eval.Evaluate(attValue); }
574   }
575
576   r *= lunit;
577
578   new G4Orb(name,r);
579}
580
581void G4GDMLReadSolids::ParaRead(const xercesc::DOMElement* const paraElement)
582{
583   G4String name;
584   G4double lunit = 1.0;
585   G4double aunit = 1.0;
586   G4double x = 0.0;
587   G4double y = 0.0;
588   G4double z = 0.0;
589   G4double alpha = 0.0;
590   G4double theta = 0.0;
591   G4double phi = 0.0;
592
593   const xercesc::DOMNamedNodeMap* const attributes
594         = paraElement->getAttributes();
595   XMLSize_t attributeCount = attributes->getLength();
596
597   for (XMLSize_t attribute_index=0;
598        attribute_index<attributeCount; attribute_index++)
599   {
600      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
601
602      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
603        { continue; }
604
605      const xercesc::DOMAttr* const attribute
606            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
607      if (!attribute)
608      {
609        G4Exception("G4GDMLReadSolids::ParaRead()",
610                    "InvalidRead", FatalException, "No attribute found!");
611        return;
612      }
613      const G4String attName = Transcode(attribute->getName());
614      const G4String attValue = Transcode(attribute->getValue());
615
616      if (attName=="name") { name = GenerateName(attValue); } else
617      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
618      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
619      if (attName=="x") { x = eval.Evaluate(attValue); } else
620      if (attName=="y") { y = eval.Evaluate(attValue); } else
621      if (attName=="z") { z = eval.Evaluate(attValue); } else
622      if (attName=="alpha") { alpha = eval.Evaluate(attValue); } else
623      if (attName=="theta") { theta = eval.Evaluate(attValue); } else
624      if (attName=="phi") { phi = eval.Evaluate(attValue); }
625   }
626
627   x *= 0.5*lunit;
628   y *= 0.5*lunit;
629   z *= 0.5*lunit;
630   alpha *= aunit;
631   theta *= aunit;
632   phi *= aunit;
633
634   new G4Para(name,x,y,z,alpha,theta,phi);
635}
636
637void G4GDMLReadSolids::
638ParaboloidRead(const xercesc::DOMElement* const paraElement)
639{
640   G4String name;
641   G4double lunit = 1.0;
642   G4double aunit = 1.0;
643   G4double rlo = 0.0;
644   G4double rhi = 0.0;
645   G4double dz = 0.0;
646 
647   const xercesc::DOMNamedNodeMap* const attributes
648         = paraElement->getAttributes();
649   XMLSize_t attributeCount = attributes->getLength();
650
651   for (XMLSize_t attribute_index=0;
652        attribute_index<attributeCount; attribute_index++)
653   {
654      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
655
656      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
657        { continue; }
658
659      const xercesc::DOMAttr* const attribute
660            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
661      if (!attribute)
662      {
663        G4Exception("G4GDMLReadSolids::ParaboloidRead()",
664                    "InvalidRead", FatalException, "No attribute found!");
665        return;
666      }
667      const G4String attName = Transcode(attribute->getName());
668      const G4String attValue = Transcode(attribute->getValue());
669
670      if (attName=="name")  { name = GenerateName(attValue); } else
671      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
672      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
673      if (attName=="rlo")   { rlo =  eval.Evaluate(attValue); } else
674      if (attName=="rhi")   { rhi = eval.Evaluate(attValue); } else
675      if (attName=="dz")    { dz = eval.Evaluate(attValue); } 
676   }     
677
678   rlo *= 1.*lunit;
679   rhi *= 1.*lunit;
680   dz *= 1.*lunit;
681 
682   new G4Paraboloid(name,dz,rlo,rhi);
683}
684
685void G4GDMLReadSolids::
686PolyconeRead(const xercesc::DOMElement* const polyconeElement) 
687{
688   G4String name;
689   G4double lunit = 1.0;
690   G4double aunit = 1.0;
691   G4double startphi = 0.0;
692   G4double deltaphi = 0.0;
693
694   const xercesc::DOMNamedNodeMap* const attributes
695         = polyconeElement->getAttributes();
696   XMLSize_t attributeCount = attributes->getLength();
697
698   for (XMLSize_t attribute_index=0;
699        attribute_index<attributeCount; attribute_index++)
700   {
701      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
702
703      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
704        { continue; }
705
706      const xercesc::DOMAttr* const attribute
707            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
708      if (!attribute)
709      {
710        G4Exception("G4GDMLReadSolids::PolyconeRead()",
711                    "InvalidRead", FatalException, "No attribute found!");
712        return;
713      }
714      const G4String attName = Transcode(attribute->getName());
715      const G4String attValue = Transcode(attribute->getValue());
716
717      if (attName=="name") { name = GenerateName(attValue); } else
718      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
719      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
720      if (attName=="startphi") { startphi = eval.Evaluate(attValue); }else
721      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
722   }
723
724   startphi *= aunit;
725   deltaphi *= aunit;
726
727   std::vector<zplaneType> zplaneList;
728
729   for (xercesc::DOMNode* iter = polyconeElement->getFirstChild();
730        iter != 0; iter = iter->getNextSibling())
731   {
732      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
733
734      const xercesc::DOMElement* const child
735            = dynamic_cast<xercesc::DOMElement*>(iter);
736      if (!child)
737      {
738        G4Exception("G4GDMLReadSolids::PolyconeRead()",
739                    "InvalidRead", FatalException, "No child found!");
740        return;
741      }
742      const G4String tag = Transcode(child->getTagName());
743
744      if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
745   }
746
747   G4int numZPlanes = zplaneList.size();
748
749   G4double* rmin_array = new G4double[numZPlanes];
750   G4double* rmax_array = new G4double[numZPlanes];
751   G4double* z_array    = new G4double[numZPlanes];
752
753   for (G4int i=0; i<numZPlanes; i++)
754   { 
755      rmin_array[i] = zplaneList[i].rmin*lunit;
756      rmax_array[i] = zplaneList[i].rmax*lunit;
757      z_array[i]    = zplaneList[i].z*lunit;
758   }
759
760   new G4Polycone(name,startphi,deltaphi,numZPlanes,
761                  z_array,rmin_array,rmax_array);
762}
763
764void G4GDMLReadSolids::
765PolyhedraRead(const xercesc::DOMElement* const polyhedraElement)
766{
767   G4String name;
768   G4double lunit = 1.0;
769   G4double aunit = 1.0;
770   G4double startphi = 0.0;
771   G4double deltaphi = 0.0;
772   G4int numsides = 0;
773
774   const xercesc::DOMNamedNodeMap* const attributes
775         = polyhedraElement->getAttributes();
776   XMLSize_t attributeCount = attributes->getLength();
777
778   for (XMLSize_t attribute_index=0;
779        attribute_index<attributeCount; attribute_index++)
780   {
781      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
782
783      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
784        { continue; }
785
786      const xercesc::DOMAttr* const attribute
787            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
788      if (!attribute)
789      {
790        G4Exception("G4GDMLReadSolids::PolyhedraRead()",
791                    "InvalidRead", FatalException, "No attribute found!");
792        return;
793      }
794      const G4String attName = Transcode(attribute->getName());
795      const G4String attValue = Transcode(attribute->getValue());
796
797      if (attName=="name") { name = GenerateName(attValue); } else
798      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
799      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
800      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
801      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
802      if (attName=="numsides") { numsides = eval.EvaluateInteger(attValue); }
803   }
804
805   startphi *= aunit;
806   deltaphi *= aunit;
807
808   std::vector<zplaneType> zplaneList;
809
810   for (xercesc::DOMNode* iter = polyhedraElement->getFirstChild();
811        iter != 0; iter = iter->getNextSibling())
812   {
813      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)  { continue; }
814
815      const xercesc::DOMElement* const child
816            = dynamic_cast<xercesc::DOMElement*>(iter);
817      if (!child)
818      {
819        G4Exception("G4GDMLReadSolids::PolyhedraRead()",
820                    "InvalidRead", FatalException, "No child found!");
821        return;
822      }
823      const G4String tag = Transcode(child->getTagName());
824
825      if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
826   }
827
828   G4int numZPlanes = zplaneList.size();
829
830   G4double* rmin_array = new G4double[numZPlanes];
831   G4double* rmax_array = new G4double[numZPlanes];
832   G4double* z_array = new G4double[numZPlanes];
833
834   for (G4int i=0; i<numZPlanes; i++)
835   { 
836      rmin_array[i] = zplaneList[i].rmin*lunit;
837      rmax_array[i] = zplaneList[i].rmax*lunit;
838      z_array[i] = zplaneList[i].z*lunit;
839   }
840
841   new G4Polyhedra(name,startphi,deltaphi,numsides,numZPlanes,
842                   z_array,rmin_array,rmax_array);
843
844}
845
846G4QuadrangularFacet* G4GDMLReadSolids::
847QuadrangularRead(const xercesc::DOMElement* const quadrangularElement)
848{
849   G4ThreeVector vertex1;
850   G4ThreeVector vertex2;
851   G4ThreeVector vertex3;
852   G4ThreeVector vertex4;
853   G4FacetVertexType type = ABSOLUTE;
854
855   const xercesc::DOMNamedNodeMap* const attributes
856         = quadrangularElement->getAttributes();
857   XMLSize_t attributeCount = attributes->getLength();
858
859   for (XMLSize_t attribute_index=0;
860        attribute_index<attributeCount; attribute_index++)
861   {
862      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
863
864      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
865        { continue; }
866
867      const xercesc::DOMAttr* const attribute
868            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
869      if (!attribute)
870      {
871        G4Exception("G4GDMLReadSolids::QuadrangularRead()",
872                    "InvalidRead", FatalException, "No attribute found!");
873        return 0;
874      }
875      const G4String attName = Transcode(attribute->getName());
876      const G4String attValue = Transcode(attribute->getValue());
877
878      if (attName=="vertex1")
879        { vertex1 = GetPosition(GenerateName(attValue)); } else
880      if (attName=="vertex2")
881        { vertex2 = GetPosition(GenerateName(attValue)); } else
882      if (attName=="vertex3")
883        { vertex3 = GetPosition(GenerateName(attValue)); } else
884      if (attName=="vertex4")
885        { vertex4 = GetPosition(GenerateName(attValue)); } else
886      if (attName=="type")
887        { if (attValue=="RELATIVE") { type = RELATIVE; } }
888   }
889
890   return new G4QuadrangularFacet(vertex1,vertex2,vertex3,vertex4,type);
891}
892
893void G4GDMLReadSolids::
894ReflectedSolidRead(const xercesc::DOMElement* const reflectedSolidElement)
895{
896   G4String name;
897   G4double lunit = 1.0;
898   G4double aunit = 1.0;
899   G4String solid;
900   G4ThreeVector scale(1.0,1.0,1.0);
901   G4ThreeVector rotation;
902   G4ThreeVector position;
903
904   const xercesc::DOMNamedNodeMap* const attributes
905         = reflectedSolidElement->getAttributes();
906   XMLSize_t attributeCount = attributes->getLength();
907
908   for (XMLSize_t attribute_index=0;
909        attribute_index<attributeCount; attribute_index++)
910   {
911      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
912
913      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
914        { continue; }
915
916      const xercesc::DOMAttr* const attribute
917            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
918      if (!attribute)
919      {
920        G4Exception("G4GDMLReadSolids::ReflectedSolidRead()",
921                    "InvalidRead", FatalException, "No attribute found!");
922        return;
923      }
924      const G4String attName = Transcode(attribute->getName());
925      const G4String attValue = Transcode(attribute->getValue());
926
927      if (attName=="name") { name = GenerateName(attValue); } else
928      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
929      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
930      if (attName=="solid") { solid = GenerateName(attValue); } else
931      if (attName=="sx") { scale.setX(eval.Evaluate(attValue)); } else
932      if (attName=="sy") { scale.setY(eval.Evaluate(attValue)); } else
933      if (attName=="sz") { scale.setZ(eval.Evaluate(attValue)); } else
934      if (attName=="rx") { rotation.setX(eval.Evaluate(attValue)); } else
935      if (attName=="ry") { rotation.setY(eval.Evaluate(attValue)); } else
936      if (attName=="rz") { rotation.setZ(eval.Evaluate(attValue)); } else
937      if (attName=="dx") { position.setX(eval.Evaluate(attValue)); } else
938      if (attName=="dy") { position.setY(eval.Evaluate(attValue)); } else
939      if (attName=="dz") { position.setZ(eval.Evaluate(attValue)); }
940   }
941
942   rotation *= aunit;
943   position *= lunit;
944
945   G4Transform3D transform(GetRotationMatrix(rotation),position);
946   transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
947         
948   new G4ReflectedSolid(name,GetSolid(solid),transform);
949}
950
951G4ExtrudedSolid::ZSection G4GDMLReadSolids::
952SectionRead(const xercesc::DOMElement* const sectionElement,G4double lunit) 
953{
954   G4double zPosition = 0.0;
955   G4TwoVector Offset;
956   G4double scalingFactor = 1.0;
957
958   const xercesc::DOMNamedNodeMap* const attributes
959         = sectionElement->getAttributes();
960   XMLSize_t attributeCount = attributes->getLength();
961
962   for (XMLSize_t attribute_index=0;
963        attribute_index<attributeCount; attribute_index++)
964   {
965      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
966
967      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
968        { continue; }
969
970      const xercesc::DOMAttr* const attribute
971            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
972      if (!attribute)
973      {
974        G4Exception("G4GDMLReadSolids::SectionRead()",
975                    "InvalidRead", FatalException, "No attribute found!");
976        return G4ExtrudedSolid::ZSection(zPosition,Offset,scalingFactor);
977      }
978      const G4String attName = Transcode(attribute->getName());
979      const G4String attValue = Transcode(attribute->getValue());
980
981      if (attName=="zPosition")
982        { zPosition = eval.Evaluate(attValue)*lunit; } else
983      if (attName=="xOffset")
984        { Offset.setX(eval.Evaluate(attValue)*lunit); } else
985      if (attName=="yOffset")
986        { Offset.setY(eval.Evaluate(attValue)*lunit); } else
987      if (attName=="scalingFactor")
988        { scalingFactor = eval.Evaluate(attValue); }
989   }
990
991   return G4ExtrudedSolid::ZSection(zPosition,Offset,scalingFactor);
992}
993
994void G4GDMLReadSolids::
995SphereRead(const xercesc::DOMElement* const sphereElement)
996{
997   G4String name;
998   G4double lunit = 1.0;
999   G4double aunit = 1.0;
1000   G4double rmin = 0.0;
1001   G4double rmax = 0.0;
1002   G4double startphi = 0.0;
1003   G4double deltaphi = 0.0;
1004   G4double starttheta = 0.0;
1005   G4double deltatheta = 0.0;
1006
1007   const xercesc::DOMNamedNodeMap* const attributes
1008         = sphereElement->getAttributes();
1009   XMLSize_t attributeCount = attributes->getLength();
1010
1011   for (XMLSize_t attribute_index=0;
1012        attribute_index<attributeCount; attribute_index++)
1013   {
1014      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1015
1016      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1017        { continue; }
1018
1019      const xercesc::DOMAttr* const attribute
1020            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1021      if (!attribute)
1022      {
1023        G4Exception("G4GDMLReadSolids::SphereRead()",
1024                    "InvalidRead", FatalException, "No attribute found!");
1025        return;
1026      }
1027      const G4String attName = Transcode(attribute->getName());
1028      const G4String attValue = Transcode(attribute->getValue());
1029
1030      if (attName=="name") { name = GenerateName(attValue); } else
1031      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1032      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1033      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1034      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1035      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1036      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
1037      if (attName=="starttheta") { starttheta = eval.Evaluate(attValue); } else
1038      if (attName=="deltatheta") { deltatheta = eval.Evaluate(attValue); }
1039   }
1040
1041   rmin *= lunit;
1042   rmax *= lunit;
1043   startphi *= aunit;
1044   deltaphi *= aunit;
1045   starttheta *= aunit;
1046   deltatheta *= aunit;
1047
1048   new G4Sphere(name,rmin,rmax,startphi,deltaphi,starttheta,deltatheta);
1049}
1050
1051void G4GDMLReadSolids::
1052TessellatedRead(const xercesc::DOMElement* const tessellatedElement)
1053{
1054   G4String name;
1055
1056   const xercesc::DOMNamedNodeMap* const attributes
1057         = tessellatedElement->getAttributes();
1058   XMLSize_t attributeCount = attributes->getLength();
1059
1060   for (XMLSize_t attribute_index=0;
1061        attribute_index<attributeCount; attribute_index++)
1062   {
1063      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1064
1065      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1066        { continue; }
1067
1068      const xercesc::DOMAttr* const attribute
1069            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1070      if (!attribute)
1071      {
1072        G4Exception("G4GDMLReadSolids::TessellatedRead()",
1073                    "InvalidRead", FatalException, "No attribute found!");
1074        return;
1075      }
1076      const G4String attName = Transcode(attribute->getName());
1077      const G4String attValue = Transcode(attribute->getValue());
1078
1079      if (attName=="name") { name = GenerateName(attValue); }
1080   }
1081   
1082   G4TessellatedSolid *tessellated = new G4TessellatedSolid(name);
1083
1084   for (xercesc::DOMNode* iter = tessellatedElement->getFirstChild();
1085        iter != 0; iter = iter->getNextSibling())
1086   {
1087      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
1088
1089      const xercesc::DOMElement* const child
1090            = dynamic_cast<xercesc::DOMElement*>(iter);
1091      if (!child)
1092      {
1093        G4Exception("G4GDMLReadSolids::TessellatedRead()",
1094                    "InvalidRead", FatalException, "No child found!");
1095        return;
1096      }
1097      const G4String tag = Transcode(child->getTagName());
1098
1099      if (tag=="triangular")
1100        { tessellated->AddFacet(TriangularRead(child)); } else
1101      if (tag=="quadrangular")
1102        { tessellated->AddFacet(QuadrangularRead(child)); }
1103   }
1104
1105   tessellated->SetSolidClosed(true);
1106}
1107
1108void G4GDMLReadSolids::TetRead(const xercesc::DOMElement* const tetElement)
1109{
1110   G4String name;
1111   G4ThreeVector vertex1;
1112   G4ThreeVector vertex2;
1113   G4ThreeVector vertex3;
1114   G4ThreeVector vertex4;
1115   
1116   const xercesc::DOMNamedNodeMap* const attributes
1117         = tetElement->getAttributes();
1118   XMLSize_t attributeCount = attributes->getLength();
1119
1120   for (XMLSize_t attribute_index=0;
1121        attribute_index<attributeCount;attribute_index++)
1122   {
1123      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1124
1125      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1126        { continue; }
1127
1128      const xercesc::DOMAttr* const attribute
1129            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1130      if (!attribute)
1131      {
1132        G4Exception("G4GDMLReadSolids::TetRead()",
1133                    "InvalidRead", FatalException, "No attribute found!");
1134        return;
1135      }
1136      const G4String attName = Transcode(attribute->getName());
1137      const G4String attValue = Transcode(attribute->getValue());
1138
1139      if (attName=="name")
1140        { name = GenerateName(attValue); } else
1141      if (attName=="vertex1")
1142        { vertex1 = GetPosition(GenerateName(attValue)); } else
1143      if (attName=="vertex2")
1144        { vertex2 = GetPosition(GenerateName(attValue)); } else
1145      if (attName=="vertex3")
1146        { vertex3 = GetPosition(GenerateName(attValue)); } else
1147      if (attName=="vertex4")
1148        { vertex4 = GetPosition(GenerateName(attValue)); }
1149   }
1150
1151   new G4Tet(name,vertex1,vertex2,vertex3,vertex4);
1152}
1153
1154void G4GDMLReadSolids::TorusRead(const xercesc::DOMElement* const torusElement)
1155{
1156   G4String name;
1157   G4double lunit = 1.0;
1158   G4double aunit = 1.0;
1159   G4double rmin = 0.0;
1160   G4double rmax = 0.0;
1161   G4double rtor = 0.0;
1162   G4double startphi = 0.0;
1163   G4double deltaphi = 0.0;
1164
1165   const xercesc::DOMNamedNodeMap* const attributes
1166         = torusElement->getAttributes();
1167   XMLSize_t attributeCount = attributes->getLength();
1168
1169   for (XMLSize_t attribute_index=0;
1170        attribute_index<attributeCount; attribute_index++)
1171   {
1172      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1173
1174      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1175        { continue; }
1176
1177      const xercesc::DOMAttr* const attribute
1178            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1179      if (!attribute)
1180      {
1181        G4Exception("G4GDMLReadSolids::TorusRead()",
1182                    "InvalidRead", FatalException, "No attribute found!");
1183        return;
1184      }
1185      const G4String attName = Transcode(attribute->getName());
1186      const G4String attValue = Transcode(attribute->getValue());
1187
1188      if (attName=="name") { name = GenerateName(attValue); } else
1189      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1190      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1191      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1192      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1193      if (attName=="rtor") { rtor = eval.Evaluate(attValue); } else
1194      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1195      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1196   }
1197
1198   rmin *= lunit;
1199   rmax *= lunit;
1200   rtor *= lunit;
1201   startphi *= aunit;
1202   deltaphi *= aunit;
1203
1204   new G4Torus(name,rmin,rmax,rtor,startphi,deltaphi);
1205}
1206
1207void G4GDMLReadSolids::
1208GenTrapRead(const xercesc::DOMElement* const gtrapElement)
1209{
1210   G4String name;
1211   G4double lunit = 1.0;
1212   G4double dz =0.0;
1213   G4double v1x=0.0, v1y=0.0, v2x=0.0, v2y=0.0, v3x=0.0, v3y=0.0,
1214            v4x=0.0, v4y=0.0, v5x=0.0, v5y=0.0, v6x=0.0, v6y=0.0,
1215            v7x=0.0, v7y=0.0, v8x=0.0, v8y=0.0;
1216
1217   const xercesc::DOMNamedNodeMap* const attributes
1218         = gtrapElement->getAttributes();
1219   XMLSize_t attributeCount = attributes->getLength();
1220
1221   for (XMLSize_t attribute_index=0;
1222        attribute_index<attributeCount; attribute_index++)
1223   {
1224      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1225
1226      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1227        { continue; }
1228
1229      const xercesc::DOMAttr* const attribute
1230            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1231      if (!attribute)
1232      {
1233        G4Exception("G4GDMLReadSolids::GenTrapRead()",
1234                    "InvalidRead", FatalException, "No attribute found!");
1235        return;
1236      }
1237      const G4String attName = Transcode(attribute->getName());
1238      const G4String attValue = Transcode(attribute->getValue());
1239
1240      if (attName=="name") { name = GenerateName(attValue); } else
1241      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1242      if (attName=="dz") { dz = eval.Evaluate(attValue); } else
1243      if (attName=="v1x") { v1x = eval.Evaluate(attValue); } else
1244      if (attName=="v1y") { v1y = eval.Evaluate(attValue); } else
1245      if (attName=="v2x") { v2x = eval.Evaluate(attValue); } else
1246      if (attName=="v2y") { v2y = eval.Evaluate(attValue); } else
1247      if (attName=="v3x") { v3x = eval.Evaluate(attValue); } else
1248      if (attName=="v3y") { v3y = eval.Evaluate(attValue); } else
1249      if (attName=="v4x") { v4x = eval.Evaluate(attValue); } else
1250      if (attName=="v4y") { v4y = eval.Evaluate(attValue); } else
1251      if (attName=="v5x") { v5x = eval.Evaluate(attValue); } else
1252      if (attName=="v5y") { v5y = eval.Evaluate(attValue); } else
1253      if (attName=="v6x") { v6x = eval.Evaluate(attValue); } else
1254      if (attName=="v6y") { v6y = eval.Evaluate(attValue); } else
1255      if (attName=="v7x") { v7x = eval.Evaluate(attValue); } else
1256      if (attName=="v7y") { v7y = eval.Evaluate(attValue); } else
1257      if (attName=="v8x") { v8x = eval.Evaluate(attValue); } else
1258      if (attName=="v8y") { v8y = eval.Evaluate(attValue); }
1259   }
1260
1261   dz *= lunit;
1262   std::vector<G4TwoVector> vertices;
1263   vertices.push_back(G4TwoVector(v1x,v1y));
1264   vertices.push_back(G4TwoVector(v2x,v2y));
1265   vertices.push_back(G4TwoVector(v3x,v3y));
1266   vertices.push_back(G4TwoVector(v4x,v4y));
1267   vertices.push_back(G4TwoVector(v5x,v5y));
1268   vertices.push_back(G4TwoVector(v6x,v6y));
1269   vertices.push_back(G4TwoVector(v7x,v7y));
1270   vertices.push_back(G4TwoVector(v8x,v8y));
1271   new G4GenericTrap(name,dz,vertices);
1272}
1273
1274void G4GDMLReadSolids::TrapRead(const xercesc::DOMElement* const trapElement)
1275{
1276   G4String name;
1277   G4double lunit = 1.0;
1278   G4double aunit = 1.0;
1279   G4double z = 0.0;
1280   G4double theta = 0.0;
1281   G4double phi = 0.0;
1282   G4double y1 = 0.0;
1283   G4double x1 = 0.0;
1284   G4double x2 = 0.0;
1285   G4double alpha1 = 0.0;
1286   G4double y2 = 0.0;
1287   G4double x3 = 0.0;
1288   G4double x4 = 0.0;
1289   G4double alpha2 = 0.0;
1290
1291   const xercesc::DOMNamedNodeMap* const attributes
1292         = trapElement->getAttributes();
1293   XMLSize_t attributeCount = attributes->getLength();
1294
1295   for (XMLSize_t attribute_index=0;
1296        attribute_index<attributeCount; attribute_index++)
1297   {
1298      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1299
1300      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1301        { continue; }
1302
1303      const xercesc::DOMAttr* const attribute
1304            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1305      if (!attribute)
1306      {
1307        G4Exception("G4GDMLReadSolids::TrapRead()",
1308                    "InvalidRead", FatalException, "No attribute found!");
1309        return;
1310      }
1311      const G4String attName = Transcode(attribute->getName());
1312      const G4String attValue = Transcode(attribute->getValue());
1313
1314      if (attName=="name") { name = GenerateName(attValue); } else
1315      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1316      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1317      if (attName=="z") { z = eval.Evaluate(attValue); } else
1318      if (attName=="theta") { theta = eval.Evaluate(attValue); } else
1319      if (attName=="phi") { phi = eval.Evaluate(attValue); } else
1320      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1321      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1322      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1323      if (attName=="alpha1") { alpha1 = eval.Evaluate(attValue); } else
1324      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1325      if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1326      if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1327      if (attName=="alpha2") { alpha2 = eval.Evaluate(attValue); }
1328   }
1329
1330   z *= 0.5*lunit;
1331   theta *= aunit;
1332   phi *= aunit;
1333   y1 *= 0.5*lunit;
1334   x1 *= 0.5*lunit;
1335   x2 *= 0.5*lunit;
1336   alpha1 *= aunit;
1337   y2 *= 0.5*lunit;
1338   x3 *= 0.5*lunit;
1339   x4 *= 0.5*lunit;
1340   alpha2 *= aunit;
1341
1342   new G4Trap(name,z,theta,phi,y1,x1,x2,alpha1,y2,x3,x4,alpha2);
1343}
1344
1345void G4GDMLReadSolids::TrdRead(const xercesc::DOMElement* const trdElement)
1346{
1347   G4String name;
1348   G4double lunit = 1.0;
1349   G4double aunit = 1.0;
1350   G4double x1 = 0.0;
1351   G4double x2 = 0.0;
1352   G4double y1 = 0.0;
1353   G4double y2 = 0.0;
1354   G4double z = 0.0;
1355
1356   const xercesc::DOMNamedNodeMap* const attributes = trdElement->getAttributes();
1357   XMLSize_t attributeCount = attributes->getLength();
1358
1359   for (XMLSize_t attribute_index=0;
1360        attribute_index<attributeCount; attribute_index++)
1361   {
1362      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1363
1364      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1365        { continue; }
1366
1367      const xercesc::DOMAttr* const attribute
1368            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1369      if (!attribute)
1370      {
1371        G4Exception("G4GDMLReadSolids::TrdRead()",
1372                    "InvalidRead", FatalException, "No attribute found!");
1373        return;
1374      }
1375      const G4String attName = Transcode(attribute->getName());
1376      const G4String attValue = Transcode(attribute->getValue());
1377
1378      if (attName=="name") { name = GenerateName(attValue); } else
1379      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1380      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1381      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1382      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1383      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1384      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1385      if (attName=="z") { z = eval.Evaluate(attValue); }
1386   }
1387
1388   x1 *= 0.5*lunit;
1389   x2 *= 0.5*lunit;
1390   y1 *= 0.5*lunit;
1391   y2 *= 0.5*lunit;
1392   z *= 0.5*lunit;
1393
1394   new G4Trd(name,x1,x2,y1,y2,z);
1395}
1396
1397G4TriangularFacet* G4GDMLReadSolids::
1398TriangularRead(const xercesc::DOMElement* const triangularElement)
1399{
1400   G4ThreeVector vertex1;
1401   G4ThreeVector vertex2;
1402   G4ThreeVector vertex3;
1403   G4FacetVertexType type = ABSOLUTE;
1404
1405   const xercesc::DOMNamedNodeMap* const attributes
1406         = triangularElement->getAttributes();
1407   XMLSize_t attributeCount = attributes->getLength();
1408
1409   for (XMLSize_t attribute_index=0;
1410        attribute_index<attributeCount; attribute_index++)
1411   {
1412      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1413
1414      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1415        { continue; }
1416
1417      const xercesc::DOMAttr* const attribute
1418            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1419      if (!attribute)
1420      {
1421        G4Exception("G4GDMLReadSolids::TriangularRead()",
1422                    "InvalidRead", FatalException, "No attribute found!");
1423        return 0;
1424      }
1425      const G4String attName = Transcode(attribute->getName());
1426      const G4String attValue = Transcode(attribute->getValue());
1427
1428      if (attName=="vertex1")
1429        { vertex1 = GetPosition(GenerateName(attValue)); } else
1430      if (attName=="vertex2")
1431        { vertex2 = GetPosition(GenerateName(attValue)); } else
1432      if (attName=="vertex3")
1433        { vertex3 = GetPosition(GenerateName(attValue)); } else
1434      if (attName=="type")
1435        { if (attValue=="RELATIVE") { type = RELATIVE; } }
1436   }
1437
1438   return new G4TriangularFacet(vertex1,vertex2,vertex3,type);
1439}
1440
1441void G4GDMLReadSolids::TubeRead(const xercesc::DOMElement* const tubeElement)
1442{
1443   G4String name;
1444   G4double lunit = 1.0;
1445   G4double aunit = 1.0;
1446   G4double rmin = 0.0;
1447   G4double rmax = 0.0;
1448   G4double z = 0.0;
1449   G4double startphi = 0.0;
1450   G4double deltaphi = 0.0;
1451
1452   const xercesc::DOMNamedNodeMap* const attributes
1453         = tubeElement->getAttributes();
1454   XMLSize_t attributeCount = attributes->getLength();
1455
1456   for (XMLSize_t attribute_index=0;
1457        attribute_index<attributeCount; attribute_index++)
1458   {
1459      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1460
1461      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1462        { continue; }
1463
1464      const xercesc::DOMAttr* const attribute
1465            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1466      if (!attribute)
1467      {
1468        G4Exception("G4GDMLReadSolids::TubeRead()",
1469                    "InvalidRead", FatalException, "No attribute found!");
1470        return;
1471      }
1472      const G4String attName = Transcode(attribute->getName());
1473      const G4String attValue = Transcode(attribute->getValue());
1474
1475      if (attName=="name") { name = GenerateName(attValue); } else
1476      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1477      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1478      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1479      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1480      if (attName=="z") { z = eval.Evaluate(attValue); } else
1481      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1482      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1483   }
1484
1485   rmin *= lunit;
1486   rmax *= lunit;
1487   z *= 0.5*lunit;
1488   startphi *= aunit;
1489   deltaphi *= aunit;
1490
1491   new G4Tubs(name,rmin,rmax,z,startphi,deltaphi);
1492}
1493
1494void G4GDMLReadSolids::
1495TwistedboxRead(const xercesc::DOMElement* const twistedboxElement)
1496{
1497   G4String name;
1498   G4double lunit = 1.0;
1499   G4double aunit = 1.0;
1500   G4double PhiTwist = 0.0;
1501   G4double x = 0.0;
1502   G4double y = 0.0;
1503   G4double z = 0.0;
1504
1505   const xercesc::DOMNamedNodeMap* const attributes
1506         = twistedboxElement->getAttributes();
1507   XMLSize_t attributeCount = attributes->getLength();
1508
1509   for (XMLSize_t attribute_index=0;
1510        attribute_index<attributeCount; attribute_index++)
1511   {
1512      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1513
1514      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1515        { continue; }
1516
1517      const xercesc::DOMAttr* const attribute
1518            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1519      if (!attribute)
1520      {
1521        G4Exception("G4GDMLReadSolids::TwistedboxRead()",
1522                    "InvalidRead", FatalException, "No attribute found!");
1523        return;
1524      }
1525      const G4String attName = Transcode(attribute->getName());
1526      const G4String attValue = Transcode(attribute->getValue());
1527
1528      if (attName=="name") { name = GenerateName(attValue); } else
1529      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1530      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1531      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1532      if (attName=="x") { x = eval.Evaluate(attValue); } else
1533      if (attName=="y") { y = eval.Evaluate(attValue); } else
1534      if (attName=="z") { z = eval.Evaluate(attValue); }
1535   }
1536
1537   PhiTwist *= aunit;
1538   x *= 0.5*lunit;
1539   y *= 0.5*lunit;
1540   z *= 0.5*lunit;
1541
1542   new G4TwistedBox(name,PhiTwist,x,y,z);
1543}
1544
1545void G4GDMLReadSolids::
1546TwistedtrapRead(const xercesc::DOMElement* const twistedtrapElement)
1547{
1548   G4String name;
1549   G4double lunit = 1.0;
1550   G4double aunit = 1.0;
1551   G4double PhiTwist = 0.0;
1552   G4double z = 0.0;
1553   G4double Theta = 0.0;
1554   G4double Phi = 0.0;
1555   G4double y1 = 0.0;
1556   G4double x1 = 0.0;
1557   G4double x2 = 0.0;
1558   G4double y2 = 0.0;
1559   G4double x3 = 0.0;
1560   G4double x4 = 0.0;
1561   G4double Alph = 0.0;
1562
1563   const xercesc::DOMNamedNodeMap* const attributes
1564         = twistedtrapElement->getAttributes();
1565   XMLSize_t attributeCount = attributes->getLength();
1566
1567   for (XMLSize_t attribute_index=0;
1568        attribute_index<attributeCount; attribute_index++)
1569   {
1570      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1571
1572      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1573        { continue; }
1574
1575      const xercesc::DOMAttr* const attribute
1576            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1577      if (!attribute)
1578      {
1579        G4Exception("G4GDMLReadSolids::TwistedtrapRead()",
1580                    "InvalidRead", FatalException, "No attribute found!");
1581        return;
1582      }
1583      const G4String attName = Transcode(attribute->getName());
1584      const G4String attValue = Transcode(attribute->getValue());
1585
1586      if (attName=="name") { name = GenerateName(attValue); } else
1587      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1588      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1589      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1590      if (attName=="z") { z = eval.Evaluate(attValue); } else
1591      if (attName=="Theta") { Theta = eval.Evaluate(attValue); } else
1592      if (attName=="Phi") { Phi = eval.Evaluate(attValue); } else
1593      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1594      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1595      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1596      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1597      if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1598      if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1599      if (attName=="Alph") { Alph = eval.Evaluate(attValue); }
1600   }
1601
1602
1603   PhiTwist *= aunit;
1604   z *= 0.5*lunit;
1605   Theta *= aunit;
1606   Phi *= aunit;
1607   Alph *= aunit;
1608   y1 *= 0.5*lunit;
1609   x1 *= 0.5*lunit;
1610   x2 *= 0.5*lunit;
1611   y2 *= 0.5*lunit;
1612   x3 *= 0.5*lunit;
1613   x4 *= 0.5*lunit;
1614
1615   new G4TwistedTrap(name,PhiTwist,z,Theta,Phi,y1,x1,x2,y2,x3,x4,Alph);
1616}
1617
1618void G4GDMLReadSolids::
1619TwistedtrdRead(const xercesc::DOMElement* const twistedtrdElement)
1620{
1621   G4String name;
1622   G4double lunit = 1.0;
1623   G4double aunit = 1.0;
1624   G4double x1 = 0.0;
1625   G4double x2 = 0.0;
1626   G4double y1 = 0.0;
1627   G4double y2 = 0.0;
1628   G4double z = 0.0;
1629   G4double PhiTwist = 0.0;
1630
1631   const xercesc::DOMNamedNodeMap* const attributes
1632         = twistedtrdElement->getAttributes();
1633   XMLSize_t attributeCount = attributes->getLength();
1634
1635   for (XMLSize_t attribute_index=0;
1636        attribute_index<attributeCount; attribute_index++)
1637   {
1638      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1639
1640      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1641        { continue; }
1642
1643      const xercesc::DOMAttr* const attribute
1644            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1645      if (!attribute)
1646      {
1647        G4Exception("G4GDMLReadSolids::TwistedtrdRead()",
1648                    "InvalidRead", FatalException, "No attribute found!");
1649        return;
1650      }
1651      const G4String attName = Transcode(attribute->getName());
1652      const G4String attValue = Transcode(attribute->getValue());
1653
1654      if (attName=="name") { name = GenerateName(attValue); } else
1655      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1656      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1657      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1658      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1659      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1660      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1661      if (attName=="z") { z = eval.Evaluate(attValue); } else
1662      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); }
1663   }
1664
1665   x1 *= 0.5*lunit;
1666   x2 *= 0.5*lunit;
1667   y1 *= 0.5*lunit;
1668   y2 *= 0.5*lunit;
1669   z *= 0.5*lunit;
1670   PhiTwist *= aunit;
1671
1672   new G4TwistedTrd(name,x1,x2,y1,y2,z,PhiTwist);
1673}
1674
1675void G4GDMLReadSolids::
1676TwistedtubsRead(const xercesc::DOMElement* const twistedtubsElement)
1677{
1678   G4String name;
1679   G4double lunit = 1.0;
1680   G4double aunit = 1.0;
1681   G4double twistedangle = 0.0;
1682   G4double endinnerrad = 0.0;
1683   G4double endouterrad = 0.0;
1684   G4double zlen = 0.0;
1685   G4double phi = 0.0;
1686
1687   const xercesc::DOMNamedNodeMap* const attributes
1688         = twistedtubsElement->getAttributes();
1689   XMLSize_t attributeCount = attributes->getLength();
1690
1691   for (XMLSize_t attribute_index=0;
1692        attribute_index<attributeCount; attribute_index++)
1693   {
1694      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1695
1696      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1697        { continue; }
1698
1699      const xercesc::DOMAttr* const attribute
1700            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1701      if (!attribute)
1702      {
1703        G4Exception("G4GDMLReadSolids::TwistedtubsRead()",
1704                    "InvalidRead", FatalException, "No attribute found!");
1705        return;
1706      }
1707      const G4String attName = Transcode(attribute->getName());
1708      const G4String attValue = Transcode(attribute->getValue());
1709
1710      if (attName=="name") { name = GenerateName(attValue); } else
1711      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1712      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1713      if (attName=="twistedangle") { twistedangle=eval.Evaluate(attValue); } else
1714      if (attName=="endinnerrad")  { endinnerrad=eval.Evaluate(attValue);  } else
1715      if (attName=="endouterrad")  { endouterrad=eval.Evaluate(attValue);  } else
1716      if (attName=="zlen") { zlen = eval.Evaluate(attValue); } else
1717      if (attName=="phi") { phi = eval.Evaluate(attValue); }
1718   }
1719
1720   twistedangle *= aunit;
1721   endinnerrad *= lunit;
1722   endouterrad *= lunit;
1723   zlen *= 0.5*lunit;
1724   phi *= aunit;
1725
1726   new G4TwistedTubs(name,twistedangle,endinnerrad,endouterrad,zlen,phi);
1727}
1728
1729G4TwoVector G4GDMLReadSolids::
1730TwoDimVertexRead(const xercesc::DOMElement* const element, G4double lunit)
1731{
1732   G4TwoVector vec;
1733   
1734   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1735   XMLSize_t attributeCount = attributes->getLength();
1736
1737   for (XMLSize_t attribute_index=0;
1738        attribute_index<attributeCount; attribute_index++)
1739   {
1740      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1741
1742      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1743        { continue; }
1744
1745      const xercesc::DOMAttr* const attribute
1746            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1747      if (!attribute)
1748      {
1749        G4Exception("G4GDMLReadSolids::TwoDimVertexRead()",
1750                    "InvalidRead", FatalException, "No attribute found!");
1751        return vec;
1752      }
1753      const G4String attName = Transcode(attribute->getName());
1754      const G4String attValue = Transcode(attribute->getValue());
1755
1756      if (attName=="x") { vec.setX(eval.Evaluate(attValue)*lunit); } else
1757      if (attName=="y") { vec.setY(eval.Evaluate(attValue)*lunit); }
1758   }
1759
1760   return vec;
1761}
1762
1763G4GDMLReadSolids::zplaneType G4GDMLReadSolids::
1764ZplaneRead(const xercesc::DOMElement* const zplaneElement)
1765{
1766   zplaneType zplane = {0.,0.,0.};
1767
1768   const xercesc::DOMNamedNodeMap* const attributes
1769         = zplaneElement->getAttributes();
1770   XMLSize_t attributeCount = attributes->getLength();
1771
1772   for (XMLSize_t attribute_index=0;
1773        attribute_index<attributeCount; attribute_index++)
1774   {
1775      xercesc::DOMNode* node = attributes->item(attribute_index);
1776
1777      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
1778
1779      const xercesc::DOMAttr* const attribute
1780            = dynamic_cast<xercesc::DOMAttr*>(node);   
1781      if (!attribute)
1782      {
1783        G4Exception("G4GDMLReadSolids::ZplaneRead()",
1784                    "InvalidRead", FatalException, "No attribute found!");
1785        return zplane;
1786      }
1787      const G4String attName = Transcode(attribute->getName());
1788      const G4String attValue = Transcode(attribute->getValue());
1789
1790      if (attName=="rmin") { zplane.rmin = eval.Evaluate(attValue); } else
1791      if (attName=="rmax") { zplane.rmax = eval.Evaluate(attValue); } else
1792      if (attName=="z") { zplane.z = eval.Evaluate(attValue); }
1793   }
1794
1795   return zplane;
1796}
1797
1798void G4GDMLReadSolids::
1799OpticalSurfaceRead(const xercesc::DOMElement* const opticalsurfaceElement)
1800{
1801   G4String name;
1802   G4String smodel;
1803   G4String sfinish;
1804   G4String stype;
1805   G4double value = 0.0;
1806
1807   const xercesc::DOMNamedNodeMap* const attributes
1808         = opticalsurfaceElement->getAttributes();
1809   XMLSize_t attributeCount = attributes->getLength();
1810
1811   for (XMLSize_t attribute_index=0;
1812        attribute_index<attributeCount; attribute_index++)
1813   {
1814      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1815
1816      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1817        { continue; }
1818
1819      const xercesc::DOMAttr* const attribute
1820            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1821      if (!attribute)
1822      {
1823        G4Exception("G4GDMLReadSolids::OpticalSurfaceRead()",
1824                    "InvalidRead", FatalException, "No attribute found!");
1825        return;
1826      }
1827      const G4String attName = Transcode(attribute->getName());
1828      const G4String attValue = Transcode(attribute->getValue());
1829
1830      if (attName=="name") { name = GenerateName(attValue); } else
1831      if (attName=="model") { smodel = attValue; } else
1832      if (attName=="finish") { sfinish = attValue; } else
1833      if (attName=="type") { stype = attValue; } else
1834      if (attName=="value") { value = eval.Evaluate(attValue); }
1835   }
1836
1837   G4OpticalSurfaceModel model; 
1838   G4OpticalSurfaceFinish finish;
1839   G4SurfaceType type;   
1840   
1841   if ((smodel=="glisur") || (smodel=="0")) { model = glisur; } else
1842   if ((smodel=="unified") || (smodel=="1")) { model = unified; }
1843   else { model = LUT; }
1844
1845   if ((sfinish=="polished") || (sfinish=="0"))
1846      { finish = polished; } else
1847   if ((sfinish=="polishedfrontpainted") || (sfinish=="1"))
1848      { finish = polishedfrontpainted; } else
1849   if ((sfinish=="polishedbackpainted") || (sfinish=="2"))
1850      { finish = polishedbackpainted; } else
1851   if ((sfinish=="ground") || (sfinish=="3"))
1852      { finish = ground; } else
1853   if ((sfinish=="groundfrontpainted") || (sfinish=="4"))
1854      { finish = groundfrontpainted; } else
1855   if ((sfinish=="groundbackpainted") || (sfinish=="5"))
1856      { finish = groundbackpainted; } else
1857   if ((sfinish=="polishedlumirrorair") || (sfinish=="6"))
1858      { finish = polishedlumirrorair; } else
1859   if ((sfinish=="polishedlumirrorglue") || (sfinish=="7"))
1860      { finish = polishedlumirrorglue; } else
1861   if ((sfinish=="polishedair") || (sfinish=="8"))
1862      { finish = polishedair; } else
1863   if ((sfinish=="polishedteflonair") || (sfinish=="9"))
1864      { finish = polishedteflonair; } else
1865   if ((sfinish=="polishedtioair") || (sfinish=="10"))
1866      { finish = polishedtioair; } else
1867   if ((sfinish=="polishedtyvekair") || (sfinish=="11"))
1868      { finish = polishedtyvekair; } else
1869   if ((sfinish=="polishedvm2000air") || (sfinish=="12"))
1870      { finish = polishedvm2000air; } else
1871   if ((sfinish=="polishedvm2000glue") || (sfinish=="13"))
1872      { finish = polishedvm2000glue; } else
1873   if ((sfinish=="etchedlumirrorair") || (sfinish=="14"))
1874      { finish = etchedlumirrorair; } else
1875   if ((sfinish=="etchedlumirrorglue") || (sfinish=="15"))
1876      { finish = etchedlumirrorglue; } else
1877   if ((sfinish=="etchedair") || (sfinish=="16"))
1878      { finish = etchedair; } else
1879   if ((sfinish=="etchedteflonair") || (sfinish=="17"))
1880      { finish = etchedteflonair; } else
1881   if ((sfinish=="etchedtioair") || (sfinish=="18"))
1882      { finish = etchedtioair; } else
1883   if ((sfinish=="etchedtyvekair") || (sfinish=="19"))
1884      { finish = etchedtyvekair; } else
1885   if ((sfinish=="etchedvm2000air") || (sfinish=="20"))
1886      { finish = etchedvm2000air; } else
1887   if ((sfinish=="etchedvm2000glue") || (sfinish=="21"))
1888      { finish = etchedvm2000glue; } else
1889   if ((sfinish=="groundlumirrorair") || (sfinish=="22"))
1890      { finish = groundlumirrorair; } else
1891   if ((sfinish=="groundlumirrorglue") || (sfinish=="23"))
1892      { finish = groundlumirrorglue; } else
1893   if ((sfinish=="groundair") || (sfinish=="24"))
1894      { finish = groundair; } else
1895   if ((sfinish=="groundteflonair") || (sfinish=="25"))
1896      { finish = groundteflonair; } else
1897   if ((sfinish=="groundtioair") || (sfinish=="26"))
1898      { finish = groundtioair; } else
1899   if ((sfinish=="groundtyvekair") || (sfinish=="27"))
1900      { finish = groundtyvekair; } else
1901   if ((sfinish=="groundvm2000air") || (sfinish=="28"))
1902      { finish = groundvm2000air; }
1903   else { finish = groundvm2000glue; }
1904
1905   if ((stype=="dielectric_metal") || (stype=="0"))
1906      { type = dielectric_metal; } else
1907   if ((stype=="dielectric_dielectric") || (stype=="1"))
1908      { type = dielectric_dielectric; } else
1909   if ((stype=="dielectric_LUT") || (stype=="2"))
1910      { type = dielectric_LUT; } else
1911   if ((stype=="firsov") || (stype=="3"))
1912      { type = firsov; }
1913   else { type = x_ray; }
1914
1915   new G4OpticalSurface(name,model,finish,type,value);
1916}
1917
1918void G4GDMLReadSolids::SolidsRead(const xercesc::DOMElement* const solidsElement)
1919{
1920   G4cout << "G4GDML: Reading solids..." << G4endl;
1921
1922   for (xercesc::DOMNode* iter = solidsElement->getFirstChild();
1923        iter != 0; iter = iter->getNextSibling())
1924   {
1925      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)  { continue; }
1926
1927      const xercesc::DOMElement* const child
1928            = dynamic_cast<xercesc::DOMElement*>(iter);
1929      if (!child)
1930      {
1931        G4Exception("G4GDMLReadSolids::SolidsRead()",
1932                    "InvalidRead", FatalException, "No child found!");
1933        return;
1934      }
1935      const G4String tag = Transcode(child->getTagName());
1936      if (tag=="define") { DefineRead(child);  }  else 
1937      if (tag=="box")    { BoxRead(child); } else
1938      if (tag=="cone")   { ConeRead(child); } else
1939      if (tag=="elcone") { ElconeRead(child); } else
1940      if (tag=="ellipsoid") { EllipsoidRead(child); }else
1941      if (tag=="eltube") { EltubeRead(child); } else
1942      if (tag=="xtru") { XtruRead(child); } else
1943      if (tag=="hype") { HypeRead(child); } else
1944      if (tag=="intersection") { BooleanRead(child,INTERSECTION); } else
1945      if (tag=="orb") { OrbRead(child); } else
1946      if (tag=="para") { ParaRead(child); } else
1947      if (tag=="paraboloid") { ParaboloidRead(child); } else
1948      if (tag=="polycone") { PolyconeRead(child); } else
1949      if (tag=="polyhedra") { PolyhedraRead(child); } else
1950      if (tag=="reflectedSolid") { ReflectedSolidRead(child); } else
1951      if (tag=="sphere") { SphereRead(child); } else
1952      if (tag=="subtraction") { BooleanRead(child,SUBTRACTION); } else
1953      if (tag=="tessellated") { TessellatedRead(child); } else
1954      if (tag=="tet") { TetRead(child); } else
1955      if (tag=="torus") { TorusRead(child); } else
1956      if (tag=="arb8") { GenTrapRead(child); } else
1957      if (tag=="trap") { TrapRead(child); } else
1958      if (tag=="trd") { TrdRead(child); } else
1959      if (tag=="tube") { TubeRead(child); } else
1960      if (tag=="twistedbox") { TwistedboxRead(child); } else
1961      if (tag=="twistedtrap") { TwistedtrapRead(child); } else
1962      if (tag=="twistedtrd") { TwistedtrdRead(child); } else
1963      if (tag=="twistedtubs") { TwistedtubsRead(child); } else
1964      if (tag=="union") { BooleanRead(child,UNION); } else
1965      if (tag=="opticalsurface") { OpticalSurfaceRead(child); } else
1966      if (tag=="loop") { LoopRead(child,&G4GDMLRead::SolidsRead); }
1967      else
1968      {
1969        G4String error_msg = "Unknown tag in solids: " + tag;
1970        G4Exception("G4GDMLReadSolids::SolidsRead()", "ReadError",
1971                    FatalException, error_msg);
1972      }
1973   }
1974}
1975
1976G4VSolid* G4GDMLReadSolids::GetSolid(const G4String& ref) const
1977{
1978   G4VSolid* solidPtr = G4SolidStore::GetInstance()->GetSolid(ref,false);
1979
1980   if (!solidPtr)
1981   {
1982     G4String error_msg = "Referenced solid '" + ref + "' was not found!";
1983     G4Exception("G4GDMLReadSolids::GetSolid()", "ReadError",
1984                 FatalException, error_msg);
1985   }
1986
1987   return solidPtr;
1988}
1989
1990G4SurfaceProperty* G4GDMLReadSolids::
1991GetSurfaceProperty(const G4String& ref) const
1992{
1993   const G4SurfacePropertyTable* surfaceList
1994         = G4SurfaceProperty::GetSurfacePropertyTable();
1995   const size_t surfaceCount = surfaceList->size();
1996
1997   for (size_t i=0; i<surfaceCount; i++)
1998   {
1999      if ((*surfaceList)[i]->GetName() == ref)  { return (*surfaceList)[i]; }
2000   }
2001
2002   G4String error_msg = "Referenced optical surface '" + ref + "' was not found!";
2003   G4Exception("G4GDMLReadSolids::GetSurfaceProperty()", "ReadError",
2004               FatalException, error_msg);
2005
2006   return 0;
2007}
Note: See TracBrowser for help on using the repository browser.