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