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

Last change on this file since 1302 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

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