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