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

Last change on this file since 1160 was 987, checked in by garnier, 15 years ago

fichiers manquants

File size: 56.6 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.22 2008/11/21 09:32:46 gcosmo Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// class G4GDMLReadSolids Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
35#include "G4GDMLReadSolids.hh"
36
37void G4GDMLReadSolids::
38BooleanRead(const xercesc::DOMElement* const booleanElement, const BooleanOp op)
39{
40   G4String name;
41   G4String first;
42   G4String second;
43   G4ThreeVector position(0.0,0.0,0.0);
44   G4ThreeVector rotation(0.0,0.0,0.0);
45   G4ThreeVector firstposition(0.0,0.0,0.0);
46   G4ThreeVector firstrotation(0.0,0.0,0.0);
47
48   const xercesc::DOMNamedNodeMap* const attributes
49         = booleanElement->getAttributes();
50   XMLSize_t attributeCount = attributes->getLength();
51
52   for (XMLSize_t attribute_index=0;
53        attribute_index<attributeCount; attribute_index++)
54   {
55      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
56
57      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
58        { continue; }
59
60      const xercesc::DOMAttr* const attribute
61            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
62      const G4String attName = Transcode(attribute->getName());
63      const G4String attValue = Transcode(attribute->getValue());
64
65      if (attName=="name")  { name = GenerateName(attValue); }
66   }
67
68   for (xercesc::DOMNode* iter = booleanElement->getFirstChild();
69        iter != 0;iter = iter->getNextSibling())
70   {
71      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)  { continue; }
72
73      const xercesc::DOMElement* const child
74            = dynamic_cast<xercesc::DOMElement*>(iter);
75      const G4String tag = Transcode(child->getTagName());
76
77      if (tag=="first") { first = RefRead(child); } else
78      if (tag=="second") { second = RefRead(child); } else
79      if (tag=="position") { VectorRead(child,position); } else
80      if (tag=="rotation") { VectorRead(child,rotation); } else
81      if (tag=="positionref")
82        { position = GetPosition(GenerateName(RefRead(child))); } else
83      if (tag=="rotationref")
84        { rotation = GetRotation(GenerateName(RefRead(child))); } else
85      if (tag=="firstposition") { VectorRead(child,firstposition); } else
86      if (tag=="firstrotation") { VectorRead(child,firstrotation); } else
87      if (tag=="firstpositionref")
88        { firstposition = GetPosition(GenerateName(RefRead(child))); } else
89      if (tag=="firstrotationref")
90        { firstrotation = GetRotation(GenerateName(RefRead(child))); } 
91      else
92      {
93        G4String error_msg = "Unknown tag in boolean solid: " + tag;
94        G4Exception("G4GDMLReadSolids::BooleanRead()", "ReadError",
95                    FatalException, error_msg);
96      }
97   }
98
99   G4VSolid* firstSolid = GetSolid(GenerateName(first));
100   G4VSolid* secondSolid = GetSolid(GenerateName(second));
101
102   G4Transform3D transform(GetRotationMatrix(rotation),position);
103
104   if (( (firstrotation.x()!=0.0) || (firstrotation.y()!=0.0)
105                                  || (firstrotation.z()!=0.0))
106    || ( (firstposition.x()!=0.0) || (firstposition.y()!=0.0)
107                                  || (firstposition.z()!=0.0)))
108   { 
109      G4Transform3D firsttransform(GetRotationMatrix(firstrotation),
110                                   firstposition);
111      firstSolid = new G4DisplacedSolid(GenerateName("displaced_"+first),
112                                        firstSolid, firsttransform);
113   }
114
115   if (op==UNION)
116     { new G4UnionSolid(name,firstSolid,secondSolid,transform); } else
117   if (op==SUBTRACTION)
118     { new G4SubtractionSolid(name,firstSolid,secondSolid,transform); } else
119   if (op==INTERSECTION)
120     { new G4IntersectionSolid(name,firstSolid,secondSolid,transform); }
121}
122
123void G4GDMLReadSolids::BoxRead(const xercesc::DOMElement* const boxElement)
124{
125   G4String name;
126   G4double lunit = 1.0;
127   G4double aunit = 1.0;
128   G4double x = 0.0;
129   G4double y = 0.0;
130   G4double z = 0.0;
131
132   const xercesc::DOMNamedNodeMap* const attributes
133         = boxElement->getAttributes();
134   XMLSize_t attributeCount = attributes->getLength();
135
136   for (XMLSize_t attribute_index=0;
137        attribute_index<attributeCount; attribute_index++)
138   {
139      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
140
141      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
142        { continue; }
143
144      const xercesc::DOMAttr* const attribute
145            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
146      const G4String attName = Transcode(attribute->getName());
147      const G4String attValue = Transcode(attribute->getValue());
148
149      if (attName=="name") { name = GenerateName(attValue); } else
150      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
151      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
152      if (attName=="x") { x = eval.Evaluate(attValue); } else
153      if (attName=="y") { y = eval.Evaluate(attValue); } else
154      if (attName=="z") { z = eval.Evaluate(attValue); }
155   }
156
157   x *= 0.5*lunit;
158   y *= 0.5*lunit;
159   z *= 0.5*lunit;
160
161   new G4Box(name,x,y,z);
162}
163
164void G4GDMLReadSolids::ConeRead(const xercesc::DOMElement* const coneElement)
165{
166   G4String name;
167   G4double lunit = 1.0;
168   G4double aunit = 1.0;
169   G4double rmin1 = 0.0;
170   G4double rmax1 = 0.0;
171   G4double rmin2 = 0.0;
172   G4double rmax2 = 0.0;
173   G4double z = 0.0;
174   G4double startphi = 0.0;
175   G4double deltaphi = 0.0;
176
177   const xercesc::DOMNamedNodeMap* const attributes
178         = coneElement->getAttributes();
179   XMLSize_t attributeCount = attributes->getLength();
180
181   for (XMLSize_t attribute_index=0;
182        attribute_index<attributeCount; attribute_index++)
183   {
184      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
185
186      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
187        { continue; }
188
189      const xercesc::DOMAttr* const attribute
190            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
191      const G4String attName = Transcode(attribute->getName());
192      const G4String attValue = Transcode(attribute->getValue());
193
194      if (attName=="name") { name = GenerateName(attValue); } else
195      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
196      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
197      if (attName=="rmin1") { rmin1 = eval.Evaluate(attValue); } else
198      if (attName=="rmax1") { rmax1 = eval.Evaluate(attValue); } else
199      if (attName=="rmin2") { rmin2 = eval.Evaluate(attValue); } else
200      if (attName=="rmax2") { rmax2 = eval.Evaluate(attValue); } else
201      if (attName=="z") { z = eval.Evaluate(attValue); } else
202      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
203      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
204   }
205
206   rmin1 *= lunit;
207   rmax1 *= lunit;
208   rmin2 *= lunit;
209   rmax2 *= lunit;
210   z *= 0.5*lunit;
211   startphi *= aunit;
212   deltaphi *= aunit;
213
214   new G4Cons(name,rmin1,rmax1,rmin2,rmax2,z,startphi,deltaphi);
215}
216
217void G4GDMLReadSolids::
218ElconeRead(const xercesc::DOMElement* const elconeElement)
219{
220   G4String name;
221   G4double lunit = 1.0;
222   G4double dx = 0.0;
223   G4double dy = 0.0;
224   G4double zmax = 0.0;
225   G4double zcut = 0.0;
226
227   const xercesc::DOMNamedNodeMap* const attributes
228         = elconeElement->getAttributes();
229   XMLSize_t attributeCount = attributes->getLength();
230
231   for (XMLSize_t attribute_index=0;
232        attribute_index<attributeCount; attribute_index++)
233   {
234      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
235
236      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
237        { continue; }
238
239      const xercesc::DOMAttr* const attribute
240            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
241      const G4String attName = Transcode(attribute->getName());
242      const G4String attValue = Transcode(attribute->getValue());
243
244      if (attName=="name") { name = GenerateName(attValue); } else
245      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
246      if (attName=="dx") { dx = eval.Evaluate(attValue); } else
247      if (attName=="dy") { dy = eval.Evaluate(attValue); } else
248      if (attName=="zmax") { zmax = eval.Evaluate(attValue); } else
249      if (attName=="zcut") { zcut = eval.Evaluate(attValue); }
250   }
251
252   dx *= lunit;
253   dy *= lunit;
254   zmax *= lunit;
255   zcut *= lunit;
256
257   new G4EllipticalCone(name,dx,dy,zmax,zcut);
258}
259
260void G4GDMLReadSolids::
261EllipsoidRead(const xercesc::DOMElement* const ellipsoidElement)
262{
263   G4String name;
264   G4double lunit = 1.0;
265   G4double ax = 0.0;
266   G4double by = 0.0;
267   G4double cz = 0.0;
268   G4double zcut1 = 0.0;
269   G4double zcut2 = 0.0; 
270
271   const xercesc::DOMNamedNodeMap* const attributes
272         = ellipsoidElement->getAttributes();
273   XMLSize_t attributeCount = attributes->getLength();
274
275   for (XMLSize_t attribute_index=0;
276        attribute_index<attributeCount; attribute_index++)
277   {
278      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
279
280      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
281        { continue; }
282
283      const xercesc::DOMAttr* const attribute
284            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
285      const G4String attName = Transcode(attribute->getName());
286      const G4String attValue = Transcode(attribute->getValue());
287
288      if (attName=="name") { name  = GenerateName(attValue); } else
289      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
290      if (attName=="ax") { ax = eval.Evaluate(attValue); } else
291      if (attName=="by") { by = eval.Evaluate(attValue); } else
292      if (attName=="cz") { cz = eval.Evaluate(attValue); } else
293      if (attName=="zcut1") { zcut1 = eval.Evaluate(attValue); } else
294      if (attName=="zcut2") { zcut2 = eval.Evaluate(attValue); }
295   }
296
297   ax *= lunit;
298   by *= lunit;
299   cz *= lunit;
300   zcut1 *= lunit;
301   zcut2 *= lunit; 
302
303   new G4Ellipsoid(name,ax,by,cz,zcut1,zcut2);
304}
305
306void G4GDMLReadSolids::
307EltubeRead(const xercesc::DOMElement* const eltubeElement)
308{
309   G4String name;
310   G4double lunit = 1.0;
311   G4double dx = 0.0;
312   G4double dy = 0.0;
313   G4double dz = 0.0;
314
315   const xercesc::DOMNamedNodeMap* const attributes
316         = eltubeElement->getAttributes();
317   XMLSize_t attributeCount = attributes->getLength();
318
319   for (XMLSize_t attribute_index=0;
320        attribute_index<attributeCount; attribute_index++)
321   {
322      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
323
324      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
325        { continue; }
326
327      const xercesc::DOMAttr* const attribute
328            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
329      const G4String attName = Transcode(attribute->getName());
330      const G4String attValue = Transcode(attribute->getValue());
331
332      if (attName=="name") { name = GenerateName(attValue); } else
333      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
334      if (attName=="dx") { dx = eval.Evaluate(attValue); } else
335      if (attName=="dy") { dy = eval.Evaluate(attValue); } else
336      if (attName=="dz") { dz = eval.Evaluate(attValue); }
337   }
338
339   dx *= lunit;
340   dy *= lunit;
341   dz *= lunit;
342
343   new G4EllipticalTube(name,dx,dy,dz);
344}
345
346void G4GDMLReadSolids::XtruRead(const xercesc::DOMElement* const xtruElement)
347{
348   G4String name;
349   G4double lunit = 1.0;
350
351   const xercesc::DOMNamedNodeMap* const attributes
352         = xtruElement->getAttributes();
353   XMLSize_t attributeCount = attributes->getLength();
354
355   for (XMLSize_t attribute_index=0;
356        attribute_index<attributeCount; attribute_index++)
357   {
358      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
359
360      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
361        { continue; } 
362
363      const xercesc::DOMAttr* const attribute
364            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
365      const G4String attName = Transcode(attribute->getName());
366      const G4String attValue = Transcode(attribute->getValue());
367
368      if (attName=="name") { name = GenerateName(attValue); } else
369      if (attName=="lunit") { lunit = eval.Evaluate(attValue); }
370   }
371
372   std::vector<G4TwoVector> twoDimVertexList;
373   std::vector<G4ExtrudedSolid::ZSection> sectionList;
374
375   for (xercesc::DOMNode* iter = xtruElement->getFirstChild();
376        iter != 0; iter = iter->getNextSibling())
377   {
378      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
379
380      const xercesc::DOMElement* const child
381            = dynamic_cast<xercesc::DOMElement*>(iter);
382      const G4String tag = Transcode(child->getTagName());
383
384      if (tag=="twoDimVertex")
385        { twoDimVertexList.push_back(TwoDimVertexRead(child,lunit)); } else
386      if (tag=="section")
387        { sectionList.push_back(SectionRead(child,lunit)); }
388   }
389
390   new G4ExtrudedSolid(name,twoDimVertexList,sectionList);
391}
392
393void G4GDMLReadSolids::HypeRead(const xercesc::DOMElement* const hypeElement)
394{
395   G4String name;
396   G4double lunit = 1.0;
397   G4double aunit = 1.0;
398   G4double rmin = 0.0;
399   G4double rmax = 0.0;
400   G4double inst = 0.0;
401   G4double outst = 0.0;
402   G4double z = 0.0;
403
404   const xercesc::DOMNamedNodeMap* const attributes
405         = hypeElement->getAttributes();
406   XMLSize_t attributeCount = attributes->getLength();
407
408   for (XMLSize_t attribute_index=0;
409        attribute_index<attributeCount; attribute_index++)
410   {
411      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
412
413      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
414        { continue; }
415
416      const xercesc::DOMAttr* const attribute
417            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
418      const G4String attName = Transcode(attribute->getName());
419      const G4String attValue = Transcode(attribute->getValue());
420
421      if (attName=="name") { name = GenerateName(attValue); } else
422      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
423      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
424      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
425      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
426      if (attName=="inst") { inst = eval.Evaluate(attValue); } else
427      if (attName=="outst") { outst = eval.Evaluate(attValue); } else
428      if (attName=="z") { z = eval.Evaluate(attValue); }
429   }
430
431   rmin *= lunit;
432   rmax *= lunit;
433   inst *= aunit;
434   outst *= aunit;
435   z *= 0.5*lunit;
436
437   new G4Hype(name,rmin,rmax,inst,outst,z);
438}
439
440void G4GDMLReadSolids::OrbRead(const xercesc::DOMElement* const orbElement)
441{
442   G4String name;
443   G4double lunit = 1.0;
444   G4double r = 0.0;
445
446   const xercesc::DOMNamedNodeMap* const attributes
447         = orbElement->getAttributes();
448   XMLSize_t attributeCount = attributes->getLength();
449
450   for (XMLSize_t attribute_index=0;
451        attribute_index<attributeCount; attribute_index++)
452   {
453      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
454
455      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
456        { continue; }
457
458      const xercesc::DOMAttr* const attribute
459            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
460      const G4String attName = Transcode(attribute->getName());
461      const G4String attValue = Transcode(attribute->getValue());
462
463      if (attName=="name") { name = GenerateName(attValue); } else
464      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
465      if (attName=="r") { r = eval.Evaluate(attValue); }
466   }
467
468   r *= lunit;
469
470   new G4Orb(name,r);
471}
472
473void G4GDMLReadSolids::ParaRead(const xercesc::DOMElement* const paraElement)
474{
475   G4String name;
476   G4double lunit = 1.0;
477   G4double aunit = 1.0;
478   G4double x = 0.0;
479   G4double y = 0.0;
480   G4double z = 0.0;
481   G4double alpha = 0.0;
482   G4double theta = 0.0;
483   G4double phi = 0.0;
484
485   const xercesc::DOMNamedNodeMap* const attributes
486         = paraElement->getAttributes();
487   XMLSize_t attributeCount = attributes->getLength();
488
489   for (XMLSize_t attribute_index=0;
490        attribute_index<attributeCount; attribute_index++)
491   {
492      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
493
494      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
495        { continue; }
496
497      const xercesc::DOMAttr* const attribute
498            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
499      const G4String attName = Transcode(attribute->getName());
500      const G4String attValue = Transcode(attribute->getValue());
501
502      if (attName=="name") { name = GenerateName(attValue); } else
503      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
504      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
505      if (attName=="x") { x = eval.Evaluate(attValue); } else
506      if (attName=="y") { y = eval.Evaluate(attValue); } else
507      if (attName=="z") { z = eval.Evaluate(attValue); } else
508      if (attName=="alpha") { alpha = eval.Evaluate(attValue); } else
509      if (attName=="theta") { theta = eval.Evaluate(attValue); } else
510      if (attName=="phi") { phi = eval.Evaluate(attValue); }
511   }
512
513   x *= 0.5*lunit;
514   y *= 0.5*lunit;
515   z *= 0.5*lunit;
516   alpha *= aunit;
517   theta *= aunit;
518   phi *= aunit;
519
520   new G4Para(name,x,y,z,alpha,theta,phi);
521}
522
523void G4GDMLReadSolids::
524ParaboloidRead(const xercesc::DOMElement* const paraElement)
525{
526   G4String name;
527   G4double lunit = 1.0;
528   G4double aunit = 1.0;
529   G4double rlo = 0.0;
530   G4double rhi = 0.0;
531   G4double dz = 0.0;
532 
533   const xercesc::DOMNamedNodeMap* const attributes
534         = paraElement->getAttributes();
535   XMLSize_t attributeCount = attributes->getLength();
536
537   for (XMLSize_t attribute_index=0;
538        attribute_index<attributeCount; attribute_index++)
539   {
540      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
541
542      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
543        { continue; }
544
545      const xercesc::DOMAttr* const attribute
546            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
547      const G4String attName = Transcode(attribute->getName());
548      const G4String attValue = Transcode(attribute->getValue());
549
550      if (attName=="name")  { name = GenerateName(attValue); } else
551      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
552      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
553      if (attName=="rlo")   { rlo =  eval.Evaluate(attValue); } else
554      if (attName=="rhi")   { rhi = eval.Evaluate(attValue); } else
555      if (attName=="dz")    { dz = eval.Evaluate(attValue); } 
556   }     
557
558   rlo *= 1.*lunit;
559   rhi *= 1.*lunit;
560   dz *= 1.*lunit;
561 
562   new G4Paraboloid(name,dz,rlo,rhi);
563}
564
565void G4GDMLReadSolids::
566PolyconeRead(const xercesc::DOMElement* const polyconeElement) 
567{
568   G4String name;
569   G4double lunit = 1.0;
570   G4double aunit = 1.0;
571   G4double startphi = 0.0;
572   G4double deltaphi = 0.0;
573
574   const xercesc::DOMNamedNodeMap* const attributes
575         = polyconeElement->getAttributes();
576   XMLSize_t attributeCount = attributes->getLength();
577
578   for (XMLSize_t attribute_index=0;
579        attribute_index<attributeCount; attribute_index++)
580   {
581      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
582
583      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
584        { continue; }
585
586      const xercesc::DOMAttr* const attribute
587            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
588      const G4String attName = Transcode(attribute->getName());
589      const G4String attValue = Transcode(attribute->getValue());
590
591      if (attName=="name") { name = GenerateName(attValue); } else
592      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
593      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
594      if (attName=="startphi") { startphi = eval.Evaluate(attValue); }else
595      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
596   }
597
598   startphi *= aunit;
599   deltaphi *= aunit;
600
601   std::vector<zplaneType> zplaneList;
602
603   for (xercesc::DOMNode* iter = polyconeElement->getFirstChild();
604        iter != 0; iter = iter->getNextSibling())
605   {
606      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
607
608      const xercesc::DOMElement* const child
609            = dynamic_cast<xercesc::DOMElement*>(iter);
610      const G4String tag = Transcode(child->getTagName());
611
612      if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
613   }
614
615   G4int numZPlanes = zplaneList.size();
616
617   G4double* rmin_array = new G4double[numZPlanes];
618   G4double* rmax_array = new G4double[numZPlanes];
619   G4double* z_array    = new G4double[numZPlanes];
620
621   for (G4int i=0; i<numZPlanes; i++)
622   { 
623      rmin_array[i] = zplaneList[i].rmin*lunit;
624      rmax_array[i] = zplaneList[i].rmax*lunit;
625      z_array[i]    = zplaneList[i].z*lunit;
626   }
627
628   new G4Polycone(name,startphi,deltaphi,numZPlanes,
629                  z_array,rmin_array,rmax_array);
630}
631
632void G4GDMLReadSolids::
633PolyhedraRead(const xercesc::DOMElement* const polyhedraElement)
634{
635   G4String name;
636   G4double lunit = 1.0;
637   G4double aunit = 1.0;
638   G4double startphi = 0.0;
639   G4double deltaphi = 0.0;
640   G4int numsides = 0;
641
642   const xercesc::DOMNamedNodeMap* const attributes
643         = polyhedraElement->getAttributes();
644   XMLSize_t attributeCount = attributes->getLength();
645
646   for (XMLSize_t attribute_index=0;
647        attribute_index<attributeCount; attribute_index++)
648   {
649      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
650
651      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
652        { continue; }
653
654      const xercesc::DOMAttr* const attribute
655            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
656      const G4String attName = Transcode(attribute->getName());
657      const G4String attValue = Transcode(attribute->getValue());
658
659      if (attName=="name") { name = GenerateName(attValue); } else
660      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
661      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
662      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
663      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
664      if (attName=="numsides") { numsides = eval.EvaluateInteger(attValue); }
665   }
666
667   startphi *= aunit;
668   deltaphi *= aunit;
669
670   std::vector<zplaneType> zplaneList;
671
672   for (xercesc::DOMNode* iter = polyhedraElement->getFirstChild();
673        iter != 0; iter = iter->getNextSibling())
674   {
675      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)  { continue; }
676
677      const xercesc::DOMElement* const child
678            = dynamic_cast<xercesc::DOMElement*>(iter);
679      const G4String tag = Transcode(child->getTagName());
680
681      if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
682   }
683
684   G4int numZPlanes = zplaneList.size();
685
686   G4double* rmin_array = new G4double[numZPlanes];
687   G4double* rmax_array = new G4double[numZPlanes];
688   G4double* z_array = new G4double[numZPlanes];
689
690   for (G4int i=0; i<numZPlanes; i++)
691   { 
692      rmin_array[i] = zplaneList[i].rmin*lunit;
693      rmax_array[i] = zplaneList[i].rmax*lunit;
694      z_array[i] = zplaneList[i].z*lunit;
695   }
696
697   new G4Polyhedra(name,startphi,deltaphi,numsides,numZPlanes,
698                   z_array,rmin_array,rmax_array);
699
700}
701
702G4QuadrangularFacet* G4GDMLReadSolids::
703QuadrangularRead(const xercesc::DOMElement* const quadrangularElement)
704{
705   G4ThreeVector vertex1;
706   G4ThreeVector vertex2;
707   G4ThreeVector vertex3;
708   G4ThreeVector vertex4;
709   G4FacetVertexType type = ABSOLUTE;
710
711   const xercesc::DOMNamedNodeMap* const attributes
712         = quadrangularElement->getAttributes();
713   XMLSize_t attributeCount = attributes->getLength();
714
715   for (XMLSize_t attribute_index=0;
716        attribute_index<attributeCount; attribute_index++)
717   {
718      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
719
720      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
721        { continue; }
722
723      const xercesc::DOMAttr* const attribute
724            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
725      const G4String attName = Transcode(attribute->getName());
726      const G4String attValue = Transcode(attribute->getValue());
727
728      if (attName=="vertex1")
729        { vertex1 = GetPosition(GenerateName(attValue)); } else
730      if (attName=="vertex2")
731        { vertex2 = GetPosition(GenerateName(attValue)); } else
732      if (attName=="vertex3")
733        { vertex3 = GetPosition(GenerateName(attValue)); } else
734      if (attName=="vertex4")
735        { vertex4 = GetPosition(GenerateName(attValue)); } else
736      if (attName=="type")
737        { if (attValue=="RELATIVE") { type = RELATIVE; } }
738   }
739
740   return new G4QuadrangularFacet(vertex1,vertex2,vertex3,vertex4,type);
741}
742
743void G4GDMLReadSolids::
744ReflectedSolidRead(const xercesc::DOMElement* const reflectedSolidElement)
745{
746   G4String name;
747   G4double lunit = 1.0;
748   G4double aunit = 1.0;
749   G4String solid;
750   G4ThreeVector scale(1.0,1.0,1.0);
751   G4ThreeVector rotation;
752   G4ThreeVector position;
753
754   const xercesc::DOMNamedNodeMap* const attributes
755         = reflectedSolidElement->getAttributes();
756   XMLSize_t attributeCount = attributes->getLength();
757
758   for (XMLSize_t attribute_index=0;
759        attribute_index<attributeCount; attribute_index++)
760   {
761      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
762
763      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
764        { continue; }
765
766      const xercesc::DOMAttr* const attribute
767            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
768      const G4String attName = Transcode(attribute->getName());
769      const G4String attValue = Transcode(attribute->getValue());
770
771      if (attName=="name") { name = GenerateName(attValue); } else
772      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
773      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
774      if (attName=="solid") { solid = GenerateName(attValue); } else
775      if (attName=="sx") { scale.setX(eval.Evaluate(attValue)); } else
776      if (attName=="sy") { scale.setY(eval.Evaluate(attValue)); } else
777      if (attName=="sz") { scale.setZ(eval.Evaluate(attValue)); } else
778      if (attName=="rx") { rotation.setX(eval.Evaluate(attValue)); } else
779      if (attName=="ry") { rotation.setY(eval.Evaluate(attValue)); } else
780      if (attName=="rz") { rotation.setZ(eval.Evaluate(attValue)); } else
781      if (attName=="dx") { position.setX(eval.Evaluate(attValue)); } else
782      if (attName=="dy") { position.setY(eval.Evaluate(attValue)); } else
783      if (attName=="dz") { position.setZ(eval.Evaluate(attValue)); }
784   }
785
786   rotation *= aunit;
787   position *= lunit;
788
789   G4Transform3D transform(GetRotationMatrix(rotation),position);
790   transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
791         
792   new G4ReflectedSolid(name,GetSolid(solid),transform);
793}
794
795G4ExtrudedSolid::ZSection G4GDMLReadSolids::
796SectionRead(const xercesc::DOMElement* const sectionElement,G4double lunit) 
797{
798   G4double zPosition = 0.0;
799   G4TwoVector Offset;
800   G4double scalingFactor = 1.0;
801
802   const xercesc::DOMNamedNodeMap* const attributes
803         = sectionElement->getAttributes();
804   XMLSize_t attributeCount = attributes->getLength();
805
806   for (XMLSize_t attribute_index=0;
807        attribute_index<attributeCount; attribute_index++)
808   {
809      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
810
811      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
812        { continue; }
813
814      const xercesc::DOMAttr* const attribute
815            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
816      const G4String attName = Transcode(attribute->getName());
817      const G4String attValue = Transcode(attribute->getValue());
818
819      if (attName=="zPosition")
820        { zPosition = eval.Evaluate(attValue)*lunit; } else
821      if (attName=="xOffset")
822        { Offset.setX(eval.Evaluate(attValue)*lunit); } else
823      if (attName=="yOffset")
824        { Offset.setY(eval.Evaluate(attValue)*lunit); } else
825      if (attName=="scalingFactor")
826        { scalingFactor = eval.Evaluate(attValue); }
827   }
828
829   return G4ExtrudedSolid::ZSection(zPosition,Offset,scalingFactor);
830}
831
832void G4GDMLReadSolids::
833SphereRead(const xercesc::DOMElement* const sphereElement)
834{
835   G4String name;
836   G4double lunit = 1.0;
837   G4double aunit = 1.0;
838   G4double rmin = 0.0;
839   G4double rmax = 0.0;
840   G4double startphi = 0.0;
841   G4double deltaphi = 0.0;
842   G4double starttheta = 0.0;
843   G4double deltatheta = 0.0;
844
845   const xercesc::DOMNamedNodeMap* const attributes
846         = sphereElement->getAttributes();
847   XMLSize_t attributeCount = attributes->getLength();
848
849   for (XMLSize_t attribute_index=0;
850        attribute_index<attributeCount; attribute_index++)
851   {
852      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
853
854      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
855        { continue; }
856
857      const xercesc::DOMAttr* const attribute
858            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
859      const G4String attName = Transcode(attribute->getName());
860      const G4String attValue = Transcode(attribute->getValue());
861
862      if (attName=="name") { name = GenerateName(attValue); } else
863      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
864      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
865      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
866      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
867      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
868      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); } else
869      if (attName=="starttheta") { starttheta = eval.Evaluate(attValue); } else
870      if (attName=="deltatheta") { deltatheta = eval.Evaluate(attValue); }
871   }
872
873   rmin *= lunit;
874   rmax *= lunit;
875   startphi *= aunit;
876   deltaphi *= aunit;
877   starttheta *= aunit;
878   deltatheta *= aunit;
879
880   new G4Sphere(name,rmin,rmax,startphi,deltaphi,starttheta,deltatheta);
881}
882
883void G4GDMLReadSolids::
884TessellatedRead(const xercesc::DOMElement* const tessellatedElement)
885{
886   G4String name;
887
888   const xercesc::DOMNamedNodeMap* const attributes
889         = tessellatedElement->getAttributes();
890   XMLSize_t attributeCount = attributes->getLength();
891
892   for (XMLSize_t attribute_index=0;
893        attribute_index<attributeCount; attribute_index++)
894   {
895      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
896
897      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
898        { continue; }
899
900      const xercesc::DOMAttr* const attribute
901            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
902      const G4String attName = Transcode(attribute->getName());
903      const G4String attValue = Transcode(attribute->getValue());
904
905      if (attName=="name") { name = GenerateName(attValue); }
906   }
907   
908   G4TessellatedSolid *tessellated = new G4TessellatedSolid(name);
909
910   for (xercesc::DOMNode* iter = tessellatedElement->getFirstChild();
911        iter != 0; iter = iter->getNextSibling())
912   {
913      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
914
915      const xercesc::DOMElement* const child
916            = dynamic_cast<xercesc::DOMElement*>(iter);
917      const G4String tag = Transcode(child->getTagName());
918
919      if (tag=="triangular")
920        { tessellated->AddFacet(TriangularRead(child)); } else
921      if (tag=="quadrangular")
922        { tessellated->AddFacet(QuadrangularRead(child)); }
923   }
924
925   tessellated->SetSolidClosed(true);
926}
927
928void G4GDMLReadSolids::TetRead(const xercesc::DOMElement* const tetElement)
929{
930   G4String name;
931   G4ThreeVector vertex1;
932   G4ThreeVector vertex2;
933   G4ThreeVector vertex3;
934   G4ThreeVector vertex4;
935   
936   const xercesc::DOMNamedNodeMap* const attributes
937         = tetElement->getAttributes();
938   XMLSize_t attributeCount = attributes->getLength();
939
940   for (XMLSize_t attribute_index=0;
941        attribute_index<attributeCount;attribute_index++)
942   {
943      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
944
945      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
946        { continue; }
947
948      const xercesc::DOMAttr* const attribute
949            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
950      const G4String attName = Transcode(attribute->getName());
951      const G4String attValue = Transcode(attribute->getValue());
952
953      if (attName=="name")
954        { name = GenerateName(attValue); } else
955      if (attName=="vertex1")
956        { vertex1 = GetPosition(GenerateName(attValue)); } else
957      if (attName=="vertex2")
958        { vertex2 = GetPosition(GenerateName(attValue)); } else
959      if (attName=="vertex3")
960        { vertex3 = GetPosition(GenerateName(attValue)); } else
961      if (attName=="vertex4")
962        { vertex4 = GetPosition(GenerateName(attValue)); }
963   }
964
965   new G4Tet(name,vertex1,vertex2,vertex3,vertex4);
966}
967
968void G4GDMLReadSolids::TorusRead(const xercesc::DOMElement* const torusElement)
969{
970   G4String name;
971   G4double lunit = 1.0;
972   G4double aunit = 1.0;
973   G4double rmin = 0.0;
974   G4double rmax = 0.0;
975   G4double rtor = 0.0;
976   G4double startphi = 0.0;
977   G4double deltaphi = 0.0;
978
979   const xercesc::DOMNamedNodeMap* const attributes
980         = torusElement->getAttributes();
981   XMLSize_t attributeCount = attributes->getLength();
982
983   for (XMLSize_t attribute_index=0;
984        attribute_index<attributeCount; attribute_index++)
985   {
986      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
987
988      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
989        { continue; }
990
991      const xercesc::DOMAttr* const attribute
992            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
993      const G4String attName = Transcode(attribute->getName());
994      const G4String attValue = Transcode(attribute->getValue());
995
996      if (attName=="name") { name = GenerateName(attValue); } else
997      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
998      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
999      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1000      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1001      if (attName=="rtor") { rtor = eval.Evaluate(attValue); } else
1002      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1003      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1004   }
1005
1006   rmin *= lunit;
1007   rmax *= lunit;
1008   rtor *= lunit;
1009   startphi *= aunit;
1010   deltaphi *= aunit;
1011
1012   new G4Torus(name,rmin,rmax,rtor,startphi,deltaphi);
1013}
1014
1015void G4GDMLReadSolids::TrapRead(const xercesc::DOMElement* const trapElement)
1016{
1017   G4String name;
1018   G4double lunit = 1.0;
1019   G4double aunit = 1.0;
1020   G4double z = 0.0;
1021   G4double theta = 0.0;
1022   G4double phi = 0.0;
1023   G4double y1 = 0.0;
1024   G4double x1 = 0.0;
1025   G4double x2 = 0.0;
1026   G4double alpha1 = 0.0;
1027   G4double y2 = 0.0;
1028   G4double x3 = 0.0;
1029   G4double x4 = 0.0;
1030   G4double alpha2 = 0.0;
1031
1032   const xercesc::DOMNamedNodeMap* const attributes
1033         = trapElement->getAttributes();
1034   XMLSize_t attributeCount = attributes->getLength();
1035
1036   for (XMLSize_t attribute_index=0;
1037        attribute_index<attributeCount; attribute_index++)
1038   {
1039      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1040
1041      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1042        { continue; }
1043
1044      const xercesc::DOMAttr* const attribute
1045            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1046      const G4String attName = Transcode(attribute->getName());
1047      const G4String attValue = Transcode(attribute->getValue());
1048
1049      if (attName=="name") { name = GenerateName(attValue); } else
1050      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1051      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1052      if (attName=="z") { z = eval.Evaluate(attValue); } else
1053      if (attName=="theta") { theta = eval.Evaluate(attValue); } else
1054      if (attName=="phi") { phi = eval.Evaluate(attValue); } else
1055      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1056      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1057      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1058      if (attName=="alpha1") { alpha1 = eval.Evaluate(attValue); } else
1059      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1060      if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1061      if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1062      if (attName=="alpha2") { alpha2 = eval.Evaluate(attValue); }
1063   }
1064
1065   z *= 0.5*lunit;
1066   theta *= aunit;
1067   phi *= aunit;
1068   y1 *= 0.5*lunit;
1069   x1 *= 0.5*lunit;
1070   x2 *= 0.5*lunit;
1071   alpha1 *= aunit;
1072   y2 *= 0.5*lunit;
1073   x3 *= 0.5*lunit;
1074   x4 *= 0.5*lunit;
1075   alpha2 *= aunit;
1076
1077   new G4Trap(name,z,theta,phi,y1,x1,x2,alpha1,y2,x3,x4,alpha2);
1078}
1079
1080void G4GDMLReadSolids::TrdRead(const xercesc::DOMElement* const trdElement)
1081{
1082   G4String name;
1083   G4double lunit = 1.0;
1084   G4double aunit = 1.0;
1085   G4double x1 = 0.0;
1086   G4double x2 = 0.0;
1087   G4double y1 = 0.0;
1088   G4double y2 = 0.0;
1089   G4double z = 0.0;
1090
1091   const xercesc::DOMNamedNodeMap* const attributes = trdElement->getAttributes();
1092   XMLSize_t attributeCount = attributes->getLength();
1093
1094   for (XMLSize_t attribute_index=0;
1095        attribute_index<attributeCount; attribute_index++)
1096   {
1097      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1098
1099      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1100        { continue; }
1101
1102      const xercesc::DOMAttr* const attribute
1103            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1104      const G4String attName = Transcode(attribute->getName());
1105      const G4String attValue = Transcode(attribute->getValue());
1106
1107      if (attName=="name") { name = GenerateName(attValue); } else
1108      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1109      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1110      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1111      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1112      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1113      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1114      if (attName=="z") { z = eval.Evaluate(attValue); }
1115   }
1116
1117   x1 *= 0.5*lunit;
1118   x2 *= 0.5*lunit;
1119   y1 *= 0.5*lunit;
1120   y2 *= 0.5*lunit;
1121   z *= 0.5*lunit;
1122
1123   new G4Trd(name,x1,x2,y1,y2,z);
1124}
1125
1126G4TriangularFacet* G4GDMLReadSolids::
1127TriangularRead(const xercesc::DOMElement* const triangularElement)
1128{
1129   G4ThreeVector vertex1;
1130   G4ThreeVector vertex2;
1131   G4ThreeVector vertex3;
1132   G4FacetVertexType type = ABSOLUTE;
1133
1134   const xercesc::DOMNamedNodeMap* const attributes
1135         = triangularElement->getAttributes();
1136   XMLSize_t attributeCount = attributes->getLength();
1137
1138   for (XMLSize_t attribute_index=0;
1139        attribute_index<attributeCount; attribute_index++)
1140   {
1141      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1142
1143      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1144        { continue; }
1145
1146      const xercesc::DOMAttr* const attribute
1147            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1148      const G4String attName = Transcode(attribute->getName());
1149      const G4String attValue = Transcode(attribute->getValue());
1150
1151      if (attName=="vertex1")
1152        { vertex1 = GetPosition(GenerateName(attValue)); } else
1153      if (attName=="vertex2")
1154        { vertex2 = GetPosition(GenerateName(attValue)); } else
1155      if (attName=="vertex3")
1156        { vertex3 = GetPosition(GenerateName(attValue)); } else
1157      if (attName=="type")
1158        { if (attValue=="RELATIVE") { type = RELATIVE; } }
1159   }
1160
1161   return new G4TriangularFacet(vertex1,vertex2,vertex3,type);
1162}
1163
1164void G4GDMLReadSolids::TubeRead(const xercesc::DOMElement* const tubeElement)
1165{
1166   G4String name;
1167   G4double lunit = 1.0;
1168   G4double aunit = 1.0;
1169   G4double rmin = 0.0;
1170   G4double rmax = 0.0;
1171   G4double z = 0.0;
1172   G4double startphi = 0.0;
1173   G4double deltaphi = 0.0;
1174
1175   const xercesc::DOMNamedNodeMap* const attributes
1176         = tubeElement->getAttributes();
1177   XMLSize_t attributeCount = attributes->getLength();
1178
1179   for (XMLSize_t attribute_index=0;
1180        attribute_index<attributeCount; attribute_index++)
1181   {
1182      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1183
1184      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1185        { continue; }
1186
1187      const xercesc::DOMAttr* const attribute
1188            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1189      const G4String attName = Transcode(attribute->getName());
1190      const G4String attValue = Transcode(attribute->getValue());
1191
1192      if (attName=="name") { name = GenerateName(attValue); } else
1193      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1194      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1195      if (attName=="rmin") { rmin = eval.Evaluate(attValue); } else
1196      if (attName=="rmax") { rmax = eval.Evaluate(attValue); } else
1197      if (attName=="z") { z = eval.Evaluate(attValue); } else
1198      if (attName=="startphi") { startphi = eval.Evaluate(attValue); } else
1199      if (attName=="deltaphi") { deltaphi = eval.Evaluate(attValue); }
1200   }
1201
1202   rmin *= lunit;
1203   rmax *= lunit;
1204   z *= 0.5*lunit;
1205   startphi *= aunit;
1206   deltaphi *= aunit;
1207
1208   new G4Tubs(name,rmin,rmax,z,startphi,deltaphi);
1209}
1210
1211void G4GDMLReadSolids::
1212TwistedboxRead(const xercesc::DOMElement* const twistedboxElement)
1213{
1214   G4String name;
1215   G4double lunit = 1.0;
1216   G4double aunit = 1.0;
1217   G4double PhiTwist = 0.0;
1218   G4double x = 0.0;
1219   G4double y = 0.0;
1220   G4double z = 0.0;
1221
1222   const xercesc::DOMNamedNodeMap* const attributes
1223         = twistedboxElement->getAttributes();
1224   XMLSize_t attributeCount = attributes->getLength();
1225
1226   for (XMLSize_t attribute_index=0;
1227        attribute_index<attributeCount; attribute_index++)
1228   {
1229      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1230
1231      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1232        { continue; }
1233
1234      const xercesc::DOMAttr* const attribute
1235            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1236      const G4String attName = Transcode(attribute->getName());
1237      const G4String attValue = Transcode(attribute->getValue());
1238
1239      if (attName=="name") { name = GenerateName(attValue); } else
1240      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1241      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1242      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1243      if (attName=="x") { x = eval.Evaluate(attValue); } else
1244      if (attName=="y") { y = eval.Evaluate(attValue); } else
1245      if (attName=="z") { z = eval.Evaluate(attValue); }
1246   }
1247
1248   PhiTwist *= aunit;
1249   x *= 0.5*lunit;
1250   y *= 0.5*lunit;
1251   z *= 0.5*lunit;
1252
1253   new G4TwistedBox(name,PhiTwist,x,y,z);
1254}
1255
1256void G4GDMLReadSolids::
1257TwistedtrapRead(const xercesc::DOMElement* const twistedtrapElement)
1258{
1259   G4String name;
1260   G4double lunit = 1.0;
1261   G4double aunit = 1.0;
1262   G4double PhiTwist = 0.0;
1263   G4double z = 0.0;
1264   G4double Theta = 0.0;
1265   G4double Phi = 0.0;
1266   G4double y1 = 0.0;
1267   G4double x1 = 0.0;
1268   G4double x2 = 0.0;
1269   G4double y2 = 0.0;
1270   G4double x3 = 0.0;
1271   G4double x4 = 0.0;
1272   G4double Alph = 0.0;
1273
1274   const xercesc::DOMNamedNodeMap* const attributes
1275         = twistedtrapElement->getAttributes();
1276   XMLSize_t attributeCount = attributes->getLength();
1277
1278   for (XMLSize_t attribute_index=0;
1279        attribute_index<attributeCount; attribute_index++)
1280   {
1281      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1282
1283      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1284        { continue; }
1285
1286      const xercesc::DOMAttr* const attribute
1287            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1288      const G4String attName = Transcode(attribute->getName());
1289      const G4String attValue = Transcode(attribute->getValue());
1290
1291      if (attName=="name") { name = GenerateName(attValue); } else
1292      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1293      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1294      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); } else
1295      if (attName=="z") { z = eval.Evaluate(attValue); } else
1296      if (attName=="Theta") { Theta = eval.Evaluate(attValue); } else
1297      if (attName=="Phi") { Phi = eval.Evaluate(attValue); } else
1298      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1299      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1300      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1301      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1302      if (attName=="x3") { x3 = eval.Evaluate(attValue); } else
1303      if (attName=="x4") { x4 = eval.Evaluate(attValue); } else
1304      if (attName=="Alph") { Alph = eval.Evaluate(attValue); }
1305   }
1306
1307
1308   PhiTwist *= aunit;
1309   z *= 0.5*lunit;
1310   Theta *= aunit;
1311   Phi *= aunit;
1312   Alph *= aunit;
1313   y1 *= 0.5*lunit;
1314   x1 *= 0.5*lunit;
1315   x2 *= 0.5*lunit;
1316   y2 *= 0.5*lunit;
1317   x3 *= 0.5*lunit;
1318   x4 *= 0.5*lunit;
1319
1320   new G4TwistedTrap(name,PhiTwist,z,Theta,Phi,y1,x1,x2,y2,x3,x4,Alph);
1321}
1322
1323void G4GDMLReadSolids::
1324TwistedtrdRead(const xercesc::DOMElement* const twistedtrdElement)
1325{
1326   G4String name;
1327   G4double lunit = 1.0;
1328   G4double aunit = 1.0;
1329   G4double x1 = 0.0;
1330   G4double x2 = 0.0;
1331   G4double y1 = 0.0;
1332   G4double y2 = 0.0;
1333   G4double z = 0.0;
1334   G4double PhiTwist = 0.0;
1335
1336   const xercesc::DOMNamedNodeMap* const attributes
1337         = twistedtrdElement->getAttributes();
1338   XMLSize_t attributeCount = attributes->getLength();
1339
1340   for (XMLSize_t attribute_index=0;
1341        attribute_index<attributeCount; attribute_index++)
1342   {
1343      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1344
1345      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1346        { continue; }
1347
1348      const xercesc::DOMAttr* const attribute
1349            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1350      const G4String attName = Transcode(attribute->getName());
1351      const G4String attValue = Transcode(attribute->getValue());
1352
1353      if (attName=="name") { name = GenerateName(attValue); } else
1354      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1355      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1356      if (attName=="x1") { x1 = eval.Evaluate(attValue); } else
1357      if (attName=="x2") { x2 = eval.Evaluate(attValue); } else
1358      if (attName=="y1") { y1 = eval.Evaluate(attValue); } else
1359      if (attName=="y2") { y2 = eval.Evaluate(attValue); } else
1360      if (attName=="z") { z = eval.Evaluate(attValue); } else
1361      if (attName=="PhiTwist") { PhiTwist = eval.Evaluate(attValue); }
1362   }
1363
1364   x1 *= 0.5*lunit;
1365   x2 *= 0.5*lunit;
1366   y1 *= 0.5*lunit;
1367   y2 *= 0.5*lunit;
1368   z *= 0.5*lunit;
1369   PhiTwist *= aunit;
1370
1371   new G4TwistedTrd(name,x1,x2,y1,y2,z,PhiTwist);
1372}
1373
1374void G4GDMLReadSolids::
1375TwistedtubsRead(const xercesc::DOMElement* const twistedtubsElement)
1376{
1377   G4String name;
1378   G4double lunit = 1.0;
1379   G4double aunit = 1.0;
1380   G4double twistedangle = 0.0;
1381   G4double endinnerrad = 0.0;
1382   G4double endouterrad = 0.0;
1383   G4double zlen = 0.0;
1384   G4double phi = 0.0;
1385
1386   const xercesc::DOMNamedNodeMap* const attributes
1387         = twistedtubsElement->getAttributes();
1388   XMLSize_t attributeCount = attributes->getLength();
1389
1390   for (XMLSize_t attribute_index=0;
1391        attribute_index<attributeCount; attribute_index++)
1392   {
1393      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1394
1395      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1396        { continue; }
1397
1398      const xercesc::DOMAttr* const attribute
1399            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1400      const G4String attName = Transcode(attribute->getName());
1401      const G4String attValue = Transcode(attribute->getValue());
1402
1403      if (attName=="name") { name = GenerateName(attValue); } else
1404      if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
1405      if (attName=="aunit") { aunit = eval.Evaluate(attValue); } else
1406      if (attName=="twistedangle") { twistedangle=eval.Evaluate(attValue); } else
1407      if (attName=="endinnerrad")  { endinnerrad=eval.Evaluate(attValue);  } else
1408      if (attName=="endouterrad")  { endouterrad=eval.Evaluate(attValue);  } else
1409      if (attName=="zlen") { zlen = eval.Evaluate(attValue); } else
1410      if (attName=="phi") { phi = eval.Evaluate(attValue); }
1411   }
1412
1413   twistedangle *= aunit;
1414   endinnerrad *= lunit;
1415   endouterrad *= lunit;
1416   zlen *= 0.5*lunit;
1417   phi *= aunit;
1418
1419   new G4TwistedTubs(name,twistedangle,endinnerrad,endouterrad,zlen,phi);
1420}
1421
1422G4TwoVector G4GDMLReadSolids::
1423TwoDimVertexRead(const xercesc::DOMElement* const element, G4double lunit)
1424{
1425   G4TwoVector vec;
1426   
1427   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1428   XMLSize_t attributeCount = attributes->getLength();
1429
1430   for (XMLSize_t attribute_index=0;
1431        attribute_index<attributeCount; attribute_index++)
1432   {
1433      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1434
1435      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1436        { continue; }
1437
1438      const xercesc::DOMAttr* const attribute
1439            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1440      const G4String attName = Transcode(attribute->getName());
1441      const G4String attValue = Transcode(attribute->getValue());
1442
1443      if (attName=="x") { vec.setX(eval.Evaluate(attValue)*lunit); } else
1444      if (attName=="y") { vec.setY(eval.Evaluate(attValue)*lunit); }
1445   }
1446
1447   return vec;
1448}
1449
1450G4GDMLReadSolids::zplaneType G4GDMLReadSolids::
1451ZplaneRead(const xercesc::DOMElement* const zplaneElement)
1452{
1453   zplaneType zplane;
1454
1455   const xercesc::DOMNamedNodeMap* const attributes
1456         = zplaneElement->getAttributes();
1457   XMLSize_t attributeCount = attributes->getLength();
1458
1459   for (XMLSize_t attribute_index=0;
1460        attribute_index<attributeCount; attribute_index++)
1461   {
1462      xercesc::DOMNode* node = attributes->item(attribute_index);
1463
1464      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
1465
1466      const xercesc::DOMAttr* const attribute
1467            = dynamic_cast<xercesc::DOMAttr*>(node);   
1468      const G4String attName = Transcode(attribute->getName());
1469      const G4String attValue = Transcode(attribute->getValue());
1470
1471      if (attName=="rmin") { zplane.rmin = eval.Evaluate(attValue); } else
1472      if (attName=="rmax") { zplane.rmax = eval.Evaluate(attValue); } else
1473      if (attName=="z") { zplane.z = eval.Evaluate(attValue); }
1474   }
1475
1476   return zplane;
1477}
1478
1479void G4GDMLReadSolids::
1480OpticalsurfaceRead(const xercesc::DOMElement* const opticalsurfaceElement)
1481{
1482   G4String name;
1483   G4String smodel;
1484   G4String sfinish;
1485   G4String stype;
1486   G4double value = 0.0;
1487
1488   const xercesc::DOMNamedNodeMap* const attributes
1489         = opticalsurfaceElement->getAttributes();
1490   XMLSize_t attributeCount = attributes->getLength();
1491
1492   for (XMLSize_t attribute_index=0;
1493        attribute_index<attributeCount; attribute_index++)
1494   {
1495      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1496
1497      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
1498        { continue; }
1499
1500      const xercesc::DOMAttr* const attribute
1501            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1502      const G4String attName = Transcode(attribute->getName());
1503      const G4String attValue = Transcode(attribute->getValue());
1504
1505      if (attName=="name") { name = GenerateName(attValue); } else
1506      if (attName=="model") { smodel = attValue; } else
1507      if (attName=="finish") { sfinish = attValue; } else
1508      if (attName=="type") { stype = attValue; } else
1509      if (attName=="value") { value = eval.Evaluate(attValue); }
1510   }
1511
1512   G4OpticalSurfaceModel model; 
1513   G4OpticalSurfaceFinish finish;
1514   G4SurfaceType type;   
1515   
1516   if (smodel="unified") { model = unified; } else { model = glisur; }
1517
1518   if (sfinish=="polishedfrontpainted") { finish = polishedfrontpainted; } else
1519   if (sfinish=="polishedbackpainted") { finish = polishedbackpainted; } else
1520   if (sfinish=="groundfrontpainted") { finish = groundfrontpainted; } else
1521   if (sfinish=="groundbackpainted") { finish = groundbackpainted; } else
1522   if (sfinish=="ground") { finish = ground; } else { finish = polished; }
1523
1524   if (stype=="dielectric_metal") { type = dielectric_metal; } else
1525   if (stype=="x_ray") { type = x_ray; } else
1526   if (stype=="firsov") { type = firsov; } else { type = dielectric_dielectric; }
1527
1528   new G4OpticalSurface(name,model,finish,type,value);
1529}
1530
1531void G4GDMLReadSolids::SolidsRead(const xercesc::DOMElement* const solidsElement)
1532{
1533   G4cout << "G4GDML: Reading solids..." << G4endl;
1534
1535   for (xercesc::DOMNode* iter = solidsElement->getFirstChild();
1536        iter != 0; iter = iter->getNextSibling())
1537   {
1538      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)  { continue; }
1539
1540      const xercesc::DOMElement* const child
1541            = dynamic_cast<xercesc::DOMElement*>(iter);
1542      const G4String tag = Transcode(child->getTagName());
1543
1544      if (tag=="box") { BoxRead(child); } else
1545      if (tag=="cone") { ConeRead(child); } else
1546      if (tag=="elcone") { ElconeRead(child); } else
1547      if (tag=="ellipsoid") { EllipsoidRead(child); }else
1548      if (tag=="eltube") { EltubeRead(child); } else
1549      if (tag=="xtru") { XtruRead(child); } else
1550      if (tag=="hype") { HypeRead(child); } else
1551      if (tag=="intersection") { BooleanRead(child,INTERSECTION); } else
1552      if (tag=="orb") { OrbRead(child); } else
1553      if (tag=="para") { ParaRead(child); } else
1554      if (tag=="paraboloid") { ParaboloidRead(child); } else
1555      if (tag=="polycone") { PolyconeRead(child); } else
1556      if (tag=="polyhedra") { PolyhedraRead(child); } else
1557      if (tag=="reflectedSolid") { ReflectedSolidRead(child); } else
1558      if (tag=="sphere") { SphereRead(child); } else
1559      if (tag=="subtraction") { BooleanRead(child,SUBTRACTION); } else
1560      if (tag=="tessellated") { TessellatedRead(child); } else
1561      if (tag=="tet") { TetRead(child); } else
1562      if (tag=="torus") { TorusRead(child); } else
1563      if (tag=="trap") { TrapRead(child); } else
1564      if (tag=="trd") { TrdRead(child); } else
1565      if (tag=="tube") { TubeRead(child); } else
1566      if (tag=="twistedbox") { TwistedboxRead(child); } else
1567      if (tag=="twistedtrap") { TwistedtrapRead(child); } else
1568      if (tag=="twistedtrd") { TwistedtrdRead(child); } else
1569      if (tag=="twistedtubs") { TwistedtubsRead(child); } else
1570      if (tag=="union") { BooleanRead(child,UNION); } else
1571      if (tag=="opticalsurface") { OpticalsurfaceRead(child); } else
1572      if (tag=="loop") { LoopRead(child,&G4GDMLRead::SolidsRead); }
1573      else
1574      {
1575        G4String error_msg = "Unknown tag in solids: " + tag;
1576        G4Exception("G4GDMLReadSolids::SolidsRead()", "ReadError",
1577                    FatalException, error_msg);
1578      }
1579   }
1580}
1581
1582G4VSolid* G4GDMLReadSolids::GetSolid(const G4String& ref) const
1583{
1584   G4VSolid* solidPtr = G4SolidStore::GetInstance()->GetSolid(ref,false);
1585
1586   if (!solidPtr)
1587   {
1588     G4String error_msg = "Referenced solid '" + ref + "' was not found!";
1589     G4Exception("G4GDMLReadSolids::GetSolid()", "ReadError",
1590                 FatalException, error_msg);
1591   }
1592
1593   return solidPtr;
1594}
1595
1596G4SurfaceProperty* G4GDMLReadSolids::
1597GetSurfaceProperty(const G4String& ref) const
1598{
1599   const G4SurfacePropertyTable* surfaceList
1600         = G4SurfaceProperty::GetSurfacePropertyTable();
1601   const size_t surfaceCount = surfaceList->size();
1602
1603   for (size_t i=0; i<surfaceCount; i++)
1604   {
1605      if ((*surfaceList)[i]->GetName() == ref)  { return (*surfaceList)[i]; }
1606   }
1607
1608   G4String error_msg = "Referenced optical surface '" + ref + "' was not found!";
1609   G4Exception("G4GDMLReadSolids::GetSurfaceProperty()", "ReadError",
1610               FatalException, error_msg);
1611
1612   return 0;
1613}
Note: See TracBrowser for help on using the repository browser.