Geant4  10.03.p02
 All Classes 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 96190 2016-03-29 08:07:36Z gcosmo $
27 //
28 // class G4GDMLReadDefine Implementation
29 //
30 // Original author: Zoltan Torzsok, November 2007
31 //
32 // --------------------------------------------------------------------
33 
34 #include "G4GDMLReadDefine.hh"
35 #include "G4UnitsTable.hh"
36 
38  : m(0), rows(0), cols(0)
39 {
40 }
41 
42 G4GDMLMatrix::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 
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 
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 
90 {
91  delete [] m;
92 }
93 
94 void 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 
104 G4double 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 
114 size_t G4GDMLMatrix::GetRows() const
115 {
116  return rows;
117 }
118 
119 size_t G4GDMLMatrix::GetCols() const
120 {
121  return cols;
122 }
123 
125 {
126 }
127 
129 {
130 }
131 
134 {
135  G4RotationMatrix rot;
136 
137  rot.rotateX(angles.x());
138  rot.rotateY(angles.y());
139  rot.rotateZ(angles.z());
140  rot.rectify(); // Rectify matrix from possible roundoff errors
141 
142  return rot;
143 }
144 
145 void
146 G4GDMLReadDefine::ConstantRead(const xercesc::DOMElement* const constantElement)
147 {
148  G4String name = "";
149  G4double value = 0.0;
150 
151  const xercesc::DOMNamedNodeMap* const attributes
152  = constantElement->getAttributes();
153  XMLSize_t attributeCount = attributes->getLength();
154 
155  for (XMLSize_t attribute_index=0;
156  attribute_index<attributeCount; attribute_index++)
157  {
158  xercesc::DOMNode* node = attributes->item(attribute_index);
159 
160  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
161 
162  const xercesc::DOMAttr* const attribute
163  = dynamic_cast<xercesc::DOMAttr*>(node);
164  if (!attribute)
165  {
166  G4Exception("G4GDMLRead::ConstantRead()", "InvalidRead",
167  FatalException, "No attribute found!");
168  return;
169  }
170  const G4String attName = Transcode(attribute->getName());
171  const G4String attValue = Transcode(attribute->getValue());
172 
173  if (attName=="name") { name = attValue; } else
174  if (attName=="value") { value = eval.Evaluate(attValue); }
175  }
176 
177  eval.DefineConstant(name,value);
178 }
179 
180 void
181 G4GDMLReadDefine::ExpressionRead(const xercesc::DOMElement* const expElement)
182 {
183  G4String name = "";
184  G4double value = 0.0;
185 
186  const xercesc::DOMNamedNodeMap* const attributes
187  = expElement->getAttributes();
188  XMLSize_t attributeCount = attributes->getLength();
189 
190  for (XMLSize_t attribute_index=0;
191  attribute_index<attributeCount; attribute_index++)
192  {
193  xercesc::DOMNode* node = attributes->item(attribute_index);
194 
195  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
196 
197  const xercesc::DOMAttr* const attribute
198  = dynamic_cast<xercesc::DOMAttr*>(node);
199  if (!attribute)
200  {
201  G4Exception("G4GDMLRead::ExpressionRead()", "InvalidRead",
202  FatalException, "No attribute found!");
203  return;
204  }
205  const G4String attName = Transcode(attribute->getName());
206  const G4String attValue = Transcode(attribute->getValue());
207 
208  if (attName=="name") { name = attValue; }
209  }
210 
211  const G4String expValue = Transcode(expElement->getTextContent());
212  value = eval.Evaluate(expValue);
213  eval.DefineConstant(name,value);
214 }
215 
216 void
217 G4GDMLReadDefine::MatrixRead(const xercesc::DOMElement* const matrixElement)
218 {
219  G4String name = "";
220  G4int coldim = 0;
221  G4String values = "";
222 
223  const xercesc::DOMNamedNodeMap* const attributes
224  = matrixElement->getAttributes();
225  XMLSize_t attributeCount = attributes->getLength();
226 
227  for (XMLSize_t attribute_index=0;
228  attribute_index<attributeCount; attribute_index++)
229  {
230  xercesc::DOMNode* node = attributes->item(attribute_index);
231 
232  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
233 
234  const xercesc::DOMAttr* const attribute
235  = dynamic_cast<xercesc::DOMAttr*>(node);
236  if (!attribute)
237  {
238  G4Exception("G4GDMLRead::MatrixRead()", "InvalidRead",
239  FatalException, "No attribute found!");
240  return;
241  }
242  const G4String attName = Transcode(attribute->getName());
243  const G4String attValue = Transcode(attribute->getValue());
244 
245  if (attName=="name") { name = GenerateName(attValue); } else
246  if (attName=="coldim") { coldim = eval.EvaluateInteger(attValue); } else
247  if (attName=="values") { values = attValue; }
248  }
249 
250  std::stringstream MatrixValueStream(values);
251  std::vector<G4double> valueList;
252 
253  while (!MatrixValueStream.eof())
254  {
255  G4String MatrixValue;
256  MatrixValueStream >> MatrixValue;
257  valueList.push_back(eval.Evaluate(MatrixValue));
258  }
259 
260  eval.DefineMatrix(name,coldim,valueList);
261 
262  G4GDMLMatrix matrix(valueList.size()/coldim,coldim);
263 
264  for (size_t i=0;i<valueList.size();i++)
265  {
266  matrix.Set(i/coldim,i%coldim,valueList[i]);
267  }
268 
269  matrixMap[name] = matrix;
270 }
271 
272 void
273 G4GDMLReadDefine::PositionRead(const xercesc::DOMElement* const positionElement)
274 {
275  G4String name = "";
276  G4double unit = 1.0;
277  G4ThreeVector position(0.,0.,0.);
278 
279  const xercesc::DOMNamedNodeMap* const attributes
280  = positionElement->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  if (!attribute)
293  {
294  G4Exception("G4GDMLRead::PositionRead()", "InvalidRead",
295  FatalException, "No attribute found!");
296  return;
297  }
298  const G4String attName = Transcode(attribute->getName());
299  const G4String attValue = Transcode(attribute->getValue());
300 
301  if (attName=="name") { name = GenerateName(attValue); } else
302  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
303  if (G4UnitDefinition::GetCategory(attValue)!="Length") {
304  G4Exception("G4GDMLReadDefine::PositionRead()", "InvalidRead",
305  FatalException, "Invalid unit for length!"); }
306  } else
307  if (attName=="x") { position.setX(eval.Evaluate(attValue)); } else
308  if (attName=="y") { position.setY(eval.Evaluate(attValue)); } else
309  if (attName=="z") { position.setZ(eval.Evaluate(attValue)); }
310  }
311 
312  positionMap[name] = position*unit;
313 }
314 
315 void
316 G4GDMLReadDefine::RotationRead(const xercesc::DOMElement* const rotationElement)
317 {
318  G4String name = "";
319  G4double unit = 1.0;
320  G4ThreeVector rotation(0.,0.,0.);
321 
322  const xercesc::DOMNamedNodeMap* const attributes
323  = rotationElement->getAttributes();
324  XMLSize_t attributeCount = attributes->getLength();
325 
326  for (XMLSize_t attribute_index=0;
327  attribute_index<attributeCount; attribute_index++)
328  {
329  xercesc::DOMNode* node = attributes->item(attribute_index);
330 
331  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
332 
333  const xercesc::DOMAttr* const attribute
334  = dynamic_cast<xercesc::DOMAttr*>(node);
335  if (!attribute)
336  {
337  G4Exception("G4GDMLRead::RotationRead()", "InvalidRead",
338  FatalException, "No attribute found!");
339  return;
340  }
341  const G4String attName = Transcode(attribute->getName());
342  const G4String attValue = Transcode(attribute->getValue());
343 
344  if (attName=="name") { name = GenerateName(attValue); } else
345  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
346  if (G4UnitDefinition::GetCategory(attValue)!="Angle") {
347  G4Exception("G4GDMLReadDefine::RotationRead()", "InvalidRead",
348  FatalException, "Invalid unit for angle!"); }
349  } else
350  if (attName=="x") { rotation.setX(eval.Evaluate(attValue)); } else
351  if (attName=="y") { rotation.setY(eval.Evaluate(attValue)); } else
352  if (attName=="z") { rotation.setZ(eval.Evaluate(attValue)); }
353  }
354 
355  rotationMap[name] = rotation*unit;
356 }
357 
358 void G4GDMLReadDefine::ScaleRead(const xercesc::DOMElement* const scaleElement)
359 {
360  G4String name = "";
361  G4ThreeVector scale(1.0,1.0,1.0);
362 
363  const xercesc::DOMNamedNodeMap* const attributes
364  = scaleElement->getAttributes();
365  XMLSize_t attributeCount = attributes->getLength();
366 
367  for (XMLSize_t attribute_index=0;
368  attribute_index<attributeCount; attribute_index++)
369  {
370  xercesc::DOMNode* node = attributes->item(attribute_index);
371 
372  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
373 
374  const xercesc::DOMAttr* const attribute
375  = dynamic_cast<xercesc::DOMAttr*>(node);
376  if (!attribute)
377  {
378  G4Exception("G4GDMLRead::ScaleRead()", "InvalidRead",
379  FatalException, "No attribute found!");
380  return;
381  }
382  const G4String attName = Transcode(attribute->getName());
383  const G4String attValue = Transcode(attribute->getValue());
384 
385  if (attName=="name") { name = GenerateName(attValue); } else
386  if (attName=="x") { scale.setX(eval.Evaluate(attValue)); } else
387  if (attName=="y") { scale.setY(eval.Evaluate(attValue)); } else
388  if (attName=="z") { scale.setZ(eval.Evaluate(attValue)); }
389  }
390 
391  scaleMap[name] = scale;
392 }
393 
394 void
395 G4GDMLReadDefine::VariableRead(const xercesc::DOMElement* const variableElement)
396 {
397  G4String name = "";
398  G4double value = 0.0;
399 
400  const xercesc::DOMNamedNodeMap* const attributes
401  = variableElement->getAttributes();
402  XMLSize_t attributeCount = attributes->getLength();
403 
404  for (XMLSize_t attribute_index=0;
405  attribute_index<attributeCount; attribute_index++)
406  {
407  xercesc::DOMNode* node = attributes->item(attribute_index);
408 
409  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
410 
411  const xercesc::DOMAttr* const attribute
412  = dynamic_cast<xercesc::DOMAttr*>(node);
413  if (!attribute)
414  {
415  G4Exception("G4GDMLRead::VariableRead()", "InvalidRead",
416  FatalException, "No attribute found!");
417  return;
418  }
419  const G4String attName = Transcode(attribute->getName());
420  const G4String attValue = Transcode(attribute->getValue());
421 
422  if (attName=="name") { name = attValue; } else
423  if (attName=="value") { value = eval.Evaluate(attValue); }
424  }
425 
426  eval.DefineVariable(name,value);
427 }
428 
429 void G4GDMLReadDefine::QuantityRead(const xercesc::DOMElement* const element)
430 {
431  G4String name = "";
432  G4double unit = 1.0;
433  G4double value = 0.0;
434 
435  const xercesc::DOMNamedNodeMap* const attributes
436  = element->getAttributes();
437  XMLSize_t attributeCount = attributes->getLength();
438 
439  for (XMLSize_t attribute_index=0;
440  attribute_index<attributeCount; attribute_index++)
441  {
442  xercesc::DOMNode* node = attributes->item(attribute_index);
443 
444  if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
445 
446  const xercesc::DOMAttr* const attribute
447  = dynamic_cast<xercesc::DOMAttr*>(node);
448  if (!attribute)
449  {
450  G4Exception("G4GDMLRead::QuantityRead()", "InvalidRead",
451  FatalException, "No attribute found!");
452  return;
453  }
454  const G4String attName = Transcode(attribute->getName());
455  const G4String attValue = Transcode(attribute->getValue());
456 
457  if (attName=="name") { name = attValue; } else
458  if (attName=="value") { value = eval.Evaluate(attValue); } else
459  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
460  }
461 
462  quantityMap[name] = value*unit;
463  eval.DefineConstant(name,value*unit);
464 }
465 
466 void
467 G4GDMLReadDefine::DefineRead(const xercesc::DOMElement* const defineElement)
468 {
469  G4cout << "G4GDML: Reading definitions..." << G4endl;
470 
471  for (xercesc::DOMNode* iter = defineElement->getFirstChild();
472  iter != 0;iter = iter->getNextSibling())
473  {
474  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
475 
476  const xercesc::DOMElement* const child
477  = dynamic_cast<xercesc::DOMElement*>(iter);
478  if (!child)
479  {
480  G4Exception("G4GDMLRead::DefineRead()", "InvalidRead",
481  FatalException, "No child found!");
482  return;
483  }
484  const G4String tag = Transcode(child->getTagName());
485 
486  if (tag=="constant") { ConstantRead(child); } else
487  if (tag=="matrix") { MatrixRead(child); } else
488  if (tag=="position") { PositionRead(child); } else
489  if (tag=="rotation") { RotationRead(child); } else
490  if (tag=="scale") { ScaleRead(child); } else
491  if (tag=="variable") { VariableRead(child); } else
492  if (tag=="quantity") { QuantityRead(child); } else
493  if (tag=="expression") { ExpressionRead(child); }
494  else
495  {
496  G4String error_msg = "Unknown tag in define: "+tag;
497  G4Exception("G4GDMLReadDefine::defineRead()", "ReadError",
498  FatalException, error_msg);
499  }
500  }
501 }
502 
503 void
504 G4GDMLReadDefine::VectorRead(const xercesc::DOMElement* const vectorElement,
505  G4ThreeVector& vec)
506 {
507  G4double unit = 1.0;
508 
509  const xercesc::DOMNamedNodeMap* const attributes
510  = vectorElement->getAttributes();
511  XMLSize_t attributeCount = attributes->getLength();
512 
513  for (XMLSize_t attribute_index=0;
514  attribute_index<attributeCount; attribute_index++)
515  {
516  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
517 
518  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
519  { continue; }
520 
521  const xercesc::DOMAttr* const attribute
522  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
523  if (!attribute)
524  {
525  G4Exception("G4GDMLRead::VectorRead()", "InvalidRead",
526  FatalException, "No attribute found!");
527  return;
528  }
529  const G4String attName = Transcode(attribute->getName());
530  const G4String attValue = Transcode(attribute->getValue());
531 
532  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); } else
533  if (attName=="x") { vec.setX(eval.Evaluate(attValue)); } else
534  if (attName=="y") { vec.setY(eval.Evaluate(attValue)); } else
535  if (attName=="z") { vec.setZ(eval.Evaluate(attValue)); }
536  }
537 
538  vec *= unit;
539 }
540 
541 G4String G4GDMLReadDefine::RefRead(const xercesc::DOMElement* const element)
542 {
543  G4String ref;
544 
545  const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
546  XMLSize_t attributeCount = attributes->getLength();
547 
548  for (XMLSize_t attribute_index=0;
549  attribute_index<attributeCount; attribute_index++)
550  {
551  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
552 
553  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
554  { continue; }
555 
556  const xercesc::DOMAttr* const attribute
557  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
558  if (!attribute)
559  {
560  G4Exception("G4GDMLRead::Read()", "InvalidRead",
561  FatalException, "No attribute found!");
562  return ref;
563  }
564  const G4String attName = Transcode(attribute->getName());
565  const G4String attValue = Transcode(attribute->getValue());
566 
567  if (attName=="ref") { ref = attValue; }
568  }
569 
570  return ref;
571 }
572 
574 {
575  return eval.IsVariable(ref);
576 }
577 
579 {
580  return eval.GetConstant(ref);
581 }
582 
584 {
585  return eval.GetVariable(ref);
586 }
587 
589 {
590  if (quantityMap.find(ref) == quantityMap.end())
591  {
592  G4String error_msg = "Quantity '"+ref+"' was not found!";
593  G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError",
594  FatalException, error_msg);
595  }
596  return quantityMap[ref];
597 }
598 
600 {
601  if (positionMap.find(ref) == positionMap.end())
602  {
603  G4String error_msg = "Position '"+ref+"' was not found!";
604  G4Exception("G4GDMLReadDefine::getPosition()", "ReadError",
605  FatalException, error_msg);
606  }
607  return positionMap[ref];
608 }
609 
611 {
612  if (rotationMap.find(ref) == rotationMap.end())
613  {
614  G4String error_msg = "Rotation '"+ref+"' was not found!";
615  G4Exception("G4GDMLReadDefine::getRotation()", "ReadError",
616  FatalException, error_msg);
617  }
618  return rotationMap[ref];
619 }
620 
622 {
623  if (scaleMap.find(ref) == scaleMap.end())
624  {
625  G4String error_msg = "Scale '"+ref+"' was not found!";
626  G4Exception("G4GDMLReadDefine::getScale()", "ReadError",
627  FatalException, error_msg);
628  }
629  return scaleMap[ref];
630 }
631 
633 {
634  if (matrixMap.find(ref) == matrixMap.end())
635  {
636  G4String error_msg = "Matrix '"+ref+"' was not found!";
637  G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError",
638  FatalException, error_msg);
639  }
640  return matrixMap[ref];
641 }
const XML_Char * name
Definition: expat.h:151
G4int EvaluateInteger(const G4String &)
G4double GetQuantity(const G4String &)
std::map< G4String, G4ThreeVector > rotationMap
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
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)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
virtual void DefineRead(const xercesc::DOMElement *const)
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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)
static G4double GetValueOf(const G4String &)
#define position
Definition: xmlparse.cc:622
G4GLOB_DLL std::ostream G4cout
static constexpr double m
Definition: G4SIunits.hh:129
void MatrixRead(const xercesc::DOMElement *const)
const XML_Char int const XML_Char * value
Definition: expat.h:331
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
static G4String GetCategory(const G4String &)
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)
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
#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
tuple c
Definition: test.py:13
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 &)