source: trunk/source/persistency/gdml/src/G4GDMLReadStructure.cc@ 1349

Last change on this file since 1349 was 1347, checked in by garnier, 15 years ago

geant4 tag 9.4

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