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

Last change on this file since 1102 was 987, checked in by garnier, 15 years ago

fichiers manquants

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