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

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

import all except CVS

File size: 51.1 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//
27// $Id: G4GDMLSolids.cc,v 1.23.2.1 2008/01/16 09:43:09 ztorzsok Exp $
28// GEANT4 tag $ Name:$
29//
30// class G4GDMLSolids Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLSolids.hh"
37
38void G4GDMLSolids::booleanRead(const xercesc::DOMElement* const element,const BooleanOp op) {
39
40   G4String name;
41   G4String first;
42   G4String second;
43
44   G4ThreeVector position;
45   G4ThreeVector rotation;
46
47   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
48   XMLSize_t attributeCount = attributes->getLength();
49
50   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
51
52      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
53
54      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
55
56      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
57
58      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
59      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
60
61      if (attribute_name=="name") name = attribute_value;
62   }
63
64   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
65
66      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
67
68      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
69
70      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
71
72      if (tag=="first") first = refRead(child); else
73      if (tag=="second") second = refRead(child); else
74      if (tag=="position") position = positionRead(child); else
75      if (tag=="rotation") rotation = rotationRead(child); else
76      G4Exception("GDML: Unknown tag in boolean solid: "+tag);
77   }
78
79   G4RotationMatrix rot;
80
81   rot.rotateX(rotation.x());
82   rot.rotateY(rotation.y());
83   rot.rotateZ(rotation.z());
84
85   G4Transform3D transform(rot,position);
86
87   G4VSolid* firstSolid = getSolid(GenerateName(first));
88   G4VSolid* secondSolid = getSolid(GenerateName(second));
89
90   if (op==UNION) new G4UnionSolid(GenerateName(name),firstSolid,secondSolid,transform); else
91   if (op==SUBTRACTION) new G4SubtractionSolid(GenerateName(name),firstSolid,secondSolid,transform); else
92   if (op==INTERSECTION) new G4IntersectionSolid(GenerateName(name),firstSolid,secondSolid,transform);
93}
94
95void G4GDMLSolids::boxRead(const xercesc::DOMElement* const element) {
96
97   G4String name;
98   G4String lunit("1");
99   G4String x;
100   G4String y;
101   G4String z;
102
103   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
104   XMLSize_t attributeCount = attributes->getLength();
105
106   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
107
108      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
109
110      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
111
112      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
113
114      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
115      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
116
117      if (attribute_name=="name" ) name  = attribute_value; else
118      if (attribute_name=="lunit") lunit = attribute_value; else
119      if (attribute_name=="x") x = attribute_value; else
120      if (attribute_name=="y") y = attribute_value; else
121      if (attribute_name=="z") z = attribute_value;
122   }
123
124   G4double _lunit = eval.Evaluate(lunit);
125
126   G4double _x = eval.Evaluate(x)*_lunit;
127   G4double _y = eval.Evaluate(y)*_lunit;
128   G4double _z = eval.Evaluate(z)*_lunit;
129
130   _x *= 0.5;
131   _y *= 0.5;
132   _z *= 0.5;
133
134   new G4Box(GenerateName(name),_x,_y,_z);
135}
136
137void G4GDMLSolids::coneRead(const xercesc::DOMElement* const element) {
138
139   G4String name;
140   G4String lunit("1");
141   G4String aunit("1");
142   G4String rmin1;
143   G4String rmax1;
144   G4String rmin2;
145   G4String rmax2;
146   G4String z;
147   G4String startphi;
148   G4String deltaphi;
149
150   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
151   XMLSize_t attributeCount = attributes->getLength();
152
153   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
154
155      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
156
157      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
158
159      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
160
161      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
162      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
163
164      if (attribute_name=="name") name = attribute_value; else
165      if (attribute_name=="lunit") lunit = attribute_value; else
166      if (attribute_name=="aunit") aunit = attribute_value; else
167      if (attribute_name=="rmin1") rmin1 = attribute_value; else
168      if (attribute_name=="rmax1") rmax1 = attribute_value; else
169      if (attribute_name=="rmin2") rmin2 = attribute_value; else
170      if (attribute_name=="rmax2") rmax2 = attribute_value; else
171      if (attribute_name=="z") z = attribute_value; else
172      if (attribute_name=="startphi") startphi = attribute_value; else
173      if (attribute_name=="deltaphi") deltaphi = attribute_value;
174   }
175
176   G4double _lunit = eval.Evaluate(lunit);
177   G4double _aunit = eval.Evaluate(aunit);
178
179   G4double _rmin1 = eval.Evaluate(rmin1)*_lunit;
180   G4double _rmax1 = eval.Evaluate(rmax1)*_lunit;
181   G4double _rmin2 = eval.Evaluate(rmin2)*_lunit;
182   G4double _rmax2 = eval.Evaluate(rmax2)*_lunit;
183   G4double _z = eval.Evaluate(z)*_lunit;
184   G4double _startphi = eval.Evaluate(startphi)*_aunit;
185   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;;
186
187   _z *= 0.5;
188
189   new G4Cons(GenerateName(name),_rmin1,_rmax1,_rmin2,_rmax2,_z,_startphi,_deltaphi);
190}
191
192void G4GDMLSolids::ellipsoidRead(const xercesc::DOMElement* const element) {
193
194   G4String name;
195   G4String lunit("1");
196   G4String ax;
197   G4String by;
198   G4String cz;
199   G4String zcut1;
200   G4String zcut2; 
201
202   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
203   XMLSize_t attributeCount = attributes->getLength();
204
205   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
206
207      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
208
209      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
210
211      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
212
213      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
214      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
215
216      if (attribute_name=="name") name  = attribute_value; else
217      if (attribute_name=="lunit") lunit = attribute_value; else
218      if (attribute_name=="ax") ax = attribute_value; else
219      if (attribute_name=="by") by = attribute_value; else
220      if (attribute_name=="cz") cz = attribute_value; else
221      if (attribute_name=="zcut1") zcut1 = attribute_value; else
222      if (attribute_name=="zcut2") zcut2 = attribute_value;
223   }
224
225   G4double _lunit = eval.Evaluate(lunit);
226
227   G4double _ax = eval.Evaluate(ax)*_lunit;
228   G4double _by = eval.Evaluate(by)*_lunit;
229   G4double _cz = eval.Evaluate(cz)*_lunit;
230   G4double _zcut1 = eval.Evaluate(zcut1)*_lunit;
231   G4double _zcut2 = eval.Evaluate(zcut2)*_lunit; 
232
233   new G4Ellipsoid(GenerateName(name),_ax,_by,_cz,_zcut1,_zcut2);
234}
235
236void G4GDMLSolids::eltubeRead(const xercesc::DOMElement* const element) {
237
238   G4String name;
239   G4String lunit("1");
240   G4String dx;
241   G4String dy;
242   G4String dz;
243
244   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
245   XMLSize_t attributeCount = attributes->getLength();
246
247   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
248
249      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
250
251      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
252
253      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
254
255      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
256      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
257
258      if (attribute_name=="name") name  = attribute_value; else
259      if (attribute_name=="lunit") lunit = attribute_value; else
260      if (attribute_name=="dx") dx = attribute_value; else
261      if (attribute_name=="dy") dy = attribute_value; else
262      if (attribute_name=="dz") dz = attribute_value;
263   }
264
265   G4double _lunit = eval.Evaluate(lunit);
266
267   G4double _dx = eval.Evaluate(dx)*_lunit;
268   G4double _dy = eval.Evaluate(dy)*_lunit;
269   G4double _dz = eval.Evaluate(dz)*_lunit;
270
271   new G4EllipticalTube(GenerateName(name),_dx,_dy,_dz);
272}
273
274void G4GDMLSolids::hypeRead(const xercesc::DOMElement* const element) {
275
276   G4String name;
277   G4String lunit("1");
278   G4String aunit("1");
279   G4String rmin;
280   G4String rmax;
281   G4String inst;
282   G4String outst;
283   G4String z;
284
285   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
286   XMLSize_t attributeCount = attributes->getLength();
287
288   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
289
290      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
291
292      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
293
294      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
295
296      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
297      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
298
299      if (attribute_name=="name") name  = attribute_value; else
300      if (attribute_name=="lunit") lunit = attribute_value; else
301      if (attribute_name=="aunit") aunit = attribute_value; else
302      if (attribute_name=="rmin") rmin = attribute_value; else
303      if (attribute_name=="rmax") rmax = attribute_value; else
304      if (attribute_name=="inst") inst = attribute_value; else
305      if (attribute_name=="outst") outst = attribute_value; else
306      if (attribute_name=="z") z = attribute_value;
307   }
308
309   G4double _lunit = eval.Evaluate(lunit);
310   G4double _aunit = eval.Evaluate(aunit);
311
312   G4double _rmin = eval.Evaluate(rmin)*_lunit;
313   G4double _rmax = eval.Evaluate(rmax)*_lunit;;
314   G4double _inst = eval.Evaluate(inst)*_aunit;
315   G4double _outst = eval.Evaluate(outst)*_aunit;
316   G4double _z = eval.Evaluate(z)*_lunit;
317
318   _z *= 0.5;
319
320   new G4Hype(GenerateName(name),_rmin,_rmax,_inst,_outst,_z);
321}
322
323void G4GDMLSolids::loopRead(const xercesc::DOMElement* const element) {
324
325   G4String var;
326   G4String from;
327   G4String to;
328   G4String step;
329
330   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
331   XMLSize_t attributeCount = attributes->getLength();
332
333   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
334
335      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
336
337      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
338
339      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
340
341      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
342      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
343
344      if (attribute_name=="var") var = attribute_value; else
345      if (attribute_name=="from") from = attribute_value; else
346      if (attribute_name=="to") to = attribute_value; else
347      if (attribute_name=="step") step = attribute_value;
348   }
349
350   eval.checkVariable(var);
351
352   G4int _var = eval.EvaluateInteger(var);
353   G4int _from = eval.EvaluateInteger(from);
354   G4int _to = eval.EvaluateInteger(to);
355   G4int _step = eval.EvaluateInteger(step);
356   
357   if (!from.empty()) _var = _from;
358
359   while (_var <= _to) {
360   
361      eval.setVariable(var,_var);
362      solidsRead(element);
363
364      _var += _step;
365   }
366}
367
368void G4GDMLSolids::orbRead(const xercesc::DOMElement* const element) {
369
370   G4String name;
371   G4String lunit("1");
372   G4String r;
373
374   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
375   XMLSize_t attributeCount = attributes->getLength();
376
377   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
378
379      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
380
381      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
382
383      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
384
385      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
386      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
387
388      if (attribute_name=="name") name  = attribute_value; else
389      if (attribute_name=="lunit") lunit = attribute_value; else
390      if (attribute_name=="r") r = attribute_value;
391   }
392
393   G4double _lunit = eval.Evaluate(lunit);
394   
395   G4double _r = eval.Evaluate(r)*_lunit;
396
397   new G4Orb(GenerateName(name),_r);
398}
399
400void G4GDMLSolids::paraRead(const xercesc::DOMElement* const element) {
401
402   G4String name;
403   G4String lunit("1");
404   G4String aunit("1");
405   G4String x;
406   G4String y;
407   G4String z;
408   G4String alpha;
409   G4String theta;
410   G4String phi;
411
412   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
413   XMLSize_t attributeCount = attributes->getLength();
414
415   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
416
417      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
418
419      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
420
421      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
422
423      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
424      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
425
426      if (attribute_name=="name" ) name  = attribute_value; else
427      if (attribute_name=="lunit") lunit = attribute_value; else
428      if (attribute_name=="aunit") aunit = attribute_value; else
429      if (attribute_name=="x") x = attribute_value; else
430      if (attribute_name=="y") y = attribute_value; else
431      if (attribute_name=="z") z = attribute_value; else
432      if (attribute_name=="alpha") alpha = attribute_value; else
433      if (attribute_name=="theta") theta = attribute_value; else
434      if (attribute_name=="phi") phi = attribute_value;
435   }
436
437   G4double _lunit = eval.Evaluate(lunit);
438   G4double _aunit = eval.Evaluate(aunit);
439
440   G4double _x = eval.Evaluate(x)*_lunit;
441   G4double _y = eval.Evaluate(y)*_lunit;
442   G4double _z = eval.Evaluate(z)*_lunit;
443   G4double _alpha = eval.Evaluate(alpha)*_aunit;
444   G4double _theta = eval.Evaluate(theta)*_aunit;
445   G4double _phi = eval.Evaluate(phi)*_aunit;
446
447   _x *= 0.5;
448   _y *= 0.5;
449   _z *= 0.5;
450
451   new G4Para(GenerateName(name),_x,_y,_z,_alpha,_theta,_phi);
452}
453
454void G4GDMLSolids::polyconeRead(const xercesc::DOMElement* const element) {
455
456   G4String name;
457   G4String lunit("1");
458   G4String aunit("1");
459   G4String startphi;
460   G4String deltaphi;
461
462   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
463   XMLSize_t attributeCount = attributes->getLength();
464
465   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
466
467      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
468
469      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
470
471      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
472
473      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
474      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
475
476      if (attribute_name=="name") name = attribute_value; else
477      if (attribute_name=="lunit") lunit = attribute_value; else
478      if (attribute_name=="aunit") aunit = attribute_value; else
479      if (attribute_name=="startphi") startphi = attribute_value; else
480      if (attribute_name=="deltaphi") deltaphi = attribute_value;
481   }
482
483   G4double _lunit = eval.Evaluate(lunit);
484   G4double _aunit = eval.Evaluate(aunit);
485
486   G4double _startphi = eval.Evaluate(startphi)*_aunit;
487   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;
488
489   std::vector<zplaneType> zplaneList;
490
491   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
492
493      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
494
495      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
496
497      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
498
499      if (tag=="zplane") zplaneList.push_back(zplaneRead(child,_lunit));
500   }
501
502   G4int numZPlanes = zplaneList.size();
503
504   G4double *rmin_array = new G4double[numZPlanes];
505   G4double *rmax_array = new G4double[numZPlanes];
506   G4double* z_array    = new G4double[numZPlanes];
507
508   for (G4int i=0;i<numZPlanes;i++) {
509   
510      rmin_array[i] = zplaneList[i].rmin;
511      rmax_array[i] = zplaneList[i].rmax;
512      z_array[i]    = zplaneList[i].z;
513   }
514
515   new G4Polycone(GenerateName(name),_startphi,_deltaphi,numZPlanes,z_array,rmin_array,rmax_array);
516}
517
518void G4GDMLSolids::polyhedraRead(const xercesc::DOMElement* const element) {
519
520   G4String name;
521   G4String lunit("1");
522   G4String aunit("1");
523   G4String startphi;
524   G4String deltaphi;
525   G4String numsides;
526
527   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
528   XMLSize_t attributeCount = attributes->getLength();
529
530   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
531
532      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
533
534      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
535
536      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
537
538      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
539      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
540
541      if (attribute_name=="name") name     = attribute_value; else
542      if (attribute_name=="lunit") lunit    = attribute_value; else
543      if (attribute_name=="aunit" ) aunit    = attribute_value; else
544      if (attribute_name=="startphi") startphi = attribute_value; else
545      if (attribute_name=="deltaphi") deltaphi = attribute_value; else
546      if (attribute_name=="numsides") numsides = attribute_value;
547   }
548
549   G4double _lunit = eval.Evaluate(lunit);
550   G4double _aunit = eval.Evaluate(aunit);
551
552   G4double _startphi = eval.Evaluate(startphi)*_aunit;
553   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;
554
555   G4int _numsides = eval.EvaluateInteger(numsides);
556
557   std::vector<zplaneType> zplaneList;
558
559   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
560
561      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
562
563      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
564
565      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
566
567      if (tag=="zplane") zplaneList.push_back(zplaneRead(child,_lunit));
568   }
569
570   G4int numZPlanes = zplaneList.size();
571
572   G4double *rmin_array = new G4double[numZPlanes];
573   G4double *rmax_array = new G4double[numZPlanes];
574   G4double* z_array    = new G4double[numZPlanes];
575
576   for (G4int i=0;i<numZPlanes;i++) {
577   
578      rmin_array[i] = zplaneList[i].rmin;
579      rmax_array[i] = zplaneList[i].rmax;
580      z_array[i]    = zplaneList[i].z;
581   }
582
583   new G4Polyhedra(GenerateName(name),_startphi,_deltaphi,_numsides,numZPlanes,z_array,rmin_array,rmax_array);
584}
585
586G4QuadrangularFacet* G4GDMLSolids::quadrangularRead(const xercesc::DOMElement* const element) {
587
588   G4String v1;
589   G4String v2;
590   G4String v3;
591   G4String v4;
592   G4String type;
593
594   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
595   XMLSize_t attributeCount = attributes->getLength();
596
597   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
598
599      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
600
601      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
602
603      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
604
605      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
606      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
607
608      if (attribute_name=="vertex1") v1 = attribute_value; else
609      if (attribute_name=="vertex2") v2 = attribute_value; else
610      if (attribute_name=="vertex3") v3 = attribute_value; else
611      if (attribute_name=="vertex4") v4 = attribute_value; else
612      if (attribute_name=="type") type = attribute_value;
613   }
614
615   const G4ThreeVector* ptr1 = getPosition(GenerateName(v1));
616   const G4ThreeVector* ptr2 = getPosition(GenerateName(v2));
617   const G4ThreeVector* ptr3 = getPosition(GenerateName(v3));
618   const G4ThreeVector* ptr4 = getPosition(GenerateName(v4));
619
620   return new G4QuadrangularFacet(*ptr1,*ptr2,*ptr3,*ptr4,(type=="RELATIVE")?(RELATIVE):(ABSOLUTE));
621}
622
623G4String G4GDMLSolids::refRead(const xercesc::DOMElement* const element) {
624
625   G4String ref;
626
627   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
628   XMLSize_t attributeCount = attributes->getLength();
629
630   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
631
632      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
633
634      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
635
636      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
637
638      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
639      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
640
641      if (attribute_name=="ref") ref = attribute_value;
642   }
643
644   return ref;
645}
646
647void G4GDMLSolids::reflectedSolidRead(const xercesc::DOMElement* const element) {
648
649   G4String name;
650   G4String lunit("1");
651   G4String aunit("1");
652   G4String solid;
653   G4String sx;
654   G4String sy;
655   G4String sz;
656   G4String rx;
657   G4String ry;
658   G4String rz;
659   G4String dx;
660   G4String dy;
661   G4String dz;
662
663   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
664   XMLSize_t attributeCount = attributes->getLength();
665
666   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
667
668      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
669
670      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
671
672      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
673
674      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
675      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
676
677      if (attribute_name=="name" ) name  = attribute_value; else
678      if (attribute_name=="lunit") lunit = attribute_value; else
679      if (attribute_name=="aunit") aunit = attribute_value; else
680      if (attribute_name=="solid") solid = attribute_value; else
681      if (attribute_name=="sx") sx = attribute_value; else
682      if (attribute_name=="sy") sy = attribute_value; else
683      if (attribute_name=="sz") sz = attribute_value; else
684      if (attribute_name=="rx") rx = attribute_value; else
685      if (attribute_name=="ry") ry = attribute_value; else
686      if (attribute_name=="rz") rz = attribute_value; else
687      if (attribute_name=="dx") dx = attribute_value; else
688      if (attribute_name=="dy") dy = attribute_value; else
689      if (attribute_name=="dz") dz = attribute_value;
690   }
691
692   G4double _lunit = eval.Evaluate(lunit);
693   G4double _aunit = eval.Evaluate(aunit);
694
695   G4double _sx = eval.Evaluate(sx);
696   G4double _sy = eval.Evaluate(sy);
697   G4double _sz = eval.Evaluate(sz);
698   G4double _rx = eval.Evaluate(rx)*_aunit;
699   G4double _ry = eval.Evaluate(ry)*_aunit;
700   G4double _rz = eval.Evaluate(rz)*_aunit;
701   G4double _dx = eval.Evaluate(dx)*_lunit;
702   G4double _dy = eval.Evaluate(dy)*_lunit;
703   G4double _dz = eval.Evaluate(dz)*_lunit;
704
705   G4VSolid* solidPtr = getSolid(GenerateName(solid));
706
707   G4RotationMatrix rot;
708   
709   rot.rotateX(_rx);
710   rot.rotateY(_ry);
711   rot.rotateZ(_rz);
712
713   G4ThreeVector trans(_dx,_dy,_dz);
714
715   G4Scale3D scale(_sx,_sy,_sz);
716
717   G4Transform3D transform(rot,trans);
718   transform = transform*scale;
719         
720   new G4ReflectedSolid(GenerateName(name),solidPtr,transform);
721}
722
723G4ExtrudedSolid::ZSection G4GDMLSolids::sectionRead(const xercesc::DOMElement* const element,G4double _lunit) {
724
725   G4String zPosition;
726   G4String xOffset;
727   G4String yOffset;
728   G4String scalingFactor;
729
730   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
731   XMLSize_t attributeCount = attributes->getLength();
732
733   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
734
735      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
736
737      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
738
739      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
740
741      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
742      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
743
744      if (attribute_name=="zPosition") zPosition = attribute_value; else
745      if (attribute_name=="xOffset") xOffset = attribute_value; else
746      if (attribute_name=="yOffset") yOffset = attribute_value; else
747      if (attribute_name=="scalingFactor") scalingFactor = attribute_value;
748   }
749
750   G4double _zPosition = eval.Evaluate(zPosition)*_lunit;
751   G4double _xOffset = eval.Evaluate(xOffset)*_lunit;
752   G4double _yOffset = eval.Evaluate(yOffset)*_lunit;
753   G4double _scalingFactor = eval.Evaluate(scalingFactor);
754
755   return G4ExtrudedSolid::ZSection(_zPosition,G4TwoVector(_xOffset,_yOffset),_scalingFactor);
756}
757
758void G4GDMLSolids::sphereRead(const xercesc::DOMElement* const element) {
759
760   G4String name;
761   G4String lunit("1");
762   G4String aunit("1");
763   G4String rmin;
764   G4String rmax;
765   G4String startphi;
766   G4String deltaphi;
767   G4String starttheta;
768   G4String deltatheta;
769
770   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
771   XMLSize_t attributeCount = attributes->getLength();
772
773   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
774
775      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
776
777      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
778
779      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
780
781      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
782      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
783
784      if (attribute_name=="name") name = attribute_value; else
785      if (attribute_name=="lunit") lunit = attribute_value; else
786      if (attribute_name=="aunit") aunit = attribute_value; else
787      if (attribute_name=="rmin") rmin = attribute_value; else
788      if (attribute_name=="rmax") rmax = attribute_value; else
789      if (attribute_name=="startphi") startphi = attribute_value; else
790      if (attribute_name=="deltaphi") deltaphi = attribute_value; else
791      if (attribute_name=="starttheta") starttheta = attribute_value; else
792      if (attribute_name=="deltatheta") deltatheta = attribute_value;
793   }
794
795   G4double _lunit = eval.Evaluate(lunit);
796   G4double _aunit = eval.Evaluate(aunit);
797
798   G4double _rmin = eval.Evaluate(rmin)*_lunit;
799   G4double _rmax = eval.Evaluate(rmax)*_lunit;
800   G4double _startphi = eval.Evaluate(startphi)*_aunit;
801   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;
802   G4double _starttheta = eval.Evaluate(starttheta)*_aunit;
803   G4double _deltatheta = eval.Evaluate(deltatheta)*_aunit;
804
805   new G4Sphere(GenerateName(name),_rmin,_rmax,_startphi,_deltaphi,_starttheta,_deltatheta);
806}
807
808void G4GDMLSolids::tessellatedRead(const xercesc::DOMElement* const element) {
809
810   G4String name;
811
812   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
813   XMLSize_t attributeCount = attributes->getLength();
814
815   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
816
817      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
818
819      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
820
821      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
822
823      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
824      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
825
826      if (attribute_name=="name") name = attribute_value;
827   }
828   
829   G4TessellatedSolid *tessellated = new G4TessellatedSolid(GenerateName(name));
830
831   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
832
833      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
834
835      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
836
837      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
838
839      if (tag=="triangular") tessellated->AddFacet(triangularRead(child)); else
840      if (tag=="quadrangular") tessellated->AddFacet(quadrangularRead(child));
841   }
842
843   tessellated->SetSolidClosed(true);
844}
845
846void G4GDMLSolids::tetRead(const xercesc::DOMElement* const element) {
847
848   G4String name;
849   G4String vertex1;
850   G4String vertex2;
851   G4String vertex3;
852   G4String vertex4;
853
854   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
855   XMLSize_t attributeCount = attributes->getLength();
856
857   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
858
859      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
860
861      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
862
863      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
864
865      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
866      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
867
868      if (attribute_name=="name") name = attribute_value; else
869      if (attribute_name=="vertex1") vertex1 = attribute_value; else
870      if (attribute_name=="vertex2") vertex2 = attribute_value; else
871      if (attribute_name=="vertex3") vertex3 = attribute_value; else
872      if (attribute_name=="vertex4") vertex4 = attribute_value;
873   }
874   
875   const G4ThreeVector* ptr1 = getPosition(GenerateName(vertex1));
876   const G4ThreeVector* ptr2 = getPosition(GenerateName(vertex2));
877   const G4ThreeVector* ptr3 = getPosition(GenerateName(vertex3));
878   const G4ThreeVector* ptr4 = getPosition(GenerateName(vertex4));
879
880   new G4Tet(GenerateName(name),*ptr1,*ptr2,*ptr3,*ptr4);
881}
882
883void G4GDMLSolids::torusRead(const xercesc::DOMElement* const element) {
884
885   G4String name;
886   G4String lunit("1");
887   G4String aunit("1");
888   G4String rmin;
889   G4String rmax;
890   G4String rtor;
891   G4String startphi;
892   G4String deltaphi;
893
894   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
895   XMLSize_t attributeCount = attributes->getLength();
896
897   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
898
899      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
900
901      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
902
903      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
904
905      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
906      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
907
908      if (attribute_name=="name") name = attribute_value; else
909      if (attribute_name=="lunit") lunit = attribute_value; else
910      if (attribute_name=="aunit") aunit = attribute_value; else
911      if (attribute_name=="rmin") rmin = attribute_value; else
912      if (attribute_name=="rmax") rmax = attribute_value; else
913      if (attribute_name=="rtor") rtor = attribute_value; else
914      if (attribute_name=="startphi") startphi = attribute_value; else
915      if (attribute_name=="deltaphi") deltaphi = attribute_value;
916   }
917
918   G4double _lunit = eval.Evaluate(lunit);
919   G4double _aunit = eval.Evaluate(aunit);
920
921   G4double _rmin = eval.Evaluate(rmin)*_lunit;
922   G4double _rmax = eval.Evaluate(rmax)*_lunit;
923   G4double _rtor = eval.Evaluate(rtor)*_lunit;
924   G4double _startphi = eval.Evaluate(startphi)*_aunit;
925   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;
926
927   new G4Torus(GenerateName(name),_rmin,_rmax,_rtor,_startphi,_deltaphi);
928}
929
930void G4GDMLSolids::trapRead(const xercesc::DOMElement* const element) {
931
932   G4String name;
933   G4String lunit("1");
934   G4String aunit("1");
935   G4String z;
936   G4String theta;
937   G4String phi;
938   G4String y1;
939   G4String x1;
940   G4String x2;
941   G4String alpha1;
942   G4String y2;
943   G4String x3;
944   G4String x4;
945   G4String alpha2;
946
947   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
948   XMLSize_t attributeCount = attributes->getLength();
949
950   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
951
952      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
953
954      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
955
956      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
957
958      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
959      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
960
961      if (attribute_name=="name") name = attribute_value; else
962      if (attribute_name=="lunit") lunit = attribute_value; else
963      if (attribute_name=="aunit") aunit = attribute_value; else
964      if (attribute_name=="z") z = attribute_value; else
965      if (attribute_name=="theta") theta  = attribute_value; else
966      if (attribute_name=="phi") phi = attribute_value; else
967      if (attribute_name=="y1") y1 = attribute_value; else
968      if (attribute_name=="x1") x1 = attribute_value; else
969      if (attribute_name=="x2") x2 = attribute_value; else
970      if (attribute_name=="alpha1") alpha1 = attribute_value; else
971      if (attribute_name=="y2") y2 = attribute_value; else
972      if (attribute_name=="x3") x3 = attribute_value; else
973      if (attribute_name=="x4") x4 = attribute_value; else
974      if (attribute_name=="alpha2") alpha2 = attribute_value;
975   }
976
977   G4double _lunit = eval.Evaluate(lunit);
978   G4double _aunit = eval.Evaluate(aunit);
979
980   G4double _z = eval.Evaluate(z)*_lunit;
981   G4double _theta  = eval.Evaluate(theta)*_aunit;
982   G4double _phi = eval.Evaluate(phi)*_aunit;
983   G4double _y1 = eval.Evaluate(y1)*_lunit;
984   G4double _x1 = eval.Evaluate(x1)*_lunit;
985   G4double _x2 = eval.Evaluate(x2)*_lunit;
986   G4double _alpha1 = eval.Evaluate(alpha1)*_aunit;
987   G4double _y2 = eval.Evaluate(y2)*_lunit;
988   G4double _x3 = eval.Evaluate(x3)*_lunit;
989   G4double _x4 = eval.Evaluate(x4)*_lunit;
990   G4double _alpha2 = eval.Evaluate(alpha2)*_aunit;
991
992   _z *= 0.5;
993   _y1 *= 0.5;
994   _x1 *= 0.5;
995   _x2 *= 0.5;
996   _y2 *= 0.5;
997   _x3 *= 0.5;
998   _x4 *= 0.5;
999
1000   new G4Trap(GenerateName(name),_z,_theta,_phi,_y1,_x1,_x2,_alpha1,_y2,_x3,_x4,_alpha2);
1001}
1002
1003void G4GDMLSolids::trdRead(const xercesc::DOMElement* const element) {
1004
1005   G4String name;
1006   G4String lunit("1");
1007   G4String x1;
1008   G4String x2;
1009   G4String y1;
1010   G4String y2;
1011   G4String z;
1012
1013   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1014   XMLSize_t attributeCount = attributes->getLength();
1015
1016   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1017
1018      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1019
1020      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1021
1022      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1023
1024      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1025      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1026
1027      if (attribute_name=="name") name = attribute_value; else
1028      if (attribute_name=="lunit") lunit = attribute_value; else
1029      if (attribute_name=="x1") x1 = attribute_value; else
1030      if (attribute_name=="x2") x2 = attribute_value; else
1031      if (attribute_name=="y1") y1 = attribute_value; else
1032      if (attribute_name=="y2") y2 = attribute_value; else
1033      if (attribute_name=="z") z = attribute_value;
1034   }
1035
1036   G4double _lunit = eval.Evaluate(lunit);
1037
1038   G4double _x1 = eval.Evaluate(x1)*_lunit;
1039   G4double _x2 = eval.Evaluate(x2)*_lunit;
1040   G4double _y1 = eval.Evaluate(y1)*_lunit;
1041   G4double _y2 = eval.Evaluate(y2)*_lunit;
1042   G4double _z  = eval.Evaluate(z)*_lunit;
1043
1044   _x1 *= 0.5;
1045   _x2 *= 0.5;
1046   _y1 *= 0.5;
1047   _y2 *= 0.5;
1048   _z *= 0.5;
1049
1050   new G4Trd(GenerateName(name),_x1,_x2,_y1,_y2,_z);
1051}
1052
1053G4TriangularFacet* G4GDMLSolids::triangularRead(const xercesc::DOMElement* const element) {
1054
1055   G4String v1;
1056   G4String v2;
1057   G4String v3;
1058   G4String type;
1059
1060   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1061   XMLSize_t attributeCount = attributes->getLength();
1062
1063   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1064
1065      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1066
1067      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1068
1069      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1070
1071      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1072      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1073
1074      if (attribute_name=="vertex1") v1 = attribute_value; else
1075      if (attribute_name=="vertex2") v2 = attribute_value; else
1076      if (attribute_name=="vertex3") v3 = attribute_value; else
1077      if (attribute_name=="type") type = attribute_value;
1078   }
1079
1080   const G4ThreeVector* ptr1 = getPosition(GenerateName(v1));
1081   const G4ThreeVector* ptr2 = getPosition(GenerateName(v2));
1082   const G4ThreeVector* ptr3 = getPosition(GenerateName(v3));
1083
1084   return new G4TriangularFacet(*ptr1,*ptr2,*ptr3,(type=="RELATIVE")?(RELATIVE):(ABSOLUTE));
1085}
1086
1087void G4GDMLSolids::tubeRead(const xercesc::DOMElement* const element) {
1088
1089   G4String name;
1090   G4String lunit("1");
1091   G4String aunit("1");
1092   G4String rmin;
1093   G4String rmax;
1094   G4String z;
1095   G4String startphi;
1096   G4String deltaphi;
1097
1098   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1099   XMLSize_t attributeCount = attributes->getLength();
1100
1101   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1102
1103      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1104
1105      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1106
1107      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1108
1109      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1110      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1111
1112      if (attribute_name=="name") name = attribute_value; else
1113      if (attribute_name=="lunit") lunit = attribute_value; else
1114      if (attribute_name=="aunit") aunit = attribute_value; else
1115      if (attribute_name=="rmin") rmin = attribute_value; else
1116      if (attribute_name=="rmax") rmax = attribute_value; else
1117      if (attribute_name=="z") z = attribute_value; else
1118      if (attribute_name=="startphi") startphi = attribute_value; else
1119      if (attribute_name=="deltaphi") deltaphi = attribute_value;
1120   }
1121
1122   G4double _lunit = eval.Evaluate(lunit);
1123   G4double _aunit = eval.Evaluate(aunit);
1124
1125   G4double _rmin = eval.Evaluate(rmin)*_lunit;
1126   G4double _rmax = eval.Evaluate(rmax)*_lunit;
1127   G4double _z = eval.Evaluate(z)*_lunit;
1128   G4double _startphi = eval.Evaluate(startphi)*_aunit;
1129   G4double _deltaphi = eval.Evaluate(deltaphi)*_aunit;
1130
1131   _z *= 0.5;
1132
1133   new G4Tubs(GenerateName(name),_rmin,_rmax,_z,_startphi,_deltaphi);
1134}
1135
1136G4TwoVector G4GDMLSolids::twoDimVertexRead(const xercesc::DOMElement* const element,G4double _lunit) {
1137
1138   G4String x;
1139   G4String y;
1140
1141   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1142   XMLSize_t attributeCount = attributes->getLength();
1143
1144   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1145
1146      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1147
1148      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1149
1150      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1151
1152      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1153      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1154
1155      if (attribute_name=="x") x = attribute_value; else
1156      if (attribute_name=="y") y = attribute_value;
1157   }
1158
1159   G4double _x = eval.Evaluate(x)*_lunit;
1160   G4double _y = eval.Evaluate(y)*_lunit;
1161
1162   return G4TwoVector(_x,_y);
1163}
1164
1165void G4GDMLSolids::xtruRead(const xercesc::DOMElement* const element) {
1166
1167   G4String name;
1168   G4String lunit("1");
1169
1170   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1171   XMLSize_t attributeCount = attributes->getLength();
1172
1173   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1174
1175      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1176
1177      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1178
1179      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1180
1181      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1182      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1183
1184      if (attribute_name=="name") name = attribute_value; else
1185      if (attribute_name=="lunit") lunit = attribute_value;
1186   }
1187
1188   G4double _lunit = eval.Evaluate(lunit);
1189
1190   std::vector<G4TwoVector> twoDimVertexList;
1191   std::vector<G4ExtrudedSolid::ZSection> sectionList;
1192
1193   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
1194
1195      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
1196
1197      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
1198
1199      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
1200
1201      if (tag=="twoDimVertex") twoDimVertexList.push_back(twoDimVertexRead(child,_lunit)); else
1202      if (tag=="section") sectionList.push_back(sectionRead(child,_lunit));     
1203   }
1204
1205   new G4ExtrudedSolid(GenerateName(name),twoDimVertexList,sectionList);
1206}
1207
1208G4GDMLSolids::zplaneType G4GDMLSolids::zplaneRead(const xercesc::DOMElement* const element,G4double _lunit) {
1209
1210   G4String rmin;
1211   G4String rmax;
1212   G4String z;
1213
1214   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1215   XMLSize_t attributeCount = attributes->getLength();
1216
1217   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1218
1219      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1220
1221      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1222
1223      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1224
1225      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1226      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1227
1228      if (attribute_name=="rmin") rmin = attribute_value; else
1229      if (attribute_name=="rmax") rmax = attribute_value; else
1230      if (attribute_name=="z") z = attribute_value;
1231   }
1232
1233   zplaneType zplane;
1234
1235   zplane.rmin = eval.Evaluate(rmin)*_lunit;
1236   zplane.rmax = eval.Evaluate(rmax)*_lunit;
1237   zplane.z = eval.Evaluate(z)*_lunit;
1238
1239   return zplane;
1240}
1241
1242void G4GDMLSolids::solidsRead(const xercesc::DOMElement* const element) {
1243
1244   for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
1245
1246      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
1247
1248      const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
1249
1250      const G4String tag = xercesc::XMLString::transcode(child->getTagName());
1251
1252      if (tag=="box") boxRead(child); else
1253      if (tag=="cone") coneRead(child); else
1254      if (tag=="ellipsoid") ellipsoidRead(child); else
1255      if (tag=="eltube") eltubeRead(child); else
1256      if (tag=="hype") hypeRead(child); else
1257      if (tag=="loop") loopRead(child); else
1258      if (tag=="orb") orbRead(child); else
1259      if (tag=="para") paraRead(child); else
1260      if (tag=="polycone") polyconeRead(child); else
1261      if (tag=="polyhedra") polyhedraRead(child); else
1262      if (tag=="sphere") sphereRead(child); else
1263      if (tag=="reflectedSolid") reflectedSolidRead(child); else
1264      if (tag=="tessellated") tessellatedRead(child); else
1265      if (tag=="tet") tetRead(child); else
1266      if (tag=="torus") torusRead(child); else
1267      if (tag=="trap") trapRead(child); else
1268      if (tag=="trd") trdRead(child); else
1269      if (tag=="tube") tubeRead(child); else
1270      if (tag=="xtru") xtruRead(child); else
1271      if (tag=="intersection") booleanRead(child,INTERSECTION); else
1272      if (tag=="subtraction") booleanRead(child,SUBTRACTION); else
1273      if (tag=="union") booleanRead(child,UNION); else
1274      G4Exception("GDML: Unknown tag in solids: "+tag);
1275   }
1276}
1277
1278G4ThreeVector G4GDMLSolids::positionRead(const xercesc::DOMElement* const element) {
1279
1280   G4String unit("1");
1281   G4String x;
1282   G4String y;
1283   G4String z;
1284
1285   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1286   XMLSize_t attributeCount = attributes->getLength();
1287
1288   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1289
1290      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1291
1292      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1293
1294      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1295
1296      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1297      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1298
1299      if (attribute_name=="unit") unit = attribute_value; else
1300      if (attribute_name=="x") x  = attribute_value; else
1301      if (attribute_name=="y") y = attribute_value; else
1302      if (attribute_name=="z") z = attribute_value;
1303   }
1304
1305   G4double _unit = eval.Evaluate(unit);
1306
1307   G4double _x = eval.Evaluate(x)*_unit;
1308   G4double _y = eval.Evaluate(y)*_unit;
1309   G4double _z = eval.Evaluate(z)*_unit;
1310   
1311   return G4ThreeVector(_x,_y,_z);
1312}
1313
1314G4ThreeVector G4GDMLSolids::rotationRead(const xercesc::DOMElement* const element) {
1315
1316   G4String unit("1");
1317   G4String x;
1318   G4String y;
1319   G4String z;
1320
1321   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1322   XMLSize_t attributeCount = attributes->getLength();
1323
1324   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1325
1326      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1327
1328      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1329
1330      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1331
1332      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1333      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1334
1335      if (attribute_name=="unit") unit = attribute_value; else
1336      if (attribute_name=="x") x = attribute_value; else
1337      if (attribute_name=="y") y = attribute_value; else
1338      if (attribute_name=="z") z = attribute_value;
1339   }
1340
1341   G4double _unit = eval.Evaluate(unit);
1342
1343   G4double _x = eval.Evaluate(x)*_unit;
1344   G4double _y = eval.Evaluate(y)*_unit;
1345   G4double _z = eval.Evaluate(z)*_unit;
1346   
1347   return G4ThreeVector(_x,_y,_z);
1348}
1349
1350G4ThreeVector G4GDMLSolids::scaleRead(const xercesc::DOMElement* const element) {
1351
1352   G4String x;
1353   G4String y;
1354   G4String z;
1355
1356   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
1357   XMLSize_t attributeCount = attributes->getLength();
1358
1359   for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
1360
1361      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
1362
1363      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
1364
1365      const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
1366
1367      const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
1368      const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
1369
1370      if (attribute_name=="x") x = attribute_value; else
1371      if (attribute_name=="y") y = attribute_value; else
1372      if (attribute_name=="z") z = attribute_value;
1373   }
1374
1375   G4double _x = eval.Evaluate(x);
1376   G4double _y = eval.Evaluate(y);
1377   G4double _z = eval.Evaluate(z);
1378   
1379   return G4ThreeVector(_x,_y,_z);
1380}
1381
1382G4VSolid* G4GDMLSolids::getSolid(const G4String& ref) const {
1383
1384   G4VSolid* solidPtr = G4SolidStore::GetInstance()->GetSolid(ref,false);
1385
1386   if (!solidPtr) G4Exception("GDML: Referenced solid '"+ref+"' was not found!");
1387
1388   return solidPtr;
1389}
Note: See TracBrowser for help on using the repository browser.