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