Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:55

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 void PHG4GDMLWriteParamvol::
0058 Box_dimensionsWrite(xercesc::DOMElement* parametersElement,
0059                     const G4Box* const box)
0060 {
0061    xercesc::DOMElement* box_dimensionsElement = NewElement("box_dimensions");
0062    box_dimensionsElement->
0063      setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
0064    box_dimensionsElement->
0065      setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
0066    box_dimensionsElement->
0067      setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
0068    box_dimensionsElement->
0069      setAttributeNode(NewAttribute("lunit","mm"));
0070    parametersElement->appendChild(box_dimensionsElement);
0071 }
0072 
0073 void PHG4GDMLWriteParamvol::
0074 Trd_dimensionsWrite(xercesc::DOMElement* parametersElement,
0075                     const G4Trd* const trd)
0076 {
0077    xercesc::DOMElement* trd_dimensionsElement = NewElement("trd_dimensions");
0078    trd_dimensionsElement->
0079      setAttributeNode(NewAttribute("x1",2.0*trd->GetXHalfLength1()/mm));
0080    trd_dimensionsElement->
0081      setAttributeNode(NewAttribute("x2",2.0*trd->GetXHalfLength2()/mm));
0082    trd_dimensionsElement->
0083      setAttributeNode(NewAttribute("y1",2.0*trd->GetYHalfLength1()/mm));
0084    trd_dimensionsElement->
0085      setAttributeNode(NewAttribute("y2",2.0*trd->GetYHalfLength2()/mm));
0086    trd_dimensionsElement->
0087      setAttributeNode(NewAttribute("z",2.0*trd->GetZHalfLength()/mm));
0088    trd_dimensionsElement->
0089      setAttributeNode(NewAttribute("lunit","mm"));
0090    parametersElement->appendChild(trd_dimensionsElement);
0091 }
0092 
0093 void PHG4GDMLWriteParamvol::
0094 Trap_dimensionsWrite(xercesc::DOMElement* parametersElement,
0095                      const G4Trap* const trap)
0096 {
0097    const G4ThreeVector simaxis = trap->GetSymAxis();
0098    const G4double phi = (simaxis.z() != 1.0)
0099                       ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
0100    const G4double theta = std::acos(simaxis.z());
0101    const G4double alpha1 = std::atan(trap->GetTanAlpha1());
0102    const G4double alpha2 = std::atan(trap->GetTanAlpha2());
0103 
0104    xercesc::DOMElement* trap_dimensionsElement = NewElement("trap");
0105    trap_dimensionsElement->
0106      setAttributeNode(NewAttribute("z",2.0*trap->GetZHalfLength()/mm));
0107    trap_dimensionsElement->
0108      setAttributeNode(NewAttribute("theta",theta/degree));
0109    trap_dimensionsElement->
0110      setAttributeNode(NewAttribute("phi",phi/degree));
0111    trap_dimensionsElement->
0112      setAttributeNode(NewAttribute("y1",2.0*trap->GetYHalfLength1()/mm));
0113    trap_dimensionsElement->
0114      setAttributeNode(NewAttribute("x1",2.0*trap->GetXHalfLength1()/mm));
0115    trap_dimensionsElement->
0116      setAttributeNode(NewAttribute("x2",2.0*trap->GetXHalfLength2()/mm));
0117    trap_dimensionsElement->
0118      setAttributeNode(NewAttribute("alpha1",alpha1/degree));
0119    trap_dimensionsElement->
0120      setAttributeNode(NewAttribute("y2",2.0*trap->GetYHalfLength2()/mm));
0121    trap_dimensionsElement->
0122      setAttributeNode(NewAttribute("x3",2.0*trap->GetXHalfLength3()/mm));
0123    trap_dimensionsElement->
0124      setAttributeNode(NewAttribute("x4",2.0*trap->GetXHalfLength4()/mm));
0125    trap_dimensionsElement->
0126      setAttributeNode(NewAttribute("alpha2",alpha2/degree));
0127    trap_dimensionsElement->
0128      setAttributeNode(NewAttribute("aunit","deg"));
0129    trap_dimensionsElement->
0130      setAttributeNode(NewAttribute("lunit","mm"));
0131    parametersElement->appendChild(trap_dimensionsElement);
0132 }
0133 
0134 void PHG4GDMLWriteParamvol::
0135 Tube_dimensionsWrite(xercesc::DOMElement* parametersElement,
0136                      const G4Tubs* const tube)
0137 {
0138    xercesc::DOMElement* tube_dimensionsElement = NewElement("tube_dimensions");
0139    tube_dimensionsElement->
0140      setAttributeNode(NewAttribute("InR",tube->GetInnerRadius()/mm));
0141    tube_dimensionsElement->
0142      setAttributeNode(NewAttribute("OutR",tube->GetOuterRadius()/mm));
0143    tube_dimensionsElement->
0144      setAttributeNode(NewAttribute("hz",2.0*tube->GetZHalfLength()/mm));
0145    tube_dimensionsElement->
0146      setAttributeNode(NewAttribute("StartPhi",tube->GetStartPhiAngle()/degree));
0147    tube_dimensionsElement->
0148      setAttributeNode(NewAttribute("DeltaPhi",tube->GetDeltaPhiAngle()/degree));
0149    tube_dimensionsElement->
0150      setAttributeNode(NewAttribute("aunit","deg"));
0151    tube_dimensionsElement->
0152      setAttributeNode(NewAttribute("lunit","mm"));
0153    parametersElement->appendChild(tube_dimensionsElement);
0154 }
0155 
0156 
0157 void PHG4GDMLWriteParamvol::
0158 Cone_dimensionsWrite(xercesc::DOMElement* parametersElement,
0159                      const G4Cons* const cone)
0160 {
0161    xercesc::DOMElement* cone_dimensionsElement = NewElement("cone_dimensions");
0162    cone_dimensionsElement->
0163      setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
0164    cone_dimensionsElement->
0165      setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
0166    cone_dimensionsElement->
0167      setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
0168    cone_dimensionsElement->
0169      setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
0170    cone_dimensionsElement->
0171      setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
0172    cone_dimensionsElement->
0173      setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
0174    cone_dimensionsElement->
0175      setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
0176    cone_dimensionsElement->
0177      setAttributeNode(NewAttribute("aunit","deg"));
0178    cone_dimensionsElement->
0179      setAttributeNode(NewAttribute("lunit","mm"));
0180    parametersElement->appendChild(cone_dimensionsElement);
0181 }
0182 
0183 void PHG4GDMLWriteParamvol::
0184 Sphere_dimensionsWrite(xercesc::DOMElement* parametersElement,
0185                        const G4Sphere* const sphere)
0186 {
0187    xercesc::DOMElement* sphere_dimensionsElement =
0188                         NewElement("sphere_dimensions");
0189    sphere_dimensionsElement->setAttributeNode(NewAttribute("rmin",
0190                              sphere->GetInnerRadius()/mm));
0191    sphere_dimensionsElement->setAttributeNode(NewAttribute("rmax",
0192                              sphere->GetOuterRadius()/mm));
0193    sphere_dimensionsElement->setAttributeNode(NewAttribute("startphi",
0194                              sphere->GetStartPhiAngle()/degree));
0195    sphere_dimensionsElement->setAttributeNode(NewAttribute("deltaphi",
0196                              sphere->GetDeltaPhiAngle()/degree));
0197    sphere_dimensionsElement->setAttributeNode(NewAttribute("starttheta",
0198                              sphere->GetStartThetaAngle()/degree));
0199    sphere_dimensionsElement->setAttributeNode(NewAttribute("deltatheta",
0200                              sphere->GetDeltaThetaAngle()/degree));
0201    sphere_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0202    sphere_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0203    parametersElement->appendChild(sphere_dimensionsElement);
0204 }
0205 
0206 void PHG4GDMLWriteParamvol::
0207 Orb_dimensionsWrite(xercesc::DOMElement* parametersElement,
0208                     const G4Orb* const orb)
0209 {
0210    xercesc::DOMElement* orb_dimensionsElement = NewElement("orb_dimensions");
0211    orb_dimensionsElement->setAttributeNode(NewAttribute("r",
0212                           orb->GetRadius()/mm));
0213    orb_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0214    parametersElement->appendChild(orb_dimensionsElement);
0215 }
0216 
0217 void PHG4GDMLWriteParamvol::
0218 Torus_dimensionsWrite(xercesc::DOMElement* parametersElement,
0219                       const G4Torus* const torus)
0220 {
0221    xercesc::DOMElement* torus_dimensionsElement =
0222                         NewElement("torus_dimensions");
0223    torus_dimensionsElement->
0224      setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
0225    torus_dimensionsElement->
0226      setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
0227    torus_dimensionsElement->
0228      setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
0229    torus_dimensionsElement->
0230      setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
0231    torus_dimensionsElement->
0232      setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
0233    torus_dimensionsElement->
0234      setAttributeNode(NewAttribute("aunit","deg"));
0235    torus_dimensionsElement->
0236      setAttributeNode(NewAttribute("lunit","mm"));
0237    parametersElement->appendChild(torus_dimensionsElement);
0238 }
0239 
0240 void PHG4GDMLWriteParamvol::
0241 Ellipsoid_dimensionsWrite(xercesc::DOMElement* parametersElement,
0242                       const G4Ellipsoid* const ellipsoid)
0243 {
0244    xercesc::DOMElement* ellipsoid_dimensionsElement =
0245                         NewElement("ellipsoid_dimensions");
0246    ellipsoid_dimensionsElement->
0247      setAttributeNode(NewAttribute("ax",ellipsoid->GetSemiAxisMax(0)/mm));
0248    ellipsoid_dimensionsElement->
0249      setAttributeNode(NewAttribute("by",ellipsoid->GetSemiAxisMax(1)/mm));
0250    ellipsoid_dimensionsElement->
0251      setAttributeNode(NewAttribute("cz",ellipsoid->GetSemiAxisMax(2)/mm));
0252    ellipsoid_dimensionsElement->
0253      setAttributeNode(NewAttribute("zcut1",ellipsoid->GetZBottomCut()/mm));
0254    ellipsoid_dimensionsElement->
0255      setAttributeNode(NewAttribute("zcut2",ellipsoid->GetZTopCut()/mm));
0256    ellipsoid_dimensionsElement->
0257      setAttributeNode(NewAttribute("lunit","mm"));
0258    parametersElement->appendChild(ellipsoid_dimensionsElement);
0259 }
0260 
0261 void PHG4GDMLWriteParamvol::
0262 Para_dimensionsWrite(xercesc::DOMElement* parametersElement,
0263                      const G4Para* const para)
0264 {
0265    const G4ThreeVector simaxis = para->GetSymAxis();
0266 
0267    const G4double alpha = std::atan(para->GetTanAlpha());
0268    const G4double theta = std::acos(simaxis.z());
0269    const G4double phi = (simaxis.z() != 1.0)
0270                       ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
0271 
0272    xercesc::DOMElement* para_dimensionsElement = NewElement("para_dimensions");
0273    para_dimensionsElement->
0274      setAttributeNode(NewAttribute("x",2.0*para->GetXHalfLength()/mm));
0275    para_dimensionsElement->
0276      setAttributeNode(NewAttribute("y",2.0*para->GetYHalfLength()/mm));
0277    para_dimensionsElement->
0278      setAttributeNode(NewAttribute("z",2.0*para->GetZHalfLength()/mm));
0279    para_dimensionsElement->
0280      setAttributeNode(NewAttribute("alpha",alpha/degree));
0281    para_dimensionsElement->
0282      setAttributeNode(NewAttribute("theta",theta/degree));
0283    para_dimensionsElement->
0284      setAttributeNode(NewAttribute("phi",phi/degree));
0285    para_dimensionsElement->
0286      setAttributeNode(NewAttribute("aunit","deg"));
0287    para_dimensionsElement->
0288      setAttributeNode(NewAttribute("lunit","mm"));
0289    parametersElement->appendChild(para_dimensionsElement);
0290 }
0291 
0292 void PHG4GDMLWriteParamvol::
0293 Hype_dimensionsWrite(xercesc::DOMElement* parametersElement,
0294                      const G4Hype* const hype)
0295 {
0296    xercesc::DOMElement* hype_dimensionsElement = NewElement("hype_dimensions");
0297    hype_dimensionsElement->
0298      setAttributeNode(NewAttribute("rmin",hype->GetInnerRadius()/mm));
0299    hype_dimensionsElement->
0300      setAttributeNode(NewAttribute("rmax",hype->GetOuterRadius()/mm));
0301    hype_dimensionsElement->
0302      setAttributeNode(NewAttribute("inst",hype->GetInnerStereo()/degree));
0303    hype_dimensionsElement->
0304      setAttributeNode(NewAttribute("outst",hype->GetOuterStereo()/degree));
0305    hype_dimensionsElement->
0306      setAttributeNode(NewAttribute("z",2.0*hype->GetZHalfLength()/mm));
0307    hype_dimensionsElement->
0308      setAttributeNode(NewAttribute("aunit","deg"));
0309    hype_dimensionsElement->
0310      setAttributeNode(NewAttribute("lunit","mm"));
0311    parametersElement->appendChild(hype_dimensionsElement);
0312 }
0313 
0314 void PHG4GDMLWriteParamvol::
0315 Polycone_dimensionsWrite(xercesc::DOMElement* parametersElement,
0316                          const G4Polycone* const pcone)
0317 {
0318    xercesc::DOMElement* pcone_dimensionsElement
0319      = NewElement("polycone_dimensions");
0320 
0321    pcone_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
0322                       pcone->GetOriginalParameters()->Num_z_planes));
0323    pcone_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
0324                       pcone->GetOriginalParameters()->Start_angle/degree));
0325    pcone_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
0326                       pcone->GetOriginalParameters()->Opening_angle/degree));
0327    pcone_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0328    pcone_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0329 
0330    parametersElement->appendChild(pcone_dimensionsElement);
0331    const size_t num_zplanes = pcone->GetOriginalParameters()->Num_z_planes;
0332    const G4double* z_array = pcone->GetOriginalParameters()->Z_values;
0333    const G4double* rmin_array = pcone->GetOriginalParameters()->Rmin;
0334    const G4double* rmax_array = pcone->GetOriginalParameters()->Rmax;
0335 
0336    for (size_t i=0; i<num_zplanes; i++)
0337    {
0338      ZplaneWrite(pcone_dimensionsElement,z_array[i],
0339                  rmin_array[i],rmax_array[i]);
0340    }
0341 }
0342 
0343 void PHG4GDMLWriteParamvol::
0344 Polyhedra_dimensionsWrite(xercesc::DOMElement* parametersElement,
0345                           const G4Polyhedra* const polyhedra)
0346 {
0347    xercesc::DOMElement* polyhedra_dimensionsElement
0348      = NewElement("polyhedra_dimensions");
0349 
0350    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numRZ",
0351                   polyhedra->GetOriginalParameters()->Num_z_planes));
0352    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("numSide",
0353                   polyhedra->GetOriginalParameters()->numSide));
0354    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("startPhi",
0355                   polyhedra->GetOriginalParameters()->Start_angle/degree));
0356    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("openPhi",
0357                   polyhedra->GetOriginalParameters()->Opening_angle/degree));
0358    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
0359    polyhedra_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
0360 
0361    parametersElement->appendChild(polyhedra_dimensionsElement);
0362    const size_t num_zplanes = polyhedra->GetOriginalParameters()->Num_z_planes;
0363    const G4double* z_array = polyhedra->GetOriginalParameters()->Z_values;
0364    const G4double* rmin_array = polyhedra->GetOriginalParameters()->Rmin;
0365    const G4double* rmax_array = polyhedra->GetOriginalParameters()->Rmax;
0366 
0367    for (size_t i=0; i<num_zplanes; i++)
0368    {
0369      ZplaneWrite(polyhedra_dimensionsElement,z_array[i],
0370                  rmin_array[i],rmax_array[i]);
0371    }
0372 }
0373 
0374 void PHG4GDMLWriteParamvol::
0375 ParametersWrite(xercesc::DOMElement* paramvolElement,
0376                 const G4VPhysicalVolume* const paramvol,const G4int& index)
0377 {
0378    paramvol->GetParameterisation()
0379      ->ComputeTransformation(index, const_cast<G4VPhysicalVolume*>(paramvol));
0380    G4ThreeVector Angles;
0381    G4String name = GenerateName(paramvol->GetName(),paramvol);
0382    std::stringstream os; 
0383    os.precision(15);
0384    os << index;     
0385    G4String sncopie = os.str(); 
0386 
0387    xercesc::DOMElement* parametersElement = NewElement("parameters");
0388    parametersElement->setAttributeNode(NewAttribute("number",index+1));
0389 
0390    PositionWrite(parametersElement, name+sncopie+"_pos",
0391                  paramvol->GetObjectTranslation());
0392    Angles=GetAngles(paramvol->GetObjectRotationValue());
0393    if (Angles.mag2()>DBL_EPSILON)
0394    {
0395      RotationWrite(parametersElement, name+sncopie+"_rot",
0396                    GetAngles(paramvol->GetObjectRotationValue()));
0397    }
0398    paramvolElement->appendChild(parametersElement);
0399 
0400    G4VSolid* solid = paramvol->GetLogicalVolume()->GetSolid();
0401 
0402    if (G4Box* box = dynamic_cast<G4Box*>(solid))
0403    {
0404       paramvol->GetParameterisation()->ComputeDimensions(*box,index,
0405                 const_cast<G4VPhysicalVolume*>(paramvol));
0406       Box_dimensionsWrite(parametersElement,box);
0407    } else
0408    if (G4Trd* trd = dynamic_cast<G4Trd*>(solid))
0409    {
0410       paramvol->GetParameterisation()->ComputeDimensions(*trd,index,
0411                 const_cast<G4VPhysicalVolume*>(paramvol));
0412       Trd_dimensionsWrite(parametersElement,trd);
0413    } else
0414    if (G4Trap* trap = dynamic_cast<G4Trap*>(solid))
0415    {
0416       paramvol->GetParameterisation()->ComputeDimensions(*trap,index,
0417                 const_cast<G4VPhysicalVolume*>(paramvol));
0418       Trap_dimensionsWrite(parametersElement,trap);
0419    } else
0420    if (G4Tubs* tube = dynamic_cast<G4Tubs*>(solid))
0421    {
0422       paramvol->GetParameterisation()->ComputeDimensions(*tube,index,
0423                 const_cast<G4VPhysicalVolume*>(paramvol));
0424       Tube_dimensionsWrite(parametersElement,tube);
0425    } else
0426    if (G4Cons* cone = dynamic_cast<G4Cons*>(solid))
0427    {
0428       paramvol->GetParameterisation()->ComputeDimensions(*cone,index,
0429                 const_cast<G4VPhysicalVolume*>(paramvol));
0430       Cone_dimensionsWrite(parametersElement,cone);
0431    } else
0432    if (G4Sphere* sphere = dynamic_cast<G4Sphere*>(solid))
0433    {
0434       paramvol->GetParameterisation()->ComputeDimensions(*sphere,index,
0435                 const_cast<G4VPhysicalVolume*>(paramvol));
0436       Sphere_dimensionsWrite(parametersElement,sphere);
0437    } else
0438    if (G4Orb* orb = dynamic_cast<G4Orb*>(solid))
0439    {
0440       paramvol->GetParameterisation()->ComputeDimensions(*orb,index,
0441                 const_cast<G4VPhysicalVolume*>(paramvol));
0442       Orb_dimensionsWrite(parametersElement,orb);
0443    } else
0444    if (G4Torus* torus = dynamic_cast<G4Torus*>(solid))
0445    {
0446       paramvol->GetParameterisation()->ComputeDimensions(*torus,index,
0447                 const_cast<G4VPhysicalVolume*>(paramvol));
0448       Torus_dimensionsWrite(parametersElement,torus);
0449    } else
0450    if (G4Ellipsoid* ellipsoid = dynamic_cast<G4Ellipsoid*>(solid))
0451    {
0452       paramvol->GetParameterisation()->ComputeDimensions(*ellipsoid,index,
0453                 const_cast<G4VPhysicalVolume*>(paramvol));
0454       Ellipsoid_dimensionsWrite(parametersElement,ellipsoid);
0455    } else
0456    if (G4Para* para = dynamic_cast<G4Para*>(solid))
0457    {
0458       paramvol->GetParameterisation()->ComputeDimensions(*para,index,
0459                 const_cast<G4VPhysicalVolume*>(paramvol));
0460       Para_dimensionsWrite(parametersElement,para);
0461    } else
0462    if (G4Hype* hype = dynamic_cast<G4Hype*>(solid))
0463    {
0464       paramvol->GetParameterisation()->ComputeDimensions(*hype,index,
0465                 const_cast<G4VPhysicalVolume*>(paramvol));
0466       Hype_dimensionsWrite(parametersElement,hype);
0467    }else
0468    if (G4Polycone* pcone = dynamic_cast<G4Polycone*>(solid))
0469    {
0470       paramvol->GetParameterisation()->ComputeDimensions(*pcone,index,
0471                 const_cast<G4VPhysicalVolume*>(paramvol));
0472       Polycone_dimensionsWrite(parametersElement,pcone);
0473    }else
0474    if (G4Polyhedra* polyhedra = dynamic_cast<G4Polyhedra*>(solid))
0475    {
0476       paramvol->GetParameterisation()->ComputeDimensions(*polyhedra,index,
0477                 const_cast<G4VPhysicalVolume*>(paramvol));
0478       Polyhedra_dimensionsWrite(parametersElement,polyhedra);
0479    }
0480    else
0481    {
0482      G4String error_msg = "Solid '" + solid->GetName()
0483                         + "' cannot be used in parameterised volume!";
0484      G4Exception("PHG4GDMLWriteParamvol::ParametersWrite()",
0485                  "InvalidSetup", FatalException, error_msg);
0486    }
0487 }
0488 
0489 void PHG4GDMLWriteParamvol::
0490 ParamvolWrite(xercesc::DOMElement* volumeElement,
0491               const G4VPhysicalVolume* const paramvol)
0492 {
0493    const G4String volumeref =
0494                   GenerateName(paramvol->GetLogicalVolume()->GetName(),
0495                                paramvol->GetLogicalVolume());
0496    xercesc::DOMElement* paramvolElement = NewElement("paramvol");
0497    paramvolElement->setAttributeNode(NewAttribute("ncopies",
0498                                      paramvol->GetMultiplicity()));
0499    xercesc::DOMElement* volumerefElement = NewElement("volumeref");
0500    volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
0501 
0502    xercesc::DOMElement* algorithmElement =
0503      NewElement("parameterised_position_size");
0504    paramvolElement->appendChild(volumerefElement);
0505    paramvolElement->appendChild(algorithmElement);
0506    ParamvolAlgorithmWrite(algorithmElement,paramvol);
0507    volumeElement->appendChild(paramvolElement);
0508 }
0509 
0510 void PHG4GDMLWriteParamvol::
0511 ParamvolAlgorithmWrite(xercesc::DOMElement* paramvolElement,
0512                        const G4VPhysicalVolume* const paramvol)
0513 {
0514    const G4String volumeref =
0515                   GenerateName(paramvol->GetLogicalVolume()->GetName(),
0516                                paramvol->GetLogicalVolume());
0517   
0518    const G4int parameterCount = paramvol->GetMultiplicity();
0519 
0520    for (G4int i=0; i<parameterCount; i++)
0521    {
0522      ParametersWrite(paramvolElement,paramvol,i);
0523    }
0524 }