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