Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GDMLWriteParamvol Class Reference

#include <G4GDMLWriteParamvol.hh>

Inheritance diagram for G4GDMLWriteParamvol:
Collaboration diagram for G4GDMLWriteParamvol:

Public Member Functions

virtual void ParamvolWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const)
 
virtual void ParamvolAlgorithmWrite (xercesc::DOMElement *paramvolElement, const G4VPhysicalVolume *const paramvol)
 
- Public Member Functions inherited from G4GDMLWriteSetup
virtual void SetupWrite (xercesc::DOMElement *, const G4LogicalVolume *const)
 
- Public Member Functions inherited from G4GDMLWriteSolids
virtual void AddSolid (const G4VSolid *const)
 
virtual void SolidsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteMaterials
void AddIsotope (const G4Isotope *const)
 
void AddElement (const G4Element *const)
 
void AddMaterial (const G4Material *const)
 
virtual void MaterialsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteDefine
G4ThreeVector GetAngles (const G4RotationMatrix &)
 
void ScaleWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
 
void RotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void PositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void FirstrotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void FirstpositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void AddPosition (const G4String &name, const G4ThreeVector &pos)
 
virtual void DefineWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWrite
G4Transform3D Write (const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
 
void AddModule (const G4VPhysicalVolume *const topVol)
 
void AddModule (const G4int depth)
 
void AddAuxiliary (G4GDMLAuxStructType myaux)
 
virtual void StructureWrite (xercesc::DOMElement *)=0
 
virtual G4Transform3D TraverseVolumeTree (const G4LogicalVolume *const, const G4int)=0
 
virtual void SurfacesWrite ()=0
 
virtual void ExtensionWrite (xercesc::DOMElement *)
 
virtual void UserinfoWrite (xercesc::DOMElement *)
 
virtual void AddExtension (xercesc::DOMElement *, const G4LogicalVolume *const)
 
G4String GenerateName (const G4String &, const void *const)
 

Protected Member Functions

 G4GDMLWriteParamvol ()
 
virtual ~G4GDMLWriteParamvol ()
 
void Box_dimensionsWrite (xercesc::DOMElement *, const G4Box *const)
 
void Trd_dimensionsWrite (xercesc::DOMElement *, const G4Trd *const)
 
void Trap_dimensionsWrite (xercesc::DOMElement *, const G4Trap *const)
 
void Tube_dimensionsWrite (xercesc::DOMElement *, const G4Tubs *const)
 
void Cone_dimensionsWrite (xercesc::DOMElement *, const G4Cons *const)
 
void Sphere_dimensionsWrite (xercesc::DOMElement *, const G4Sphere *const)
 
void Orb_dimensionsWrite (xercesc::DOMElement *, const G4Orb *const)
 
void Torus_dimensionsWrite (xercesc::DOMElement *, const G4Torus *const)
 
void Ellipsoid_dimensionsWrite (xercesc::DOMElement *, const G4Ellipsoid *const)
 
void Para_dimensionsWrite (xercesc::DOMElement *, const G4Para *const)
 
void Hype_dimensionsWrite (xercesc::DOMElement *, const G4Hype *const)
 
void Polycone_dimensionsWrite (xercesc::DOMElement *, const G4Polycone *const)
 
void Polyhedra_dimensionsWrite (xercesc::DOMElement *, const G4Polyhedra *const)
 
void ParametersWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const, const G4int &)
 
- Protected Member Functions inherited from G4GDMLWriteSetup
 G4GDMLWriteSetup ()
 
virtual ~G4GDMLWriteSetup ()
 
- Protected Member Functions inherited from G4GDMLWriteSolids
 G4GDMLWriteSolids ()
 
virtual ~G4GDMLWriteSolids ()
 
void MultiUnionWrite (xercesc::DOMElement *solElement, const G4MultiUnion *const)
 
void BooleanWrite (xercesc::DOMElement *, const G4BooleanSolid *const)
 
void ScaledWrite (xercesc::DOMElement *, const G4ScaledSolid *const)
 
void BoxWrite (xercesc::DOMElement *, const G4Box *const)
 
void ConeWrite (xercesc::DOMElement *, const G4Cons *const)
 
void ElconeWrite (xercesc::DOMElement *, const G4EllipticalCone *const)
 
void EllipsoidWrite (xercesc::DOMElement *, const G4Ellipsoid *const)
 
void EltubeWrite (xercesc::DOMElement *, const G4EllipticalTube *const)
 
void XtruWrite (xercesc::DOMElement *, const G4ExtrudedSolid *const)
 
void HypeWrite (xercesc::DOMElement *, const G4Hype *const)
 
void OrbWrite (xercesc::DOMElement *, const G4Orb *const)
 
void ParaWrite (xercesc::DOMElement *, const G4Para *const)
 
void ParaboloidWrite (xercesc::DOMElement *, const G4Paraboloid *const)
 
void PolyconeWrite (xercesc::DOMElement *, const G4Polycone *const)
 
void GenericPolyconeWrite (xercesc::DOMElement *, const G4GenericPolycone *const)
 
void PolyhedraWrite (xercesc::DOMElement *, const G4Polyhedra *const)
 
void SphereWrite (xercesc::DOMElement *, const G4Sphere *const)
 
void TessellatedWrite (xercesc::DOMElement *, const G4TessellatedSolid *const)
 
void TetWrite (xercesc::DOMElement *, const G4Tet *const)
 
void TorusWrite (xercesc::DOMElement *, const G4Torus *const)
 
void GenTrapWrite (xercesc::DOMElement *, const G4GenericTrap *const)
 
void TrapWrite (xercesc::DOMElement *, const G4Trap *const)
 
void TrdWrite (xercesc::DOMElement *, const G4Trd *const)
 
void TubeWrite (xercesc::DOMElement *, const G4Tubs *const)
 
void CutTubeWrite (xercesc::DOMElement *, const G4CutTubs *const)
 
void TwistedboxWrite (xercesc::DOMElement *, const G4TwistedBox *const)
 
void TwistedtrapWrite (xercesc::DOMElement *, const G4TwistedTrap *const)
 
void TwistedtrdWrite (xercesc::DOMElement *, const G4TwistedTrd *const)
 
void TwistedtubsWrite (xercesc::DOMElement *, const G4TwistedTubs *const)
 
void ZplaneWrite (xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
 
void RZPointWrite (xercesc::DOMElement *, const G4double &, const G4double &)
 
void OpticalSurfaceWrite (xercesc::DOMElement *, const G4OpticalSurface *const)
 
- Protected Member Functions inherited from G4GDMLWriteMaterials
 G4GDMLWriteMaterials ()
 
virtual ~G4GDMLWriteMaterials ()
 
void AtomWrite (xercesc::DOMElement *, const G4double &)
 
void DWrite (xercesc::DOMElement *, const G4double &)
 
void PWrite (xercesc::DOMElement *, const G4double &)
 
void TWrite (xercesc::DOMElement *, const G4double &)
 
void MEEWrite (xercesc::DOMElement *, const G4double &)
 
void IsotopeWrite (const G4Isotope *const)
 
void ElementWrite (const G4Element *const)
 
void MaterialWrite (const G4Material *const)
 
void PropertyWrite (xercesc::DOMElement *, const G4Material *const)
 
void PropertyVectorWrite (const G4String &, const G4PhysicsOrderedFreeVector *const)
 
- Protected Member Functions inherited from G4GDMLWriteDefine
 G4GDMLWriteDefine ()
 
virtual ~G4GDMLWriteDefine ()
 
void Scale_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Rotation_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Position_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
- Protected Member Functions inherited from G4GDMLWrite
 G4GDMLWrite ()
 
virtual ~G4GDMLWrite ()
 
VolumeMapType & VolumeMap ()
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4String &)
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4double &)
 
xercesc::DOMElement * NewElement (const G4String &)
 
G4String Modularize (const G4VPhysicalVolume *const topvol, const G4int depth)
 
void AddAuxInfo (G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
 
G4bool FileExists (const G4String &) const
 
PhysVolumeMapType & PvolumeMap ()
 
DepthMapType & DepthMap ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GDMLWrite
static void SetAddPointerToName (G4bool)
 
- Protected Attributes inherited from G4GDMLWriteSolids
std::vector< const G4VSolid * > solidList
 
xercesc::DOMElement * solidsElement
 
- Protected Attributes inherited from G4GDMLWriteMaterials
std::vector< const G4Isotope * > isotopeList
 
std::vector< const G4Element * > elementList
 
std::vector< const G4Material * > materialList
 
xercesc::DOMElement * materialsElement
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc
 
xercesc::DOMElement * extElement
 
xercesc::DOMElement * userinfoElement
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 
- Static Protected Attributes inherited from G4GDMLWriteSolids
static const G4int maxTransforms = 8
 
- Static Protected Attributes inherited from G4GDMLWriteDefine
static const G4double kRelativePrecision = DBL_EPSILON
 
static const G4double kAngularPrecision = DBL_EPSILON
 
static const G4double kLinearPrecision = DBL_EPSILON
 
- Static Protected Attributes inherited from G4GDMLWrite
static G4bool addPointerToName = true
 

Detailed Description

Definition at line 60 of file G4GDMLWriteParamvol.hh.

Constructor & Destructor Documentation

G4GDMLWriteParamvol::G4GDMLWriteParamvol ( )
protected

Definition at line 58 of file G4GDMLWriteParamvol.cc.

59 {
60 }
G4GDMLWriteParamvol::~G4GDMLWriteParamvol ( )
protectedvirtual

Definition at line 63 of file G4GDMLWriteParamvol.cc.

64 {
65 }

Member Function Documentation

void G4GDMLWriteParamvol::Box_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Box * const  box 
)
protected

Definition at line 68 of file G4GDMLWriteParamvol.cc.

70 {
71  xercesc::DOMElement* box_dimensionsElement = NewElement("box_dimensions");
72  box_dimensionsElement->
73  setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
74  box_dimensionsElement->
75  setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
76  box_dimensionsElement->
77  setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
78  box_dimensionsElement->
79  setAttributeNode(NewAttribute("lunit","mm"));
80  parametersElement->appendChild(box_dimensionsElement);
81 }
G4double GetXHalfLength() const
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
G4double GetYHalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Cone_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Cons * const  cone 
)
protected

Definition at line 168 of file G4GDMLWriteParamvol.cc.

170 {
171  xercesc::DOMElement* cone_dimensionsElement = NewElement("cone_dimensions");
172  cone_dimensionsElement->
173  setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
174  cone_dimensionsElement->
175  setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
176  cone_dimensionsElement->
177  setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
178  cone_dimensionsElement->
179  setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
180  cone_dimensionsElement->
181  setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
182  cone_dimensionsElement->
183  setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
184  cone_dimensionsElement->
185  setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
186  cone_dimensionsElement->
187  setAttributeNode(NewAttribute("aunit","deg"));
188  cone_dimensionsElement->
189  setAttributeNode(NewAttribute("lunit","mm"));
190  parametersElement->appendChild(cone_dimensionsElement);
191 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetOuterRadiusMinusZ() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetStartPhiAngle() const
G4double GetInnerRadiusPlusZ() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetInnerRadiusMinusZ() const
G4double GetOuterRadiusPlusZ() const
G4double GetZHalfLength() const
G4double GetDeltaPhiAngle() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Ellipsoid_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Ellipsoid * const  ellipsoid 
)
protected

Definition at line 251 of file G4GDMLWriteParamvol.cc.

253 {
254  xercesc::DOMElement* ellipsoid_dimensionsElement =
255  NewElement("ellipsoid_dimensions");
256  ellipsoid_dimensionsElement->
257  setAttributeNode(NewAttribute("ax",ellipsoid->GetSemiAxisMax(0)/mm));
258  ellipsoid_dimensionsElement->
259  setAttributeNode(NewAttribute("by",ellipsoid->GetSemiAxisMax(1)/mm));
260  ellipsoid_dimensionsElement->
261  setAttributeNode(NewAttribute("cz",ellipsoid->GetSemiAxisMax(2)/mm));
262  ellipsoid_dimensionsElement->
263  setAttributeNode(NewAttribute("zcut1",ellipsoid->GetZBottomCut()/mm));
264  ellipsoid_dimensionsElement->
265  setAttributeNode(NewAttribute("zcut2",ellipsoid->GetZTopCut()/mm));
266  ellipsoid_dimensionsElement->
267  setAttributeNode(NewAttribute("lunit","mm"));
268  parametersElement->appendChild(ellipsoid_dimensionsElement);
269 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetZTopCut() const
G4double GetZBottomCut() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetSemiAxisMax(G4int i) const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Hype_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Hype * const  hype 
)
protected

Definition at line 303 of file G4GDMLWriteParamvol.cc.

305 {
306  xercesc::DOMElement* hype_dimensionsElement = NewElement("hype_dimensions");
307  hype_dimensionsElement->
308  setAttributeNode(NewAttribute("rmin",hype->GetInnerRadius()/mm));
309  hype_dimensionsElement->
310  setAttributeNode(NewAttribute("rmax",hype->GetOuterRadius()/mm));
311  hype_dimensionsElement->
312  setAttributeNode(NewAttribute("inst",hype->GetInnerStereo()/degree));
313  hype_dimensionsElement->
314  setAttributeNode(NewAttribute("outst",hype->GetOuterStereo()/degree));
315  hype_dimensionsElement->
316  setAttributeNode(NewAttribute("z",2.0*hype->GetZHalfLength()/mm));
317  hype_dimensionsElement->
318  setAttributeNode(NewAttribute("aunit","deg"));
319  hype_dimensionsElement->
320  setAttributeNode(NewAttribute("lunit","mm"));
321  parametersElement->appendChild(hype_dimensionsElement);
322 }
G4double GetOuterStereo() const
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetInnerStereo() const
G4double GetInnerRadius() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetOuterRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZHalfLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Orb_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Orb * const  orb 
)
protected

Definition at line 217 of file G4GDMLWriteParamvol.cc.

219 {
220  xercesc::DOMElement* orb_dimensionsElement = NewElement("orb_dimensions");
221  orb_dimensionsElement->setAttributeNode(NewAttribute("r",
222  orb->GetRadius()/mm));
223  orb_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
224  parametersElement->appendChild(orb_dimensionsElement);
225 }
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Para_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Para * const  para 
)
protected

Definition at line 272 of file G4GDMLWriteParamvol.cc.

274 {
275  const G4ThreeVector simaxis = para->GetSymAxis();
276 
277  const G4double alpha = std::atan(para->GetTanAlpha());
278  const G4double theta = std::acos(simaxis.z());
279  const G4double phi = (simaxis.z() != 1.0)
280  ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
281 
282  xercesc::DOMElement* para_dimensionsElement = NewElement("para_dimensions");
283  para_dimensionsElement->
284  setAttributeNode(NewAttribute("x",2.0*para->GetXHalfLength()/mm));
285  para_dimensionsElement->
286  setAttributeNode(NewAttribute("y",2.0*para->GetYHalfLength()/mm));
287  para_dimensionsElement->
288  setAttributeNode(NewAttribute("z",2.0*para->GetZHalfLength()/mm));
289  para_dimensionsElement->
290  setAttributeNode(NewAttribute("alpha",alpha/degree));
291  para_dimensionsElement->
292  setAttributeNode(NewAttribute("theta",theta/degree));
293  para_dimensionsElement->
294  setAttributeNode(NewAttribute("phi",phi/degree));
295  para_dimensionsElement->
296  setAttributeNode(NewAttribute("aunit","deg"));
297  para_dimensionsElement->
298  setAttributeNode(NewAttribute("lunit","mm"));
299  parametersElement->appendChild(para_dimensionsElement);
300 }
G4ThreeVector GetSymAxis() const
static constexpr double mm
Definition: G4SIunits.hh:115
double x() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
double z() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetXHalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetTanAlpha() const
double y() const
G4double GetZHalfLength() const
double G4double
Definition: G4Types.hh:76
static const G4double alpha
G4double GetYHalfLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::ParametersWrite ( xercesc::DOMElement *  paramvolElement,
const G4VPhysicalVolume * const  paramvol,
const G4int index 
)
protected

Definition at line 385 of file G4GDMLWriteParamvol.cc.

387 {
388  paramvol->GetParameterisation()
389  ->ComputeTransformation(index, const_cast<G4VPhysicalVolume*>(paramvol));
390  G4ThreeVector Angles;
391  G4String name = GenerateName(paramvol->GetName(),paramvol);
392  std::stringstream os;
393  os.precision(15);
394  os << index;
395  G4String sncopie = os.str();
396 
397  xercesc::DOMElement* parametersElement = NewElement("parameters");
398  parametersElement->setAttributeNode(NewAttribute("number",index+1));
399 
400  PositionWrite(parametersElement, name+sncopie+"_pos",
401  paramvol->GetObjectTranslation());
402  Angles=GetAngles(paramvol->GetObjectRotationValue());
403  if (Angles.mag2()>DBL_EPSILON)
404  {
405  RotationWrite(parametersElement, name+sncopie+"_rot",
406  GetAngles(paramvol->GetObjectRotationValue()));
407  }
408  paramvolElement->appendChild(parametersElement);
409 
410  G4VSolid* solid = paramvol->GetLogicalVolume()->GetSolid();
411 
412  if (G4Box* box = dynamic_cast<G4Box*>(solid))
413  {
414  paramvol->GetParameterisation()->ComputeDimensions(*box,index,
415  const_cast<G4VPhysicalVolume*>(paramvol));
416  Box_dimensionsWrite(parametersElement,box);
417  } else
418  if (G4Trd* trd = dynamic_cast<G4Trd*>(solid))
419  {
420  paramvol->GetParameterisation()->ComputeDimensions(*trd,index,
421  const_cast<G4VPhysicalVolume*>(paramvol));
422  Trd_dimensionsWrite(parametersElement,trd);
423  } else
424  if (G4Trap* trap = dynamic_cast<G4Trap*>(solid))
425  {
426  paramvol->GetParameterisation()->ComputeDimensions(*trap,index,
427  const_cast<G4VPhysicalVolume*>(paramvol));
428  Trap_dimensionsWrite(parametersElement,trap);
429  } else
430  if (G4Tubs* tube = dynamic_cast<G4Tubs*>(solid))
431  {
432  paramvol->GetParameterisation()->ComputeDimensions(*tube,index,
433  const_cast<G4VPhysicalVolume*>(paramvol));
434  Tube_dimensionsWrite(parametersElement,tube);
435  } else
436  if (G4Cons* cone = dynamic_cast<G4Cons*>(solid))
437  {
438  paramvol->GetParameterisation()->ComputeDimensions(*cone,index,
439  const_cast<G4VPhysicalVolume*>(paramvol));
440  Cone_dimensionsWrite(parametersElement,cone);
441  } else
442  if (G4Sphere* sphere = dynamic_cast<G4Sphere*>(solid))
443  {
444  paramvol->GetParameterisation()->ComputeDimensions(*sphere,index,
445  const_cast<G4VPhysicalVolume*>(paramvol));
446  Sphere_dimensionsWrite(parametersElement,sphere);
447  } else
448  if (G4Orb* orb = dynamic_cast<G4Orb*>(solid))
449  {
450  paramvol->GetParameterisation()->ComputeDimensions(*orb,index,
451  const_cast<G4VPhysicalVolume*>(paramvol));
452  Orb_dimensionsWrite(parametersElement,orb);
453  } else
454  if (G4Torus* torus = dynamic_cast<G4Torus*>(solid))
455  {
456  paramvol->GetParameterisation()->ComputeDimensions(*torus,index,
457  const_cast<G4VPhysicalVolume*>(paramvol));
458  Torus_dimensionsWrite(parametersElement,torus);
459  } else
460  if (G4Ellipsoid* ellipsoid = dynamic_cast<G4Ellipsoid*>(solid))
461  {
462  paramvol->GetParameterisation()->ComputeDimensions(*ellipsoid,index,
463  const_cast<G4VPhysicalVolume*>(paramvol));
464  Ellipsoid_dimensionsWrite(parametersElement,ellipsoid);
465  } else
466  if (G4Para* para = dynamic_cast<G4Para*>(solid))
467  {
468  paramvol->GetParameterisation()->ComputeDimensions(*para,index,
469  const_cast<G4VPhysicalVolume*>(paramvol));
470  Para_dimensionsWrite(parametersElement,para);
471  } else
472  if (G4Hype* hype = dynamic_cast<G4Hype*>(solid))
473  {
474  paramvol->GetParameterisation()->ComputeDimensions(*hype,index,
475  const_cast<G4VPhysicalVolume*>(paramvol));
476  Hype_dimensionsWrite(parametersElement,hype);
477  }else
478  if (G4Polycone* pcone = dynamic_cast<G4Polycone*>(solid))
479  {
480  paramvol->GetParameterisation()->ComputeDimensions(*pcone,index,
481  const_cast<G4VPhysicalVolume*>(paramvol));
482  Polycone_dimensionsWrite(parametersElement,pcone);
483  }else
484  if (G4Polyhedra* polyhedra = dynamic_cast<G4Polyhedra*>(solid))
485  {
486  paramvol->GetParameterisation()->ComputeDimensions(*polyhedra,index,
487  const_cast<G4VPhysicalVolume*>(paramvol));
488  Polyhedra_dimensionsWrite(parametersElement,polyhedra);
489  }
490  else
491  {
492  G4String error_msg = "Solid '" + solid->GetName()
493  + "' cannot be used in parameterised volume!";
494  G4Exception("G4GDMLWriteParamvol::ParametersWrite()",
495  "InvalidSetup", FatalException, error_msg);
496  }
497 }
G4String GetName() const
void Box_dimensionsWrite(xercesc::DOMElement *, const G4Box *const)
const XML_Char * name
Definition: expat.h:151
void Polycone_dimensionsWrite(xercesc::DOMElement *, const G4Polycone *const)
Definition: G4Para.hh:77
void Orb_dimensionsWrite(xercesc::DOMElement *, const G4Orb *const)
Definition: G4Box.hh:64
void Hype_dimensionsWrite(xercesc::DOMElement *, const G4Hype *const)
Definition: G4Tubs.hh:85
G4VSolid * GetSolid() const
void Polyhedra_dimensionsWrite(xercesc::DOMElement *, const G4Polyhedra *const)
Definition: G4Trd.hh:72
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
void Tube_dimensionsWrite(xercesc::DOMElement *, const G4Tubs *const)
void Cone_dimensionsWrite(xercesc::DOMElement *, const G4Cons *const)
const G4String & GetName() const
Definition: G4Hype.hh:67
Definition: G4Cons.hh:83
virtual G4VPVParameterisation * GetParameterisation() const =0
G4RotationMatrix GetObjectRotationValue() const
#define DBL_EPSILON
Definition: templates.hh:87
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
Definition: G4Orb.hh:61
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void Ellipsoid_dimensionsWrite(xercesc::DOMElement *, const G4Ellipsoid *const)
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
void Trap_dimensionsWrite(xercesc::DOMElement *, const G4Trap *const)
void Trd_dimensionsWrite(xercesc::DOMElement *, const G4Trd *const)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4LogicalVolume * GetLogicalVolume() const
void Para_dimensionsWrite(xercesc::DOMElement *, const G4Para *const)
void Torus_dimensionsWrite(xercesc::DOMElement *, const G4Torus *const)
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4ThreeVector GetAngles(const G4RotationMatrix &)
void Sphere_dimensionsWrite(xercesc::DOMElement *, const G4Sphere *const)
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
G4ThreeVector GetObjectTranslation() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::ParamvolAlgorithmWrite ( xercesc::DOMElement *  paramvolElement,
const G4VPhysicalVolume *const  paramvol 
)
virtual

Definition at line 521 of file G4GDMLWriteParamvol.cc.

523 {
524  const G4String volumeref =
525  GenerateName(paramvol->GetLogicalVolume()->GetName(),
526  paramvol->GetLogicalVolume());
527 
528  const G4int parameterCount = paramvol->GetMultiplicity();
529 
530  for (G4int i=0; i<parameterCount; i++)
531  {
532  ParametersWrite(paramvolElement,paramvol,i);
533  }
534 }
int G4int
Definition: G4Types.hh:78
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
void ParametersWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const, const G4int &)
const G4String & GetName() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::ParamvolWrite ( xercesc::DOMElement *  volumeElement,
const G4VPhysicalVolume * const  paramvol 
)
virtual

Definition at line 500 of file G4GDMLWriteParamvol.cc.

502 {
503  const G4String volumeref =
504  GenerateName(paramvol->GetLogicalVolume()->GetName(),
505  paramvol->GetLogicalVolume());
506  xercesc::DOMElement* paramvolElement = NewElement("paramvol");
507  paramvolElement->setAttributeNode(NewAttribute("ncopies",
508  paramvol->GetMultiplicity()));
509  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
510  volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
511 
512  xercesc::DOMElement* algorithmElement =
513  NewElement("parameterised_position_size");
514  paramvolElement->appendChild(volumerefElement);
515  paramvolElement->appendChild(algorithmElement);
516  ParamvolAlgorithmWrite(algorithmElement,paramvol);
517  volumeElement->appendChild(paramvolElement);
518 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
virtual void ParamvolAlgorithmWrite(xercesc::DOMElement *paramvolElement, const G4VPhysicalVolume *const paramvol)
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetMultiplicity() const
const G4String & GetName() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Polycone_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Polycone * const  pcone 
)
protected

Definition at line 325 of file G4GDMLWriteParamvol.cc.

327 {
328  xercesc::DOMElement* pcone_dimensionsElement
329  = NewElement("polycone_dimensions");
330 
331  pcone_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
333  pcone_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
335  pcone_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
337  pcone_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
338  pcone_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
339 
340  parametersElement->appendChild(pcone_dimensionsElement);
341  const size_t num_zplanes = pcone->GetOriginalParameters()->Num_z_planes;
342  const G4double* z_array = pcone->GetOriginalParameters()->Z_values;
343  const G4double* rmin_array = pcone->GetOriginalParameters()->Rmin;
344  const G4double* rmax_array = pcone->GetOriginalParameters()->Rmax;
345 
346  for (size_t i=0; i<num_zplanes; i++)
347  {
348  ZplaneWrite(pcone_dimensionsElement,z_array[i],
349  rmin_array[i],rmax_array[i]);
350  }
351 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
G4PolyconeHistorical * GetOriginalParameters() const
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Polyhedra_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Polyhedra * const  polyhedra 
)
protected

Definition at line 354 of file G4GDMLWriteParamvol.cc.

356 {
357  xercesc::DOMElement* polyhedra_dimensionsElement
358  = NewElement("polyhedra_dimensions");
359 
360  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
361  polyhedra->GetOriginalParameters()->Num_z_planes));
362  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numSide",
363  polyhedra->GetOriginalParameters()->numSide));
364  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
365  polyhedra->GetOriginalParameters()->Start_angle/degree));
366  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
368  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
369  polyhedra_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
370 
371  parametersElement->appendChild(polyhedra_dimensionsElement);
372  const size_t num_zplanes = polyhedra->GetOriginalParameters()->Num_z_planes;
373  const G4double* z_array = polyhedra->GetOriginalParameters()->Z_values;
374  const G4double* rmin_array = polyhedra->GetOriginalParameters()->Rmin;
375  const G4double* rmax_array = polyhedra->GetOriginalParameters()->Rmax;
376 
377  for (size_t i=0; i<num_zplanes; i++)
378  {
379  ZplaneWrite(polyhedra_dimensionsElement,z_array[i],
380  rmin_array[i],rmax_array[i]);
381  }
382 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4PolyhedraHistorical * GetOriginalParameters() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Sphere_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Sphere * const  sphere 
)
protected

Definition at line 194 of file G4GDMLWriteParamvol.cc.

196 {
197  xercesc::DOMElement* sphere_dimensionsElement =
198  NewElement("sphere_dimensions");
199  sphere_dimensionsElement->setAttributeNode(NewAttribute("rmin",
200  sphere->GetInnerRadius()/mm));
201  sphere_dimensionsElement->setAttributeNode(NewAttribute("rmax",
202  sphere->GetOuterRadius()/mm));
203  sphere_dimensionsElement->setAttributeNode(NewAttribute("startphi",
204  sphere->GetStartPhiAngle()/degree));
205  sphere_dimensionsElement->setAttributeNode(NewAttribute("deltaphi",
206  sphere->GetDeltaPhiAngle()/degree));
207  sphere_dimensionsElement->setAttributeNode(NewAttribute("starttheta",
208  sphere->GetStartThetaAngle()/degree));
209  sphere_dimensionsElement->setAttributeNode(NewAttribute("deltatheta",
210  sphere->GetDeltaThetaAngle()/degree));
211  sphere_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
212  sphere_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
213  parametersElement->appendChild(sphere_dimensionsElement);
214 }
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDeltaPhiAngle() const
G4double GetStartThetaAngle() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetInnerRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetDeltaThetaAngle() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Torus_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Torus * const  torus 
)
protected

Definition at line 228 of file G4GDMLWriteParamvol.cc.

230 {
231  xercesc::DOMElement* torus_dimensionsElement =
232  NewElement("torus_dimensions");
233  torus_dimensionsElement->
234  setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
235  torus_dimensionsElement->
236  setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
237  torus_dimensionsElement->
238  setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
239  torus_dimensionsElement->
240  setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
241  torus_dimensionsElement->
242  setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
243  torus_dimensionsElement->
244  setAttributeNode(NewAttribute("aunit","deg"));
245  torus_dimensionsElement->
246  setAttributeNode(NewAttribute("lunit","mm"));
247  parametersElement->appendChild(torus_dimensionsElement);
248 }
G4double GetSPhi() const
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetRmax() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetRtor() const
G4double GetRmin() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetDPhi() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Trap_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Trap * const  trap 
)
protected

Definition at line 104 of file G4GDMLWriteParamvol.cc.

106 {
107  const G4ThreeVector simaxis = trap->GetSymAxis();
108  const G4double phi = (simaxis.z() != 1.0)
109  ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
110  const G4double theta = std::acos(simaxis.z());
111  const G4double alpha1 = std::atan(trap->GetTanAlpha1());
112  const G4double alpha2 = std::atan(trap->GetTanAlpha2());
113 
114  xercesc::DOMElement* trap_dimensionsElement = NewElement("trap");
115  trap_dimensionsElement->
116  setAttributeNode(NewAttribute("z",2.0*trap->GetZHalfLength()/mm));
117  trap_dimensionsElement->
118  setAttributeNode(NewAttribute("theta",theta/degree));
119  trap_dimensionsElement->
120  setAttributeNode(NewAttribute("phi",phi/degree));
121  trap_dimensionsElement->
122  setAttributeNode(NewAttribute("y1",2.0*trap->GetYHalfLength1()/mm));
123  trap_dimensionsElement->
124  setAttributeNode(NewAttribute("x1",2.0*trap->GetXHalfLength1()/mm));
125  trap_dimensionsElement->
126  setAttributeNode(NewAttribute("x2",2.0*trap->GetXHalfLength2()/mm));
127  trap_dimensionsElement->
128  setAttributeNode(NewAttribute("alpha1",alpha1/degree));
129  trap_dimensionsElement->
130  setAttributeNode(NewAttribute("y2",2.0*trap->GetYHalfLength2()/mm));
131  trap_dimensionsElement->
132  setAttributeNode(NewAttribute("x3",2.0*trap->GetXHalfLength3()/mm));
133  trap_dimensionsElement->
134  setAttributeNode(NewAttribute("x4",2.0*trap->GetXHalfLength4()/mm));
135  trap_dimensionsElement->
136  setAttributeNode(NewAttribute("alpha2",alpha2/degree));
137  trap_dimensionsElement->
138  setAttributeNode(NewAttribute("aunit","deg"));
139  trap_dimensionsElement->
140  setAttributeNode(NewAttribute("lunit","mm"));
141  parametersElement->appendChild(trap_dimensionsElement);
142 }
G4double GetXHalfLength4() const
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetYHalfLength2() const
double x() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetXHalfLength2() const
double z() const
G4double GetTanAlpha2() const
G4double GetXHalfLength1() const
G4double GetXHalfLength3() const
static constexpr double degree
Definition: G4SIunits.hh:144
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
double y() const
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4double GetTanAlpha1() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Trd_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Trd * const  trd 
)
protected

Definition at line 84 of file G4GDMLWriteParamvol.cc.

86 {
87  xercesc::DOMElement* trd_dimensionsElement = NewElement("trd_dimensions");
88  trd_dimensionsElement->
89  setAttributeNode(NewAttribute("x1",2.0*trd->GetXHalfLength1()/mm));
90  trd_dimensionsElement->
91  setAttributeNode(NewAttribute("x2",2.0*trd->GetXHalfLength2()/mm));
92  trd_dimensionsElement->
93  setAttributeNode(NewAttribute("y1",2.0*trd->GetYHalfLength1()/mm));
94  trd_dimensionsElement->
95  setAttributeNode(NewAttribute("y2",2.0*trd->GetYHalfLength2()/mm));
96  trd_dimensionsElement->
97  setAttributeNode(NewAttribute("z",2.0*trd->GetZHalfLength()/mm));
98  trd_dimensionsElement->
99  setAttributeNode(NewAttribute("lunit","mm"));
100  parametersElement->appendChild(trd_dimensionsElement);
101 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetXHalfLength1() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteParamvol::Tube_dimensionsWrite ( xercesc::DOMElement *  parametersElement,
const G4Tubs * const  tube 
)
protected

Definition at line 145 of file G4GDMLWriteParamvol.cc.

147 {
148  xercesc::DOMElement* tube_dimensionsElement = NewElement("tube_dimensions");
149  tube_dimensionsElement->
150  setAttributeNode(NewAttribute("InR",tube->GetInnerRadius()/mm));
151  tube_dimensionsElement->
152  setAttributeNode(NewAttribute("OutR",tube->GetOuterRadius()/mm));
153  tube_dimensionsElement->
154  setAttributeNode(NewAttribute("hz",2.0*tube->GetZHalfLength()/mm));
155  tube_dimensionsElement->
156  setAttributeNode(NewAttribute("StartPhi",tube->GetStartPhiAngle()/degree));
157  tube_dimensionsElement->
158  setAttributeNode(NewAttribute("DeltaPhi",tube->GetDeltaPhiAngle()/degree));
159  tube_dimensionsElement->
160  setAttributeNode(NewAttribute("aunit","deg"));
161  tube_dimensionsElement->
162  setAttributeNode(NewAttribute("lunit","mm"));
163  parametersElement->appendChild(tube_dimensionsElement);
164 }
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDeltaPhiAngle() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetStartPhiAngle() const
G4double GetInnerRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZHalfLength() const
G4double GetOuterRadius() const

Here is the call graph for this function:

Here is the caller graph for this function:


The documentation for this class was generated from the following files: