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

Last change on this file since 1330 was 1315, checked in by garnier, 15 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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