Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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$
27 //
28 // class G4GDMLReadParamvol Implementation
29 //
30 // History:
31 // - Created. Zoltan Torzsok, November 2007
32 // -------------------------------------------------------------------------
33 
34 #include "G4GDMLReadParamvol.hh"
35 
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 Para_dimensionsRead( const xercesc::DOMElement* const element,
457 {
458  G4double lunit = 1.0;
459  G4double aunit = 1.0;
460 
461  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
462  XMLSize_t attributeCount = attributes->getLength();
463 
464  for (XMLSize_t attribute_index=0;
465  attribute_index<attributeCount; attribute_index++)
466  {
467  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
468 
469  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
470  { continue; }
471 
472  const xercesc::DOMAttr* const attribute
473  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
474  if (!attribute)
475  {
476  G4Exception("G4GDMLReadParamvol::Para_dimensionsRead()",
477  "InvalidRead", FatalException, "No attribute found!");
478  return;
479  }
480  const G4String attName = Transcode(attribute->getName());
481  const G4String attValue = Transcode(attribute->getValue());
482 
483  if (attName=="lunit")
484  { lunit = eval.Evaluate(attValue); } else
485  if (attName=="aunit")
486  { aunit = eval.Evaluate(attValue); } else
487  if (attName=="x")
488  { parameter.dimension[0] = eval.Evaluate(attValue); } else
489  if (attName=="y")
490  { parameter.dimension[1] = eval.Evaluate(attValue); } else
491  if (attName=="z")
492  { parameter.dimension[2] = eval.Evaluate(attValue); } else
493  if (attName=="alpha")
494  { parameter.dimension[3] = eval.Evaluate(attValue); } else
495  if (attName=="theta")
496  { parameter.dimension[4] = eval.Evaluate(attValue); } else
497  if (attName=="phi")
498  { parameter.dimension[5] = eval.Evaluate(attValue); }
499  }
500 
501  parameter.dimension[0] = 0.5*lunit;
502  parameter.dimension[1] = 0.5*lunit;
503  parameter.dimension[2] = 0.5*lunit;
504  parameter.dimension[3] = aunit;
505  parameter.dimension[4] = aunit;
506  parameter.dimension[5] = aunit;
507 }
508 
510 Hype_dimensionsRead( const xercesc::DOMElement* const element,
512 {
513  G4double lunit = 1.0;
514  G4double aunit = 1.0;
515 
516  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
517  XMLSize_t attributeCount = attributes->getLength();
518 
519  for (XMLSize_t attribute_index=0;
520  attribute_index<attributeCount; attribute_index++)
521  {
522  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
523 
524  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
525  { continue; }
526 
527  const xercesc::DOMAttr* const attribute
528  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
529  if (!attribute)
530  {
531  G4Exception("G4GDMLReadParamvol::Hype_dimensionsRead()",
532  "InvalidRead", FatalException, "No attribute found!");
533  return;
534  }
535  const G4String attName = Transcode(attribute->getName());
536  const G4String attValue = Transcode(attribute->getValue());
537 
538  if (attName=="lunit")
539  { lunit = eval.Evaluate(attValue); } else
540  if (attName=="aunit")
541  { aunit = eval.Evaluate(attValue); } else
542  if (attName=="rmin")
543  { parameter.dimension[0] = eval.Evaluate(attValue); } else
544  if (attName=="rmax")
545  { parameter.dimension[1] = eval.Evaluate(attValue); } else
546  if (attName=="inst")
547  { parameter.dimension[2] = eval.Evaluate(attValue); } else
548  if (attName=="outst")
549  { parameter.dimension[3] = eval.Evaluate(attValue); } else
550  if (attName=="z")
551  { parameter.dimension[4] = eval.Evaluate(attValue); }
552  }
553 
554  parameter.dimension[0] = lunit;
555  parameter.dimension[1] = lunit;
556  parameter.dimension[2] = aunit;
557  parameter.dimension[3] = aunit;
558  parameter.dimension[4] = 0.5*lunit;
559 }
560 
562 ParametersRead(const xercesc::DOMElement* const element) {
563 
564  G4ThreeVector rotation(0.0,0.0,0.0);
565  G4ThreeVector position(0.0,0.0,0.0);
566 
568 
569  for (xercesc::DOMNode* iter = element->getFirstChild();
570  iter != 0; iter = iter->getNextSibling())
571  {
572  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
573 
574  const xercesc::DOMElement* const child
575  = dynamic_cast<xercesc::DOMElement*>(iter);
576  if (!child)
577  {
578  G4Exception("G4GDMLReadParamvol::ParametersRead()",
579  "InvalidRead", FatalException, "No child found!");
580  return;
581  }
582  const G4String tag = Transcode(child->getTagName());
583  if (tag=="rotation") { VectorRead(child,rotation); } else
584  if (tag=="position") { VectorRead(child,position); } else
585  if (tag=="positionref")
586  { position = GetPosition(GenerateName(RefRead(child))); } else
587  if (tag=="rotationref")
588  { rotation = GetRotation(GenerateName(RefRead(child))); } else
589  if (tag=="box_dimensions") { Box_dimensionsRead(child,parameter); } else
590  if (tag=="trd_dimensions") { Trd_dimensionsRead(child,parameter); } else
591  if (tag=="trap_dimensions") { Trap_dimensionsRead(child,parameter); } else
592  if (tag=="tube_dimensions") { Tube_dimensionsRead(child,parameter); } else
593  if (tag=="cone_dimensions") { Cone_dimensionsRead(child,parameter); } else
594  if (tag=="sphere_dimensions") { Cone_dimensionsRead(child,parameter); } else
595  if (tag=="orb_dimensions") { Cone_dimensionsRead(child,parameter); } else
596  if (tag=="torus_dimensions") { Cone_dimensionsRead(child,parameter); } else
597  if (tag=="para_dimensions") { Cone_dimensionsRead(child,parameter); } else
598  if (tag=="hype_dimensions") { Hype_dimensionsRead(child,parameter); }
599  else
600  {
601  G4String error_msg = "Unknown tag in parameters: " + tag;
602  G4Exception("G4GDMLReadParamvol::ParametersRead()", "ReadError",
603  FatalException, error_msg);
604  }
605  }
606 
607  parameter.pRot = new G4RotationMatrix();
608 
609  parameter.pRot->rotateX(rotation.x());
610  parameter.pRot->rotateY(rotation.y());
611  parameter.pRot->rotateZ(rotation.z());
612 
613  parameter.position = position;
614 
615  parameterisation->AddParameter(parameter);
616 }
617 
619 ParameterisedRead(const xercesc::DOMElement* const element)
620 {
621  for (xercesc::DOMNode* iter = element->getFirstChild();
622  iter != 0; iter = iter->getNextSibling())
623  {
624  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
625 
626  const xercesc::DOMElement* const child
627  = dynamic_cast<xercesc::DOMElement*>(iter);
628  if (!child)
629  {
630  G4Exception("G4GDMLReadParamvol::ParameterisedRead()",
631  "InvalidRead", FatalException, "No child found!");
632  return;
633  }
634  const G4String tag = Transcode(child->getTagName());
635 
636  if (tag=="parameters")
637  {
638  const xercesc::DOMNamedNodeMap* const attributes
639  = element->getAttributes();
640  XMLSize_t attributeCount = attributes->getLength();
641  for (XMLSize_t attribute_index=0;
642  attribute_index<attributeCount; attribute_index++)
643  {
644  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
645 
646  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
647  { continue; }
648 
649  const xercesc::DOMAttr* const attribute
650  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
651  if (!attribute)
652  {
653  G4Exception("G4GDMLReadParamvol::ParameterisedRead()",
654  "InvalidRead", FatalException, "No attribute found!");
655  return;
656  }
657  const G4String attName = Transcode(attribute->getName());
658  const G4String attValue = Transcode(attribute->getValue());
659 
660  if (attName=="number") { eval.Evaluate(attValue); }
661  }
662  ParametersRead(child);
663  }
664  else
665  {
666  if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); }
667  }
668  }
669 }
670 
672 Paramvol_contentRead(const xercesc::DOMElement* const element)
673 {
674  for (xercesc::DOMNode* iter = element->getFirstChild();
675  iter != 0; iter = iter->getNextSibling())
676  {
677  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
678 
679  const xercesc::DOMElement* const child
680  = dynamic_cast<xercesc::DOMElement*>(iter);
681  if (!child)
682  {
683  G4Exception("G4GDMLReadParamvol::Paramvol_contentRead()", "InvalidRead",
684  FatalException, "No child found!");
685  return;
686  }
687  const G4String tag = Transcode(child->getTagName());
688  if (tag=="parameterised_position_size") { ParameterisedRead(child); }else
689  if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); }
690  }
691 }
692 
694 ParamvolRead(const xercesc::DOMElement* const element, G4LogicalVolume* mother)
695 {
696  G4String volumeref;
697 
699 
700  for (xercesc::DOMNode* iter = element->getFirstChild();
701  iter != 0; iter = iter->getNextSibling())
702  {
703  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
704 
705  const xercesc::DOMElement* const child
706  = dynamic_cast<xercesc::DOMElement*>(iter);
707  if (!child)
708  {
709  G4Exception("G4GDMLReadParamvol::ParamvolRead()", "InvalidRead",
710  FatalException, "No child found!");
711  return;
712  }
713  const G4String tag = Transcode(child->getTagName());
714 
715  if (tag=="volumeref") { volumeref = RefRead(child); }
716 
717  }
718 
719  Paramvol_contentRead(element);
720 
721  G4LogicalVolume* logvol = GetVolume(GenerateName(volumeref));
722 
723  if (parameterisation->GetSize()==0)
724  {
725  G4Exception("G4GDMLReadParamvol::ParamvolRead()",
726  "ReadError", FatalException,
727  "No parameters are defined in parameterised volume!");
728  }
729  G4String pv_name = logvol->GetName() + "_param";
730  new G4PVParameterised(pv_name, logvol, mother, kUndefined,
732 }