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

Last change on this file since 1287 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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