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

Last change on this file since 1342 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 17.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.25 2010/02/17 18:06:25 gcosmo Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-01 $
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 ((rows<=0) || (cols<=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      const G4String attName = Transcode(attribute->getName());
164      const G4String attValue = Transcode(attribute->getValue());
165
166      if (attName=="name")  { name = attValue; }  else
167      if (attName=="value") { value = eval.Evaluate(attValue); }
168   }
169
170   eval.DefineConstant(name,value);
171}
172
173void
174G4GDMLReadDefine::ExpressionRead(const xercesc::DOMElement* const expElement)
175{
176   G4String name  = "";
177   G4double value = 0.0;
178
179   const xercesc::DOMNamedNodeMap* const attributes
180         = expElement->getAttributes();
181   XMLSize_t attributeCount = attributes->getLength();
182
183   for (XMLSize_t attribute_index=0;
184        attribute_index<attributeCount; attribute_index++)
185   {
186      xercesc::DOMNode* node = attributes->item(attribute_index);
187
188      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
189
190      const xercesc::DOMAttr* const attribute
191            = dynamic_cast<xercesc::DOMAttr*>(node);   
192      const G4String attName = Transcode(attribute->getName());
193      const G4String attValue = Transcode(attribute->getValue());
194
195      if (attName=="name")  { name = attValue; }
196   }
197
198   const G4String expValue = Transcode(expElement->getTextContent());
199   value = eval.Evaluate(expValue);
200   eval.DefineConstant(name,value);
201}
202
203void
204G4GDMLReadDefine::MatrixRead(const xercesc::DOMElement* const matrixElement) 
205{
206   G4String name = "";
207   G4int coldim  = 0;
208   G4String values = "";
209
210   const xercesc::DOMNamedNodeMap* const attributes
211         = matrixElement->getAttributes();
212   XMLSize_t attributeCount = attributes->getLength();
213
214   for (XMLSize_t attribute_index=0;
215        attribute_index<attributeCount; attribute_index++)
216   {
217      xercesc::DOMNode* node = attributes->item(attribute_index);
218
219      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
220
221      const xercesc::DOMAttr* const attribute
222            = dynamic_cast<xercesc::DOMAttr*>(node);   
223      const G4String attName = Transcode(attribute->getName());
224      const G4String attValue = Transcode(attribute->getValue());
225
226      if (attName=="name")   { name  = GenerateName(attValue); } else
227      if (attName=="coldim") { coldim = eval.EvaluateInteger(attValue); } else
228      if (attName=="values") { values = attValue; }
229   }
230
231   std::stringstream MatrixValueStream(values);
232   std::vector<G4double> valueList;
233
234   while (!MatrixValueStream.eof())
235   {
236      G4String MatrixValue;
237      MatrixValueStream >> MatrixValue;
238      valueList.push_back(eval.Evaluate(MatrixValue));
239   }
240
241   eval.DefineMatrix(name,coldim,valueList);
242
243   G4GDMLMatrix matrix(valueList.size()/coldim,coldim);
244
245   for (size_t i=0;i<valueList.size();i++)
246   {
247      matrix.Set(i/coldim,i%coldim,valueList[i]);
248   }
249
250   matrixMap[name] = matrix;
251}
252
253void
254G4GDMLReadDefine::PositionRead(const xercesc::DOMElement* const positionElement)
255{
256   G4String name = "";
257   G4double unit = 1.0;
258   G4ThreeVector position(0.,0.,0.);
259
260   const xercesc::DOMNamedNodeMap* const attributes
261         = positionElement->getAttributes();
262   XMLSize_t attributeCount = attributes->getLength();
263
264   for (XMLSize_t attribute_index=0;
265        attribute_index<attributeCount; attribute_index++)
266   {
267      xercesc::DOMNode* node = attributes->item(attribute_index);
268
269      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
270
271      const xercesc::DOMAttr* const attribute
272            = dynamic_cast<xercesc::DOMAttr*>(node);   
273      const G4String attName = Transcode(attribute->getName());
274      const G4String attValue = Transcode(attribute->getValue());
275
276      if (attName=="name") { name = GenerateName(attValue); }  else
277      if (attName=="unit") { unit = eval.Evaluate(attValue); } else
278      if (attName=="x") { position.setX(eval.Evaluate(attValue)); } else
279      if (attName=="y") { position.setY(eval.Evaluate(attValue)); } else
280      if (attName=="z") { position.setZ(eval.Evaluate(attValue)); }
281   }
282
283   positionMap[name] = position*unit;
284}
285
286void
287G4GDMLReadDefine::RotationRead(const xercesc::DOMElement* const rotationElement)
288{
289   G4String name = "";
290   G4double unit = 1.0;
291   G4ThreeVector rotation(0.,0.,0.);
292
293   const xercesc::DOMNamedNodeMap* const attributes
294         = rotationElement->getAttributes();
295   XMLSize_t attributeCount = attributes->getLength();
296
297   for (XMLSize_t attribute_index=0;
298        attribute_index<attributeCount; attribute_index++)
299   {
300      xercesc::DOMNode* node = attributes->item(attribute_index);
301
302      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
303
304      const xercesc::DOMAttr* const attribute
305            = dynamic_cast<xercesc::DOMAttr*>(node);   
306      const G4String attName = Transcode(attribute->getName());
307      const G4String attValue = Transcode(attribute->getValue());
308
309      if (attName=="name") { name = GenerateName(attValue); }  else
310      if (attName=="unit") { unit = eval.Evaluate(attValue); } else
311      if (attName=="x") { rotation.setX(eval.Evaluate(attValue)); } else
312      if (attName=="y") { rotation.setY(eval.Evaluate(attValue)); } else
313      if (attName=="z") { rotation.setZ(eval.Evaluate(attValue)); }
314   }
315
316   rotationMap[name] = rotation*unit;
317}
318
319void G4GDMLReadDefine::ScaleRead(const xercesc::DOMElement* const scaleElement)
320{
321   G4String name = "";
322   G4ThreeVector scale(1.0,1.0,1.0);
323
324   const xercesc::DOMNamedNodeMap* const attributes
325         = scaleElement->getAttributes();
326   XMLSize_t attributeCount = attributes->getLength();
327
328   for (XMLSize_t attribute_index=0;
329        attribute_index<attributeCount; attribute_index++)
330   {
331      xercesc::DOMNode* node = attributes->item(attribute_index);
332
333      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
334
335      const xercesc::DOMAttr* const attribute
336            = dynamic_cast<xercesc::DOMAttr*>(node);   
337      const G4String attName = Transcode(attribute->getName());
338      const G4String attValue = Transcode(attribute->getValue());
339
340      if (attName=="name") { name = GenerateName(attValue); }    else
341      if (attName=="x") { scale.setX(eval.Evaluate(attValue)); } else
342      if (attName=="y") { scale.setY(eval.Evaluate(attValue)); } else
343      if (attName=="z") { scale.setZ(eval.Evaluate(attValue)); }
344   }
345
346   scaleMap[name] = scale;
347}
348
349void
350G4GDMLReadDefine::VariableRead(const xercesc::DOMElement* const variableElement)
351{
352   G4String name  = "";
353   G4double value = 0.0;
354
355   const xercesc::DOMNamedNodeMap* const attributes
356         = variableElement->getAttributes();
357   XMLSize_t attributeCount = attributes->getLength();
358
359   for (XMLSize_t attribute_index=0;
360        attribute_index<attributeCount; attribute_index++)
361   {
362      xercesc::DOMNode* node = attributes->item(attribute_index);
363
364      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
365
366      const xercesc::DOMAttr* const attribute
367            = dynamic_cast<xercesc::DOMAttr*>(node);   
368      const G4String attName = Transcode(attribute->getName());
369      const G4String attValue = Transcode(attribute->getValue());
370
371      if (attName=="name")  { name = attValue; } else
372      if (attName=="value") { value = eval.Evaluate(attValue); }
373   }
374
375   eval.DefineVariable(name,value);
376}
377
378void G4GDMLReadDefine::QuantityRead(const xercesc::DOMElement* const element)
379{
380   G4String name = "";
381   G4double unit = 1.0;
382   G4double value = 0.0;
383
384   const xercesc::DOMNamedNodeMap* const attributes
385         = element->getAttributes();
386   XMLSize_t attributeCount = attributes->getLength();
387
388   for (XMLSize_t attribute_index=0;
389        attribute_index<attributeCount; attribute_index++)
390   {
391      xercesc::DOMNode* node = attributes->item(attribute_index);
392
393      if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
394
395      const xercesc::DOMAttr* const attribute
396            = dynamic_cast<xercesc::DOMAttr*>(node);   
397      const G4String attName = Transcode(attribute->getName());
398      const G4String attValue = Transcode(attribute->getValue());
399
400      if (attName=="name") { name = attValue; } else
401      if (attName=="value") { value = eval.Evaluate(attValue); } else
402      if (attName=="unit") { unit = eval.Evaluate(attValue); }
403   }
404
405   quantityMap[name] = value*unit;
406   eval.DefineConstant(name,value*unit);
407}
408
409void
410G4GDMLReadDefine::DefineRead(const xercesc::DOMElement* const defineElement)
411{
412   G4cout << "G4GDML: Reading definitions..." << G4endl;
413
414   for (xercesc::DOMNode* iter = defineElement->getFirstChild();
415        iter != 0;iter = iter->getNextSibling())
416   {
417      if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
418
419      const xercesc::DOMElement* const child
420            = dynamic_cast<xercesc::DOMElement*>(iter);
421      const G4String tag = Transcode(child->getTagName());
422
423      if (tag=="constant") { ConstantRead(child); } else
424      if (tag=="matrix")   { MatrixRead(child); }   else
425      if (tag=="position") { PositionRead(child); } else
426      if (tag=="rotation") { RotationRead(child); } else
427      if (tag=="scale")    { ScaleRead(child); }    else
428      if (tag=="variable") { VariableRead(child); } else
429      if (tag=="quantity") { QuantityRead(child); } else
430      if (tag=="expression") { ExpressionRead(child); }
431      else
432      {
433        G4String error_msg = "Unknown tag in define: "+tag;
434        G4Exception("G4GDMLReadDefine::defineRead()", "ReadError",
435                    FatalException, error_msg);
436      }
437   }
438}
439
440void
441G4GDMLReadDefine::VectorRead(const xercesc::DOMElement* const vectorElement,
442                             G4ThreeVector& vec)
443{
444   G4double unit = 1.0;
445
446   const xercesc::DOMNamedNodeMap* const attributes
447         = vectorElement->getAttributes();
448   XMLSize_t attributeCount = attributes->getLength();
449
450   for (XMLSize_t attribute_index=0;
451        attribute_index<attributeCount; attribute_index++)
452   {
453      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
454
455      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
456      { continue; }
457
458      const xercesc::DOMAttr* const attribute
459            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
460      const G4String attName = Transcode(attribute->getName());
461      const G4String attValue = Transcode(attribute->getValue());
462
463      if (attName=="unit") { unit = eval.Evaluate(attValue); } else
464      if (attName=="x") { vec.setX(eval.Evaluate(attValue)); } else
465      if (attName=="y") { vec.setY(eval.Evaluate(attValue)); } else
466      if (attName=="z") { vec.setZ(eval.Evaluate(attValue)); }
467   }
468
469   vec *= unit;
470}
471
472G4String G4GDMLReadDefine::RefRead(const xercesc::DOMElement* const element)
473{
474   G4String ref;
475
476   const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
477   XMLSize_t attributeCount = attributes->getLength();
478
479   for (XMLSize_t attribute_index=0;
480        attribute_index<attributeCount; attribute_index++)
481   {
482      xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
483
484      if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
485      { continue; }
486
487      const xercesc::DOMAttr* const attribute
488            = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
489      const G4String attName = Transcode(attribute->getName());
490      const G4String attValue = Transcode(attribute->getValue());
491
492      if (attName=="ref") { ref = attValue; }
493   }
494
495   return ref;
496}
497
498G4bool G4GDMLReadDefine::IsValidID(const G4String& ref) const
499{
500   return eval.IsVariable(ref);
501}
502
503G4double G4GDMLReadDefine::GetConstant(const G4String& ref)
504{
505   return eval.GetConstant(ref);
506}
507
508G4double G4GDMLReadDefine::GetVariable(const G4String& ref)
509{
510   return eval.GetVariable(ref);
511}
512
513G4double G4GDMLReadDefine::GetQuantity(const G4String& ref)
514{
515   if (quantityMap.find(ref) == quantityMap.end())
516   {
517     G4String error_msg = "Quantity '"+ref+"' was not found!";
518     G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError",
519                 FatalException, error_msg);
520   }
521   return quantityMap[ref];
522}
523
524G4ThreeVector G4GDMLReadDefine::GetPosition(const G4String& ref)
525{
526   if (positionMap.find(ref) == positionMap.end())
527   {
528     G4String error_msg = "Position '"+ref+"' was not found!";
529     G4Exception("G4GDMLReadDefine::getPosition()", "ReadError",
530                 FatalException, error_msg);
531   }
532   return positionMap[ref];
533}
534
535G4ThreeVector G4GDMLReadDefine::GetRotation(const G4String& ref)
536{
537   if (rotationMap.find(ref) == rotationMap.end())
538   {
539     G4String error_msg = "Rotation '"+ref+"' was not found!";
540     G4Exception("G4GDMLReadDefine::getRotation()", "ReadError",
541                 FatalException, error_msg);
542   }
543   return rotationMap[ref];
544}
545
546G4ThreeVector G4GDMLReadDefine::GetScale(const G4String& ref)
547{
548   if (scaleMap.find(ref) == scaleMap.end())
549   {
550     G4String error_msg = "Scale '"+ref+"' was not found!";
551     G4Exception("G4GDMLReadDefine::getScale()", "ReadError",
552                 FatalException, error_msg);
553   }
554   return scaleMap[ref];
555}
556
557G4GDMLMatrix G4GDMLReadDefine::GetMatrix(const G4String& ref)
558{
559   if (matrixMap.find(ref) == matrixMap.end())
560   {
561     G4String error_msg = "Matrix '"+ref+"' was not found!";
562     G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError",
563                 FatalException, error_msg);
564   }
565   return matrixMap[ref];
566}
Note: See TracBrowser for help on using the repository browser.