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