Geant4  10.02
G4GDMLReadStructure.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: G4GDMLReadStructure.cc 91689 2015-07-31 09:58:27Z gcosmo $
27 //
28 // class G4GDMLReadStructure Implementation
29 //
30 // Original author: Zoltan Torzsok, November 2007
31 //
32 // --------------------------------------------------------------------
33 
34 #include "G4GDMLReadStructure.hh"
35 
36 #include "G4UnitsTable.hh"
37 #include "G4LogicalVolume.hh"
38 #include "G4VPhysicalVolume.hh"
39 #include "G4PVPlacement.hh"
40 #include "G4LogicalVolumeStore.hh"
41 #include "G4PhysicalVolumeStore.hh"
42 #include "G4AssemblyVolume.hh"
43 #include "G4ReflectionFactory.hh"
44 #include "G4PVDivisionFactory.hh"
46 #include "G4LogicalSkinSurface.hh"
47 #include "G4VisAttributes.hh"
48 
50  : G4GDMLReadParamvol(), pMotherLogical(0), strip(false)
51 {
52 }
53 
55 {
56 }
57 
58 
60 BorderSurfaceRead(const xercesc::DOMElement* const bordersurfaceElement)
61 {
62  G4String name;
63  G4VPhysicalVolume* pv1 = 0;
64  G4VPhysicalVolume* pv2 = 0;
65  G4SurfaceProperty* prop = 0;
66  G4int index = 0;
67 
68  const xercesc::DOMNamedNodeMap* const attributes
69  = bordersurfaceElement->getAttributes();
70  XMLSize_t attributeCount = attributes->getLength();
71 
72  for (XMLSize_t attribute_index=0;
73  attribute_index<attributeCount; attribute_index++)
74  {
75  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
76 
77  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
78  { continue; }
79 
80  const xercesc::DOMAttr* const attribute
81  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
82  if (!attribute)
83  {
84  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
85  "InvalidRead", FatalException, "No attribute found!");
86  return;
87  }
88  const G4String attName = Transcode(attribute->getName());
89  const G4String attValue = Transcode(attribute->getValue());
90 
91  if (attName=="name")
92  { name = GenerateName(attValue); } else
93  if (attName=="surfaceproperty")
94  { prop = GetSurfaceProperty(GenerateName(attValue)); }
95  }
96 
97  for (xercesc::DOMNode* iter = bordersurfaceElement->getFirstChild();
98  iter != 0; iter = iter->getNextSibling())
99  {
100  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
101 
102  const xercesc::DOMElement* const child
103  = dynamic_cast<xercesc::DOMElement*>(iter);
104  if (!child)
105  {
106  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
107  "InvalidRead", FatalException, "No child found!");
108  return;
109  }
110  const G4String tag = Transcode(child->getTagName());
111 
112  if (tag != "physvolref") { continue; }
113 
114  if (index==0)
115  { pv1 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
116  if (index==1)
117  { pv2 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
118  break;
119  }
120 
121  new G4LogicalBorderSurface(Strip(name),pv1,pv2,prop);
122 }
123 
125 DivisionvolRead(const xercesc::DOMElement* const divisionvolElement)
126 {
127  G4String name;
128  G4double unit = 1.0;
129  G4double width = 0.0;
130  G4double offset = 0.0;
131  G4int number = 0;
132  EAxis axis = kUndefined;
133  G4LogicalVolume* logvol = 0;
134 
135  const xercesc::DOMNamedNodeMap* const attributes
136  = divisionvolElement->getAttributes();
137  XMLSize_t attributeCount = attributes->getLength();
138 
139  for (XMLSize_t attribute_index=0;
140  attribute_index<attributeCount; attribute_index++)
141  {
142  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
143 
144  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
145  { continue; }
146 
147  const xercesc::DOMAttr* const attribute
148  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
149  if (!attribute)
150  {
151  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
152  "InvalidRead", FatalException, "No attribute found!");
153  return;
154  }
155  const G4String attName = Transcode(attribute->getName());
156  const G4String attValue = Transcode(attribute->getValue());
157 
158  if (attName=="name") { name = attValue; } else
159  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); } else
160  if (attName=="width") { width = eval.Evaluate(attValue); } else
161  if (attName=="offset") { offset = eval.Evaluate(attValue); } else
162  if (attName=="number") { number = eval.EvaluateInteger(attValue); } else
163  if (attName=="axis")
164  {
165  if (attValue=="kXAxis") { axis = kXAxis; } else
166  if (attValue=="kYAxis") { axis = kYAxis; } else
167  if (attValue=="kZAxis") { axis = kZAxis; } else
168  if (attValue=="kRho") { axis = kRho; } else
169  if (attValue=="kPhi") { axis = kPhi; }
170  }
171  }
172 
173  width *= unit;
174  offset *= unit;
175 
176  for (xercesc::DOMNode* iter = divisionvolElement->getFirstChild();
177  iter != 0;iter = iter->getNextSibling())
178  {
179  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
180 
181  const xercesc::DOMElement* const child
182  = dynamic_cast<xercesc::DOMElement*>(iter);
183  if (!child)
184  {
185  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
186  "InvalidRead", FatalException, "No child found!");
187  return;
188  }
189  const G4String tag = Transcode(child->getTagName());
190 
191  if (tag=="volumeref") { logvol = GetVolume(GenerateName(RefRead(child))); }
192  }
193 
194  if (!logvol) { return; }
195 
198 
199  G4String pv_name = logvol->GetName() + "_div";
200  if ((number != 0) && (width == 0.0))
201  {
203  ->Divide(pv_name,logvol,pMotherLogical,axis,number,offset);
204  }
205  else if ((number == 0) && (width != 0.0))
206  {
208  ->Divide(pv_name,logvol,pMotherLogical,axis,width,offset);
209  }
210  else
211  {
213  ->Divide(pv_name,logvol,pMotherLogical,axis,number,width,offset);
214  }
215 
216  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
217  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
218 }
219 
221 FileRead(const xercesc::DOMElement* const fileElement)
222 {
223  G4String name;
224  G4String volname;
225 
226  const xercesc::DOMNamedNodeMap* const attributes
227  = fileElement->getAttributes();
228  XMLSize_t attributeCount = attributes->getLength();
229 
230  for (XMLSize_t attribute_index=0;
231  attribute_index<attributeCount; attribute_index++)
232  {
233  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
234 
235  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
236  { continue; }
237 
238  const xercesc::DOMAttr* const attribute
239  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
240  if (!attribute)
241  {
242  G4Exception("G4GDMLReadStructure::FileRead()",
243  "InvalidRead", FatalException, "No attribute found!");
244  return 0;
245  }
246  const G4String attName = Transcode(attribute->getName());
247  const G4String attValue = Transcode(attribute->getValue());
248 
249  if (attName=="name") { name = attValue; } else
250  if (attName=="volname") { volname = attValue; }
251  }
252 
253  const G4bool isModule = true;
254  G4GDMLReadStructure structure;
255  structure.Read(name,validate,isModule);
256 
257  // Register existing auxiliar information defined in child module
258  //
259  const G4GDMLAuxMapType* aux = structure.GetAuxMap();
260  if (!aux->empty())
261  {
262  G4GDMLAuxMapType::const_iterator pos;
263  for (pos = aux->begin(); pos != aux->end(); ++pos)
264  {
265  auxMap.insert(std::make_pair(pos->first,pos->second));
266  }
267  }
268 
269  // Return volume structure from child module
270  //
271  if (volname.empty())
272  {
273  return structure.GetVolume(structure.GetSetup("Default"));
274  }
275  else
276  {
277  return structure.GetVolume(structure.GenerateName(volname));
278  }
279 }
280 
282 PhysvolRead(const xercesc::DOMElement* const physvolElement,
283  G4AssemblyVolume* pAssembly)
284 {
285  G4String name;
286  G4LogicalVolume* logvol = 0;
287  G4AssemblyVolume* assembly = 0;
288  G4ThreeVector position(0.0,0.0,0.0);
289  G4ThreeVector rotation(0.0,0.0,0.0);
290  G4ThreeVector scale(1.0,1.0,1.0);
291  G4int copynumber = 0;
292 
293  const xercesc::DOMNamedNodeMap* const attributes
294  = physvolElement->getAttributes();
295  XMLSize_t attributeCount = attributes->getLength();
296 
297  for (XMLSize_t attribute_index=0;
298  attribute_index<attributeCount; attribute_index++)
299  {
300  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
301 
302  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
303  { continue; }
304 
305  const xercesc::DOMAttr* const attribute
306  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
307  if (!attribute)
308  {
309  G4Exception("G4GDMLReadStructure::PhysvolRead()",
310  "InvalidRead", FatalException, "No attribute found!");
311  return;
312  }
313  const G4String attName = Transcode(attribute->getName());
314  const G4String attValue = Transcode(attribute->getValue());
315 
316  if (attName=="name") { name = attValue; }
317  if (attName=="copynumber") { copynumber = eval.EvaluateInteger(attValue); }
318  }
319 
320  for (xercesc::DOMNode* iter = physvolElement->getFirstChild();
321  iter != 0; iter = iter->getNextSibling())
322  {
323  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
324 
325  const xercesc::DOMElement* const child
326  = dynamic_cast<xercesc::DOMElement*>(iter);
327  if (!child)
328  {
329  G4Exception("G4GDMLReadStructure::PhysvolRead()",
330  "InvalidRead", FatalException, "No child found!");
331  return;
332  }
333  const G4String tag = Transcode(child->getTagName());
334 
335  if (tag=="volumeref")
336  {
337  const G4String& child_name = GenerateName(RefRead(child));
338  assembly = GetAssembly(child_name);
339  if (!assembly) { logvol = GetVolume(child_name); }
340  }
341  else if (tag=="file")
342  { logvol = FileRead(child); }
343  else if (tag=="position")
344  { VectorRead(child,position); }
345  else if (tag=="rotation")
346  { VectorRead(child,rotation); }
347  else if (tag=="scale")
348  { VectorRead(child,scale); }
349  else if (tag=="positionref")
350  { position = GetPosition(GenerateName(RefRead(child))); }
351  else if (tag=="rotationref")
352  { rotation = GetRotation(GenerateName(RefRead(child))); }
353  else if (tag=="scaleref")
354  { scale = GetScale(GenerateName(RefRead(child))); }
355  else
356  {
357  G4String error_msg = "Unknown tag in physvol: " + tag;
358  G4Exception("G4GDMLReadStructure::PhysvolRead()", "ReadError",
359  FatalException, error_msg);
360  return;
361  }
362  }
363 
364  G4Transform3D transform(GetRotationMatrix(rotation).inverse(),position);
365  transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
366 
367  if (pAssembly) // Fill assembly structure
368  {
369  if (!logvol) { return; }
370  pAssembly->AddPlacedVolume(logvol, transform);
371  }
372  else // Generate physical volume tree or do assembly imprint
373  {
374  if (assembly)
375  {
376  assembly->MakeImprint(pMotherLogical, transform, 0, check);
377  }
378  else
379  {
380  if (!logvol) { return; }
381  G4String pv_name = logvol->GetName() + "_PV";
383  ->Place(transform,pv_name,logvol,pMotherLogical,false,copynumber,check);
384 
385  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
386  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
387  }
388  }
389 }
390 
392 ReplicavolRead(const xercesc::DOMElement* const replicavolElement, G4int number)
393 {
394  G4LogicalVolume* logvol = 0;
395  for (xercesc::DOMNode* iter = replicavolElement->getFirstChild();
396  iter != 0; iter = iter->getNextSibling())
397  {
398  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
399 
400  const xercesc::DOMElement* const child
401  = dynamic_cast<xercesc::DOMElement*>(iter);
402  if (!child)
403  {
404  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
405  "InvalidRead", FatalException, "No child found!");
406  return;
407  }
408  const G4String tag = Transcode(child->getTagName());
409 
410  if (tag=="volumeref")
411  {
412  logvol = GetVolume(GenerateName(RefRead(child)));
413  }
414  else if (tag=="replicate_along_axis")
415  {
416  if (logvol) { ReplicaRead(child,logvol,number); }
417  }
418  else
419  {
420  G4String error_msg = "Unknown tag in ReplicavolRead: " + tag;
421  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
422  "ReadError", FatalException, error_msg);
423  }
424  }
425 }
426 
428 ReplicaRead(const xercesc::DOMElement* const replicaElement,
429  G4LogicalVolume* logvol, G4int number)
430 {
431  G4double width = 0.0;
432  G4double offset = 0.0;
433  G4ThreeVector position(0.0,0.0,0.0);
434  G4ThreeVector rotation(0.0,0.0,0.0);
435  EAxis axis = kUndefined;
436  G4String name;
437 
438  for (xercesc::DOMNode* iter = replicaElement->getFirstChild();
439  iter != 0; iter = iter->getNextSibling())
440  {
441  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
442 
443  const xercesc::DOMElement* const child
444  = dynamic_cast<xercesc::DOMElement*>(iter);
445  if (!child)
446  {
447  G4Exception("G4GDMLReadStructure::ReplicaRead()",
448  "InvalidRead", FatalException, "No child found!");
449  return;
450  }
451  const G4String tag = Transcode(child->getTagName());
452 
453  if (tag=="position")
454  { VectorRead(child,position); } else
455  if (tag=="rotation")
456  { VectorRead(child,rotation); } else
457  if (tag=="positionref")
458  { position = GetPosition(GenerateName(RefRead(child))); } else
459  if (tag=="rotationref")
460  { rotation = GetRotation(GenerateName(RefRead(child))); } else
461  if (tag=="direction")
462  { axis=AxisRead(child); } else
463  if (tag=="width")
464  { width=QuantityRead(child); } else
465  if (tag=="offset")
466  { offset=QuantityRead(child); }
467  else
468  {
469  G4String error_msg = "Unknown tag in ReplicaRead: " + tag;
470  G4Exception("G4GDMLReadStructure::ReplicaRead()", "ReadError",
471  FatalException, error_msg);
472  }
473  }
474 
475  G4String pv_name = logvol->GetName() + "_PV";
477  ->Replicate(pv_name,logvol,pMotherLogical,axis,number,width,offset);
478 
479  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
480  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
481 
482 }
483 
485 AxisRead(const xercesc::DOMElement* const axisElement)
486 {
487  EAxis axis = kUndefined;
488 
489  const xercesc::DOMNamedNodeMap* const attributes
490  = axisElement->getAttributes();
491  XMLSize_t attributeCount = attributes->getLength();
492 
493  for (XMLSize_t attribute_index=0;
494  attribute_index<attributeCount; attribute_index++)
495  {
496  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
497 
498  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
499  { continue; }
500 
501  const xercesc::DOMAttr* const attribute
502  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
503  if (!attribute)
504  {
505  G4Exception("G4GDMLReadStructure::AxisRead()",
506  "InvalidRead", FatalException, "No attribute found!");
507  return axis;
508  }
509  const G4String attName = Transcode(attribute->getName());
510  const G4String attValue = Transcode(attribute->getValue());
511  if (attName=="x")
512  { if( eval.Evaluate(attValue)==1.) {axis=kXAxis;} }
513  else if (attName=="y")
514  { if( eval.Evaluate(attValue)==1.) {axis=kYAxis;} }
515  else if (attName=="z")
516  { if( eval.Evaluate(attValue)==1.) {axis=kZAxis;} }
517  else if (attName=="rho")
518  { if( eval.Evaluate(attValue)==1.) {axis=kRho;} }
519  else if (attName=="phi")
520  { if( eval.Evaluate(attValue)==1.) {axis=kPhi;} }
521  }
522 
523  return axis;
524 }
525 
527 QuantityRead(const xercesc::DOMElement* const readElement)
528 {
529  G4double value = 0.0;
530  G4double unit = 0.0;
531  const xercesc::DOMNamedNodeMap* const attributes
532  = readElement->getAttributes();
533  XMLSize_t attributeCount = attributes->getLength();
534 
535  for (XMLSize_t attribute_index=0;
536  attribute_index<attributeCount; attribute_index++)
537  {
538  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
539 
540  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
541  { continue; }
542  const xercesc::DOMAttr* const attribute
543  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
544  if (!attribute)
545  {
546  G4Exception("G4GDMLReadStructure::QuantityRead()",
547  "InvalidRead", FatalException, "No attribute found!");
548  return value;
549  }
550  const G4String attName = Transcode(attribute->getName());
551  const G4String attValue = Transcode(attribute->getValue());
552 
553  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); } else
554  if (attName=="value"){ value= eval.Evaluate(attValue); }
555  }
556 
557  return value*unit;
558 }
559 
561 VolumeRead(const xercesc::DOMElement* const volumeElement)
562 {
563  G4VSolid* solidPtr = 0;
564  G4Material* materialPtr = 0;
565  G4GDMLAuxListType auxList;
566 
567  XMLCh *name_attr = xercesc::XMLString::transcode("name");
568  const G4String name = Transcode(volumeElement->getAttribute(name_attr));
569  xercesc::XMLString::release(&name_attr);
570 
571  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
572  iter != 0; iter = iter->getNextSibling())
573  {
574  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
575 
576  const xercesc::DOMElement* const child
577  = dynamic_cast<xercesc::DOMElement*>(iter);
578  if (!child)
579  {
580  G4Exception("G4GDMLReadStructure::VolumeRead()",
581  "InvalidRead", FatalException, "No child found!");
582  return;
583  }
584  const G4String tag = Transcode(child->getTagName());
585 
586  if (tag=="auxiliary")
587  { auxList.push_back(AuxiliaryRead(child)); } else
588  if (tag=="materialref")
589  { materialPtr = GetMaterial(GenerateName(RefRead(child),true)); } else
590  if (tag=="solidref")
591  { solidPtr = GetSolid(GenerateName(RefRead(child))); }
592  }
593 
594  pMotherLogical = new G4LogicalVolume(solidPtr,materialPtr,
595  GenerateName(name),0,0,0);
596 
597  if (!auxList.empty()) { auxMap[pMotherLogical] = auxList; }
598 
599  Volume_contentRead(volumeElement);
600 }
601 
603 AssemblyRead(const xercesc::DOMElement* const assemblyElement)
604 {
605  XMLCh *name_attr = xercesc::XMLString::transcode("name");
606  const G4String name = Transcode(assemblyElement->getAttribute(name_attr));
607  xercesc::XMLString::release(&name_attr);
608 
609  G4AssemblyVolume* pAssembly = new G4AssemblyVolume();
610  assemblyMap.insert(std::make_pair(GenerateName(name), pAssembly));
611 
612  for (xercesc::DOMNode* iter = assemblyElement->getFirstChild();
613  iter != 0; iter = iter->getNextSibling())
614  {
615  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
616  const xercesc::DOMElement* const child
617  = dynamic_cast<xercesc::DOMElement*>(iter);
618  if (!child)
619  {
620  G4Exception("G4GDMLReadStructure::AssemblyRead()",
621  "InvalidRead", FatalException, "No child found!");
622  return;
623  }
624  const G4String tag = Transcode(child->getTagName());
625 
626  if (tag=="physvol")
627  {
628  PhysvolRead(child, pAssembly);
629  }
630  else
631  {
632  G4cout << "Unsupported GDML tag '" << tag
633  << "' for Geant4 assembly structure !" << G4endl;
634  }
635  }
636 }
637 
639 SkinSurfaceRead(const xercesc::DOMElement* const skinsurfaceElement)
640 {
641  G4String name;
642  G4LogicalVolume* logvol = 0;
643  G4SurfaceProperty* prop = 0;
644 
645  const xercesc::DOMNamedNodeMap* const attributes
646  = skinsurfaceElement->getAttributes();
647  XMLSize_t attributeCount = attributes->getLength();
648 
649  for (XMLSize_t attribute_index=0;
650  attribute_index<attributeCount; attribute_index++)
651  {
652  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
653 
654  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
655  { continue; }
656 
657  const xercesc::DOMAttr* const attribute
658  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
659  if (!attribute)
660  {
661  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
662  "InvalidRead", FatalException, "No attribute found!");
663  return;
664  }
665  const G4String attName = Transcode(attribute->getName());
666  const G4String attValue = Transcode(attribute->getValue());
667 
668  if (attName=="name")
669  { name = GenerateName(attValue); } else
670  if (attName=="surfaceproperty")
671  { prop = GetSurfaceProperty(GenerateName(attValue)); }
672  }
673 
674  for (xercesc::DOMNode* iter = skinsurfaceElement->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("G4GDMLReadStructure::SkinsurfaceRead()",
684  "InvalidRead", FatalException, "No child found!");
685  return;
686  }
687  const G4String tag = Transcode(child->getTagName());
688 
689  if (tag=="volumeref")
690  {
691  logvol = GetVolume(GenerateName(RefRead(child)));
692  }
693  else
694  {
695  G4String error_msg = "Unknown tag in skinsurface: " + tag;
696  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
697  FatalException, error_msg);
698  }
699  }
700 
701  new G4LogicalSkinSurface(Strip(name),logvol,prop);
702 }
703 
705 Volume_contentRead(const xercesc::DOMElement* const volumeElement)
706 {
707  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
708  iter != 0; iter = iter->getNextSibling())
709  {
710  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
711 
712  const xercesc::DOMElement* const child
713  = dynamic_cast<xercesc::DOMElement*>(iter);
714  if (!child)
715  {
716  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
717  "InvalidRead", FatalException, "No child found!");
718  return;
719  }
720  const G4String tag = Transcode(child->getTagName());
721 
722  if ((tag=="auxiliary") || (tag=="materialref") || (tag=="solidref"))
723  {
724  // These are already processed in VolumeRead()
725  }
726  else if (tag=="paramvol")
727  {
729  }
730  else if (tag=="physvol")
731  {
732  PhysvolRead(child);
733  }
734  else if (tag=="replicavol")
735  {
736  G4int number = 1;
737  const xercesc::DOMNamedNodeMap* const attributes
738  = child->getAttributes();
739  XMLSize_t attributeCount = attributes->getLength();
740  for (XMLSize_t attribute_index=0;
741  attribute_index<attributeCount; attribute_index++)
742  {
743  xercesc::DOMNode* attribute_node
744  = attributes->item(attribute_index);
745  if (attribute_node->getNodeType()!=xercesc::DOMNode::ATTRIBUTE_NODE)
746  {
747  continue;
748  }
749  const xercesc::DOMAttr* const attribute
750  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
751  if (!attribute)
752  {
753  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
754  "InvalidRead", FatalException, "No attribute found!");
755  return;
756  }
757  const G4String attName = Transcode(attribute->getName());
758  const G4String attValue = Transcode(attribute->getValue());
759  if (attName=="number")
760  {
761  number = eval.EvaluateInteger(attValue);
762  }
763  }
764  ReplicavolRead(child,number);
765  }
766  else if (tag=="divisionvol")
767  {
768  DivisionvolRead(child);
769  }
770  else if (tag=="loop")
771  {
773  }
774  else
775  {
776  G4cout << "Treating unknown GDML tag in volume '" << tag
777  << "' as GDML extension..." << G4endl;
778  }
779  }
780 }
781 
783 StructureRead(const xercesc::DOMElement* const structureElement)
784 {
785  G4cout << "G4GDML: Reading structure..." << G4endl;
786 
787  for (xercesc::DOMNode* iter = structureElement->getFirstChild();
788  iter != 0; iter = iter->getNextSibling())
789  {
790  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
791 
792  const xercesc::DOMElement* const child
793  = dynamic_cast<xercesc::DOMElement*>(iter);
794  if (!child)
795  {
796  G4Exception("G4GDMLReadStructure::StructureRead()",
797  "InvalidRead", FatalException, "No child found!");
798  return;
799  }
800  const G4String tag = Transcode(child->getTagName());
801 
802  if (tag=="bordersurface") { BorderSurfaceRead(child); } else
803  if (tag=="skinsurface") { SkinSurfaceRead(child); } else
804  if (tag=="volume") { VolumeRead(child); } else
805  if (tag=="assembly") { AssemblyRead(child); } else
806  if (tag=="loop") { LoopRead(child,&G4GDMLRead::StructureRead); }
807  else
808  {
809  G4String error_msg = "Unknown tag in structure: " + tag;
810  G4Exception("G4GDMLReadStructure::StructureRead()",
811  "ReadError", FatalException, error_msg);
812  }
813  }
814 }
815 
817 GetPhysvol(const G4String& ref) const
818 {
819  G4VPhysicalVolume* physvolPtr =
821 
822  if (!physvolPtr)
823  {
824  G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
825  G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
826  FatalException, error_msg);
827  }
828 
829  return physvolPtr;
830 }
831 
833 GetVolume(const G4String& ref) const
834 {
835  G4LogicalVolume *volumePtr
837 
838  if (!volumePtr)
839  {
840  G4String error_msg = "Referenced volume '" + ref + "' was not found!";
841  G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError",
842  FatalException, error_msg);
843  }
844 
845  return volumePtr;
846 }
847 
849 GetAssembly(const G4String& ref) const
850 {
851  G4GDMLAssemblyMapType::const_iterator pos = assemblyMap.find(ref);
852  if (pos != assemblyMap.end()) { return pos->second; }
853  return 0;
854 }
855 
858 {
859  G4GDMLAuxMapType::const_iterator pos = auxMap.find(logvol);
860  if (pos != auxMap.end()) { return pos->second; }
861  else { return G4GDMLAuxListType(); }
862 }
863 
865 GetWorldVolume(const G4String& setupName)
866 {
867  G4String sname = GetSetup(setupName);
868  if (sname == "") { return 0; }
869 
870  G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
872 
873  G4VPhysicalVolume* pvWorld = 0;
874 
875  if(setuptoPV[setupName])
876  {
877  pvWorld = setuptoPV[setupName];
878  }
879  else
880  {
881  pvWorld = new G4PVPlacement(0,G4ThreeVector(0,0,0),volume,
882  volume->GetName()+"_PV",0,0,0);
883  setuptoPV[setupName] = pvWorld;
884  }
885  return pvWorld;
886 }
887 
889 {
890  eval.Clear();
891  setuptoPV.clear();
892 }
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
Definition: geomdefs.hh:54
G4AssemblyVolume * GetAssembly(const G4String &) const
G4LogicalVolume * FileRead(const xercesc::DOMElement *const)
G4int EvaluateInteger(const G4String &)
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
CLHEP::Hep3Vector G4ThreeVector
const G4GDMLAuxMapType * GetAuxMap() const
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
G4int first(char) const
G4String name
Definition: TRTMaterials.hh:40
virtual void Volume_contentRead(const xercesc::DOMElement *const)=0
Definition: xmlparse.cc:187
G4bool dostrip
Definition: G4GDMLRead.hh:158
G4double QuantityRead(const xercesc::DOMElement *const readElement)
void DivisionvolRead(const xercesc::DOMElement *const)
#define width
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
G4ThreeVector GetScale(const G4String &)
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
void ReplicavolRead(const xercesc::DOMElement *const, G4int number)
G4PhysicalVolumesPair Replicate(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofReplicas, G4double width, G4double offset=0)
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:81
int G4int
Definition: G4Types.hh:78
G4VPhysicalVolume * GetPhysvol(const G4String &) const
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4ReflectionFactory * Instance()
G4VSolid * GetSolid(const G4String &) const
virtual void StructureRead(const xercesc::DOMElement *const)=0
static G4PhysicalVolumeStore * GetInstance()
G4String RefRead(const xercesc::DOMElement *const)
static G4PVDivisionFactory * GetInstance()
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
G4GLOB_DLL std::ostream G4cout
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
G4GDMLAssemblyMapType assemblyMap
G4LogicalVolume * GetVolume(const G4String &) const
bool G4bool
Definition: G4Types.hh:79
void ReplicaRead(const xercesc::DOMElement *const replicaElement, G4LogicalVolume *logvol, G4int number)
static G4LogicalVolumeStore * GetInstance()
HepGeom::Transform3D G4Transform3D
HepGeom::Scale3D G4Scale3D
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
std::map< std::string, G4VPhysicalVolume * > setuptoPV
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
void BorderSurfaceRead(const xercesc::DOMElement *const)
void AddPlacedVolume(G4LogicalVolume *pPlacedVolume, G4ThreeVector &translation, G4RotationMatrix *rotation)
G4GDMLAuxMapType auxMap
EAxis
Definition: geomdefs.hh:54
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:97
G4VPhysicalVolume * GetWorldVolume(const G4String &)
EAxis AxisRead(const xercesc::DOMElement *const axisElement)
void SkinSurfaceRead(const xercesc::DOMElement *const)
G4bool validate
Definition: G4GDMLRead.hh:156
G4GDMLAuxStructType AuxiliaryRead(const xercesc::DOMElement *const auxElem)
Definition: G4GDMLRead.cc:262
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
static const G4VisAttributes Invisible
virtual void Volume_contentRead(const xercesc::DOMElement *const)
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector GetRotation(const G4String &)
G4GDMLAuxListType GetVolumeAuxiliaryInformation(G4LogicalVolume *) const
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
std::map< G4LogicalVolume *, G4GDMLAuxListType > G4GDMLAuxMapType
Definition: geomdefs.hh:54
virtual void StructureRead(const xercesc::DOMElement *const)
void AssemblyRead(const xercesc::DOMElement *const)
G4VPhysicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
G4bool check
Definition: G4GDMLRead.hh:157
void SetVisAttributes(const G4VisAttributes *pVA)
static const G4double pos
G4String GetSetup(const G4String &)
void Read(const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
Definition: G4GDMLRead.cc:361
G4double Evaluate(const G4String &)
void PhysvolRead(const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
G4ThreeVector GetPosition(const G4String &)
virtual void VolumeRead(const xercesc::DOMElement *const)