source: trunk/source/persistency/gdml/src/G4GDMLReadDefine.cc@ 1193

Last change on this file since 1193 was 987, checked in by garnier, 17 years ago

fichiers manquants

File size: 15.4 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: G4GDMLReadDefine.cc,v 1.20 2008/07/16 15:46:34 gcosmo Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// class G4GDMLReadDefine Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
35#include "G4GDMLReadDefine.hh"
36
37G4GDMLMatrix::G4GDMLMatrix()
38{
39 rows = 0;
40 cols = 0;
41 m = 0;
42}
43
44G4GDMLMatrix::G4GDMLMatrix(size_t rows0,size_t cols0)
45{
46 rows = rows0;
47 cols = cols0;
48 m = new G4double[rows*cols];
49}
50
51G4GDMLMatrix::~G4GDMLMatrix()
52{
53 if (m) { delete [] m; }
54}
55
56void G4GDMLMatrix::Set(size_t r,size_t c,G4double a)
57{
58 if (r>=rows || c>=cols)
59 {
60 G4Exception("G4GDMLMatrix::set()", "InvalidSetup",
61 FatalException, "Index out of range!");
62 }
63 m[cols*r+c] = a;
64}
65
66G4double G4GDMLMatrix::Get(size_t r,size_t c) const
67{
68 if (r>=rows || c>=cols)
69 {
70 G4Exception("G4GDMLMatrix::get()", "InvalidSetup",
71 FatalException, "Index out of range!");
72 }
73 return m[cols*r+c];
74}
75
76size_t G4GDMLMatrix::GetRows() const
77{
78 return rows;
79}
80
81size_t G4GDMLMatrix::GetCols() const
82{
83 return cols;
84}
85
86G4RotationMatrix
87G4GDMLReadDefine::GetRotationMatrix(const G4ThreeVector& angles)
88{
89 G4RotationMatrix rot;
90
91 rot.rotateX(angles.x());
92 rot.rotateY(angles.y());
93 rot.rotateZ(angles.z());
94
95 return rot;
96}
97
98void
99G4GDMLReadDefine::ConstantRead(const xercesc::DOMElement* const constantElement)
100{
101 G4String name = "";
102 G4double value = 0.0;
103
104 const xercesc::DOMNamedNodeMap* const attributes
105 = constantElement->getAttributes();
106 XMLSize_t attributeCount = attributes->getLength();
107
108 for (XMLSize_t attribute_index=0;
109 attribute_index<attributeCount; attribute_index++)
110 {
111 xercesc::DOMNode* node = attributes->item(attribute_index);
112
113 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
114
115 const xercesc::DOMAttr* const attribute
116 = dynamic_cast<xercesc::DOMAttr*>(node);
117 const G4String attName = Transcode(attribute->getName());
118 const G4String attValue = Transcode(attribute->getValue());
119
120 if (attName=="name") { name = attValue; } else
121 if (attName=="value") { value = eval.Evaluate(attValue); }
122 }
123
124 eval.DefineConstant(name,value);
125}
126
127void
128G4GDMLReadDefine::MatrixRead(const xercesc::DOMElement* const matrixElement)
129{
130 G4String name = "";
131 G4int coldim = 0;
132 G4String values = "";
133
134 const xercesc::DOMNamedNodeMap* const attributes
135 = matrixElement->getAttributes();
136 XMLSize_t attributeCount = attributes->getLength();
137
138 for (XMLSize_t attribute_index=0;
139 attribute_index<attributeCount; attribute_index++)
140 {
141 xercesc::DOMNode* node = attributes->item(attribute_index);
142
143 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
144
145 const xercesc::DOMAttr* const attribute
146 = dynamic_cast<xercesc::DOMAttr*>(node);
147 const G4String attName = Transcode(attribute->getName());
148 const G4String attValue = Transcode(attribute->getValue());
149
150 if (attName=="name") { name = GenerateName(attValue); } else
151 if (attName=="coldim") { coldim = eval.EvaluateInteger(attValue); } else
152 if (attName=="values") { values = attValue; }
153 }
154
155 std::stringstream MatrixValueStream(values);
156 std::vector<G4double> valueList;
157
158 while (!MatrixValueStream.eof())
159 {
160 G4String MatrixValue;
161 MatrixValueStream >> MatrixValue;
162 valueList.push_back(eval.Evaluate(MatrixValue));
163 }
164
165 eval.DefineMatrix(name,coldim,valueList);
166
167 G4GDMLMatrix matrix(valueList.size()/coldim,coldim);
168
169 for (size_t i=0;i<valueList.size();i++)
170 {
171 matrix.Set(i/coldim,i%coldim,valueList[i]);
172 }
173
174 matrixMap[name] = matrix;
175}
176
177void
178G4GDMLReadDefine::PositionRead(const xercesc::DOMElement* const positionElement)
179{
180 G4String name = "";
181 G4double unit = 1.0;
182 G4ThreeVector position(0.,0.,0.);
183
184 const xercesc::DOMNamedNodeMap* const attributes
185 = positionElement->getAttributes();
186 XMLSize_t attributeCount = attributes->getLength();
187
188 for (XMLSize_t attribute_index=0;
189 attribute_index<attributeCount; attribute_index++)
190 {
191 xercesc::DOMNode* node = attributes->item(attribute_index);
192
193 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
194
195 const xercesc::DOMAttr* const attribute
196 = dynamic_cast<xercesc::DOMAttr*>(node);
197 const G4String attName = Transcode(attribute->getName());
198 const G4String attValue = Transcode(attribute->getValue());
199
200 if (attName=="name") { name = GenerateName(attValue); } else
201 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
202 if (attName=="x") { position.setX(eval.Evaluate(attValue)); } else
203 if (attName=="y") { position.setY(eval.Evaluate(attValue)); } else
204 if (attName=="z") { position.setZ(eval.Evaluate(attValue)); }
205 }
206
207 positionMap[name] = position*unit;
208}
209
210void
211G4GDMLReadDefine::RotationRead(const xercesc::DOMElement* const rotationElement)
212{
213 G4String name = "";
214 G4double unit = 1.0;
215 G4ThreeVector rotation(0.,0.,0.);
216
217 const xercesc::DOMNamedNodeMap* const attributes
218 = rotationElement->getAttributes();
219 XMLSize_t attributeCount = attributes->getLength();
220
221 for (XMLSize_t attribute_index=0;
222 attribute_index<attributeCount; attribute_index++)
223 {
224 xercesc::DOMNode* node = attributes->item(attribute_index);
225
226 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
227
228 const xercesc::DOMAttr* const attribute
229 = dynamic_cast<xercesc::DOMAttr*>(node);
230 const G4String attName = Transcode(attribute->getName());
231 const G4String attValue = Transcode(attribute->getValue());
232
233 if (attName=="name") { name = GenerateName(attValue); } else
234 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
235 if (attName=="x") { rotation.setX(eval.Evaluate(attValue)); } else
236 if (attName=="y") { rotation.setY(eval.Evaluate(attValue)); } else
237 if (attName=="z") { rotation.setZ(eval.Evaluate(attValue)); }
238 }
239
240 rotationMap[name] = rotation*unit;
241}
242
243void G4GDMLReadDefine::ScaleRead(const xercesc::DOMElement* const scaleElement)
244{
245 G4String name = "";
246 G4ThreeVector scale(1.0,1.0,1.0);
247
248 const xercesc::DOMNamedNodeMap* const attributes
249 = scaleElement->getAttributes();
250 XMLSize_t attributeCount = attributes->getLength();
251
252 for (XMLSize_t attribute_index=0;
253 attribute_index<attributeCount; attribute_index++)
254 {
255 xercesc::DOMNode* node = attributes->item(attribute_index);
256
257 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
258
259 const xercesc::DOMAttr* const attribute
260 = dynamic_cast<xercesc::DOMAttr*>(node);
261 const G4String attName = Transcode(attribute->getName());
262 const G4String attValue = Transcode(attribute->getValue());
263
264 if (attName=="name") { name = GenerateName(attValue); } else
265 if (attName=="x") { scale.setX(eval.Evaluate(attValue)); } else
266 if (attName=="y") { scale.setY(eval.Evaluate(attValue)); } else
267 if (attName=="z") { scale.setZ(eval.Evaluate(attValue)); }
268 }
269
270 scaleMap[name] = scale;
271}
272
273void
274G4GDMLReadDefine::VariableRead(const xercesc::DOMElement* const variableElement)
275{
276 G4String name = "";
277 G4double value = 0.0;
278
279 const xercesc::DOMNamedNodeMap* const attributes
280 = variableElement->getAttributes();
281 XMLSize_t attributeCount = attributes->getLength();
282
283 for (XMLSize_t attribute_index=0;
284 attribute_index<attributeCount; attribute_index++)
285 {
286 xercesc::DOMNode* node = attributes->item(attribute_index);
287
288 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
289
290 const xercesc::DOMAttr* const attribute
291 = dynamic_cast<xercesc::DOMAttr*>(node);
292 const G4String attName = Transcode(attribute->getName());
293 const G4String attValue = Transcode(attribute->getValue());
294
295 if (attName=="name") { name = attValue; } else
296 if (attName=="value") { value = eval.Evaluate(attValue); }
297 }
298
299 eval.DefineVariable(name,value);
300}
301
302void G4GDMLReadDefine::QuantityRead(const xercesc::DOMElement* const element)
303{
304 G4String name = "";
305 G4double unit = 1.0;
306 G4double value = 0.0;
307
308 const xercesc::DOMNamedNodeMap* const attributes
309 = element->getAttributes();
310 XMLSize_t attributeCount = attributes->getLength();
311
312 for (XMLSize_t attribute_index=0;
313 attribute_index<attributeCount; attribute_index++)
314 {
315 xercesc::DOMNode* node = attributes->item(attribute_index);
316
317 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
318
319 const xercesc::DOMAttr* const attribute
320 = dynamic_cast<xercesc::DOMAttr*>(node);
321 const G4String attName = Transcode(attribute->getName());
322 const G4String attValue = Transcode(attribute->getValue());
323
324 if (attName=="name") { name = attValue; } else
325 if (attName=="value") { value = eval.Evaluate(attValue); } else
326 if (attName=="unit") { unit = eval.Evaluate(attValue); }
327 }
328
329 quantityMap[name] = value*unit;
330}
331
332void
333G4GDMLReadDefine::DefineRead(const xercesc::DOMElement* const defineElement)
334{
335 G4cout << "G4GDML: Reading definitions..." << G4endl;
336
337 for (xercesc::DOMNode* iter = defineElement->getFirstChild();
338 iter != 0;iter = iter->getNextSibling())
339 {
340 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
341
342 const xercesc::DOMElement* const child
343 = dynamic_cast<xercesc::DOMElement*>(iter);
344 const G4String tag = Transcode(child->getTagName());
345
346 if (tag=="constant") { ConstantRead(child); } else
347 if (tag=="matrix") { MatrixRead(child); } else
348 if (tag=="position") { PositionRead(child); } else
349 if (tag=="rotation") { RotationRead(child); } else
350 if (tag=="scale") { ScaleRead(child); } else
351 if (tag=="variable") { VariableRead(child); } else
352 if (tag=="quantity") { QuantityRead(child); }
353 else
354 {
355 G4String error_msg = "Unknown tag in define: "+tag;
356 G4Exception("G4GDMLReadDefine::defineRead()", "ReadError",
357 FatalException, error_msg);
358 }
359 }
360}
361
362void
363G4GDMLReadDefine::VectorRead(const xercesc::DOMElement* const vectorElement,
364 G4ThreeVector& vec)
365{
366 G4double unit = 1.0;
367
368 const xercesc::DOMNamedNodeMap* const attributes
369 = vectorElement->getAttributes();
370 XMLSize_t attributeCount = attributes->getLength();
371
372 for (XMLSize_t attribute_index=0;
373 attribute_index<attributeCount; attribute_index++)
374 {
375 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
376
377 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
378 { continue; }
379
380 const xercesc::DOMAttr* const attribute
381 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
382 const G4String attName = Transcode(attribute->getName());
383 const G4String attValue = Transcode(attribute->getValue());
384
385 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
386 if (attName=="x") { vec.setX(eval.Evaluate(attValue)); } else
387 if (attName=="y") { vec.setY(eval.Evaluate(attValue)); } else
388 if (attName=="z") { vec.setZ(eval.Evaluate(attValue)); }
389 }
390
391 vec *= unit;
392}
393
394G4String G4GDMLReadDefine::RefRead(const xercesc::DOMElement* const element)
395{
396 G4String ref;
397
398 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
399 XMLSize_t attributeCount = attributes->getLength();
400
401 for (XMLSize_t attribute_index=0;
402 attribute_index<attributeCount; attribute_index++)
403 {
404 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
405
406 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
407 { continue; }
408
409 const xercesc::DOMAttr* const attribute
410 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
411 const G4String attName = Transcode(attribute->getName());
412 const G4String attValue = Transcode(attribute->getValue());
413
414 if (attName=="ref") { ref = attValue; }
415 }
416
417 return ref;
418}
419
420G4double G4GDMLReadDefine::GetConstant(const G4String& ref)
421{
422 return eval.GetConstant(ref);
423}
424
425G4double G4GDMLReadDefine::GetVariable(const G4String& ref)
426{
427 return eval.GetVariable(ref);
428}
429
430G4double G4GDMLReadDefine::GetQuantity(const G4String& ref)
431{
432 if (quantityMap.find(ref) == quantityMap.end())
433 {
434 G4String error_msg = "Quantity '"+ref+"' was not found!";
435 G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError",
436 FatalException, error_msg);
437 }
438 return quantityMap[ref];
439}
440
441G4ThreeVector G4GDMLReadDefine::GetPosition(const G4String& ref)
442{
443 if (positionMap.find(ref) == positionMap.end())
444 {
445 G4String error_msg = "Position '"+ref+"' was not found!";
446 G4Exception("G4GDMLReadDefine::getPosition()", "ReadError",
447 FatalException, error_msg);
448 }
449 return positionMap[ref];
450}
451
452G4ThreeVector G4GDMLReadDefine::GetRotation(const G4String& ref)
453{
454 if (rotationMap.find(ref) == rotationMap.end())
455 {
456 G4String error_msg = "Rotation '"+ref+"' was not found!";
457 G4Exception("G4GDMLReadDefine::getRotation()", "ReadError",
458 FatalException, error_msg);
459 }
460 return rotationMap[ref];
461}
462
463G4ThreeVector G4GDMLReadDefine::GetScale(const G4String& ref)
464{
465 if (scaleMap.find(ref) == scaleMap.end())
466 {
467 G4String error_msg = "Scale '"+ref+"' was not found!";
468 G4Exception("G4GDMLReadDefine::getScale()", "ReadError",
469 FatalException, error_msg);
470 }
471 return scaleMap[ref];
472}
473
474G4GDMLMatrix G4GDMLReadDefine::GetMatrix(const G4String& ref)
475{
476 if (matrixMap.find(ref) == matrixMap.end())
477 {
478 G4String error_msg = "Matrix '"+ref+"' was not found!";
479 G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError",
480 FatalException, error_msg);
481 }
482 return matrixMap[ref];
483}
Note: See TracBrowser for help on using the repository browser.