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

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

geant4 tag 9.4

File size: 19.3 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.26 2010/10/14 16:19:40 gcosmo Exp $
27// GEANT4 tag $Name: gdml-V09-03-09 $
28//
29// class G4GDMLReadDefine Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
35#include "G4GDMLReadDefine.hh"
36
37G4GDMLMatrix::G4GDMLMatrix()
38 : m(0), rows(0), cols(0)
39{
40}
41
42G4GDMLMatrix::G4GDMLMatrix(size_t rows0, size_t cols0)
43{
44 if ((rows0==0) || (cols0==0))
45 {
46 G4Exception("G4GDMLMatrix::G4GDMLMatrix(r,c)", "InvalidSetup",
47 FatalException, "Zero indeces as arguments!?");
48 }
49 rows = rows0;
50 cols = cols0;
51 m = new G4double[rows*cols];
52}
53
54G4GDMLMatrix::G4GDMLMatrix(const G4GDMLMatrix& rhs)
55 : m(0), rows(0), cols(0)
56{
57 if (rhs.m)
58 {
59 rows = rhs.rows;
60 cols = rhs.cols;
61 m = new G4double[rows*cols];
62 for (size_t i=0; i<rows*cols; i++) { m[i] = rhs.m[i]; }
63 }
64}
65
66G4GDMLMatrix& G4GDMLMatrix::operator=(const G4GDMLMatrix& rhs)
67{
68 // Check assignment to self
69 //
70 if (this == &rhs) { return *this; }
71
72 // Copy data
73 //
74 rows = rhs.rows;
75 cols = rhs.cols;
76 if (rhs.m)
77 {
78 m = new G4double[rows*cols];
79 for (size_t i=0; i<rows*cols; i++) { m[i] = rhs.m[i]; }
80 }
81 else
82 {
83 m = 0;
84 }
85
86 return *this;
87}
88
89G4GDMLMatrix::~G4GDMLMatrix()
90{
91 delete [] m;
92}
93
94void G4GDMLMatrix::Set(size_t r,size_t c,G4double a)
95{
96 if (r>=rows || c>=cols)
97 {
98 G4Exception("G4GDMLMatrix::set()", "InvalidSetup",
99 FatalException, "Index out of range!");
100 }
101 m[cols*r+c] = a;
102}
103
104G4double G4GDMLMatrix::Get(size_t r,size_t c) const
105{
106 if (r>=rows || c>=cols)
107 {
108 G4Exception("G4GDMLMatrix::get()", "InvalidSetup",
109 FatalException, "Index out of range!");
110 }
111 return m[cols*r+c];
112}
113
114size_t G4GDMLMatrix::GetRows() const
115{
116 return rows;
117}
118
119size_t G4GDMLMatrix::GetCols() const
120{
121 return cols;
122}
123
124G4GDMLReadDefine::G4GDMLReadDefine() : G4GDMLRead()
125{
126}
127
128G4GDMLReadDefine::~G4GDMLReadDefine()
129{
130}
131
132G4RotationMatrix
133G4GDMLReadDefine::GetRotationMatrix(const G4ThreeVector& angles)
134{
135 G4RotationMatrix rot;
136
137 rot.rotateX(angles.x());
138 rot.rotateY(angles.y());
139 rot.rotateZ(angles.z());
140
141 return rot;
142}
143
144void
145G4GDMLReadDefine::ConstantRead(const xercesc::DOMElement* const constantElement)
146{
147 G4String name = "";
148 G4double value = 0.0;
149
150 const xercesc::DOMNamedNodeMap* const attributes
151 = constantElement->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* node = attributes->item(attribute_index);
158
159 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
160
161 const xercesc::DOMAttr* const attribute
162 = dynamic_cast<xercesc::DOMAttr*>(node);
163 if (!attribute)
164 {
165 G4Exception("G4GDMLRead::ConstantRead()", "InvalidRead",
166 FatalException, "No attribute found!");
167 return;
168 }
169 const G4String attName = Transcode(attribute->getName());
170 const G4String attValue = Transcode(attribute->getValue());
171
172 if (attName=="name") { name = attValue; } else
173 if (attName=="value") { value = eval.Evaluate(attValue); }
174 }
175
176 eval.DefineConstant(name,value);
177}
178
179void
180G4GDMLReadDefine::ExpressionRead(const xercesc::DOMElement* const expElement)
181{
182 G4String name = "";
183 G4double value = 0.0;
184
185 const xercesc::DOMNamedNodeMap* const attributes
186 = expElement->getAttributes();
187 XMLSize_t attributeCount = attributes->getLength();
188
189 for (XMLSize_t attribute_index=0;
190 attribute_index<attributeCount; attribute_index++)
191 {
192 xercesc::DOMNode* node = attributes->item(attribute_index);
193
194 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
195
196 const xercesc::DOMAttr* const attribute
197 = dynamic_cast<xercesc::DOMAttr*>(node);
198 if (!attribute)
199 {
200 G4Exception("G4GDMLRead::ExpressionRead()", "InvalidRead",
201 FatalException, "No attribute found!");
202 return;
203 }
204 const G4String attName = Transcode(attribute->getName());
205 const G4String attValue = Transcode(attribute->getValue());
206
207 if (attName=="name") { name = attValue; }
208 }
209
210 const G4String expValue = Transcode(expElement->getTextContent());
211 value = eval.Evaluate(expValue);
212 eval.DefineConstant(name,value);
213}
214
215void
216G4GDMLReadDefine::MatrixRead(const xercesc::DOMElement* const matrixElement)
217{
218 G4String name = "";
219 G4int coldim = 0;
220 G4String values = "";
221
222 const xercesc::DOMNamedNodeMap* const attributes
223 = matrixElement->getAttributes();
224 XMLSize_t attributeCount = attributes->getLength();
225
226 for (XMLSize_t attribute_index=0;
227 attribute_index<attributeCount; attribute_index++)
228 {
229 xercesc::DOMNode* node = attributes->item(attribute_index);
230
231 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
232
233 const xercesc::DOMAttr* const attribute
234 = dynamic_cast<xercesc::DOMAttr*>(node);
235 if (!attribute)
236 {
237 G4Exception("G4GDMLRead::MatrixRead()", "InvalidRead",
238 FatalException, "No attribute found!");
239 return;
240 }
241 const G4String attName = Transcode(attribute->getName());
242 const G4String attValue = Transcode(attribute->getValue());
243
244 if (attName=="name") { name = GenerateName(attValue); } else
245 if (attName=="coldim") { coldim = eval.EvaluateInteger(attValue); } else
246 if (attName=="values") { values = attValue; }
247 }
248
249 std::stringstream MatrixValueStream(values);
250 std::vector<G4double> valueList;
251
252 while (!MatrixValueStream.eof())
253 {
254 G4String MatrixValue;
255 MatrixValueStream >> MatrixValue;
256 valueList.push_back(eval.Evaluate(MatrixValue));
257 }
258
259 eval.DefineMatrix(name,coldim,valueList);
260
261 G4GDMLMatrix matrix(valueList.size()/coldim,coldim);
262
263 for (size_t i=0;i<valueList.size();i++)
264 {
265 matrix.Set(i/coldim,i%coldim,valueList[i]);
266 }
267
268 matrixMap[name] = matrix;
269}
270
271void
272G4GDMLReadDefine::PositionRead(const xercesc::DOMElement* const positionElement)
273{
274 G4String name = "";
275 G4double unit = 1.0;
276 G4ThreeVector position(0.,0.,0.);
277
278 const xercesc::DOMNamedNodeMap* const attributes
279 = positionElement->getAttributes();
280 XMLSize_t attributeCount = attributes->getLength();
281
282 for (XMLSize_t attribute_index=0;
283 attribute_index<attributeCount; attribute_index++)
284 {
285 xercesc::DOMNode* node = attributes->item(attribute_index);
286
287 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
288
289 const xercesc::DOMAttr* const attribute
290 = dynamic_cast<xercesc::DOMAttr*>(node);
291 if (!attribute)
292 {
293 G4Exception("G4GDMLRead::PositionRead()", "InvalidRead",
294 FatalException, "No attribute found!");
295 return;
296 }
297 const G4String attName = Transcode(attribute->getName());
298 const G4String attValue = Transcode(attribute->getValue());
299
300 if (attName=="name") { name = GenerateName(attValue); } else
301 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
302 if (attName=="x") { position.setX(eval.Evaluate(attValue)); } else
303 if (attName=="y") { position.setY(eval.Evaluate(attValue)); } else
304 if (attName=="z") { position.setZ(eval.Evaluate(attValue)); }
305 }
306
307 positionMap[name] = position*unit;
308}
309
310void
311G4GDMLReadDefine::RotationRead(const xercesc::DOMElement* const rotationElement)
312{
313 G4String name = "";
314 G4double unit = 1.0;
315 G4ThreeVector rotation(0.,0.,0.);
316
317 const xercesc::DOMNamedNodeMap* const attributes
318 = rotationElement->getAttributes();
319 XMLSize_t attributeCount = attributes->getLength();
320
321 for (XMLSize_t attribute_index=0;
322 attribute_index<attributeCount; attribute_index++)
323 {
324 xercesc::DOMNode* node = attributes->item(attribute_index);
325
326 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
327
328 const xercesc::DOMAttr* const attribute
329 = dynamic_cast<xercesc::DOMAttr*>(node);
330 if (!attribute)
331 {
332 G4Exception("G4GDMLRead::RotationRead()", "InvalidRead",
333 FatalException, "No attribute found!");
334 return;
335 }
336 const G4String attName = Transcode(attribute->getName());
337 const G4String attValue = Transcode(attribute->getValue());
338
339 if (attName=="name") { name = GenerateName(attValue); } else
340 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
341 if (attName=="x") { rotation.setX(eval.Evaluate(attValue)); } else
342 if (attName=="y") { rotation.setY(eval.Evaluate(attValue)); } else
343 if (attName=="z") { rotation.setZ(eval.Evaluate(attValue)); }
344 }
345
346 rotationMap[name] = rotation*unit;
347}
348
349void G4GDMLReadDefine::ScaleRead(const xercesc::DOMElement* const scaleElement)
350{
351 G4String name = "";
352 G4ThreeVector scale(1.0,1.0,1.0);
353
354 const xercesc::DOMNamedNodeMap* const attributes
355 = scaleElement->getAttributes();
356 XMLSize_t attributeCount = attributes->getLength();
357
358 for (XMLSize_t attribute_index=0;
359 attribute_index<attributeCount; attribute_index++)
360 {
361 xercesc::DOMNode* node = attributes->item(attribute_index);
362
363 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
364
365 const xercesc::DOMAttr* const attribute
366 = dynamic_cast<xercesc::DOMAttr*>(node);
367 if (!attribute)
368 {
369 G4Exception("G4GDMLRead::ScaleRead()", "InvalidRead",
370 FatalException, "No attribute found!");
371 return;
372 }
373 const G4String attName = Transcode(attribute->getName());
374 const G4String attValue = Transcode(attribute->getValue());
375
376 if (attName=="name") { name = GenerateName(attValue); } else
377 if (attName=="x") { scale.setX(eval.Evaluate(attValue)); } else
378 if (attName=="y") { scale.setY(eval.Evaluate(attValue)); } else
379 if (attName=="z") { scale.setZ(eval.Evaluate(attValue)); }
380 }
381
382 scaleMap[name] = scale;
383}
384
385void
386G4GDMLReadDefine::VariableRead(const xercesc::DOMElement* const variableElement)
387{
388 G4String name = "";
389 G4double value = 0.0;
390
391 const xercesc::DOMNamedNodeMap* const attributes
392 = variableElement->getAttributes();
393 XMLSize_t attributeCount = attributes->getLength();
394
395 for (XMLSize_t attribute_index=0;
396 attribute_index<attributeCount; attribute_index++)
397 {
398 xercesc::DOMNode* node = attributes->item(attribute_index);
399
400 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
401
402 const xercesc::DOMAttr* const attribute
403 = dynamic_cast<xercesc::DOMAttr*>(node);
404 if (!attribute)
405 {
406 G4Exception("G4GDMLRead::VariableRead()", "InvalidRead",
407 FatalException, "No attribute found!");
408 return;
409 }
410 const G4String attName = Transcode(attribute->getName());
411 const G4String attValue = Transcode(attribute->getValue());
412
413 if (attName=="name") { name = attValue; } else
414 if (attName=="value") { value = eval.Evaluate(attValue); }
415 }
416
417 eval.DefineVariable(name,value);
418}
419
420void G4GDMLReadDefine::QuantityRead(const xercesc::DOMElement* const element)
421{
422 G4String name = "";
423 G4double unit = 1.0;
424 G4double value = 0.0;
425
426 const xercesc::DOMNamedNodeMap* const attributes
427 = element->getAttributes();
428 XMLSize_t attributeCount = attributes->getLength();
429
430 for (XMLSize_t attribute_index=0;
431 attribute_index<attributeCount; attribute_index++)
432 {
433 xercesc::DOMNode* node = attributes->item(attribute_index);
434
435 if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
436
437 const xercesc::DOMAttr* const attribute
438 = dynamic_cast<xercesc::DOMAttr*>(node);
439 if (!attribute)
440 {
441 G4Exception("G4GDMLRead::QuantityRead()", "InvalidRead",
442 FatalException, "No attribute found!");
443 return;
444 }
445 const G4String attName = Transcode(attribute->getName());
446 const G4String attValue = Transcode(attribute->getValue());
447
448 if (attName=="name") { name = attValue; } else
449 if (attName=="value") { value = eval.Evaluate(attValue); } else
450 if (attName=="unit") { unit = eval.Evaluate(attValue); }
451 }
452
453 quantityMap[name] = value*unit;
454 eval.DefineConstant(name,value*unit);
455}
456
457void
458G4GDMLReadDefine::DefineRead(const xercesc::DOMElement* const defineElement)
459{
460 G4cout << "G4GDML: Reading definitions..." << G4endl;
461
462 for (xercesc::DOMNode* iter = defineElement->getFirstChild();
463 iter != 0;iter = iter->getNextSibling())
464 {
465 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
466
467 const xercesc::DOMElement* const child
468 = dynamic_cast<xercesc::DOMElement*>(iter);
469 if (!child)
470 {
471 G4Exception("G4GDMLRead::DefineRead()", "InvalidRead",
472 FatalException, "No child found!");
473 return;
474 }
475 const G4String tag = Transcode(child->getTagName());
476
477 if (tag=="constant") { ConstantRead(child); } else
478 if (tag=="matrix") { MatrixRead(child); } else
479 if (tag=="position") { PositionRead(child); } else
480 if (tag=="rotation") { RotationRead(child); } else
481 if (tag=="scale") { ScaleRead(child); } else
482 if (tag=="variable") { VariableRead(child); } else
483 if (tag=="quantity") { QuantityRead(child); } else
484 if (tag=="expression") { ExpressionRead(child); }
485 else
486 {
487 G4String error_msg = "Unknown tag in define: "+tag;
488 G4Exception("G4GDMLReadDefine::defineRead()", "ReadError",
489 FatalException, error_msg);
490 }
491 }
492}
493
494void
495G4GDMLReadDefine::VectorRead(const xercesc::DOMElement* const vectorElement,
496 G4ThreeVector& vec)
497{
498 G4double unit = 1.0;
499
500 const xercesc::DOMNamedNodeMap* const attributes
501 = vectorElement->getAttributes();
502 XMLSize_t attributeCount = attributes->getLength();
503
504 for (XMLSize_t attribute_index=0;
505 attribute_index<attributeCount; attribute_index++)
506 {
507 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
508
509 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
510 { continue; }
511
512 const xercesc::DOMAttr* const attribute
513 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
514 if (!attribute)
515 {
516 G4Exception("G4GDMLRead::VectorRead()", "InvalidRead",
517 FatalException, "No attribute found!");
518 return;
519 }
520 const G4String attName = Transcode(attribute->getName());
521 const G4String attValue = Transcode(attribute->getValue());
522
523 if (attName=="unit") { unit = eval.Evaluate(attValue); } else
524 if (attName=="x") { vec.setX(eval.Evaluate(attValue)); } else
525 if (attName=="y") { vec.setY(eval.Evaluate(attValue)); } else
526 if (attName=="z") { vec.setZ(eval.Evaluate(attValue)); }
527 }
528
529 vec *= unit;
530}
531
532G4String G4GDMLReadDefine::RefRead(const xercesc::DOMElement* const element)
533{
534 G4String ref;
535
536 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
537 XMLSize_t attributeCount = attributes->getLength();
538
539 for (XMLSize_t attribute_index=0;
540 attribute_index<attributeCount; attribute_index++)
541 {
542 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
543
544 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
545 { continue; }
546
547 const xercesc::DOMAttr* const attribute
548 = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
549 if (!attribute)
550 {
551 G4Exception("G4GDMLRead::Read()", "InvalidRead",
552 FatalException, "No attribute found!");
553 return ref;
554 }
555 const G4String attName = Transcode(attribute->getName());
556 const G4String attValue = Transcode(attribute->getValue());
557
558 if (attName=="ref") { ref = attValue; }
559 }
560
561 return ref;
562}
563
564G4bool G4GDMLReadDefine::IsValidID(const G4String& ref) const
565{
566 return eval.IsVariable(ref);
567}
568
569G4double G4GDMLReadDefine::GetConstant(const G4String& ref)
570{
571 return eval.GetConstant(ref);
572}
573
574G4double G4GDMLReadDefine::GetVariable(const G4String& ref)
575{
576 return eval.GetVariable(ref);
577}
578
579G4double G4GDMLReadDefine::GetQuantity(const G4String& ref)
580{
581 if (quantityMap.find(ref) == quantityMap.end())
582 {
583 G4String error_msg = "Quantity '"+ref+"' was not found!";
584 G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError",
585 FatalException, error_msg);
586 }
587 return quantityMap[ref];
588}
589
590G4ThreeVector G4GDMLReadDefine::GetPosition(const G4String& ref)
591{
592 if (positionMap.find(ref) == positionMap.end())
593 {
594 G4String error_msg = "Position '"+ref+"' was not found!";
595 G4Exception("G4GDMLReadDefine::getPosition()", "ReadError",
596 FatalException, error_msg);
597 }
598 return positionMap[ref];
599}
600
601G4ThreeVector G4GDMLReadDefine::GetRotation(const G4String& ref)
602{
603 if (rotationMap.find(ref) == rotationMap.end())
604 {
605 G4String error_msg = "Rotation '"+ref+"' was not found!";
606 G4Exception("G4GDMLReadDefine::getRotation()", "ReadError",
607 FatalException, error_msg);
608 }
609 return rotationMap[ref];
610}
611
612G4ThreeVector G4GDMLReadDefine::GetScale(const G4String& ref)
613{
614 if (scaleMap.find(ref) == scaleMap.end())
615 {
616 G4String error_msg = "Scale '"+ref+"' was not found!";
617 G4Exception("G4GDMLReadDefine::getScale()", "ReadError",
618 FatalException, error_msg);
619 }
620 return scaleMap[ref];
621}
622
623G4GDMLMatrix G4GDMLReadDefine::GetMatrix(const G4String& ref)
624{
625 if (matrixMap.find(ref) == matrixMap.end())
626 {
627 G4String error_msg = "Matrix '"+ref+"' was not found!";
628 G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError",
629 FatalException, error_msg);
630 }
631 return matrixMap[ref];
632}
Note: See TracBrowser for help on using the repository browser.