Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:19:14

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // $Id: PHG4GDMLWriteParamvol.cc 77913 2013-11-29 10:59:07Z gcosmo $
0028 //
0029 // class G4GDMLParamVol Implementation
0030 //
0031 // Original author: Zoltan Torzsok, November 2007
0032 //
0033 // --------------------------------------------------------------------
0034 
0035 #include "PHG4GDMLWriteParamvol.hh"
0036 #include "PHG4GDMLWriteSolids.hh"
0037 
0038 #include <Geant4/G4SystemOfUnits.hh>
0039 #include <Geant4/G4Box.hh>
0040 #include <Geant4/G4Trd.hh>
0041 #include <Geant4/G4Trap.hh>
0042 #include <Geant4/G4Tubs.hh>
0043 #include <Geant4/G4Cons.hh>
0044 #include <Geant4/G4Sphere.hh>
0045 #include <Geant4/G4Orb.hh>
0046 #include <Geant4/G4Torus.hh>
0047 #include <Geant4/G4Ellipsoid.hh>
0048 #include <Geant4/G4Para.hh>
0049 #include <Geant4/G4Hype.hh>
0050 #include <Geant4/G4Polycone.hh>
0051 #include <Geant4/G4Polyhedra.hh>
0052 #include <Geant4/G4LogicalVolume.hh>
0053 #include <Geant4/G4VPhysicalVolume.hh>
0054 #include <Geant4/G4PVParameterised.hh>
0055 #include <Geant4/G4VPVParameterisation.hh>
0056 
0057 PHG4GDMLWriteParamvol::
0058 PHG4GDMLWriteParamvol() : PHG4GDMLWriteSetup()
0059 {
0060 }
0061 
0062 PHG4GDMLWriteParamvol::
0063 ~PHG4GDMLWriteParamvol()
0064 {
0065 }
0066 
0067 void PHG4GDMLWriteParamvol::
0068 Box_dimensionsWrite(xercesc::DOMElement* parametersElement,
0069                     const G4Box* const box)
0070 {
0071    xercesc::DOMElement* box_dimensionsElement = NewElement("box_dimensions");
0072    box_dimensionsElement->
0073      setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
0074    box_dimensionsElement->
0075      setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
0076    box_dimensionsElement->
0077      setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
0078    box_dimensionsElement->
0079      setAttributeNode(NewAttribute("lunit","mm"));
0080    parametersElement->appendChild(box_dimensionsElement);
0081 }
0082 
0083 void PHG4GDMLWriteParamvol::
0084 Trd_dimensionsWrite(xercesc::DOMElement* parametersElement,
0085                     const G4Trd* const trd)
0086 {
0087    xercesc::DOMElement* trd_dimensionsElement = NewElement("trd_dimensions");
0088    trd_dimensionsElement->
0089      setAttributeNode(NewAttribute("x1",2.0*trd->GetXHalfLength1()/mm));
0090    trd_dimensionsElement->
0091      setAttributeNode(NewAttribute("x2",2.0*trd->GetXHalfLength2()/mm));
0092    trd_dimensionsElement->
0093      setAttributeNode(NewAttribute("y1",2.0*trd->GetYHalfLength1()/mm));
0094    trd_dimensionsElement->
0095      setAttributeNode(NewAttribute("y2",2.0*trd->GetYHalfLength2()/mm));
0096    trd_dimensionsElement->
0097      setAttributeNode(NewAttribute("z",2.0*trd->GetZHalfLength()/mm));
0098    trd_dimensionsElement->
0099      setAttributeNode(NewAttribute("lunit","mm"));
0100    parametersElement->appendChild(trd_dimensionsElement);
0101 }
0102 
0103 void PHG4GDMLWriteParamvol::
0104 Trap_dimensionsWrite(xercesc::DOMElement* parametersElement,
0105                      const G4Trap* const trap)
0106 {
0107    const G4ThreeVector simaxis = trap->GetSymAxis();
0108    const G4double phi = (simaxis.z() != 1.0)
0109                       ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
0110    const G4double theta = std::acos(simaxis.z());
0111    const G4double alpha1 = std::atan(trap->GetTanAlpha1());
0112    const G4double alpha2 = std::atan(trap->GetTanAlpha2());
0113 
0114    xercesc::DOMElement* trap_dimensionsElement = NewElement("trap");
0115    trap_dimensionsElement->
0116      setAttributeNode(NewAttribute("z",2.0*trap->GetZHalfLength()/mm));
0117    trap_dimensionsElement->
0118      setAttributeNode(NewAttribute("theta",theta/degree));
0119    trap_dimensionsElement->
0120      setAttributeNode(NewAttribute("phi",phi/degree));
0121    trap_dimensionsElement->
0122      setAttributeNode(NewAttribute("y1",2.0*trap->GetYHalfLength1()/mm));
0123    trap_dimensionsElement->
0124      setAttributeNode(NewAttribute("x1",2.0*trap->GetXHalfLength1()/mm));
0125    trap_dimensionsElement->
0126      setAttributeNode(NewAttribute("x2",2.0*trap->GetXHalfLength2()/mm));
0127    trap_dimensionsElement->
0128      setAttributeNode(NewAttribute("alpha1",alpha1/degree));
0129    trap_dimensionsElement->
0130      setAttributeNode(NewAttribute("y2",2.0*trap->GetYHalfLength2()/mm));
0131    trap_dimensionsElement->
0132      setAttributeNode(NewAttribute("x3",2.0*trap->GetXHalfLength3()/mm));
0133    trap_dimensionsElement->
0134      setAttributeNode(NewAttribute("x4",2.0*trap->GetXHalfLength4()/mm));
0135    trap_dimensionsElement->
0136      setAttributeNode(NewAttribute("alpha2",alpha2/degree));
0137    trap_dimensionsElement->
0138      setAttributeNode(NewAttribute("aunit","deg"));
0139    trap_dimensionsElement->
0140      setAttributeNode(NewAttribute("lunit","mm"));
0141    parametersElement->appendChild(trap_dimensionsElement);
0142 }
0143 
0144 void PHG4GDMLWriteParamvol::
0145 Tube_dimensionsWrite(xercesc::DOMElement* parametersElement,
0146                      const G4Tubs* const tube)
0147 {
0148    xercesc::DOMElement* tube_dimensionsElement = NewElement("tube_dimensions");
0149    tube_dimensionsElement->
0150      setAttributeNode(NewAttribute("InR",tube->GetInnerRadius()/mm));
0151    tube_dimensionsElement->
0152      setAttributeNode(NewAttribute("OutR",tube->GetOuterRadius()/mm));
0153    tube_dimensionsElement->
0154      setAttributeNode(NewAttribute("hz",2.0*tube->GetZHalfLength()/mm));
0155    tube_dimensionsElement->
0156      setAttributeNode(NewAttribute("StartPhi",tube->GetStartPhiAngle()/degree));
0157    tube_dimensionsElement->
0158      setAttributeNode(NewAttribute("DeltaPhi",tube->GetDeltaPhiAngle()/degree));
0159    tube_dimensionsElement->
0160      setAttributeNode(NewAttribute("aunit","deg"));
0161    tube_dimensionsElement->
0162      setAttributeNode(NewAttribute("lunit","mm"));
0163    parametersElement->appendChild(tube_dimensionsElement);
0164 }
0165 
0166 
0167 void PHG4GDMLWriteParamvol::
0168 Cone_dimensionsWrite(xercesc::DOMElement* parametersElement,
0169                      const G4Cons* const cone)
0170 {
0171    xercesc::DOMElement* cone_dimensionsElement = NewElement("cone_dimensions");
0172    cone_dimensionsElement->
0173      setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
0174    cone_dimensionsElement->
0175      setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
0176    cone_dimensionsElement->
0177      setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
0178    cone_dimensionsElement->
0179      setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
0180    cone_dimensionsElement->
0181      setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
0182    cone_dimensionsElement->
0183      setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
0184    cone_dimensionsElement->
0185      setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
0186    cone_dimensionsElement->
0187      setAttributeNode(NewAttribute("aunit","deg"));
0188    cone_dimensionsElement->
0189      setAttributeNode(NewAttribute("lunit","mm"));
0190    parametersElement->appendChild(cone_dimensionsElement);
0191 }
0192 
0193 void PHG4GDMLWriteParamvol::
0194 Sphere_dimensionsWrite(xercesc::DOMElement* parametersElement,
0195                        const G4Sphere* const sphere)
0196 {
0197    xercesc::DOMElement* sphere_dimensionsElement =
0198                         NewElement("sphere_dimensions");
0199    sphere_dimensionsElement->setAttributeNode(NewAttribute("rmin",
0200                              sphere->GetInnerRadius()/mm));
0201    sphere_dimensionsElement->setAttributeNode(NewAttribute("rmax",
0202                              sphere->GetOuterRadius()/mm));
0203    sphere_dimensionsElement->setAttributeNode(NewAttribute("startphi",
0204                              sphere->GetStartPhiAngle()/degree));
0205    sphere_dimensionsElement->setAttributeNode(NewAttribute("deltaphi",
0206                              sphere->GetDeltaPhiAngle()/degree));
0207    sphere_dimensionsElement->setAttributeNode(NewAttribute("starttheta",
0208                              sphere->GetStartThetaAngle()/degree));
0209    sphere_dimensionsElement->setAttributeNode(NewAttribute("deltatheta",
0210                              sphere->GetDeltaThetaAngle()/degree));
0211    sphere_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0212    sphere_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0213    parametersElement->appendChild(sphere_dimensionsElement);
0214 }
0215 
0216 void PHG4GDMLWriteParamvol::
0217 Orb_dimensionsWrite(xercesc::DOMElement* parametersElement,
0218                     const G4Orb* const orb)
0219 {
0220    xercesc::DOMElement* orb_dimensionsElement = NewElement("orb_dimensions");
0221    orb_dimensionsElement->setAttributeNode(NewAttribute("r",
0222                           orb->GetRadius()/mm));
0223    orb_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0224    parametersElement->appendChild(orb_dimensionsElement);
0225 }
0226 
0227 void PHG4GDMLWriteParamvol::
0228 Torus_dimensionsWrite(xercesc::DOMElement* parametersElement,
0229                       const G4Torus* const torus)
0230 {
0231    xercesc::DOMElement* torus_dimensionsElement =
0232                         NewElement("torus_dimensions");
0233    torus_dimensionsElement->
0234      setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
0235    torus_dimensionsElement->
0236      setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
0237    torus_dimensionsElement->
0238      setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
0239    torus_dimensionsElement->
0240      setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
0241    torus_dimensionsElement->
0242      setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
0243    torus_dimensionsElement->
0244      setAttributeNode(NewAttribute("aunit","deg"));
0245    torus_dimensionsElement->
0246      setAttributeNode(NewAttribute("lunit","mm"));
0247    parametersElement->appendChild(torus_dimensionsElement);
0248 }
0249 
0250 void PHG4GDMLWriteParamvol::
0251 Ellipsoid_dimensionsWrite(xercesc::DOMElement* parametersElement,
0252                       const G4Ellipsoid* const ellipsoid)
0253 {
0254    xercesc::DOMElement* ellipsoid_dimensionsElement =
0255                         NewElement("ellipsoid_dimensions");
0256    ellipsoid_dimensionsElement->
0257      setAttributeNode(NewAttribute("ax",ellipsoid->GetSemiAxisMax(0)/mm));
0258    ellipsoid_dimensionsElement->
0259      setAttributeNode(NewAttribute("by",ellipsoid->GetSemiAxisMax(1)/mm));
0260    ellipsoid_dimensionsElement->
0261      setAttributeNode(NewAttribute("cz",ellipsoid->GetSemiAxisMax(2)/mm));
0262    ellipsoid_dimensionsElement->
0263      setAttributeNode(NewAttribute("zcut1",ellipsoid->GetZBottomCut()/mm));
0264    ellipsoid_dimensionsElement->
0265      setAttributeNode(NewAttribute("zcut2",ellipsoid->GetZTopCut()/mm));
0266    ellipsoid_dimensionsElement->
0267      setAttributeNode(NewAttribute("lunit","mm"));
0268    parametersElement->appendChild(ellipsoid_dimensionsElement);
0269 }
0270 
0271 void PHG4GDMLWriteParamvol::
0272 Para_dimensionsWrite(xercesc::DOMElement* parametersElement,
0273                      const G4Para* const para)
0274 {
0275    const G4ThreeVector simaxis = para->GetSymAxis();
0276 
0277    const G4double alpha = std::atan(para->GetTanAlpha());
0278    const G4double theta = std::acos(simaxis.z());
0279    const G4double phi = (simaxis.z() != 1.0)
0280                       ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
0281 
0282    xercesc::DOMElement* para_dimensionsElement = NewElement("para_dimensions");
0283    para_dimensionsElement->
0284      setAttributeNode(NewAttribute("x",2.0*para->GetXHalfLength()/mm));
0285    para_dimensionsElement->
0286      setAttributeNode(NewAttribute("y",2.0*para->GetYHalfLength()/mm));
0287    para_dimensionsElement->
0288      setAttributeNode(NewAttribute("z",2.0*para->GetZHalfLength()/mm));
0289    para_dimensionsElement->
0290      setAttributeNode(NewAttribute("alpha",alpha/degree));
0291    para_dimensionsElement->
0292      setAttributeNode(NewAttribute("theta",theta/degree));
0293    para_dimensionsElement->
0294      setAttributeNode(NewAttribute("phi",phi/degree));
0295    para_dimensionsElement->
0296      setAttributeNode(NewAttribute("aunit","deg"));
0297    para_dimensionsElement->
0298      setAttributeNode(NewAttribute("lunit","mm"));
0299    parametersElement->appendChild(para_dimensionsElement);
0300 }
0301 
0302 void PHG4GDMLWriteParamvol::
0303 Hype_dimensionsWrite(xercesc::DOMElement* parametersElement,
0304                      const G4Hype* const hype)
0305 {
0306    xercesc::DOMElement* hype_dimensionsElement = NewElement("hype_dimensions");
0307    hype_dimensionsElement->
0308      setAttributeNode(NewAttribute("rmin",hype->GetInnerRadius()/mm));
0309    hype_dimensionsElement->
0310      setAttributeNode(NewAttribute("rmax",hype->GetOuterRadius()/mm));
0311    hype_dimensionsElement->
0312      setAttributeNode(NewAttribute("inst",hype->GetInnerStereo()/degree));
0313    hype_dimensionsElement->
0314      setAttributeNode(NewAttribute("outst",hype->GetOuterStereo()/degree));
0315    hype_dimensionsElement->
0316      setAttributeNode(NewAttribute("z",2.0*hype->GetZHalfLength()/mm));
0317    hype_dimensionsElement->
0318      setAttributeNode(NewAttribute("aunit","deg"));
0319    hype_dimensionsElement->
0320      setAttributeNode(NewAttribute("lunit","mm"));
0321    parametersElement->appendChild(hype_dimensionsElement);
0322 }
0323 
0324 void PHG4GDMLWriteParamvol::
0325 Polycone_dimensionsWrite(xercesc::DOMElement* parametersElement,
0326                          const G4Polycone* const pcone)
0327 {
0328    xercesc::DOMElement* pcone_dimensionsElement
0329      = NewElement("polycone_dimensions");
0330 
0331    pcone_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
0332                       pcone->GetOriginalParameters()->Num_z_planes));
0333    pcone_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
0334                       pcone->GetOriginalParameters()->Start_angle/degree));
0335    pcone_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
0336                       pcone->GetOriginalParameters()->Opening_angle/degree));
0337    pcone_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0338    pcone_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0339 
0340    parametersElement->appendChild(pcone_dimensionsElement);
0341    const size_t num_zplanes = pcone->GetOriginalParameters()->Num_z_planes;
0342    const G4double* z_array = pcone->GetOriginalParameters()->Z_values;
0343    const G4double* rmin_array = pcone->GetOriginalParameters()->Rmin;
0344    const G4double* rmax_array = pcone->GetOriginalParameters()->Rmax;
0345 
0346    for (size_t i=0; i<num_zplanes; i++)
0347    {
0348      ZplaneWrite(pcone_dimensionsElement,z_array[i],
0349                  rmin_array[i],rmax_array[i]);
0350    }
0351 }
0352 
0353 void PHG4GDMLWriteParamvol::
0354 Polyhedra_dimensionsWrite(xercesc::DOMElement* parametersElement,
0355                           const G4Polyhedra* const polyhedra)
0356 {
0357    xercesc::DOMElement* polyhedra_dimensionsElement
0358      = NewElement("polyhedra_dimensions");
0359 
0360    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
0361                   polyhedra->GetOriginalParameters()->Num_z_planes));
0362    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numSide",
0363                   polyhedra->GetOriginalParameters()->numSide));
0364    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
0365                   polyhedra->GetOriginalParameters()->Start_angle/degree));
0366    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
0367                   polyhedra->GetOriginalParameters()->Opening_angle/degree));
0368    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0369    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0370 
0371    parametersElement->appendChild(polyhedra_dimensionsElement);
0372    const size_t num_zplanes = polyhedra->GetOriginalParameters()->Num_z_planes;
0373    const G4double* z_array = polyhedra->GetOriginalParameters()->Z_values;
0374    const G4double* rmin_array = polyhedra->GetOriginalParameters()->Rmin;
0375    const G4double* rmax_array = polyhedra->GetOriginalParameters()->Rmax;
0376 
0377    for (size_t i=0; i<num_zplanes; i++)
0378    {
0379      ZplaneWrite(polyhedra_dimensionsElement,z_array[i],
0380                  rmin_array[i],rmax_array[i]);
0381    }
0382 }
0383 
0384 void PHG4GDMLWriteParamvol::
0385 ParametersWrite(xercesc::DOMElement* paramvolElement,
0386                 const G4VPhysicalVolume* const paramvol,const G4int& index)
0387 {
0388    paramvol->GetParameterisation()
0389      ->ComputeTransformation(index, const_cast<G4VPhysicalVolume*>(paramvol));
0390    G4ThreeVector Angles;
0391    G4String name = GenerateName(paramvol->GetName(),paramvol);
0392    std::stringstream os; 
0393    os.precision(15);
0394    os << index;     
0395    G4String sncopie = os.str(); 
0396 
0397    xercesc::DOMElement* parametersElement = NewElement("parameters");
0398    parametersElement->setAttributeNode(NewAttribute("number",index+1));
0399 
0400    PositionWrite(parametersElement, name+sncopie+"_pos",
0401                  paramvol->GetObjectTranslation());
0402    Angles=GetAngles(paramvol->GetObjectRotationValue());
0403    if (Angles.mag2()>DBL_EPSILON)
0404    {
0405      RotationWrite(parametersElement, name+sncopie+"_rot",
0406                    GetAngles(paramvol->GetObjectRotationValue()));
0407    }
0408    paramvolElement->appendChild(parametersElement);
0409 
0410    G4VSolid* solid = paramvol->GetLogicalVolume()->GetSolid();
0411 
0412    if (G4Box* box = dynamic_cast<G4Box*>(solid))
0413    {
0414       paramvol->GetParameterisation()->ComputeDimensions(*box,index,
0415                 const_cast<G4VPhysicalVolume*>(paramvol));
0416       Box_dimensionsWrite(parametersElement,box);
0417    } else
0418    if (G4Trd* trd = dynamic_cast<G4Trd*>(solid))
0419    {
0420       paramvol->GetParameterisation()->ComputeDimensions(*trd,index,
0421                 const_cast<G4VPhysicalVolume*>(paramvol));
0422       Trd_dimensionsWrite(parametersElement,trd);
0423    } else
0424    if (G4Trap* trap = dynamic_cast<G4Trap*>(solid))
0425    {
0426       paramvol->GetParameterisation()->ComputeDimensions(*trap,index,
0427                 const_cast<G4VPhysicalVolume*>(paramvol));
0428       Trap_dimensionsWrite(parametersElement,trap);
0429    } else
0430    if (G4Tubs* tube = dynamic_cast<G4Tubs*>(solid))
0431    {
0432       paramvol->GetParameterisation()->ComputeDimensions(*tube,index,
0433                 const_cast<G4VPhysicalVolume*>(paramvol));
0434       Tube_dimensionsWrite(parametersElement,tube);
0435    } else
0436    if (G4Cons* cone = dynamic_cast<G4Cons*>(solid))
0437    {
0438       paramvol->GetParameterisation()->ComputeDimensions(*cone,index,
0439                 const_cast<G4VPhysicalVolume*>(paramvol));
0440       Cone_dimensionsWrite(parametersElement,cone);
0441    } else
0442    if (G4Sphere* sphere = dynamic_cast<G4Sphere*>(solid))
0443    {
0444       paramvol->GetParameterisation()->ComputeDimensions(*sphere,index,
0445                 const_cast<G4VPhysicalVolume*>(paramvol));
0446       Sphere_dimensionsWrite(parametersElement,sphere);
0447    } else
0448    if (G4Orb* orb = dynamic_cast<G4Orb*>(solid))
0449    {
0450       paramvol->GetParameterisation()->ComputeDimensions(*orb,index,
0451                 const_cast<G4VPhysicalVolume*>(paramvol));
0452       Orb_dimensionsWrite(parametersElement,orb);
0453    } else
0454    if (G4Torus* torus = dynamic_cast<G4Torus*>(solid))
0455    {
0456       paramvol->GetParameterisation()->ComputeDimensions(*torus,index,
0457                 const_cast<G4VPhysicalVolume*>(paramvol));
0458       Torus_dimensionsWrite(parametersElement,torus);
0459    } else
0460    if (G4Ellipsoid* ellipsoid = dynamic_cast<G4Ellipsoid*>(solid))
0461    {
0462       paramvol->GetParameterisation()->ComputeDimensions(*ellipsoid,index,
0463                 const_cast<G4VPhysicalVolume*>(paramvol));
0464       Ellipsoid_dimensionsWrite(parametersElement,ellipsoid);
0465    } else
0466    if (G4Para* para = dynamic_cast<G4Para*>(solid))
0467    {
0468       paramvol->GetParameterisation()->ComputeDimensions(*para,index,
0469                 const_cast<G4VPhysicalVolume*>(paramvol));
0470       Para_dimensionsWrite(parametersElement,para);
0471    } else
0472    if (G4Hype* hype = dynamic_cast<G4Hype*>(solid))
0473    {
0474       paramvol->GetParameterisation()->ComputeDimensions(*hype,index,
0475                 const_cast<G4VPhysicalVolume*>(paramvol));
0476       Hype_dimensionsWrite(parametersElement,hype);
0477    }else
0478    if (G4Polycone* pcone = dynamic_cast<G4Polycone*>(solid))
0479    {
0480       paramvol->GetParameterisation()->ComputeDimensions(*pcone,index,
0481                 const_cast<G4VPhysicalVolume*>(paramvol));
0482       Polycone_dimensionsWrite(parametersElement,pcone);
0483    }else
0484    if (G4Polyhedra* polyhedra = dynamic_cast<G4Polyhedra*>(solid))
0485    {
0486       paramvol->GetParameterisation()->ComputeDimensions(*polyhedra,index,
0487                 const_cast<G4VPhysicalVolume*>(paramvol));
0488       Polyhedra_dimensionsWrite(parametersElement,polyhedra);
0489    }
0490    else
0491    {
0492      G4String error_msg = "Solid '" + solid->GetName()
0493                         + "' cannot be used in parameterised volume!";
0494      G4Exception("PHG4GDMLWriteParamvol::ParametersWrite()",
0495                  "InvalidSetup", FatalException, error_msg);
0496    }
0497 }
0498 
0499 void PHG4GDMLWriteParamvol::
0500 ParamvolWrite(xercesc::DOMElement* volumeElement,
0501               const G4VPhysicalVolume* const paramvol)
0502 {
0503    const G4String volumeref =
0504                   GenerateName(paramvol->GetLogicalVolume()->GetName(),
0505                                paramvol->GetLogicalVolume());
0506    xercesc::DOMElement* paramvolElement = NewElement("paramvol");
0507    paramvolElement->setAttributeNode(NewAttribute("ncopies",
0508                                      paramvol->GetMultiplicity()));
0509    xercesc::DOMElement* volumerefElement = NewElement("volumeref");
0510    volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
0511 
0512    xercesc::DOMElement* algorithmElement =
0513      NewElement("parameterised_position_size");
0514    paramvolElement->appendChild(volumerefElement);
0515    paramvolElement->appendChild(algorithmElement);
0516    ParamvolAlgorithmWrite(algorithmElement,paramvol);
0517    volumeElement->appendChild(paramvolElement);
0518 }
0519 
0520 void PHG4GDMLWriteParamvol::
0521 ParamvolAlgorithmWrite(xercesc::DOMElement* paramvolElement,
0522                        const G4VPhysicalVolume* const paramvol)
0523 {
0524    const G4String volumeref =
0525                   GenerateName(paramvol->GetLogicalVolume()->GetName(),
0526                                paramvol->GetLogicalVolume());
0527   
0528    const G4int parameterCount = paramvol->GetMultiplicity();
0529 
0530    for (G4int i=0; i<parameterCount; i++)
0531    {
0532      ParametersWrite(paramvolElement,paramvol,i);
0533    }
0534 }