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