Geant4_10
G4GDMLReadParamvol.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: G4GDMLReadParamvol.cc 77913 2013-11-29 10:59:07Z gcosmo $
27 //
28 // class G4GDMLReadParamvol Implementation
29 //
30 // History:
31 // - Created. Zoltan Torzsok, November 2007
32 // -------------------------------------------------------------------------
33 
34 #include "G4GDMLReadParamvol.hh"
35 #include "G4GDMLReadSolids.hh"
36 #include "G4LogicalVolume.hh"
37 #include "G4PVParameterised.hh"
38 #include "G4PVPlacement.hh"
39 #include "G4VPhysicalVolume.hh"
40 
42  : G4GDMLReadSetup(), parameterisation(0)
43 {
44 }
45 
47 {
48 }
49 
51 Box_dimensionsRead( const xercesc::DOMElement* const element,
53 {
54  G4double lunit = 1.0;
55 
56  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
57  XMLSize_t attributeCount = attributes->getLength();
58 
59  for (XMLSize_t attribute_index=0;
60  attribute_index<attributeCount; attribute_index++)
61  {
62  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
63 
64  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
65  { continue; }
66 
67  const xercesc::DOMAttr* const attribute
68  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
69  if (!attribute)
70  {
71  G4Exception("G4GDMLReadParamvol::Box_dimensionsRead()",
72  "InvalidRead", FatalException, "No attribute found!");
73  return;
74  }
75  const G4String attName = Transcode(attribute->getName());
76  const G4String attValue = Transcode(attribute->getValue());
77 
78  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
79  if (attName=="x") { parameter.dimension[0] = eval.Evaluate(attValue); } else
80  if (attName=="y") { parameter.dimension[1] = eval.Evaluate(attValue); } else
81  if (attName=="z") { parameter.dimension[2] = eval.Evaluate(attValue); }
82  }
83 
84  parameter.dimension[0] *= 0.5*lunit;
85  parameter.dimension[1] *= 0.5*lunit;
86  parameter.dimension[2] *= 0.5*lunit;
87 }
88 
90 Trd_dimensionsRead( const xercesc::DOMElement* const element,
92 {
93  G4double lunit = 1.0;
94 
95  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
96  XMLSize_t attributeCount = attributes->getLength();
97 
98  for (XMLSize_t attribute_index=0;
99  attribute_index<attributeCount; attribute_index++)
100  {
101  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
102 
103  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
104  { continue; }
105 
106  const xercesc::DOMAttr* const attribute
107  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
108  if (!attribute)
109  {
110  G4Exception("G4GDMLReadParamvol::Trd_dimensionsRead()",
111  "InvalidRead", FatalException, "No attribute found!");
112  return;
113  }
114  const G4String attName = Transcode(attribute->getName());
115  const G4String attValue = Transcode(attribute->getValue());
116 
117  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
118  if (attName=="x1") { parameter.dimension[0]=eval.Evaluate(attValue); } else
119  if (attName=="x2") { parameter.dimension[1]=eval.Evaluate(attValue); } else
120  if (attName=="y1") { parameter.dimension[2]=eval.Evaluate(attValue); } else
121  if (attName=="y2") { parameter.dimension[3]=eval.Evaluate(attValue); } else
122  if (attName=="z") { parameter.dimension[4]=eval.Evaluate(attValue); }
123  }
124 
125  parameter.dimension[0] *= 0.5*lunit;
126  parameter.dimension[1] *= 0.5*lunit;
127  parameter.dimension[2] *= 0.5*lunit;
128  parameter.dimension[3] *= 0.5*lunit;
129  parameter.dimension[4] *= 0.5*lunit;
130 }
131 
133 Trap_dimensionsRead( const xercesc::DOMElement* const element,
135 {
136  G4double lunit = 1.0;
137  G4double aunit = 1.0;
138 
139  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
140  XMLSize_t attributeCount = attributes->getLength();
141 
142  for (XMLSize_t attribute_index=0;
143  attribute_index<attributeCount; attribute_index++)
144  {
145  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
146 
147  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
148  { continue; }
149 
150  const xercesc::DOMAttr* const attribute
151  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
152  if (!attribute)
153  {
154  G4Exception("G4GDMLReadParamvol::Trap_dimensionsRead()",
155  "InvalidRead", FatalException, "No attribute found!");
156  return;
157  }
158  const G4String attName = Transcode(attribute->getName());
159  const G4String attValue = Transcode(attribute->getValue());
160 
161  if (attName=="lunit")
162  { lunit = eval.Evaluate(attValue); } else
163  if (attName=="aunit")
164  { aunit = eval.Evaluate(attValue); } else
165  if (attName=="z")
166  { parameter.dimension[0] = eval.Evaluate(attValue); } else
167  if (attName=="theta")
168  { parameter.dimension[1] = eval.Evaluate(attValue); } else
169  if (attName=="phi")
170  { parameter.dimension[2] = eval.Evaluate(attValue); } else
171  if (attName=="y1")
172  { parameter.dimension[3] = eval.Evaluate(attValue); } else
173  if (attName=="x1")
174  { parameter.dimension[4] = eval.Evaluate(attValue); } else
175  if (attName=="x2")
176  { parameter.dimension[5] = eval.Evaluate(attValue); } else
177  if (attName=="alpha1")
178  { parameter.dimension[6] = eval.Evaluate(attValue); } else
179  if (attName=="y2")
180  { parameter.dimension[7] = eval.Evaluate(attValue); } else
181  if (attName=="x3")
182  { parameter.dimension[8] = eval.Evaluate(attValue); } else
183  if (attName=="x4")
184  { parameter.dimension[9] = eval.Evaluate(attValue); } else
185  if (attName=="alpha2")
186  { parameter.dimension[10] = eval.Evaluate(attValue); }
187  }
188 
189  parameter.dimension[0] *= 0.5*lunit;
190  parameter.dimension[1] *= aunit;
191  parameter.dimension[2] *= aunit;
192  parameter.dimension[3] *= 0.5*lunit;
193  parameter.dimension[4] *= 0.5*lunit;
194  parameter.dimension[5] *= 0.5*lunit;
195  parameter.dimension[6] *= aunit;
196  parameter.dimension[7] *= 0.5*lunit;
197  parameter.dimension[8] *= 0.5*lunit;
198  parameter.dimension[9] *= 0.5*lunit;
199  parameter.dimension[10] *= aunit;
200 }
201 
203 Tube_dimensionsRead( const xercesc::DOMElement* const element,
205 {
206  G4double lunit = 1.0;
207  G4double aunit = 1.0;
208 
209  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
210  XMLSize_t attributeCount = attributes->getLength();
211 
212  for (XMLSize_t attribute_index=0;
213  attribute_index<attributeCount; attribute_index++)
214  {
215  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
216 
217  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
218  { continue; }
219 
220  const xercesc::DOMAttr* const attribute
221  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
222  if (!attribute)
223  {
224  G4Exception("G4GDMLReadParamvol::Tube_dimensionsRead()",
225  "InvalidRead", FatalException, "No attribute found!");
226  return;
227  }
228  const G4String attName = Transcode(attribute->getName());
229  const G4String attValue = Transcode(attribute->getValue());
230 
231  if (attName=="lunit")
232  { lunit = eval.Evaluate(attValue); } else
233  if (attName=="aunit")
234  { aunit = eval.Evaluate(attValue); } else
235  if (attName=="InR")
236  { parameter.dimension[0] = eval.Evaluate(attValue); } else
237  if (attName=="OutR")
238  { parameter.dimension[1] = eval.Evaluate(attValue); } else
239  if (attName=="hz")
240  { parameter.dimension[2] = eval.Evaluate(attValue); } else
241  if (attName=="StartPhi")
242  { parameter.dimension[3] = eval.Evaluate(attValue); } else
243  if (attName=="DeltaPhi")
244  { parameter.dimension[4] = eval.Evaluate(attValue); }
245  }
246 
247  parameter.dimension[0] *= lunit;
248  parameter.dimension[1] *= lunit;
249  parameter.dimension[2] *= 0.5*lunit;
250  parameter.dimension[3] *= aunit;
251  parameter.dimension[4] *= aunit;
252 }
253 
255 Cone_dimensionsRead( const xercesc::DOMElement* const element,
257 {
258  G4double lunit = 1.0;
259  G4double aunit = 1.0;
260 
261  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
262  XMLSize_t attributeCount = attributes->getLength();
263 
264  for (XMLSize_t attribute_index=0;
265  attribute_index<attributeCount; attribute_index++)
266  {
267  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
268 
269  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
270  { continue; }
271 
272  const xercesc::DOMAttr* const attribute
273  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
274  if (!attribute)
275  {
276  G4Exception("G4GDMLReadParamvol::Cone_dimensionsRead()",
277  "InvalidRead", FatalException, "No attribute found!");
278  return;
279  }
280  const G4String attName = Transcode(attribute->getName());
281  const G4String attValue = Transcode(attribute->getValue());
282 
283  if (attName=="lunit")
284  { lunit = eval.Evaluate(attValue); } else
285  if (attName=="aunit")
286  { aunit = eval.Evaluate(attValue); } else
287  if (attName=="rmin1")
288  { parameter.dimension[0] = eval.Evaluate(attValue); } else
289  if (attName=="rmax1")
290  { parameter.dimension[1] = eval.Evaluate(attValue); } else
291  if (attName=="rmin2")
292  { parameter.dimension[2] = eval.Evaluate(attValue); } else
293  if (attName=="rmax2")
294  { parameter.dimension[3] = eval.Evaluate(attValue); } else
295  if (attName=="z")
296  { parameter.dimension[4] = eval.Evaluate(attValue); } else
297  if (attName=="startphi")
298  { parameter.dimension[5] = eval.Evaluate(attValue); } else
299  if (attName=="deltaphi")
300  { parameter.dimension[6] = eval.Evaluate(attValue); }
301  }
302 
303  parameter.dimension[0] *= lunit;
304  parameter.dimension[1] *= lunit;
305  parameter.dimension[2] *= lunit;
306  parameter.dimension[3] *= lunit;
307  parameter.dimension[4] *= 0.5*lunit;
308  parameter.dimension[5] *= aunit;
309  parameter.dimension[6] *= aunit;
310 }
311 
313 Sphere_dimensionsRead( const xercesc::DOMElement* const element,
315 {
316  G4double lunit = 1.0;
317  G4double aunit = 1.0;
318 
319  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
320  XMLSize_t attributeCount = attributes->getLength();
321 
322  for (XMLSize_t attribute_index=0;
323  attribute_index<attributeCount; attribute_index++)
324  {
325  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
326 
327  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
328  { continue; }
329 
330  const xercesc::DOMAttr* const attribute
331  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
332  if (!attribute)
333  {
334  G4Exception("G4GDMLReadParamvol::Sphere_dimensionsRead()",
335  "InvalidRead", FatalException, "No attribute found!");
336  return;
337  }
338  const G4String attName = Transcode(attribute->getName());
339  const G4String attValue = Transcode(attribute->getValue());
340 
341  if (attName=="lunit")
342  { lunit = eval.Evaluate(attValue); } else
343  if (attName=="aunit")
344  { aunit = eval.Evaluate(attValue); } else
345  if (attName=="rmin")
346  { parameter.dimension[0] = eval.Evaluate(attValue); } else
347  if (attName=="rmax")
348  { parameter.dimension[1] = eval.Evaluate(attValue); } else
349  if (attName=="startphi")
350  { parameter.dimension[2] = eval.Evaluate(attValue); } else
351  if (attName=="deltaphi")
352  { parameter.dimension[3] = eval.Evaluate(attValue); } else
353  if (attName=="starttheta")
354  { parameter.dimension[4] = eval.Evaluate(attValue); } else
355  if (attName=="deltatheta")
356  { parameter.dimension[5] = eval.Evaluate(attValue); }
357  }
358 
359  parameter.dimension[0] *= lunit;
360  parameter.dimension[1] *= lunit;
361  parameter.dimension[2] *= aunit;
362  parameter.dimension[3] *= aunit;
363  parameter.dimension[4] *= aunit;
364  parameter.dimension[5] *= aunit;
365 }
366 
368 Orb_dimensionsRead( const xercesc::DOMElement* const element,
370 {
371  G4double lunit = 1.0;
372 
373  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
374  XMLSize_t attributeCount = attributes->getLength();
375 
376  for (XMLSize_t attribute_index=0;
377  attribute_index<attributeCount; attribute_index++)
378  {
379  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
380 
381  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
382  { continue; }
383 
384  const xercesc::DOMAttr* const attribute
385  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
386  if (!attribute)
387  {
388  G4Exception("G4GDMLReadParamvol::Orb_dimensionsRead()",
389  "InvalidRead", FatalException, "No attribute found!");
390  return;
391  }
392  const G4String attName = Transcode(attribute->getName());
393  const G4String attValue = Transcode(attribute->getValue());
394 
395  if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else
396  if (attName=="r") { parameter.dimension[0] = eval.Evaluate(attValue); }
397  }
398 
399  parameter.dimension[0] *= lunit;
400 }
401 
403 Torus_dimensionsRead( const xercesc::DOMElement* const element,
405 {
406  G4double lunit = 1.0;
407  G4double aunit = 1.0;
408 
409  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
410  XMLSize_t attributeCount = attributes->getLength();
411 
412  for (XMLSize_t attribute_index=0;
413  attribute_index<attributeCount; attribute_index++)
414  {
415  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
416 
417  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
418  { continue; }
419 
420  const xercesc::DOMAttr* const attribute
421  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
422  if (!attribute)
423  {
424  G4Exception("G4GDMLReadParamvol::Torus_dimensionsRead()",
425  "InvalidRead", FatalException, "No attribute found!");
426  return;
427  }
428  const G4String attName = Transcode(attribute->getName());
429  const G4String attValue = Transcode(attribute->getValue());
430 
431  if (attName=="lunit")
432  { lunit = eval.Evaluate(attValue); } else
433  if (attName=="aunit")
434  { aunit = eval.Evaluate(attValue); } else
435  if (attName=="rmin")
436  { parameter.dimension[0] = eval.Evaluate(attValue); } else
437  if (attName=="rmax")
438  { parameter.dimension[1] = eval.Evaluate(attValue); } else
439  if (attName=="rtor")
440  { parameter.dimension[2] = eval.Evaluate(attValue); } else
441  if (attName=="startphi")
442  { parameter.dimension[3] = eval.Evaluate(attValue); } else
443  if (attName=="deltaphi")
444  { parameter.dimension[4] = eval.Evaluate(attValue); }
445  }
446 
447  parameter.dimension[0] *= lunit;
448  parameter.dimension[1] *= lunit;
449  parameter.dimension[2] *= lunit;
450  parameter.dimension[3] *= aunit;
451  parameter.dimension[4] *= aunit;
452 }
453 
455 Ellipsoid_dimensionsRead( const xercesc::DOMElement* const element,
457 {
458  G4double lunit = 1.0;
459  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
460  XMLSize_t attributeCount = attributes->getLength();
461 
462  for (XMLSize_t attribute_index=0;
463  attribute_index<attributeCount; attribute_index++)
464  {
465  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
466 
467  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
468  { continue; }
469 
470  const xercesc::DOMAttr* const attribute
471  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
472  if (!attribute)
473  {
474  G4Exception("G4GDMLReadParamvol::Ellipsoid_dimensionsRead()",
475  "InvalidRead", FatalException, "No attribute found!");
476  return;
477  }
478  const G4String attName = Transcode(attribute->getName());
479  const G4String attValue = Transcode(attribute->getValue());
480 
481  if (attName=="lunit")
482  { lunit = eval.Evaluate(attValue); } else
483  if (attName=="ax")
484  { parameter.dimension[0] = eval.Evaluate(attValue); } else
485  if (attName=="by")
486  { parameter.dimension[1] = eval.Evaluate(attValue); } else
487  if (attName=="cz")
488  { parameter.dimension[2] = eval.Evaluate(attValue); } else
489  if (attName=="zcut1")
490  { parameter.dimension[3] = eval.Evaluate(attValue); } else
491  if (attName=="zcut2")
492  { parameter.dimension[4] = eval.Evaluate(attValue); }
493  }
494 
495  parameter.dimension[0] *= lunit;
496  parameter.dimension[1] *= lunit;
497  parameter.dimension[2] *= lunit;
498  parameter.dimension[3] *= lunit;
499  parameter.dimension[4] *= lunit;
500 }
501 
503 Para_dimensionsRead( const xercesc::DOMElement* const element,
505 {
506  G4double lunit = 1.0;
507  G4double aunit = 1.0;
508 
509  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
510  XMLSize_t attributeCount = attributes->getLength();
511 
512  for (XMLSize_t attribute_index=0;
513  attribute_index<attributeCount; attribute_index++)
514  {
515  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
516 
517  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
518  { continue; }
519 
520  const xercesc::DOMAttr* const attribute
521  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
522  if (!attribute)
523  {
524  G4Exception("G4GDMLReadParamvol::Para_dimensionsRead()",
525  "InvalidRead", FatalException, "No attribute found!");
526  return;
527  }
528  const G4String attName = Transcode(attribute->getName());
529  const G4String attValue = Transcode(attribute->getValue());
530 
531  if (attName=="lunit")
532  { lunit = eval.Evaluate(attValue); } else
533  if (attName=="aunit")
534  { aunit = eval.Evaluate(attValue); } else
535  if (attName=="x")
536  { parameter.dimension[0] = eval.Evaluate(attValue); } else
537  if (attName=="y")
538  { parameter.dimension[1] = eval.Evaluate(attValue); } else
539  if (attName=="z")
540  { parameter.dimension[2] = eval.Evaluate(attValue); } else
541  if (attName=="alpha")
542  { parameter.dimension[3] = eval.Evaluate(attValue); } else
543  if (attName=="theta")
544  { parameter.dimension[4] = eval.Evaluate(attValue); } else
545  if (attName=="phi")
546  { parameter.dimension[5] = eval.Evaluate(attValue); }
547  }
548 
549  parameter.dimension[0] = 0.5*lunit;
550  parameter.dimension[1] = 0.5*lunit;
551  parameter.dimension[2] = 0.5*lunit;
552  parameter.dimension[3] = aunit;
553  parameter.dimension[4] = aunit;
554  parameter.dimension[5] = aunit;
555 }
556 
558 Hype_dimensionsRead( const xercesc::DOMElement* const element,
560 {
561  G4double lunit = 1.0;
562  G4double aunit = 1.0;
563 
564  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
565  XMLSize_t attributeCount = attributes->getLength();
566 
567  for (XMLSize_t attribute_index=0;
568  attribute_index<attributeCount; attribute_index++)
569  {
570  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
571 
572  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
573  { continue; }
574 
575  const xercesc::DOMAttr* const attribute
576  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
577  if (!attribute)
578  {
579  G4Exception("G4GDMLReadParamvol::Hype_dimensionsRead()",
580  "InvalidRead", FatalException, "No attribute found!");
581  return;
582  }
583  const G4String attName = Transcode(attribute->getName());
584  const G4String attValue = Transcode(attribute->getValue());
585 
586  if (attName=="lunit")
587  { lunit = eval.Evaluate(attValue); } else
588  if (attName=="aunit")
589  { aunit = eval.Evaluate(attValue); } else
590  if (attName=="rmin")
591  { parameter.dimension[0] = eval.Evaluate(attValue); } else
592  if (attName=="rmax")
593  { parameter.dimension[1] = eval.Evaluate(attValue); } else
594  if (attName=="inst")
595  { parameter.dimension[2] = eval.Evaluate(attValue); } else
596  if (attName=="outst")
597  { parameter.dimension[3] = eval.Evaluate(attValue); } else
598  if (attName=="z")
599  { parameter.dimension[4] = eval.Evaluate(attValue); }
600  }
601 
602  parameter.dimension[0] = lunit;
603  parameter.dimension[1] = lunit;
604  parameter.dimension[2] = aunit;
605  parameter.dimension[3] = aunit;
606  parameter.dimension[4] = 0.5*lunit;
607 }
608 
610 Polycone_dimensionsRead( const xercesc::DOMElement* const element,
612 {
613  G4double lunit = 1.0;
614  G4double aunit = 1.0;
615 
616  std::vector<zplaneType> zplaneList;
617 
618  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
619  XMLSize_t attributeCount = attributes->getLength();
620 
621  for (XMLSize_t attribute_index=0;
622  attribute_index<attributeCount; attribute_index++)
623  {
624  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
625 
626  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
627  { continue; }
628 
629  const xercesc::DOMAttr* const attribute
630  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
631  if (!attribute)
632  {
633  G4Exception("G4GDMLReadParamvol::Polycone_dimensionsRead()",
634  "InvalidRead", FatalException, "No attribute found!");
635  return;
636  }
637  const G4String attName = Transcode(attribute->getName());
638  const G4String attValue = Transcode(attribute->getValue());
639 
640  if (attName=="lunit")
641  { lunit = eval.Evaluate(attValue); } else
642  if (attName=="aunit")
643  { aunit = eval.Evaluate(attValue); } else
644  if (attName=="startPhi")
645  { parameter.dimension[0] = eval.Evaluate(attValue); } else
646  if (attName=="openPhi")
647  { parameter.dimension[1] = eval.Evaluate(attValue); } else
648  if (attName=="numRZ")
649  { parameter.dimension[2] = eval.Evaluate(attValue); }
650  }
651 
652  parameter.dimension[0] *= aunit;
653  parameter.dimension[1] *= aunit;
654 
655  for (xercesc::DOMNode* iter = element->getFirstChild();
656  iter != 0; iter = iter->getNextSibling())
657  {
658  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
659 
660  const xercesc::DOMElement* const child
661  = dynamic_cast<xercesc::DOMElement*>(iter);
662  if (!child)
663  {
664  G4Exception("G4GDMLReadParamVol::Polycone_dimensionsRead()",
665  "InvalidRead", FatalException, "No child found!");
666  return;
667  }
668  const G4String tag = Transcode(child->getTagName());
669 
670  if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
671  }
672 
673  G4int numZPlanes = zplaneList.size();
674  for (G4int i=0; i<numZPlanes; i++)
675  {
676  parameter.dimension[3+i*3] = zplaneList[i].rmin*lunit;
677  parameter.dimension[4+i*3] = zplaneList[i].rmax*lunit;
678  parameter.dimension[5+i*3] = zplaneList[i].z*lunit;
679  }
680 }
681 
683 Polyhedra_dimensionsRead( const xercesc::DOMElement* const element,
685 {
686  G4double lunit = 1.0;
687  G4double aunit = 1.0;
688 
689  std::vector<zplaneType> zplaneList;
690 
691  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
692  XMLSize_t attributeCount = attributes->getLength();
693 
694  for (XMLSize_t attribute_index=0;
695  attribute_index<attributeCount; attribute_index++)
696  {
697  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
698 
699  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
700  { continue; }
701 
702  const xercesc::DOMAttr* const attribute
703  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
704  if (!attribute)
705  {
706  G4Exception("G4GDMLReadParamvol::Polycone_dimensionsRead()",
707  "InvalidRead", FatalException, "No attribute found!");
708  return;
709  }
710  const G4String attName = Transcode(attribute->getName());
711  const G4String attValue = Transcode(attribute->getValue());
712 
713  if (attName=="lunit")
714  { lunit = eval.Evaluate(attValue); } else
715  if (attName=="aunit")
716  { aunit = eval.Evaluate(attValue); } else
717  if (attName=="startPhi")
718  { parameter.dimension[0] = eval.Evaluate(attValue); } else
719  if (attName=="openPhi")
720  { parameter.dimension[1] = eval.Evaluate(attValue); } else
721  if (attName=="numRZ")
722  { parameter.dimension[2] = eval.Evaluate(attValue); } else
723  if (attName=="numSide")
724  { parameter.dimension[3] = eval.Evaluate(attValue); }
725  }
726 
727  parameter.dimension[0] *= aunit;
728  parameter.dimension[1] *= aunit;
729 
730  for (xercesc::DOMNode* iter = element->getFirstChild();
731  iter != 0; iter = iter->getNextSibling())
732  {
733  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
734 
735  const xercesc::DOMElement* const child
736  = dynamic_cast<xercesc::DOMElement*>(iter);
737  if (!child)
738  {
739  G4Exception("G4GDMLReadParamvo::PolyhedraRead()",
740  "InvalidRead", FatalException, "No child found!");
741  return;
742  }
743  const G4String tag = Transcode(child->getTagName());
744 
745  if (tag=="zplane") { zplaneList.push_back(ZplaneRead(child)); }
746  }
747 
748  G4int numZPlanes = zplaneList.size();
749  for (G4int i=0; i<numZPlanes; i++)
750  {
751  parameter.dimension[4+i*3] = zplaneList[i].rmin*lunit;
752  parameter.dimension[5+i*3] = zplaneList[i].rmax*lunit;
753  parameter.dimension[6+i*3] = zplaneList[i].z*lunit;
754  }
755 }
756 
758 ParametersRead(const xercesc::DOMElement* const element)
759 {
760  G4ThreeVector rotation(0.0,0.0,0.0);
761  G4ThreeVector position(0.0,0.0,0.0);
762 
764 
765  for (xercesc::DOMNode* iter = element->getFirstChild();
766  iter != 0; iter = iter->getNextSibling())
767  {
768  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
769 
770  const xercesc::DOMElement* const child
771  = dynamic_cast<xercesc::DOMElement*>(iter);
772  if (!child)
773  {
774  G4Exception("G4GDMLReadParamvol::ParametersRead()",
775  "InvalidRead", FatalException, "No child found!");
776  return;
777  }
778  const G4String tag = Transcode(child->getTagName());
779  if (tag=="rotation") { VectorRead(child,rotation); } else
780  if (tag=="position") { VectorRead(child,position); } else
781  if (tag=="positionref")
782  { position = GetPosition(GenerateName(RefRead(child))); } else
783  if (tag=="rotationref")
784  { rotation = GetRotation(GenerateName(RefRead(child))); } else
785  if (tag=="box_dimensions") { Box_dimensionsRead(child,parameter); } else
786  if (tag=="trd_dimensions") { Trd_dimensionsRead(child,parameter); } else
787  if (tag=="trap_dimensions") { Trap_dimensionsRead(child,parameter); } else
788  if (tag=="tube_dimensions") { Tube_dimensionsRead(child,parameter); } else
789  if (tag=="cone_dimensions") { Cone_dimensionsRead(child,parameter); } else
790  if (tag=="sphere_dimensions") { Sphere_dimensionsRead(child,parameter); } else
791  if (tag=="orb_dimensions") { Orb_dimensionsRead(child,parameter); } else
792  if (tag=="torus_dimensions") { Torus_dimensionsRead(child,parameter); } else
793  if (tag=="ellipsoid_dimensions") { Ellipsoid_dimensionsRead(child,parameter); } else
794  if (tag=="para_dimensions") { Para_dimensionsRead(child,parameter); } else
795  if (tag=="polycone_dimensions") { Polycone_dimensionsRead(child,parameter); } else
796  if (tag=="polyhedra_dimensions") { Polyhedra_dimensionsRead(child,parameter); } else
797  if (tag=="hype_dimensions") { Hype_dimensionsRead(child,parameter); }
798  else
799  {
800  G4String error_msg = "Unknown tag in parameters: " + tag;
801  G4Exception("G4GDMLReadParamvol::ParametersRead()", "ReadError",
802  FatalException, error_msg);
803  }
804  }
805 
806  parameter.pRot = new G4RotationMatrix();
807 
808  parameter.pRot->rotateX(rotation.x());
809  parameter.pRot->rotateY(rotation.y());
810  parameter.pRot->rotateZ(rotation.z());
811 
812  parameter.position = position;
813 
814  parameterisation->AddParameter(parameter);
815 }
816 
818 ParameterisedRead(const xercesc::DOMElement* const element)
819 {
820  for (xercesc::DOMNode* iter = element->getFirstChild();
821  iter != 0; iter = iter->getNextSibling())
822  {
823  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
824 
825  const xercesc::DOMElement* const child
826  = dynamic_cast<xercesc::DOMElement*>(iter);
827  if (!child)
828  {
829  G4Exception("G4GDMLReadParamvol::ParameterisedRead()",
830  "InvalidRead", FatalException, "No child found!");
831  return;
832  }
833  const G4String tag = Transcode(child->getTagName());
834 
835  if (tag=="parameters")
836  {
837  const xercesc::DOMNamedNodeMap* const attributes
838  = element->getAttributes();
839  XMLSize_t attributeCount = attributes->getLength();
840  for (XMLSize_t attribute_index=0;
841  attribute_index<attributeCount; attribute_index++)
842  {
843  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
844 
845  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
846  { continue; }
847 
848  const xercesc::DOMAttr* const attribute
849  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
850  if (!attribute)
851  {
852  G4Exception("G4GDMLReadParamvol::ParameterisedRead()",
853  "InvalidRead", FatalException, "No attribute found!");
854  return;
855  }
856  const G4String attName = Transcode(attribute->getName());
857  const G4String attValue = Transcode(attribute->getValue());
858 
859  if (attName=="number") { eval.Evaluate(attValue); }
860  }
861  ParametersRead(child);
862  }
863  else
864  {
865  if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); }
866  }
867  }
868 }
869 
871 Paramvol_contentRead(const xercesc::DOMElement* const element)
872 {
873  for (xercesc::DOMNode* iter = element->getFirstChild();
874  iter != 0; iter = iter->getNextSibling())
875  {
876  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
877 
878  const xercesc::DOMElement* const child
879  = dynamic_cast<xercesc::DOMElement*>(iter);
880  if (!child)
881  {
882  G4Exception("G4GDMLReadParamvol::Paramvol_contentRead()", "InvalidRead",
883  FatalException, "No child found!");
884  return;
885  }
886  const G4String tag = Transcode(child->getTagName());
887  if (tag=="parameterised_position_size") { ParameterisedRead(child); }else
888  if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); }
889  }
890 }
891 
893 ParamvolRead(const xercesc::DOMElement* const element, G4LogicalVolume* mother)
894 {
895  G4String volumeref;
896 
898  for (xercesc::DOMNode* iter = element->getFirstChild();
899  iter != 0; iter = iter->getNextSibling())
900  {
901  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
902 
903  const xercesc::DOMElement* const child
904  = dynamic_cast<xercesc::DOMElement*>(iter);
905  if (!child)
906  {
907  G4Exception("G4GDMLReadParamvol::ParamvolRead()", "InvalidRead",
908  FatalException, "No child found!");
909  return;
910  }
911  const G4String tag = Transcode(child->getTagName());
912 
913  if (tag=="volumeref") { volumeref = RefRead(child); }
914 
915  }
916 
917  Paramvol_contentRead(element);
918 
919  G4LogicalVolume* logvol = GetVolume(GenerateName(volumeref));
920 
921  if (parameterisation->GetSize()==0)
922  {
923  G4Exception("G4GDMLReadParamvol::ParamvolRead()",
924  "ReadError", FatalException,
925  "No parameters are defined in parameterised volume!");
926  }
927  G4String pv_name = logvol->GetName() + "_param";
928  new G4PVParameterised(pv_name, logvol, mother, kUndefined,
930 }
zplaneType ZplaneRead(const xercesc::DOMElement *const)
void Torus_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Trap_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Orb_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
G4String GetName() const
HepRotation & rotateX(double delta)
Definition: Rotation.cc:66
double x() const
CLHEP::HepRotation G4RotationMatrix
void AddParameter(const PARAMETER &)
Definition: xmlparse.cc:179
virtual void Paramvol_contentRead(const xercesc::DOMElement *const)
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
void Ellipsoid_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Trd_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
int G4int
Definition: G4Types.hh:78
void Polycone_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Polyhedra_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
double z() const
G4String RefRead(const xercesc::DOMElement *const)
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:179
void Box_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Sphere_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
G4GDMLParameterisation * parameterisation
void Tube_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Cone_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void Hype_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void Paramvol_contentRead(const xercesc::DOMElement *const)=0
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
void Para_dimensionsRead(const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
void ParameterisedRead(const xercesc::DOMElement *const)
void ParametersRead(const xercesc::DOMElement *const)
int position
Definition: filter.cc:7
double y() const
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
virtual G4LogicalVolume * GetVolume(const G4String &) const =0
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
G4ThreeVector GetRotation(const G4String &)
double G4double
Definition: G4Types.hh:76
G4bool check
Definition: G4GDMLRead.hh:146
G4double Evaluate(const G4String &)
G4ThreeVector GetPosition(const G4String &)