G4GDMLReadDefine.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 // class G4GDMLReadDefine Implementation
00029 //
00030 // Original author: Zoltan Torzsok, November 2007
00031 //
00032 // --------------------------------------------------------------------
00033 
00034 #include "G4GDMLReadDefine.hh"
00035 
00036 G4GDMLMatrix::G4GDMLMatrix()
00037   : m(0), rows(0), cols(0)
00038 {
00039 }
00040 
00041 G4GDMLMatrix::G4GDMLMatrix(size_t rows0, size_t cols0)
00042 {   
00043    if ((rows0==0) || (cols0==0))
00044    {
00045      G4Exception("G4GDMLMatrix::G4GDMLMatrix(r,c)", "InvalidSetup",
00046                  FatalException, "Zero indeces as arguments!?");
00047    }
00048    rows = rows0;
00049    cols = cols0;
00050    m = new G4double[rows*cols];
00051 }
00052 
00053 G4GDMLMatrix::G4GDMLMatrix(const G4GDMLMatrix& rhs)
00054   : m(0), rows(0), cols(0)
00055 {
00056    if (rhs.m)
00057    {
00058      rows = rhs.rows;
00059      cols = rhs.cols;
00060      m = new G4double[rows*cols];
00061      for (size_t i=0; i<rows*cols; i++)  { m[i] = rhs.m[i]; }
00062    }
00063 }
00064 
00065 G4GDMLMatrix& G4GDMLMatrix::operator=(const G4GDMLMatrix& rhs)
00066 {
00067    // Check assignment to self
00068    //
00069    if (this == &rhs)  { return *this; }
00070 
00071    // Copy data
00072    //
00073    rows = rhs.rows;
00074    cols = rhs.cols;
00075    if (rhs.m)
00076    {
00077      m = new G4double[rows*cols];
00078      for (size_t i=0; i<rows*cols; i++)  { m[i] = rhs.m[i]; }
00079    }
00080    else
00081    {
00082      m = 0;
00083    }
00084 
00085    return *this;
00086 }
00087 
00088 G4GDMLMatrix::~G4GDMLMatrix()
00089 {
00090    delete [] m;
00091 }
00092 
00093 void G4GDMLMatrix::Set(size_t r,size_t c,G4double a)
00094 {   
00095    if (r>=rows || c>=cols)
00096    {
00097      G4Exception("G4GDMLMatrix::set()", "InvalidSetup",
00098                  FatalException, "Index out of range!");
00099    }
00100    m[cols*r+c] = a;
00101 }
00102 
00103 G4double G4GDMLMatrix::Get(size_t r,size_t c) const
00104 {   
00105    if (r>=rows || c>=cols)
00106    {
00107      G4Exception("G4GDMLMatrix::get()", "InvalidSetup",
00108                  FatalException, "Index out of range!");
00109    }
00110    return m[cols*r+c];
00111 }
00112 
00113 size_t G4GDMLMatrix::GetRows() const
00114 {
00115    return rows;
00116 }
00117 
00118 size_t G4GDMLMatrix::GetCols() const
00119 {
00120    return cols;
00121 }
00122 
00123 G4GDMLReadDefine::G4GDMLReadDefine() : G4GDMLRead()
00124 {
00125 }
00126 
00127 G4GDMLReadDefine::~G4GDMLReadDefine()
00128 {
00129 }
00130 
00131 G4RotationMatrix
00132 G4GDMLReadDefine::GetRotationMatrix(const G4ThreeVector& angles)
00133 {
00134    G4RotationMatrix rot;
00135 
00136    rot.rotateX(angles.x());
00137    rot.rotateY(angles.y());
00138    rot.rotateZ(angles.z());
00139 
00140    return rot;
00141 }
00142 
00143 void
00144 G4GDMLReadDefine::ConstantRead(const xercesc::DOMElement* const constantElement)
00145 {
00146    G4String name  = "";
00147    G4double value = 0.0;
00148 
00149    const xercesc::DOMNamedNodeMap* const attributes
00150          = constantElement->getAttributes();
00151    XMLSize_t attributeCount = attributes->getLength();
00152 
00153    for (XMLSize_t attribute_index=0;
00154         attribute_index<attributeCount; attribute_index++)
00155    {
00156       xercesc::DOMNode* node = attributes->item(attribute_index);
00157 
00158       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00159 
00160       const xercesc::DOMAttr* const attribute
00161             = dynamic_cast<xercesc::DOMAttr*>(node);   
00162       if (!attribute)
00163       {
00164         G4Exception("G4GDMLRead::ConstantRead()", "InvalidRead",
00165                     FatalException, "No attribute found!");
00166         return;
00167       }
00168       const G4String attName = Transcode(attribute->getName());
00169       const G4String attValue = Transcode(attribute->getValue());
00170 
00171       if (attName=="name")  { name = attValue; }  else
00172       if (attName=="value") { value = eval.Evaluate(attValue); }
00173    }
00174 
00175    eval.DefineConstant(name,value);
00176 }
00177 
00178 void
00179 G4GDMLReadDefine::ExpressionRead(const xercesc::DOMElement* const expElement)
00180 {
00181    G4String name  = "";
00182    G4double value = 0.0;
00183 
00184    const xercesc::DOMNamedNodeMap* const attributes
00185          = expElement->getAttributes();
00186    XMLSize_t attributeCount = attributes->getLength();
00187 
00188    for (XMLSize_t attribute_index=0;
00189         attribute_index<attributeCount; attribute_index++)
00190    {
00191       xercesc::DOMNode* node = attributes->item(attribute_index);
00192 
00193       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00194 
00195       const xercesc::DOMAttr* const attribute
00196             = dynamic_cast<xercesc::DOMAttr*>(node);   
00197       if (!attribute)
00198       {
00199         G4Exception("G4GDMLRead::ExpressionRead()", "InvalidRead",
00200                     FatalException, "No attribute found!");
00201         return;
00202       }
00203       const G4String attName = Transcode(attribute->getName());
00204       const G4String attValue = Transcode(attribute->getValue());
00205 
00206       if (attName=="name")  { name = attValue; }
00207    }
00208 
00209    const G4String expValue = Transcode(expElement->getTextContent());
00210    value = eval.Evaluate(expValue);
00211    eval.DefineConstant(name,value);
00212 }
00213 
00214 void
00215 G4GDMLReadDefine::MatrixRead(const xercesc::DOMElement* const matrixElement) 
00216 {
00217    G4String name = "";
00218    G4int coldim  = 0;
00219    G4String values = "";
00220 
00221    const xercesc::DOMNamedNodeMap* const attributes
00222          = matrixElement->getAttributes();
00223    XMLSize_t attributeCount = attributes->getLength();
00224 
00225    for (XMLSize_t attribute_index=0;
00226         attribute_index<attributeCount; attribute_index++)
00227    {
00228       xercesc::DOMNode* node = attributes->item(attribute_index);
00229 
00230       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00231 
00232       const xercesc::DOMAttr* const attribute
00233             = dynamic_cast<xercesc::DOMAttr*>(node);   
00234       if (!attribute)
00235       {
00236         G4Exception("G4GDMLRead::MatrixRead()", "InvalidRead",
00237                     FatalException, "No attribute found!");
00238         return;
00239       }
00240       const G4String attName = Transcode(attribute->getName());
00241       const G4String attValue = Transcode(attribute->getValue());
00242 
00243       if (attName=="name")   { name  = GenerateName(attValue); } else
00244       if (attName=="coldim") { coldim = eval.EvaluateInteger(attValue); } else
00245       if (attName=="values") { values = attValue; }
00246    }
00247 
00248    std::stringstream MatrixValueStream(values);
00249    std::vector<G4double> valueList;
00250 
00251    while (!MatrixValueStream.eof())
00252    {
00253       G4String MatrixValue;
00254       MatrixValueStream >> MatrixValue;
00255       valueList.push_back(eval.Evaluate(MatrixValue));
00256    }
00257 
00258    eval.DefineMatrix(name,coldim,valueList);
00259 
00260    G4GDMLMatrix matrix(valueList.size()/coldim,coldim);
00261 
00262    for (size_t i=0;i<valueList.size();i++)
00263    {
00264       matrix.Set(i/coldim,i%coldim,valueList[i]);
00265    }
00266 
00267    matrixMap[name] = matrix;
00268 }
00269 
00270 void
00271 G4GDMLReadDefine::PositionRead(const xercesc::DOMElement* const positionElement)
00272 {
00273    G4String name = "";
00274    G4double unit = 1.0;
00275    G4ThreeVector position(0.,0.,0.);
00276 
00277    const xercesc::DOMNamedNodeMap* const attributes
00278          = positionElement->getAttributes();
00279    XMLSize_t attributeCount = attributes->getLength();
00280 
00281    for (XMLSize_t attribute_index=0;
00282         attribute_index<attributeCount; attribute_index++)
00283    {
00284       xercesc::DOMNode* node = attributes->item(attribute_index);
00285 
00286       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00287 
00288       const xercesc::DOMAttr* const attribute
00289             = dynamic_cast<xercesc::DOMAttr*>(node);   
00290       if (!attribute)
00291       {
00292         G4Exception("G4GDMLRead::PositionRead()", "InvalidRead",
00293                     FatalException, "No attribute found!");
00294         return;
00295       }
00296       const G4String attName = Transcode(attribute->getName());
00297       const G4String attValue = Transcode(attribute->getValue());
00298 
00299       if (attName=="name") { name = GenerateName(attValue); }  else
00300       if (attName=="unit") { unit = eval.Evaluate(attValue); } else
00301       if (attName=="x") { position.setX(eval.Evaluate(attValue)); } else
00302       if (attName=="y") { position.setY(eval.Evaluate(attValue)); } else
00303       if (attName=="z") { position.setZ(eval.Evaluate(attValue)); }
00304    }
00305 
00306    positionMap[name] = position*unit;
00307 }
00308 
00309 void
00310 G4GDMLReadDefine::RotationRead(const xercesc::DOMElement* const rotationElement)
00311 {
00312    G4String name = "";
00313    G4double unit = 1.0;
00314    G4ThreeVector rotation(0.,0.,0.);
00315 
00316    const xercesc::DOMNamedNodeMap* const attributes
00317          = rotationElement->getAttributes();
00318    XMLSize_t attributeCount = attributes->getLength();
00319 
00320    for (XMLSize_t attribute_index=0;
00321         attribute_index<attributeCount; attribute_index++)
00322    {
00323       xercesc::DOMNode* node = attributes->item(attribute_index);
00324 
00325       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00326 
00327       const xercesc::DOMAttr* const attribute
00328             = dynamic_cast<xercesc::DOMAttr*>(node);   
00329       if (!attribute)
00330       {
00331         G4Exception("G4GDMLRead::RotationRead()", "InvalidRead",
00332                     FatalException, "No attribute found!");
00333         return;
00334       }
00335       const G4String attName = Transcode(attribute->getName());
00336       const G4String attValue = Transcode(attribute->getValue());
00337 
00338       if (attName=="name") { name = GenerateName(attValue); }  else
00339       if (attName=="unit") { unit = eval.Evaluate(attValue); } else
00340       if (attName=="x") { rotation.setX(eval.Evaluate(attValue)); } else
00341       if (attName=="y") { rotation.setY(eval.Evaluate(attValue)); } else
00342       if (attName=="z") { rotation.setZ(eval.Evaluate(attValue)); }
00343    }
00344 
00345    rotationMap[name] = rotation*unit;
00346 }
00347 
00348 void G4GDMLReadDefine::ScaleRead(const xercesc::DOMElement* const scaleElement)
00349 {
00350    G4String name = "";
00351    G4ThreeVector scale(1.0,1.0,1.0);
00352 
00353    const xercesc::DOMNamedNodeMap* const attributes
00354          = scaleElement->getAttributes();
00355    XMLSize_t attributeCount = attributes->getLength();
00356 
00357    for (XMLSize_t attribute_index=0;
00358         attribute_index<attributeCount; attribute_index++)
00359    {
00360       xercesc::DOMNode* node = attributes->item(attribute_index);
00361 
00362       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00363 
00364       const xercesc::DOMAttr* const attribute
00365             = dynamic_cast<xercesc::DOMAttr*>(node);   
00366       if (!attribute)
00367       {
00368         G4Exception("G4GDMLRead::ScaleRead()", "InvalidRead",
00369                     FatalException, "No attribute found!");
00370         return;
00371       }
00372       const G4String attName = Transcode(attribute->getName());
00373       const G4String attValue = Transcode(attribute->getValue());
00374 
00375       if (attName=="name") { name = GenerateName(attValue); }    else
00376       if (attName=="x") { scale.setX(eval.Evaluate(attValue)); } else
00377       if (attName=="y") { scale.setY(eval.Evaluate(attValue)); } else
00378       if (attName=="z") { scale.setZ(eval.Evaluate(attValue)); }
00379    }
00380 
00381    scaleMap[name] = scale;
00382 }
00383 
00384 void
00385 G4GDMLReadDefine::VariableRead(const xercesc::DOMElement* const variableElement)
00386 {
00387    G4String name  = "";
00388    G4double value = 0.0;
00389 
00390    const xercesc::DOMNamedNodeMap* const attributes
00391          = variableElement->getAttributes();
00392    XMLSize_t attributeCount = attributes->getLength();
00393 
00394    for (XMLSize_t attribute_index=0;
00395         attribute_index<attributeCount; attribute_index++)
00396    {
00397       xercesc::DOMNode* node = attributes->item(attribute_index);
00398 
00399       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00400 
00401       const xercesc::DOMAttr* const attribute
00402             = dynamic_cast<xercesc::DOMAttr*>(node);   
00403       if (!attribute)
00404       {
00405         G4Exception("G4GDMLRead::VariableRead()", "InvalidRead",
00406                     FatalException, "No attribute found!");
00407         return;
00408       }
00409       const G4String attName = Transcode(attribute->getName());
00410       const G4String attValue = Transcode(attribute->getValue());
00411 
00412       if (attName=="name")  { name = attValue; } else
00413       if (attName=="value") { value = eval.Evaluate(attValue); }
00414    }
00415 
00416    eval.DefineVariable(name,value);
00417 }
00418 
00419 void G4GDMLReadDefine::QuantityRead(const xercesc::DOMElement* const element)
00420 {
00421    G4String name = "";
00422    G4double unit = 1.0;
00423    G4double value = 0.0;
00424 
00425    const xercesc::DOMNamedNodeMap* const attributes
00426          = element->getAttributes();
00427    XMLSize_t attributeCount = attributes->getLength();
00428 
00429    for (XMLSize_t attribute_index=0;
00430         attribute_index<attributeCount; attribute_index++)
00431    {
00432       xercesc::DOMNode* node = attributes->item(attribute_index);
00433 
00434       if (node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) { continue; }
00435 
00436       const xercesc::DOMAttr* const attribute
00437             = dynamic_cast<xercesc::DOMAttr*>(node);   
00438       if (!attribute)
00439       {
00440         G4Exception("G4GDMLRead::QuantityRead()", "InvalidRead",
00441                     FatalException, "No attribute found!");
00442         return;
00443       }
00444       const G4String attName = Transcode(attribute->getName());
00445       const G4String attValue = Transcode(attribute->getValue());
00446 
00447       if (attName=="name") { name = attValue; } else
00448       if (attName=="value") { value = eval.Evaluate(attValue); } else
00449       if (attName=="unit") { unit = eval.Evaluate(attValue); }
00450    }
00451 
00452    quantityMap[name] = value*unit;
00453    eval.DefineConstant(name,value*unit);
00454 }
00455 
00456 void
00457 G4GDMLReadDefine::DefineRead(const xercesc::DOMElement* const defineElement)
00458 {
00459    G4cout << "G4GDML: Reading definitions..." << G4endl;
00460 
00461    for (xercesc::DOMNode* iter = defineElement->getFirstChild();
00462         iter != 0;iter = iter->getNextSibling())
00463    {
00464       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00465 
00466       const xercesc::DOMElement* const child
00467             = dynamic_cast<xercesc::DOMElement*>(iter);
00468       if (!child)
00469       {
00470         G4Exception("G4GDMLRead::DefineRead()", "InvalidRead",
00471                     FatalException, "No child found!");
00472         return;
00473       }
00474       const G4String tag = Transcode(child->getTagName());
00475 
00476       if (tag=="constant") { ConstantRead(child); } else
00477       if (tag=="matrix")   { MatrixRead(child); }   else
00478       if (tag=="position") { PositionRead(child); } else
00479       if (tag=="rotation") { RotationRead(child); } else
00480       if (tag=="scale")    { ScaleRead(child); }    else
00481       if (tag=="variable") { VariableRead(child); } else
00482       if (tag=="quantity") { QuantityRead(child); } else
00483       if (tag=="expression") { ExpressionRead(child); }
00484       else
00485       {
00486         G4String error_msg = "Unknown tag in define: "+tag;
00487         G4Exception("G4GDMLReadDefine::defineRead()", "ReadError",
00488                     FatalException, error_msg);
00489       }
00490    }
00491 }
00492 
00493 void
00494 G4GDMLReadDefine::VectorRead(const xercesc::DOMElement* const vectorElement,
00495                              G4ThreeVector& vec)
00496 {
00497    G4double unit = 1.0;
00498 
00499    const xercesc::DOMNamedNodeMap* const attributes
00500          = vectorElement->getAttributes();
00501    XMLSize_t attributeCount = attributes->getLength();
00502 
00503    for (XMLSize_t attribute_index=0;
00504         attribute_index<attributeCount; attribute_index++)
00505    {
00506       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00507 
00508       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00509       { continue; }
00510 
00511       const xercesc::DOMAttr* const attribute
00512             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00513       if (!attribute)
00514       {
00515         G4Exception("G4GDMLRead::VectorRead()", "InvalidRead",
00516                     FatalException, "No attribute found!");
00517         return;
00518       }
00519       const G4String attName = Transcode(attribute->getName());
00520       const G4String attValue = Transcode(attribute->getValue());
00521 
00522       if (attName=="unit") { unit = eval.Evaluate(attValue); } else
00523       if (attName=="x") { vec.setX(eval.Evaluate(attValue)); } else
00524       if (attName=="y") { vec.setY(eval.Evaluate(attValue)); } else
00525       if (attName=="z") { vec.setZ(eval.Evaluate(attValue)); }
00526    }
00527 
00528    vec *= unit;
00529 }
00530 
00531 G4String G4GDMLReadDefine::RefRead(const xercesc::DOMElement* const element)
00532 {
00533    G4String ref;
00534 
00535    const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
00536    XMLSize_t attributeCount = attributes->getLength();
00537 
00538    for (XMLSize_t attribute_index=0;
00539         attribute_index<attributeCount; attribute_index++)
00540    {
00541       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00542 
00543       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00544       { continue; }
00545 
00546       const xercesc::DOMAttr* const attribute
00547             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00548       if (!attribute)
00549       {
00550         G4Exception("G4GDMLRead::Read()", "InvalidRead",
00551                     FatalException, "No attribute found!");
00552         return ref;
00553       }
00554       const G4String attName = Transcode(attribute->getName());
00555       const G4String attValue = Transcode(attribute->getValue());
00556 
00557       if (attName=="ref") { ref = attValue; }
00558    }
00559 
00560    return ref;
00561 }
00562 
00563 G4bool G4GDMLReadDefine::IsValidID(const G4String& ref) const
00564 {
00565    return eval.IsVariable(ref);
00566 }
00567 
00568 G4double G4GDMLReadDefine::GetConstant(const G4String& ref)
00569 {
00570    return eval.GetConstant(ref);
00571 }
00572 
00573 G4double G4GDMLReadDefine::GetVariable(const G4String& ref)
00574 {
00575    return eval.GetVariable(ref);
00576 }
00577 
00578 G4double G4GDMLReadDefine::GetQuantity(const G4String& ref)
00579 {
00580    if (quantityMap.find(ref) == quantityMap.end())
00581    {
00582      G4String error_msg = "Quantity '"+ref+"' was not found!";
00583      G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError",
00584                  FatalException, error_msg);
00585    }
00586    return quantityMap[ref];
00587 }
00588 
00589 G4ThreeVector G4GDMLReadDefine::GetPosition(const G4String& ref)
00590 {
00591    if (positionMap.find(ref) == positionMap.end())
00592    {
00593      G4String error_msg = "Position '"+ref+"' was not found!";
00594      G4Exception("G4GDMLReadDefine::getPosition()", "ReadError",
00595                  FatalException, error_msg);
00596    }
00597    return positionMap[ref];
00598 }
00599 
00600 G4ThreeVector G4GDMLReadDefine::GetRotation(const G4String& ref)
00601 {
00602    if (rotationMap.find(ref) == rotationMap.end())
00603    {
00604      G4String error_msg = "Rotation '"+ref+"' was not found!";
00605      G4Exception("G4GDMLReadDefine::getRotation()", "ReadError",
00606                  FatalException, error_msg);
00607    }
00608    return rotationMap[ref];
00609 }
00610 
00611 G4ThreeVector G4GDMLReadDefine::GetScale(const G4String& ref)
00612 {
00613    if (scaleMap.find(ref) == scaleMap.end())
00614    {
00615      G4String error_msg = "Scale '"+ref+"' was not found!";
00616      G4Exception("G4GDMLReadDefine::getScale()", "ReadError",
00617                  FatalException, error_msg);
00618    }
00619    return scaleMap[ref];
00620 }
00621 
00622 G4GDMLMatrix G4GDMLReadDefine::GetMatrix(const G4String& ref)
00623 {
00624    if (matrixMap.find(ref) == matrixMap.end())
00625    {
00626      G4String error_msg = "Matrix '"+ref+"' was not found!";
00627      G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError",
00628                  FatalException, error_msg);
00629    }
00630    return matrixMap[ref];
00631 }

Generated on Mon May 27 17:48:19 2013 for Geant4 by  doxygen 1.4.7