source: trunk/source/persistency/gdml/src/G4GDMLStructure.cc@ 1191

Last change on this file since 1191 was 818, checked in by garnier, 17 years ago

import all except CVS

File size: 18.5 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//
27// $Id: G4GDMLStructure.cc,v 1.25.2.1 2008/01/16 09:44:28 ztorzsok Exp $
28// GEANT4 tag $ Name:$
29//
30// class G4GDMLStructure Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLStructure.hh"
37
38EAxis G4GDMLStructure::directionRead(const xercesc::DOMElement* const element) {
39
40 G4String x;
41 G4String y;
42 G4String z;
43
44 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
45 XMLSize_t attributeCount = attributes->getLength();
46
47 for (XMLSize_t attribute_index=0;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) continue;
52
53 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
54
55 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
56 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
57
58 if (attribute_name=="x") x = attribute_value; else
59 if (attribute_name=="y") y = attribute_value; else
60 if (attribute_name=="z") z = attribute_value;
61 }
62
63 G4double _x = eval.Evaluate(x);
64 G4double _y = eval.Evaluate(y);
65 G4double _z = eval.Evaluate(z);
66
67 if (_x == 1.0 && _y == 0.0 && _z == 0.0) return kXAxis; else
68 if (_x == 0.0 && _y == 1.0 && _z == 0.0) return kYAxis; else
69 if (_x == 0.0 && _y == 0.0 && _z == 1.0) return kZAxis;
70
71 G4Exception("GDML: Only directions along axes are supported!");
72
73 return kZAxis;
74}
75
76void G4GDMLStructure::divisionvolRead(const xercesc::DOMElement* const element,G4LogicalVolume* pMother) {
77
78 G4String unit("1");
79 G4String axis;
80 G4String width;
81 G4String offset;
82 G4String number;
83 G4String volumeref;
84
85 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
86 XMLSize_t attributeCount = attributes->getLength();
87
88 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
89
90 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
91
92 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
93
94 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
95
96 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
97 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
98
99 if (attribute_name=="unit") unit = attribute_value; else
100 if (attribute_name=="axis") axis = attribute_value; else
101 if (attribute_name=="width") width = attribute_value; else
102 if (attribute_name=="offset") offset = attribute_value; else
103 if (attribute_name=="number") number = attribute_value;
104 }
105
106 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
107
108 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
109
110 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
111
112 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
113
114 if (tag=="volumeref") volumeref = refRead(child);
115 }
116
117 G4LogicalVolume* pLogical = getVolume(GenerateName(volumeref));
118
119 G4double _unit = eval.Evaluate(unit);
120
121 G4double _width = eval.Evaluate(width)*_unit;
122 G4double _offset = eval.Evaluate(offset)*_unit;
123 G4double _number = eval.Evaluate(number);
124
125 EAxis _axis = kZAxis;
126
127 if (axis=="kXAxis") _axis = kXAxis; else
128 if (axis=="kYAxis") _axis = kYAxis; else
129 if (axis=="kZAxis") _axis = kZAxis; else
130 if (axis=="kRho") _axis = kRho; else
131 if (axis=="kPhi") _axis = kPhi;
132
133 new G4PVDivision("",pLogical,pMother,_axis,(G4int)_number,_width,_offset);
134}
135
136G4LogicalVolume* G4GDMLStructure::fileRead(const xercesc::DOMElement* const element) {
137
138 G4String name;
139 G4String volname;
140
141 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
142 XMLSize_t attributeCount = attributes->getLength();
143
144 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
145
146 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
147
148 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
149
150 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
151
152 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
153 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
154
155 if (attribute_name=="name") name = attribute_value; else
156 if (attribute_name=="volname") volname = attribute_value;
157 }
158
159 G4GDMLStructure structure; // We create a new structure with a new evaluator
160
161 structure.Parse(name);
162
163 if (volname.empty()) return structure.getVolume(structure.getSetup("Default"));
164
165 return structure.getVolume(structure.GenerateName(volname));
166}
167
168void G4GDMLStructure::loopRead(const xercesc::DOMElement* const element) {
169
170 G4String var;
171 G4String from;
172 G4String to;
173 G4String step;
174
175 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
176 XMLSize_t attributeCount = attributes->getLength();
177
178 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
179
180 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
181
182 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
183
184 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
185
186 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
187 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
188
189 if (attribute_name=="var") var = attribute_value; else
190 if (attribute_name=="from") from = attribute_value; else
191 if (attribute_name=="to") to = attribute_value; else
192 if (attribute_name=="step") step = attribute_value;
193 }
194
195 eval.checkVariable(var);
196
197 G4int _var = eval.EvaluateInteger(var );
198 G4int _from = eval.EvaluateInteger(from);
199 G4int _to = eval.EvaluateInteger(to );
200 G4int _step = eval.EvaluateInteger(step);
201
202 if (!from.empty()) _var = _from;
203
204 while (_var <= _to) {
205
206 eval.setVariable(var,_var);
207 structureRead(element);
208
209 _var += _step;
210 }
211}
212
213void G4GDMLStructure::paramvolRead(const xercesc::DOMElement* const element,G4LogicalVolume *pMotherLogical) {
214
215 pMotherLogical = 0;
216
217 G4String volumeref;
218 G4String parameterised_position_size;
219
220 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
221
222 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
223
224 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
225
226 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
227
228 if (tag=="parameterised_position_size") { ; } else
229 if (tag=="volumeref") volumeref = refRead(child);
230 }
231}
232
233void G4GDMLStructure::physvolRead(const xercesc::DOMElement* const element,G4LogicalVolume *mother) {
234
235 G4String volumeref;
236 G4String positionref;
237 G4String rotationref;
238 G4String scaleref;
239
240 G4LogicalVolume* logvol = 0;
241
242 G4ThreeVector position;
243 G4ThreeVector rotation;
244 G4ThreeVector scale(1.0,1.0,1.0);
245
246 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
247
248 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
249
250 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
251
252 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
253
254 if (tag=="file") logvol = fileRead(child); else
255 if (tag=="volumeref") volumeref = refRead(child); else
256 if (tag=="position") position = positionRead(child); else
257 if (tag=="rotation") rotation = rotationRead(child); else
258 if (tag=="scale") scale = scaleRead(child); else
259 if (tag=="positionref") positionref = refRead(child); else
260 if (tag=="rotationref") rotationref = refRead(child); else
261 if (tag=="scaleref") scaleref = refRead(child);
262 }
263
264 if (!volumeref.empty()) logvol = getVolume(GenerateName(volumeref));
265
266 if (!positionref.empty()) position = *getPosition(GenerateName(positionref));
267 if (!rotationref.empty()) rotation = *getRotation(GenerateName(rotationref));
268 if (!scaleref.empty()) scale = *getScale(GenerateName(scaleref));
269
270 G4RotationMatrix Rot;
271
272 Rot.rotateX(rotation.x());
273 Rot.rotateY(rotation.y());
274 Rot.rotateZ(rotation.z());
275
276 G4Transform3D transform(Rot.inverse(),position);
277 transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
278
279 G4ReflectionFactory::Instance()->Place(transform,"",logvol,mother,false,0);
280}
281
282G4double G4GDMLStructure::quantityRead(const xercesc::DOMElement* const element) {
283
284 G4String value;
285 G4String unit("1");
286
287 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
288 XMLSize_t attributeCount = attributes->getLength();
289
290 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
291
292 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
293
294 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
295
296 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
297
298 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
299 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
300
301 if (attribute_name=="value") value = attribute_value;
302 if (attribute_name=="unit") unit = attribute_value;
303 }
304
305 return eval.Evaluate(value)*eval.Evaluate(unit);
306}
307
308G4String G4GDMLStructure::refRead(const xercesc::DOMElement* const element) {
309
310 G4String ref;
311
312 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
313 XMLSize_t attributeCount = attributes->getLength();
314
315 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
316
317 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
318
319 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
320
321 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
322
323 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
324 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
325
326 if (attribute_name=="ref") ref = attribute_value;
327 }
328
329 return ref;
330}
331
332void G4GDMLStructure::replicate_along_axisRead(const xercesc::DOMElement* const element,G4double& _width,G4double& _offset,EAxis& _axis) {
333
334 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
335
336 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
337
338 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
339
340 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
341
342 if (tag=="width") _width = quantityRead(child); else
343 if (tag=="offset") _offset = quantityRead(child); else
344 if (tag=="direction") _axis = directionRead(child);
345 }
346}
347
348void G4GDMLStructure::replicavolRead(const xercesc::DOMElement* const element,G4LogicalVolume* pMother) {
349
350 G4String volumeref;
351 G4String numb;
352
353 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
354 XMLSize_t attributeCount = attributes->getLength();
355
356 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
357
358 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
359
360 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
361
362 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
363
364 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
365 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
366
367 if (attribute_name=="numb") numb = attribute_value;
368 }
369
370 G4double _numb = eval.Evaluate(numb);
371 G4double _width = 0.0;
372 G4double _offset = 0.0;
373 EAxis _axis;
374
375 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
376
377 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
378
379 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
380
381 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
382
383 if (tag=="volumeref") volumeref = refRead(child); else
384 if (tag=="replicate_along_axis") replicate_along_axisRead(child,_width,_offset,_axis);
385 }
386
387 G4LogicalVolume* pLogical = getVolume(GenerateName(volumeref));
388
389 new G4PVReplica("",pLogical,pMother,_axis,(G4int)_numb,_width,_offset);
390}
391
392void G4GDMLStructure::volumeRead(const xercesc::DOMElement* const element) {
393
394 G4String name;
395 G4String solidref;
396 G4String materialref;
397
398 XMLCh *name_attr = xercesc::XMLString::transcode("name");
399 name = xercesc::XMLString::transcode(element->getAttribute(name_attr));
400 xercesc::XMLString::release(&name_attr);
401
402 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
403
404 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
405
406 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
407
408 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
409
410 if (tag=="materialref") materialref = refRead(child); else
411 if (tag=="solidref") solidref = refRead(child);
412 }
413
414 G4Material* materialPtr = getMaterial(GenerateName(materialref));
415 G4VSolid* solidPtr = getSolid(GenerateName(solidref));
416
417 volume_contentRead(element,new G4LogicalVolume(solidPtr,materialPtr,GenerateName(name),0,0,0));
418}
419
420void G4GDMLStructure::volume_contentRead(const xercesc::DOMElement* const element,G4LogicalVolume* volumePtr) {
421
422 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
423
424 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
425
426 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
427
428 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
429
430 if (tag=="loop") volume_loopRead(child,volumePtr); else
431 if (tag=="paramvol") paramvolRead(child,volumePtr); else
432 if (tag=="physvol") physvolRead(child,volumePtr); else
433 if (tag=="replicavol") replicavolRead(child,volumePtr); else
434 if (tag=="divisionvol") divisionvolRead(child,volumePtr);
435 }
436}
437
438void G4GDMLStructure::volume_loopRead(const xercesc::DOMElement* const element,G4LogicalVolume* volumePtr) {
439
440 G4String var;
441 G4String from;
442 G4String to;
443 G4String step;
444
445 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
446 XMLSize_t attributeCount = attributes->getLength();
447
448 for (XMLSize_t attribute_index=0;attribute_index<attributeCount;attribute_index++) {
449
450 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
451
452 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) continue;
453
454 const xercesc::DOMAttr* const attribute = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
455
456 const G4String attribute_name = xercesc::XMLString::transcode(attribute->getName());
457 const G4String attribute_value = xercesc::XMLString::transcode(attribute->getValue());
458
459 if (attribute_name=="var") var = attribute_value; else
460 if (attribute_name=="from") from = attribute_value; else
461 if (attribute_name=="to") to = attribute_value; else
462 if (attribute_name=="step") step = attribute_value;
463 }
464
465 eval.checkVariable(var);
466
467 G4int _var = eval.EvaluateInteger(var);
468 G4int _from = eval.EvaluateInteger(from);
469 G4int _to = eval.EvaluateInteger(to);
470 G4int _step = eval.EvaluateInteger(step);
471
472 if (!from.empty()) _var = _from;
473
474 while (_var <= _to) {
475
476 eval.setVariable(var,_var);
477 volume_contentRead(element,volumePtr);
478
479 _var += _step;
480 }
481}
482
483void G4GDMLStructure::structureRead(const xercesc::DOMElement* const element) {
484
485 for (xercesc::DOMNode* iter = element->getFirstChild();iter != 0;iter = iter->getNextSibling()) {
486
487 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) continue;
488
489 const xercesc::DOMElement* const child = dynamic_cast<xercesc::DOMElement*>(iter);
490
491 const G4String tag = xercesc::XMLString::transcode(child->getTagName());
492
493 if (tag=="loop") loopRead(child); else
494 if (tag=="volume") volumeRead(child); else
495 G4Exception("GDML: Unknown tag in structure: "+tag);
496 }
497}
498
499G4LogicalVolume* G4GDMLStructure::getVolume(const G4String& ref) const {
500
501 G4LogicalVolume *volumePtr = G4LogicalVolumeStore::GetInstance()->GetVolume(ref,false);
502
503 if (!volumePtr) G4Exception("GDML: Referenced volume '"+ref+"' was not found!");
504
505 return volumePtr;
506}
Note: See TracBrowser for help on using the repository browser.