Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 101687 2016-11-21 09:43:18Z 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  G4String unitname;
139 
140  for (XMLSize_t attribute_index=0;
141  attribute_index<attributeCount; attribute_index++)
142  {
143  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
144 
145  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
146  { continue; }
147 
148  const xercesc::DOMAttr* const attribute
149  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
150  if (!attribute)
151  {
152  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
153  "InvalidRead", FatalException, "No attribute found!");
154  return;
155  }
156  const G4String attName = Transcode(attribute->getName());
157  const G4String attValue = Transcode(attribute->getValue());
158 
159  if (attName=="name") { name = attValue; } else
160  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
161  unitname = G4UnitDefinition::GetCategory(attValue);
162  } else
163  if (attName=="width") { width = eval.Evaluate(attValue); } else
164  if (attName=="offset") { offset = eval.Evaluate(attValue); } else
165  if (attName=="number") { number = eval.EvaluateInteger(attValue); } else
166  if (attName=="axis")
167  {
168  if (attValue=="kXAxis") { axis = kXAxis; } else
169  if (attValue=="kYAxis") { axis = kYAxis; } else
170  if (attValue=="kZAxis") { axis = kZAxis; } else
171  if (attValue=="kRho") { axis = kRho; } else
172  if (attValue=="kPhi") { axis = kPhi; }
173  }
174  }
175 
176  if ( ((axis == kXAxis || axis == kYAxis || axis == kZAxis) && unitname!="Length") ||
177  ((axis == kRho || axis == kPhi) && unitname!="Angle")) {
178  G4Exception("G4GDMLReadStructure::DivisionvolRead()", "InvalidRead",
179  FatalException, "Invalid unit!"); }
180 
181  width *= unit;
182  offset *= unit;
183 
184  for (xercesc::DOMNode* iter = divisionvolElement->getFirstChild();
185  iter != 0;iter = iter->getNextSibling())
186  {
187  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
188 
189  const xercesc::DOMElement* const child
190  = dynamic_cast<xercesc::DOMElement*>(iter);
191  if (!child)
192  {
193  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
194  "InvalidRead", FatalException, "No child found!");
195  return;
196  }
197  const G4String tag = Transcode(child->getTagName());
198 
199  if (tag=="volumeref") { logvol = GetVolume(GenerateName(RefRead(child))); }
200  }
201 
202  if (!logvol) { return; }
203 
206 
207  G4String pv_name = logvol->GetName() + "_div";
208  if ((number != 0) && (width == 0.0))
209  {
211  ->Divide(pv_name,logvol,pMotherLogical,axis,number,offset);
212  }
213  else if ((number == 0) && (width != 0.0))
214  {
216  ->Divide(pv_name,logvol,pMotherLogical,axis,width,offset);
217  }
218  else
219  {
221  ->Divide(pv_name,logvol,pMotherLogical,axis,number,width,offset);
222  }
223 
224  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
225  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
226 }
227 
229 FileRead(const xercesc::DOMElement* const fileElement)
230 {
231  G4String name;
232  G4String volname;
233 
234  const xercesc::DOMNamedNodeMap* const attributes
235  = fileElement->getAttributes();
236  XMLSize_t attributeCount = attributes->getLength();
237 
238  for (XMLSize_t attribute_index=0;
239  attribute_index<attributeCount; attribute_index++)
240  {
241  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
242 
243  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
244  { continue; }
245 
246  const xercesc::DOMAttr* const attribute
247  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
248  if (!attribute)
249  {
250  G4Exception("G4GDMLReadStructure::FileRead()",
251  "InvalidRead", FatalException, "No attribute found!");
252  return 0;
253  }
254  const G4String attName = Transcode(attribute->getName());
255  const G4String attValue = Transcode(attribute->getValue());
256 
257  if (attName=="name") { name = attValue; } else
258  if (attName=="volname") { volname = attValue; }
259  }
260 
261  const G4bool isModule = true;
262  G4GDMLReadStructure structure;
263  structure.Read(name,validate,isModule);
264 
265  // Register existing auxiliar information defined in child module
266  //
267  const G4GDMLAuxMapType* aux = structure.GetAuxMap();
268  if (!aux->empty())
269  {
270  G4GDMLAuxMapType::const_iterator pos;
271  for (pos = aux->begin(); pos != aux->end(); ++pos)
272  {
273  auxMap.insert(std::make_pair(pos->first,pos->second));
274  }
275  }
276 
277  // Return volume structure from child module
278  //
279  if (volname.empty())
280  {
281  return structure.GetVolume(structure.GetSetup("Default"));
282  }
283  else
284  {
285  return structure.GetVolume(structure.GenerateName(volname));
286  }
287 }
288 
290 PhysvolRead(const xercesc::DOMElement* const physvolElement,
291  G4AssemblyVolume* pAssembly)
292 {
293  G4String name;
294  G4LogicalVolume* logvol = 0;
295  G4AssemblyVolume* assembly = 0;
296  G4ThreeVector position(0.0,0.0,0.0);
297  G4ThreeVector rotation(0.0,0.0,0.0);
298  G4ThreeVector scale(1.0,1.0,1.0);
299  G4int copynumber = 0;
300 
301  const xercesc::DOMNamedNodeMap* const attributes
302  = physvolElement->getAttributes();
303  XMLSize_t attributeCount = attributes->getLength();
304 
305  for (XMLSize_t attribute_index=0;
306  attribute_index<attributeCount; attribute_index++)
307  {
308  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
309 
310  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
311  { continue; }
312 
313  const xercesc::DOMAttr* const attribute
314  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
315  if (!attribute)
316  {
317  G4Exception("G4GDMLReadStructure::PhysvolRead()",
318  "InvalidRead", FatalException, "No attribute found!");
319  return;
320  }
321  const G4String attName = Transcode(attribute->getName());
322  const G4String attValue = Transcode(attribute->getValue());
323 
324  if (attName=="name") { name = attValue; }
325  if (attName=="copynumber") { copynumber = eval.EvaluateInteger(attValue); }
326  }
327 
328  for (xercesc::DOMNode* iter = physvolElement->getFirstChild();
329  iter != 0; iter = iter->getNextSibling())
330  {
331  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
332 
333  const xercesc::DOMElement* const child
334  = dynamic_cast<xercesc::DOMElement*>(iter);
335  if (!child)
336  {
337  G4Exception("G4GDMLReadStructure::PhysvolRead()",
338  "InvalidRead", FatalException, "No child found!");
339  return;
340  }
341  const G4String tag = Transcode(child->getTagName());
342 
343  if (tag=="volumeref")
344  {
345  const G4String& child_name = GenerateName(RefRead(child));
346  assembly = GetAssembly(child_name);
347  if (!assembly) { logvol = GetVolume(child_name); }
348  }
349  else if (tag=="file")
350  { logvol = FileRead(child); }
351  else if (tag=="position")
352  { VectorRead(child,position); }
353  else if (tag=="rotation")
354  { VectorRead(child,rotation); }
355  else if (tag=="scale")
356  { VectorRead(child,scale); }
357  else if (tag=="positionref")
358  { position = GetPosition(GenerateName(RefRead(child))); }
359  else if (tag=="rotationref")
360  { rotation = GetRotation(GenerateName(RefRead(child))); }
361  else if (tag=="scaleref")
362  { scale = GetScale(GenerateName(RefRead(child))); }
363  else
364  {
365  G4String error_msg = "Unknown tag in physvol: " + tag;
366  G4Exception("G4GDMLReadStructure::PhysvolRead()", "ReadError",
367  FatalException, error_msg);
368  return;
369  }
370  }
371 
372  G4Transform3D transform(GetRotationMatrix(rotation).inverse(),position);
373  transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
374 
375  if (pAssembly) // Fill assembly structure
376  {
377  if (!logvol) { return; }
378  pAssembly->AddPlacedVolume(logvol, transform);
379  }
380  else // Generate physical volume tree or do assembly imprint
381  {
382  if (assembly)
383  {
384  assembly->MakeImprint(pMotherLogical, transform, 0, check);
385  }
386  else
387  {
388  if (!logvol) { return; }
389  G4String pv_name = logvol->GetName() + "_PV";
391  ->Place(transform,pv_name,logvol,pMotherLogical,false,copynumber,check);
392 
393  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
394  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
395  }
396  }
397 }
398 
400 ReplicavolRead(const xercesc::DOMElement* const replicavolElement, G4int number)
401 {
402  G4LogicalVolume* logvol = 0;
403  for (xercesc::DOMNode* iter = replicavolElement->getFirstChild();
404  iter != 0; iter = iter->getNextSibling())
405  {
406  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
407 
408  const xercesc::DOMElement* const child
409  = dynamic_cast<xercesc::DOMElement*>(iter);
410  if (!child)
411  {
412  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
413  "InvalidRead", FatalException, "No child found!");
414  return;
415  }
416  const G4String tag = Transcode(child->getTagName());
417 
418  if (tag=="volumeref")
419  {
420  logvol = GetVolume(GenerateName(RefRead(child)));
421  }
422  else if (tag=="replicate_along_axis")
423  {
424  if (logvol) { ReplicaRead(child,logvol,number); }
425  }
426  else
427  {
428  G4String error_msg = "Unknown tag in ReplicavolRead: " + tag;
429  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
430  "ReadError", FatalException, error_msg);
431  }
432  }
433 }
434 
436 ReplicaRead(const xercesc::DOMElement* const replicaElement,
437  G4LogicalVolume* logvol, G4int number)
438 {
439  G4double width = 0.0;
440  G4double offset = 0.0;
441  G4ThreeVector position(0.0,0.0,0.0);
442  G4ThreeVector rotation(0.0,0.0,0.0);
443  EAxis axis = kUndefined;
444  G4String name;
445 
446  for (xercesc::DOMNode* iter = replicaElement->getFirstChild();
447  iter != 0; iter = iter->getNextSibling())
448  {
449  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
450 
451  const xercesc::DOMElement* const child
452  = dynamic_cast<xercesc::DOMElement*>(iter);
453  if (!child)
454  {
455  G4Exception("G4GDMLReadStructure::ReplicaRead()",
456  "InvalidRead", FatalException, "No child found!");
457  return;
458  }
459  const G4String tag = Transcode(child->getTagName());
460 
461  if (tag=="position")
462  { VectorRead(child,position); } else
463  if (tag=="rotation")
464  { VectorRead(child,rotation); } else
465  if (tag=="positionref")
466  { position = GetPosition(GenerateName(RefRead(child))); } else
467  if (tag=="rotationref")
468  { rotation = GetRotation(GenerateName(RefRead(child))); } else
469  if (tag=="direction")
470  { axis=AxisRead(child); } else
471  if (tag=="width")
472  { width=QuantityRead(child); } else
473  if (tag=="offset")
474  { offset=QuantityRead(child); }
475  else
476  {
477  G4String error_msg = "Unknown tag in ReplicaRead: " + tag;
478  G4Exception("G4GDMLReadStructure::ReplicaRead()", "ReadError",
479  FatalException, error_msg);
480  }
481  }
482 
483  G4String pv_name = logvol->GetName() + "_PV";
485  ->Replicate(pv_name,logvol,pMotherLogical,axis,number,width,offset);
486 
487  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
488  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
489 
490 }
491 
493 AxisRead(const xercesc::DOMElement* const axisElement)
494 {
495  EAxis axis = kUndefined;
496 
497  const xercesc::DOMNamedNodeMap* const attributes
498  = axisElement->getAttributes();
499  XMLSize_t attributeCount = attributes->getLength();
500 
501  for (XMLSize_t attribute_index=0;
502  attribute_index<attributeCount; attribute_index++)
503  {
504  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
505 
506  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
507  { continue; }
508 
509  const xercesc::DOMAttr* const attribute
510  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
511  if (!attribute)
512  {
513  G4Exception("G4GDMLReadStructure::AxisRead()",
514  "InvalidRead", FatalException, "No attribute found!");
515  return axis;
516  }
517  const G4String attName = Transcode(attribute->getName());
518  const G4String attValue = Transcode(attribute->getValue());
519  if (attName=="x")
520  { if( eval.Evaluate(attValue)==1.) {axis=kXAxis;} }
521  else if (attName=="y")
522  { if( eval.Evaluate(attValue)==1.) {axis=kYAxis;} }
523  else if (attName=="z")
524  { if( eval.Evaluate(attValue)==1.) {axis=kZAxis;} }
525  else if (attName=="rho")
526  { if( eval.Evaluate(attValue)==1.) {axis=kRho;} }
527  else if (attName=="phi")
528  { if( eval.Evaluate(attValue)==1.) {axis=kPhi;} }
529  }
530 
531  return axis;
532 }
533 
535 QuantityRead(const xercesc::DOMElement* const readElement)
536 {
537  G4double value = 0.0;
538  G4double unit = 0.0;
539  const xercesc::DOMNamedNodeMap* const attributes
540  = readElement->getAttributes();
541  XMLSize_t attributeCount = attributes->getLength();
542 
543  for (XMLSize_t attribute_index=0;
544  attribute_index<attributeCount; attribute_index++)
545  {
546  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
547 
548  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
549  { continue; }
550  const xercesc::DOMAttr* const attribute
551  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
552  if (!attribute)
553  {
554  G4Exception("G4GDMLReadStructure::QuantityRead()",
555  "InvalidRead", FatalException, "No attribute found!");
556  return value;
557  }
558  const G4String attName = Transcode(attribute->getName());
559  const G4String attValue = Transcode(attribute->getValue());
560 
561  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
562  if (G4UnitDefinition::GetCategory(attValue)!="Length" && G4UnitDefinition::GetCategory(attValue)!="Angle") {
563  G4Exception("G4GDMLReadStructure::QuantityRead()", "InvalidRead",
564  FatalException, "Invalid unit for lenght or angle (width, offset)!"); }
565  } else
566  if (attName=="value"){ value= eval.Evaluate(attValue); }
567  }
568 
569  return value*unit;
570 }
571 
573 VolumeRead(const xercesc::DOMElement* const volumeElement)
574 {
575  G4VSolid* solidPtr = 0;
576  G4Material* materialPtr = 0;
577  G4GDMLAuxListType auxList;
578 
579  XMLCh *name_attr = xercesc::XMLString::transcode("name");
580  const G4String name = Transcode(volumeElement->getAttribute(name_attr));
581  xercesc::XMLString::release(&name_attr);
582 
583  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
584  iter != 0; iter = iter->getNextSibling())
585  {
586  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
587 
588  const xercesc::DOMElement* const child
589  = dynamic_cast<xercesc::DOMElement*>(iter);
590  if (!child)
591  {
592  G4Exception("G4GDMLReadStructure::VolumeRead()",
593  "InvalidRead", FatalException, "No child found!");
594  return;
595  }
596  const G4String tag = Transcode(child->getTagName());
597 
598  if (tag=="auxiliary")
599  { auxList.push_back(AuxiliaryRead(child)); } else
600  if (tag=="materialref")
601  { materialPtr = GetMaterial(GenerateName(RefRead(child),true)); } else
602  if (tag=="solidref")
603  { solidPtr = GetSolid(GenerateName(RefRead(child))); }
604  }
605 
606  pMotherLogical = new G4LogicalVolume(solidPtr,materialPtr,
607  GenerateName(name),0,0,0);
608 
609  if (!auxList.empty()) { auxMap[pMotherLogical] = auxList; }
610 
611  Volume_contentRead(volumeElement);
612 }
613 
615 AssemblyRead(const xercesc::DOMElement* const assemblyElement)
616 {
617  XMLCh *name_attr = xercesc::XMLString::transcode("name");
618  const G4String name = Transcode(assemblyElement->getAttribute(name_attr));
619  xercesc::XMLString::release(&name_attr);
620 
621  G4AssemblyVolume* pAssembly = new G4AssemblyVolume();
622  assemblyMap.insert(std::make_pair(GenerateName(name), pAssembly));
623 
624  for (xercesc::DOMNode* iter = assemblyElement->getFirstChild();
625  iter != 0; iter = iter->getNextSibling())
626  {
627  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
628  const xercesc::DOMElement* const child
629  = dynamic_cast<xercesc::DOMElement*>(iter);
630  if (!child)
631  {
632  G4Exception("G4GDMLReadStructure::AssemblyRead()",
633  "InvalidRead", FatalException, "No child found!");
634  return;
635  }
636  const G4String tag = Transcode(child->getTagName());
637 
638  if (tag=="physvol")
639  {
640  PhysvolRead(child, pAssembly);
641  }
642  else
643  {
644  G4cout << "Unsupported GDML tag '" << tag
645  << "' for Geant4 assembly structure !" << G4endl;
646  }
647  }
648 }
649 
651 SkinSurfaceRead(const xercesc::DOMElement* const skinsurfaceElement)
652 {
653  G4String name;
654  G4LogicalVolume* logvol = 0;
655  G4SurfaceProperty* prop = 0;
656 
657  const xercesc::DOMNamedNodeMap* const attributes
658  = skinsurfaceElement->getAttributes();
659  XMLSize_t attributeCount = attributes->getLength();
660 
661  for (XMLSize_t attribute_index=0;
662  attribute_index<attributeCount; attribute_index++)
663  {
664  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
665 
666  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
667  { continue; }
668 
669  const xercesc::DOMAttr* const attribute
670  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
671  if (!attribute)
672  {
673  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
674  "InvalidRead", FatalException, "No attribute found!");
675  return;
676  }
677  const G4String attName = Transcode(attribute->getName());
678  const G4String attValue = Transcode(attribute->getValue());
679 
680  if (attName=="name")
681  { name = GenerateName(attValue); } else
682  if (attName=="surfaceproperty")
683  { prop = GetSurfaceProperty(GenerateName(attValue)); }
684  }
685 
686  for (xercesc::DOMNode* iter = skinsurfaceElement->getFirstChild();
687  iter != 0; iter = iter->getNextSibling())
688  {
689  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
690 
691  const xercesc::DOMElement* const child
692  = dynamic_cast<xercesc::DOMElement*>(iter);
693  if (!child)
694  {
695  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
696  "InvalidRead", FatalException, "No child found!");
697  return;
698  }
699  const G4String tag = Transcode(child->getTagName());
700 
701  if (tag=="volumeref")
702  {
703  logvol = GetVolume(GenerateName(RefRead(child)));
704  }
705  else
706  {
707  G4String error_msg = "Unknown tag in skinsurface: " + tag;
708  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
709  FatalException, error_msg);
710  }
711  }
712 
713  new G4LogicalSkinSurface(Strip(name),logvol,prop);
714 }
715 
717 Volume_contentRead(const xercesc::DOMElement* const volumeElement)
718 {
719  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
720  iter != 0; iter = iter->getNextSibling())
721  {
722  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
723 
724  const xercesc::DOMElement* const child
725  = dynamic_cast<xercesc::DOMElement*>(iter);
726  if (!child)
727  {
728  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
729  "InvalidRead", FatalException, "No child found!");
730  return;
731  }
732  const G4String tag = Transcode(child->getTagName());
733 
734  if ((tag=="auxiliary") || (tag=="materialref") || (tag=="solidref"))
735  {
736  // These are already processed in VolumeRead()
737  }
738  else if (tag=="paramvol")
739  {
741  }
742  else if (tag=="physvol")
743  {
744  PhysvolRead(child);
745  }
746  else if (tag=="replicavol")
747  {
748  G4int number = 1;
749  const xercesc::DOMNamedNodeMap* const attributes
750  = child->getAttributes();
751  XMLSize_t attributeCount = attributes->getLength();
752  for (XMLSize_t attribute_index=0;
753  attribute_index<attributeCount; attribute_index++)
754  {
755  xercesc::DOMNode* attribute_node
756  = attributes->item(attribute_index);
757  if (attribute_node->getNodeType()!=xercesc::DOMNode::ATTRIBUTE_NODE)
758  {
759  continue;
760  }
761  const xercesc::DOMAttr* const attribute
762  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
763  if (!attribute)
764  {
765  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
766  "InvalidRead", FatalException, "No attribute found!");
767  return;
768  }
769  const G4String attName = Transcode(attribute->getName());
770  const G4String attValue = Transcode(attribute->getValue());
771  if (attName=="number")
772  {
773  number = eval.EvaluateInteger(attValue);
774  }
775  }
776  ReplicavolRead(child,number);
777  }
778  else if (tag=="divisionvol")
779  {
780  DivisionvolRead(child);
781  }
782  else if (tag=="loop")
783  {
785  }
786  else
787  {
788  G4cout << "Treating unknown GDML tag in volume '" << tag
789  << "' as GDML extension..." << G4endl;
790  }
791  }
792 }
793 
795 StructureRead(const xercesc::DOMElement* const structureElement)
796 {
797  G4cout << "G4GDML: Reading structure..." << G4endl;
798 
799  for (xercesc::DOMNode* iter = structureElement->getFirstChild();
800  iter != 0; iter = iter->getNextSibling())
801  {
802  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
803 
804  const xercesc::DOMElement* const child
805  = dynamic_cast<xercesc::DOMElement*>(iter);
806  if (!child)
807  {
808  G4Exception("G4GDMLReadStructure::StructureRead()",
809  "InvalidRead", FatalException, "No child found!");
810  return;
811  }
812  const G4String tag = Transcode(child->getTagName());
813 
814  if (tag=="bordersurface") { BorderSurfaceRead(child); } else
815  if (tag=="skinsurface") { SkinSurfaceRead(child); } else
816  if (tag=="volume") { VolumeRead(child); } else
817  if (tag=="assembly") { AssemblyRead(child); } else
818  if (tag=="loop") { LoopRead(child,&G4GDMLRead::StructureRead); }
819  else
820  {
821  G4String error_msg = "Unknown tag in structure: " + tag;
822  G4Exception("G4GDMLReadStructure::StructureRead()",
823  "ReadError", FatalException, error_msg);
824  }
825  }
826 }
827 
829 GetPhysvol(const G4String& ref) const
830 {
831  G4VPhysicalVolume* physvolPtr =
833 
834  if (!physvolPtr)
835  {
836  G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
837  G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
838  FatalException, error_msg);
839  }
840 
841  return physvolPtr;
842 }
843 
845 GetVolume(const G4String& ref) const
846 {
847  G4LogicalVolume *volumePtr
849 
850  if (!volumePtr)
851  {
852  G4String error_msg = "Referenced volume '" + ref + "' was not found!";
853  G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError",
854  FatalException, error_msg);
855  }
856 
857  return volumePtr;
858 }
859 
861 GetAssembly(const G4String& ref) const
862 {
863  G4GDMLAssemblyMapType::const_iterator pos = assemblyMap.find(ref);
864  if (pos != assemblyMap.end()) { return pos->second; }
865  return 0;
866 }
867 
870 {
871  G4GDMLAuxMapType::const_iterator pos = auxMap.find(logvol);
872  if (pos != auxMap.end()) { return pos->second; }
873  else { return G4GDMLAuxListType(); }
874 }
875 
877 GetWorldVolume(const G4String& setupName)
878 {
879  G4String sname = GetSetup(setupName);
880  if (sname == "") { return 0; }
881 
882  G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
884 
885  G4VPhysicalVolume* pvWorld = 0;
886 
887  if(setuptoPV[setupName])
888  {
889  pvWorld = setuptoPV[setupName];
890  }
891  else
892  {
893  pvWorld = new G4PVPlacement(0,G4ThreeVector(0,0,0),volume,
894  volume->GetName()+"_PV",0,0,0);
895  setuptoPV[setupName] = pvWorld;
896  }
897  return pvWorld;
898 }
899 
901 {
902  eval.Clear();
903  setuptoPV.clear();
904 }
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
Definition: geomdefs.hh:54
const XML_Char * name
Definition: expat.h:151
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
double x() const
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
G4int first(char) const
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)
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
double z() 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
const XML_Char int const XML_Char * value
Definition: expat.h:331
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::Scale3D G4Scale3D
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
static G4String GetCategory(const G4String &)
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
double y() const
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
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
static const G4VisAttributes & GetInvisible()
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)