#include <G4GDMLReadParamvol.hh>
Inheritance diagram for G4GDMLReadParamvol:
Public Member Functions | |
virtual void | ParamvolRead (const xercesc::DOMElement *const, G4LogicalVolume *) |
virtual void | Paramvol_contentRead (const xercesc::DOMElement *const) |
Protected Member Functions | |
G4GDMLReadParamvol () | |
virtual | ~G4GDMLReadParamvol () |
void | Box_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Trd_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Trap_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Tube_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Cone_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Sphere_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Orb_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Torus_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Para_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | Hype_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &) |
void | ParameterisedRead (const xercesc::DOMElement *const) |
void | ParametersRead (const xercesc::DOMElement *const) |
Protected Attributes | |
G4GDMLParameterisation * | parameterisation |
Definition at line 48 of file G4GDMLReadParamvol.hh.
G4GDMLReadParamvol::G4GDMLReadParamvol | ( | ) | [protected] |
Definition at line 41 of file G4GDMLReadParamvol.cc.
00042 : G4GDMLReadSetup(), parameterisation(0) 00043 { 00044 }
G4GDMLReadParamvol::~G4GDMLReadParamvol | ( | ) | [protected, virtual] |
void G4GDMLReadParamvol::Box_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 51 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00053 { 00054 G4double lunit = 1.0; 00055 00056 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00057 XMLSize_t attributeCount = attributes->getLength(); 00058 00059 for (XMLSize_t attribute_index=0; 00060 attribute_index<attributeCount; attribute_index++) 00061 { 00062 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00063 00064 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00065 { continue; } 00066 00067 const xercesc::DOMAttr* const attribute 00068 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00069 if (!attribute) 00070 { 00071 G4Exception("G4GDMLReadParamvol::Box_dimensionsRead()", 00072 "InvalidRead", FatalException, "No attribute found!"); 00073 return; 00074 } 00075 const G4String attName = Transcode(attribute->getName()); 00076 const G4String attValue = Transcode(attribute->getValue()); 00077 00078 if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else 00079 if (attName=="x") { parameter.dimension[0] = eval.Evaluate(attValue); } else 00080 if (attName=="y") { parameter.dimension[1] = eval.Evaluate(attValue); } else 00081 if (attName=="z") { parameter.dimension[2] = eval.Evaluate(attValue); } 00082 } 00083 00084 parameter.dimension[0] *= 0.5*lunit; 00085 parameter.dimension[1] *= 0.5*lunit; 00086 parameter.dimension[2] *= 0.5*lunit; 00087 }
void G4GDMLReadParamvol::Cone_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 255 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00257 { 00258 G4double lunit = 1.0; 00259 G4double aunit = 1.0; 00260 00261 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00262 XMLSize_t attributeCount = attributes->getLength(); 00263 00264 for (XMLSize_t attribute_index=0; 00265 attribute_index<attributeCount; attribute_index++) 00266 { 00267 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00268 00269 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00270 { continue; } 00271 00272 const xercesc::DOMAttr* const attribute 00273 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00274 if (!attribute) 00275 { 00276 G4Exception("G4GDMLReadParamvol::Cone_dimensionsRead()", 00277 "InvalidRead", FatalException, "No attribute found!"); 00278 return; 00279 } 00280 const G4String attName = Transcode(attribute->getName()); 00281 const G4String attValue = Transcode(attribute->getValue()); 00282 00283 if (attName=="lunit") 00284 { lunit = eval.Evaluate(attValue); } else 00285 if (attName=="aunit") 00286 { aunit = eval.Evaluate(attValue); } else 00287 if (attName=="rmin1") 00288 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00289 if (attName=="rmax1") 00290 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00291 if (attName=="rmin2") 00292 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00293 if (attName=="rmax2") 00294 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00295 if (attName=="z") 00296 { parameter.dimension[4] = eval.Evaluate(attValue); } else 00297 if (attName=="startphi") 00298 { parameter.dimension[5] = eval.Evaluate(attValue); } else 00299 if (attName=="deltaphi") 00300 { parameter.dimension[6] = eval.Evaluate(attValue); } 00301 } 00302 00303 parameter.dimension[0] *= lunit; 00304 parameter.dimension[1] *= lunit; 00305 parameter.dimension[2] *= lunit; 00306 parameter.dimension[3] *= lunit; 00307 parameter.dimension[4] *= 0.5*lunit; 00308 parameter.dimension[5] *= aunit; 00309 parameter.dimension[6] *= aunit; 00310 }
void G4GDMLReadParamvol::Hype_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 510 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00512 { 00513 G4double lunit = 1.0; 00514 G4double aunit = 1.0; 00515 00516 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00517 XMLSize_t attributeCount = attributes->getLength(); 00518 00519 for (XMLSize_t attribute_index=0; 00520 attribute_index<attributeCount; attribute_index++) 00521 { 00522 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00523 00524 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00525 { continue; } 00526 00527 const xercesc::DOMAttr* const attribute 00528 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00529 if (!attribute) 00530 { 00531 G4Exception("G4GDMLReadParamvol::Hype_dimensionsRead()", 00532 "InvalidRead", FatalException, "No attribute found!"); 00533 return; 00534 } 00535 const G4String attName = Transcode(attribute->getName()); 00536 const G4String attValue = Transcode(attribute->getValue()); 00537 00538 if (attName=="lunit") 00539 { lunit = eval.Evaluate(attValue); } else 00540 if (attName=="aunit") 00541 { aunit = eval.Evaluate(attValue); } else 00542 if (attName=="rmin") 00543 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00544 if (attName=="rmax") 00545 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00546 if (attName=="inst") 00547 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00548 if (attName=="outst") 00549 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00550 if (attName=="z") 00551 { parameter.dimension[4] = eval.Evaluate(attValue); } 00552 } 00553 00554 parameter.dimension[0] = lunit; 00555 parameter.dimension[1] = lunit; 00556 parameter.dimension[2] = aunit; 00557 parameter.dimension[3] = aunit; 00558 parameter.dimension[4] = 0.5*lunit; 00559 }
void G4GDMLReadParamvol::Orb_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 368 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
00370 { 00371 G4double lunit = 1.0; 00372 00373 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00374 XMLSize_t attributeCount = attributes->getLength(); 00375 00376 for (XMLSize_t attribute_index=0; 00377 attribute_index<attributeCount; attribute_index++) 00378 { 00379 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00380 00381 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00382 { continue; } 00383 00384 const xercesc::DOMAttr* const attribute 00385 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00386 if (!attribute) 00387 { 00388 G4Exception("G4GDMLReadParamvol::Orb_dimensionsRead()", 00389 "InvalidRead", FatalException, "No attribute found!"); 00390 return; 00391 } 00392 const G4String attName = Transcode(attribute->getName()); 00393 const G4String attValue = Transcode(attribute->getValue()); 00394 00395 if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else 00396 if (attName=="r") { parameter.dimension[0] = eval.Evaluate(attValue); } 00397 } 00398 00399 parameter.dimension[0] *= lunit; 00400 }
void G4GDMLReadParamvol::Para_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 455 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
00457 { 00458 G4double lunit = 1.0; 00459 G4double aunit = 1.0; 00460 00461 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00462 XMLSize_t attributeCount = attributes->getLength(); 00463 00464 for (XMLSize_t attribute_index=0; 00465 attribute_index<attributeCount; attribute_index++) 00466 { 00467 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00468 00469 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00470 { continue; } 00471 00472 const xercesc::DOMAttr* const attribute 00473 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00474 if (!attribute) 00475 { 00476 G4Exception("G4GDMLReadParamvol::Para_dimensionsRead()", 00477 "InvalidRead", FatalException, "No attribute found!"); 00478 return; 00479 } 00480 const G4String attName = Transcode(attribute->getName()); 00481 const G4String attValue = Transcode(attribute->getValue()); 00482 00483 if (attName=="lunit") 00484 { lunit = eval.Evaluate(attValue); } else 00485 if (attName=="aunit") 00486 { aunit = eval.Evaluate(attValue); } else 00487 if (attName=="x") 00488 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00489 if (attName=="y") 00490 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00491 if (attName=="z") 00492 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00493 if (attName=="alpha") 00494 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00495 if (attName=="theta") 00496 { parameter.dimension[4] = eval.Evaluate(attValue); } else 00497 if (attName=="phi") 00498 { parameter.dimension[5] = eval.Evaluate(attValue); } 00499 } 00500 00501 parameter.dimension[0] = 0.5*lunit; 00502 parameter.dimension[1] = 0.5*lunit; 00503 parameter.dimension[2] = 0.5*lunit; 00504 parameter.dimension[3] = aunit; 00505 parameter.dimension[4] = aunit; 00506 parameter.dimension[5] = aunit; 00507 }
void G4GDMLReadParamvol::ParameterisedRead | ( | const xercesc::DOMElement * | const | ) | [protected] |
Definition at line 619 of file G4GDMLReadParamvol.cc.
References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::LoopRead(), ParametersRead(), G4GDMLRead::Paramvol_contentRead(), and G4GDMLRead::Transcode().
Referenced by Paramvol_contentRead().
00620 { 00621 for (xercesc::DOMNode* iter = element->getFirstChild(); 00622 iter != 0; iter = iter->getNextSibling()) 00623 { 00624 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; } 00625 00626 const xercesc::DOMElement* const child 00627 = dynamic_cast<xercesc::DOMElement*>(iter); 00628 if (!child) 00629 { 00630 G4Exception("G4GDMLReadParamvol::ParameterisedRead()", 00631 "InvalidRead", FatalException, "No child found!"); 00632 return; 00633 } 00634 const G4String tag = Transcode(child->getTagName()); 00635 00636 if (tag=="parameters") 00637 { 00638 const xercesc::DOMNamedNodeMap* const attributes 00639 = element->getAttributes(); 00640 XMLSize_t attributeCount = attributes->getLength(); 00641 for (XMLSize_t attribute_index=0; 00642 attribute_index<attributeCount; attribute_index++) 00643 { 00644 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00645 00646 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00647 { continue; } 00648 00649 const xercesc::DOMAttr* const attribute 00650 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00651 if (!attribute) 00652 { 00653 G4Exception("G4GDMLReadParamvol::ParameterisedRead()", 00654 "InvalidRead", FatalException, "No attribute found!"); 00655 return; 00656 } 00657 const G4String attName = Transcode(attribute->getName()); 00658 const G4String attValue = Transcode(attribute->getValue()); 00659 00660 if (attName=="number") { eval.Evaluate(attValue); } 00661 } 00662 ParametersRead(child); 00663 } 00664 else 00665 { 00666 if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); } 00667 } 00668 } 00669 }
void G4GDMLReadParamvol::ParametersRead | ( | const xercesc::DOMElement * | const | ) | [protected] |
Definition at line 562 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::AddParameter(), Box_dimensionsRead(), Cone_dimensionsRead(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4GDMLReadDefine::GetPosition(), G4GDMLReadDefine::GetRotation(), Hype_dimensionsRead(), parameterisation, G4GDMLParameterisation::PARAMETER::position, position, G4GDMLParameterisation::PARAMETER::pRot, G4GDMLReadDefine::RefRead(), G4GDMLRead::Transcode(), Trap_dimensionsRead(), Trd_dimensionsRead(), Tube_dimensionsRead(), and G4GDMLReadDefine::VectorRead().
Referenced by ParameterisedRead().
00562 { 00563 00564 G4ThreeVector rotation(0.0,0.0,0.0); 00565 G4ThreeVector position(0.0,0.0,0.0); 00566 00567 G4GDMLParameterisation::PARAMETER parameter; 00568 00569 for (xercesc::DOMNode* iter = element->getFirstChild(); 00570 iter != 0; iter = iter->getNextSibling()) 00571 { 00572 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; } 00573 00574 const xercesc::DOMElement* const child 00575 = dynamic_cast<xercesc::DOMElement*>(iter); 00576 if (!child) 00577 { 00578 G4Exception("G4GDMLReadParamvol::ParametersRead()", 00579 "InvalidRead", FatalException, "No child found!"); 00580 return; 00581 } 00582 const G4String tag = Transcode(child->getTagName()); 00583 if (tag=="rotation") { VectorRead(child,rotation); } else 00584 if (tag=="position") { VectorRead(child,position); } else 00585 if (tag=="positionref") 00586 { position = GetPosition(GenerateName(RefRead(child))); } else 00587 if (tag=="rotationref") 00588 { rotation = GetRotation(GenerateName(RefRead(child))); } else 00589 if (tag=="box_dimensions") { Box_dimensionsRead(child,parameter); } else 00590 if (tag=="trd_dimensions") { Trd_dimensionsRead(child,parameter); } else 00591 if (tag=="trap_dimensions") { Trap_dimensionsRead(child,parameter); } else 00592 if (tag=="tube_dimensions") { Tube_dimensionsRead(child,parameter); } else 00593 if (tag=="cone_dimensions") { Cone_dimensionsRead(child,parameter); } else 00594 if (tag=="sphere_dimensions") { Cone_dimensionsRead(child,parameter); } else 00595 if (tag=="orb_dimensions") { Cone_dimensionsRead(child,parameter); } else 00596 if (tag=="torus_dimensions") { Cone_dimensionsRead(child,parameter); } else 00597 if (tag=="para_dimensions") { Cone_dimensionsRead(child,parameter); } else 00598 if (tag=="hype_dimensions") { Hype_dimensionsRead(child,parameter); } 00599 else 00600 { 00601 G4String error_msg = "Unknown tag in parameters: " + tag; 00602 G4Exception("G4GDMLReadParamvol::ParametersRead()", "ReadError", 00603 FatalException, error_msg); 00604 } 00605 } 00606 00607 parameter.pRot = new G4RotationMatrix(); 00608 00609 parameter.pRot->rotateX(rotation.x()); 00610 parameter.pRot->rotateY(rotation.y()); 00611 parameter.pRot->rotateZ(rotation.z()); 00612 00613 parameter.position = position; 00614 00615 parameterisation->AddParameter(parameter); 00616 }
void G4GDMLReadParamvol::Paramvol_contentRead | ( | const xercesc::DOMElement * | const | ) | [virtual] |
Implements G4GDMLRead.
Definition at line 672 of file G4GDMLReadParamvol.cc.
References FatalException, G4Exception(), G4GDMLRead::LoopRead(), ParameterisedRead(), G4GDMLRead::Paramvol_contentRead(), and G4GDMLRead::Transcode().
Referenced by ParamvolRead().
00673 { 00674 for (xercesc::DOMNode* iter = element->getFirstChild(); 00675 iter != 0; iter = iter->getNextSibling()) 00676 { 00677 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; } 00678 00679 const xercesc::DOMElement* const child 00680 = dynamic_cast<xercesc::DOMElement*>(iter); 00681 if (!child) 00682 { 00683 G4Exception("G4GDMLReadParamvol::Paramvol_contentRead()", "InvalidRead", 00684 FatalException, "No child found!"); 00685 return; 00686 } 00687 const G4String tag = Transcode(child->getTagName()); 00688 if (tag=="parameterised_position_size") { ParameterisedRead(child); }else 00689 if (tag=="loop") { LoopRead(child,&G4GDMLRead::Paramvol_contentRead); } 00690 } 00691 }
void G4GDMLReadParamvol::ParamvolRead | ( | const xercesc::DOMElement * | const, | |
G4LogicalVolume * | ||||
) | [virtual] |
Definition at line 694 of file G4GDMLReadParamvol.cc.
References G4GDMLRead::check, FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4LogicalVolume::GetName(), G4GDMLParameterisation::GetSize(), G4GDMLRead::GetVolume(), kUndefined, parameterisation, Paramvol_contentRead(), G4GDMLReadDefine::RefRead(), and G4GDMLRead::Transcode().
Referenced by G4GDMLReadStructure::Volume_contentRead().
00695 { 00696 G4String volumeref; 00697 00698 parameterisation = new G4GDMLParameterisation(); 00699 00700 for (xercesc::DOMNode* iter = element->getFirstChild(); 00701 iter != 0; iter = iter->getNextSibling()) 00702 { 00703 if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; } 00704 00705 const xercesc::DOMElement* const child 00706 = dynamic_cast<xercesc::DOMElement*>(iter); 00707 if (!child) 00708 { 00709 G4Exception("G4GDMLReadParamvol::ParamvolRead()", "InvalidRead", 00710 FatalException, "No child found!"); 00711 return; 00712 } 00713 const G4String tag = Transcode(child->getTagName()); 00714 00715 if (tag=="volumeref") { volumeref = RefRead(child); } 00716 00717 } 00718 00719 Paramvol_contentRead(element); 00720 00721 G4LogicalVolume* logvol = GetVolume(GenerateName(volumeref)); 00722 00723 if (parameterisation->GetSize()==0) 00724 { 00725 G4Exception("G4GDMLReadParamvol::ParamvolRead()", 00726 "ReadError", FatalException, 00727 "No parameters are defined in parameterised volume!"); 00728 } 00729 G4String pv_name = logvol->GetName() + "_param"; 00730 new G4PVParameterised(pv_name, logvol, mother, kUndefined, 00731 parameterisation->GetSize(), parameterisation, check); 00732 }
void G4GDMLReadParamvol::Sphere_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 313 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
00315 { 00316 G4double lunit = 1.0; 00317 G4double aunit = 1.0; 00318 00319 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00320 XMLSize_t attributeCount = attributes->getLength(); 00321 00322 for (XMLSize_t attribute_index=0; 00323 attribute_index<attributeCount; attribute_index++) 00324 { 00325 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00326 00327 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00328 { continue; } 00329 00330 const xercesc::DOMAttr* const attribute 00331 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00332 if (!attribute) 00333 { 00334 G4Exception("G4GDMLReadParamvol::Sphere_dimensionsRead()", 00335 "InvalidRead", FatalException, "No attribute found!"); 00336 return; 00337 } 00338 const G4String attName = Transcode(attribute->getName()); 00339 const G4String attValue = Transcode(attribute->getValue()); 00340 00341 if (attName=="lunit") 00342 { lunit = eval.Evaluate(attValue); } else 00343 if (attName=="aunit") 00344 { aunit = eval.Evaluate(attValue); } else 00345 if (attName=="rmin") 00346 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00347 if (attName=="rmax") 00348 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00349 if (attName=="startphi") 00350 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00351 if (attName=="deltaphi") 00352 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00353 if (attName=="starttheta") 00354 { parameter.dimension[4] = eval.Evaluate(attValue); } else 00355 if (attName=="deltatheta") 00356 { parameter.dimension[5] = eval.Evaluate(attValue); } 00357 } 00358 00359 parameter.dimension[0] *= lunit; 00360 parameter.dimension[1] *= lunit; 00361 parameter.dimension[2] *= aunit; 00362 parameter.dimension[3] *= aunit; 00363 parameter.dimension[4] *= aunit; 00364 parameter.dimension[5] *= aunit; 00365 }
void G4GDMLReadParamvol::Torus_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 403 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
00405 { 00406 G4double lunit = 1.0; 00407 G4double aunit = 1.0; 00408 00409 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00410 XMLSize_t attributeCount = attributes->getLength(); 00411 00412 for (XMLSize_t attribute_index=0; 00413 attribute_index<attributeCount; attribute_index++) 00414 { 00415 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00416 00417 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00418 { continue; } 00419 00420 const xercesc::DOMAttr* const attribute 00421 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00422 if (!attribute) 00423 { 00424 G4Exception("G4GDMLReadParamvol::Torus_dimensionsRead()", 00425 "InvalidRead", FatalException, "No attribute found!"); 00426 return; 00427 } 00428 const G4String attName = Transcode(attribute->getName()); 00429 const G4String attValue = Transcode(attribute->getValue()); 00430 00431 if (attName=="lunit") 00432 { lunit = eval.Evaluate(attValue); } else 00433 if (attName=="aunit") 00434 { aunit = eval.Evaluate(attValue); } else 00435 if (attName=="rmin") 00436 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00437 if (attName=="rmax") 00438 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00439 if (attName=="rtor") 00440 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00441 if (attName=="startphi") 00442 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00443 if (attName=="deltaphi") 00444 { parameter.dimension[4] = eval.Evaluate(attValue); } 00445 } 00446 00447 parameter.dimension[0] *= lunit; 00448 parameter.dimension[1] *= lunit; 00449 parameter.dimension[2] *= lunit; 00450 parameter.dimension[3] *= aunit; 00451 parameter.dimension[4] *= aunit; 00452 }
void G4GDMLReadParamvol::Trap_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 133 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00135 { 00136 G4double lunit = 1.0; 00137 G4double aunit = 1.0; 00138 00139 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00140 XMLSize_t attributeCount = attributes->getLength(); 00141 00142 for (XMLSize_t attribute_index=0; 00143 attribute_index<attributeCount; attribute_index++) 00144 { 00145 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00146 00147 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00148 { continue; } 00149 00150 const xercesc::DOMAttr* const attribute 00151 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00152 if (!attribute) 00153 { 00154 G4Exception("G4GDMLReadParamvol::Trap_dimensionsRead()", 00155 "InvalidRead", FatalException, "No attribute found!"); 00156 return; 00157 } 00158 const G4String attName = Transcode(attribute->getName()); 00159 const G4String attValue = Transcode(attribute->getValue()); 00160 00161 if (attName=="lunit") 00162 { lunit = eval.Evaluate(attValue); } else 00163 if (attName=="aunit") 00164 { aunit = eval.Evaluate(attValue); } else 00165 if (attName=="z") 00166 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00167 if (attName=="theta") 00168 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00169 if (attName=="phi") 00170 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00171 if (attName=="y1") 00172 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00173 if (attName=="x1") 00174 { parameter.dimension[4] = eval.Evaluate(attValue); } else 00175 if (attName=="x2") 00176 { parameter.dimension[5] = eval.Evaluate(attValue); } else 00177 if (attName=="alpha1") 00178 { parameter.dimension[6] = eval.Evaluate(attValue); } else 00179 if (attName=="y2") 00180 { parameter.dimension[7] = eval.Evaluate(attValue); } else 00181 if (attName=="x3") 00182 { parameter.dimension[8] = eval.Evaluate(attValue); } else 00183 if (attName=="x4") 00184 { parameter.dimension[9] = eval.Evaluate(attValue); } else 00185 if (attName=="alpha2") 00186 { parameter.dimension[10] = eval.Evaluate(attValue); } 00187 } 00188 00189 parameter.dimension[0] *= 0.5*lunit; 00190 parameter.dimension[1] *= aunit; 00191 parameter.dimension[2] *= aunit; 00192 parameter.dimension[3] *= 0.5*lunit; 00193 parameter.dimension[4] *= 0.5*lunit; 00194 parameter.dimension[5] *= 0.5*lunit; 00195 parameter.dimension[6] *= aunit; 00196 parameter.dimension[7] *= 0.5*lunit; 00197 parameter.dimension[8] *= 0.5*lunit; 00198 parameter.dimension[9] *= 0.5*lunit; 00199 parameter.dimension[10] *= aunit; 00200 }
void G4GDMLReadParamvol::Trd_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 90 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00092 { 00093 G4double lunit = 1.0; 00094 00095 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00096 XMLSize_t attributeCount = attributes->getLength(); 00097 00098 for (XMLSize_t attribute_index=0; 00099 attribute_index<attributeCount; attribute_index++) 00100 { 00101 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00102 00103 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00104 { continue; } 00105 00106 const xercesc::DOMAttr* const attribute 00107 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00108 if (!attribute) 00109 { 00110 G4Exception("G4GDMLReadParamvol::Trd_dimensionsRead()", 00111 "InvalidRead", FatalException, "No attribute found!"); 00112 return; 00113 } 00114 const G4String attName = Transcode(attribute->getName()); 00115 const G4String attValue = Transcode(attribute->getValue()); 00116 00117 if (attName=="lunit") { lunit = eval.Evaluate(attValue); } else 00118 if (attName=="x1") { parameter.dimension[0]=eval.Evaluate(attValue); } else 00119 if (attName=="x2") { parameter.dimension[1]=eval.Evaluate(attValue); } else 00120 if (attName=="y1") { parameter.dimension[2]=eval.Evaluate(attValue); } else 00121 if (attName=="y2") { parameter.dimension[3]=eval.Evaluate(attValue); } else 00122 if (attName=="z") { parameter.dimension[4]=eval.Evaluate(attValue); } 00123 } 00124 00125 parameter.dimension[0] *= 0.5*lunit; 00126 parameter.dimension[1] *= 0.5*lunit; 00127 parameter.dimension[2] *= 0.5*lunit; 00128 parameter.dimension[3] *= 0.5*lunit; 00129 parameter.dimension[4] *= 0.5*lunit; 00130 }
void G4GDMLReadParamvol::Tube_dimensionsRead | ( | const xercesc::DOMElement * | const, | |
G4GDMLParameterisation::PARAMETER & | ||||
) | [protected] |
Definition at line 203 of file G4GDMLReadParamvol.cc.
References G4GDMLParameterisation::PARAMETER::dimension, G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().
Referenced by ParametersRead().
00205 { 00206 G4double lunit = 1.0; 00207 G4double aunit = 1.0; 00208 00209 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes(); 00210 XMLSize_t attributeCount = attributes->getLength(); 00211 00212 for (XMLSize_t attribute_index=0; 00213 attribute_index<attributeCount; attribute_index++) 00214 { 00215 xercesc::DOMNode* attribute_node = attributes->item(attribute_index); 00216 00217 if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE) 00218 { continue; } 00219 00220 const xercesc::DOMAttr* const attribute 00221 = dynamic_cast<xercesc::DOMAttr*>(attribute_node); 00222 if (!attribute) 00223 { 00224 G4Exception("G4GDMLReadParamvol::Tube_dimensionsRead()", 00225 "InvalidRead", FatalException, "No attribute found!"); 00226 return; 00227 } 00228 const G4String attName = Transcode(attribute->getName()); 00229 const G4String attValue = Transcode(attribute->getValue()); 00230 00231 if (attName=="lunit") 00232 { lunit = eval.Evaluate(attValue); } else 00233 if (attName=="aunit") 00234 { aunit = eval.Evaluate(attValue); } else 00235 if (attName=="InR") 00236 { parameter.dimension[0] = eval.Evaluate(attValue); } else 00237 if (attName=="OutR") 00238 { parameter.dimension[1] = eval.Evaluate(attValue); } else 00239 if (attName=="hz") 00240 { parameter.dimension[2] = eval.Evaluate(attValue); } else 00241 if (attName=="StartPhi") 00242 { parameter.dimension[3] = eval.Evaluate(attValue); } else 00243 if (attName=="DeltaPhi") 00244 { parameter.dimension[4] = eval.Evaluate(attValue); } 00245 } 00246 00247 parameter.dimension[0] *= lunit; 00248 parameter.dimension[1] *= lunit; 00249 parameter.dimension[2] *= 0.5*lunit; 00250 parameter.dimension[3] *= aunit; 00251 parameter.dimension[4] *= aunit; 00252 }
Definition at line 89 of file G4GDMLReadParamvol.hh.
Referenced by ParametersRead(), and ParamvolRead().