Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GDMLWriteParamvol.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4GDMLWriteParamvol.cc 77913 2013-11-29 10:59:07Z gcosmo $
28 //
29 // class G4GDMLParamVol Implementation
30 //
31 // Original author: Zoltan Torzsok, November 2007
32 //
33 // --------------------------------------------------------------------
34 
35 #include "G4GDMLWriteParamvol.hh"
36 #include "G4GDMLWriteSolids.hh"
37 
38 #include "G4SystemOfUnits.hh"
39 #include "G4Box.hh"
40 #include "G4Trd.hh"
41 #include "G4Trap.hh"
42 #include "G4Tubs.hh"
43 #include "G4Cons.hh"
44 #include "G4Sphere.hh"
45 #include "G4Orb.hh"
46 #include "G4Torus.hh"
47 #include "G4Ellipsoid.hh"
48 #include "G4Para.hh"
49 #include "G4Hype.hh"
50 #include "G4Polycone.hh"
51 #include "G4Polyhedra.hh"
52 #include "G4LogicalVolume.hh"
53 #include "G4VPhysicalVolume.hh"
54 #include "G4PVParameterised.hh"
55 #include "G4VPVParameterisation.hh"
56 
59 {
60 }
61 
64 {
65 }
66 
68 Box_dimensionsWrite(xercesc::DOMElement* parametersElement,
69  const G4Box* const box)
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 }
82 
84 Trd_dimensionsWrite(xercesc::DOMElement* parametersElement,
85  const G4Trd* const trd)
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 }
102 
104 Trap_dimensionsWrite(xercesc::DOMElement* parametersElement,
105  const G4Trap* const trap)
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 }
143 
145 Tube_dimensionsWrite(xercesc::DOMElement* parametersElement,
146  const G4Tubs* const tube)
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 }
165 
166 
168 Cone_dimensionsWrite(xercesc::DOMElement* parametersElement,
169  const G4Cons* const cone)
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 }
192 
194 Sphere_dimensionsWrite(xercesc::DOMElement* parametersElement,
195  const G4Sphere* const sphere)
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 }
215 
217 Orb_dimensionsWrite(xercesc::DOMElement* parametersElement,
218  const G4Orb* const orb)
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 }
226 
228 Torus_dimensionsWrite(xercesc::DOMElement* parametersElement,
229  const G4Torus* const torus)
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 }
249 
251 Ellipsoid_dimensionsWrite(xercesc::DOMElement* parametersElement,
252  const G4Ellipsoid* const ellipsoid)
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 }
270 
272 Para_dimensionsWrite(xercesc::DOMElement* parametersElement,
273  const G4Para* const para)
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 }
301 
303 Hype_dimensionsWrite(xercesc::DOMElement* parametersElement,
304  const G4Hype* const hype)
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 }
323 
325 Polycone_dimensionsWrite(xercesc::DOMElement* parametersElement,
326  const G4Polycone* const pcone)
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 }
352 
354 Polyhedra_dimensionsWrite(xercesc::DOMElement* parametersElement,
355  const G4Polyhedra* const polyhedra)
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 }
383 
385 ParametersWrite(xercesc::DOMElement* paramvolElement,
386  const G4VPhysicalVolume* const paramvol,const G4int& index)
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 }
498 
500 ParamvolWrite(xercesc::DOMElement* volumeElement,
501  const G4VPhysicalVolume* const paramvol)
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 }
519 
521 ParamvolAlgorithmWrite(xercesc::DOMElement* paramvolElement,
522  const G4VPhysicalVolume* const paramvol)
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 }
G4ThreeVector GetSymAxis() const
G4double GetSPhi() const
G4String GetName() const
G4double GetXHalfLength4() 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)
G4double GetXHalfLength() const
G4double GetOuterStereo() const
Definition: G4Para.hh:77
void Orb_dimensionsWrite(xercesc::DOMElement *, const G4Orb *const)
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetYHalfLength2() const
G4double GetYHalfLength1() const
G4double GetInnerStereo() const
G4double GetInnerRadius() const
double x() const
Definition: G4Box.hh:64
G4double GetZTopCut() const
G4double GetZHalfLength() const
void Hype_dimensionsWrite(xercesc::DOMElement *, const G4Hype *const)
Definition: G4Tubs.hh:85
G4double GetRmax() const
G4VSolid * GetSolid() const
G4double GetZBottomCut() const
void Polyhedra_dimensionsWrite(xercesc::DOMElement *, const G4Polyhedra *const)
G4double GetOuterRadiusMinusZ() const
virtual void ParamvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
Definition: G4Trd.hh:72
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDeltaPhiAngle() const
int G4int
Definition: G4Types.hh:78
G4double GetZHalfLength() const
G4double GetSemiAxisMax(G4int i) const
G4double GetXHalfLength2() const
void Tube_dimensionsWrite(xercesc::DOMElement *, const G4Tubs *const)
double z() const
G4double GetRtor() const
G4double GetTanAlpha2() const
G4double GetXHalfLength2() const
G4double GetRadius() const
G4double GetXHalfLength1() const
virtual void ParamvolAlgorithmWrite(xercesc::DOMElement *paramvolElement, const G4VPhysicalVolume *const paramvol)
void Cone_dimensionsWrite(xercesc::DOMElement *, const G4Cons *const)
G4double GetDeltaPhiAngle() const
G4double GetStartThetaAngle() const
const G4String & GetName() const
G4double GetRmin() const
G4double GetXHalfLength3() const
static constexpr double degree
Definition: G4SIunits.hh:144
Definition: G4Hype.hh:67
Definition: G4Cons.hh:83
virtual G4VPVParameterisation * GetParameterisation() const =0
G4RotationMatrix GetObjectRotationValue() const
#define DBL_EPSILON
Definition: templates.hh:87
G4double GetYHalfLength() const
G4double GetStartPhiAngle() const
G4double GetYHalfLength2() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetInnerRadiusPlusZ() const
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
G4double GetDPhi() const
G4double GetInnerRadius() const
G4PolyconeHistorical * GetOriginalParameters() const
Definition: G4Orb.hh:61
G4double GetXHalfLength() const
G4double GetInnerRadius() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
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
G4double GetTanAlpha() 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 &)
double y() const
void Sphere_dimensionsWrite(xercesc::DOMElement *, const G4Sphere *const)
virtual G4int GetMultiplicity() const
G4PolyhedraHistorical * GetOriginalParameters() const
G4double GetZHalfLength() const
G4double GetZHalfLength() const
G4ThreeVector GetSymAxis() const
G4double GetZHalfLength() const
void ParametersWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const, const G4int &)
G4double GetXHalfLength1() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetYHalfLength1() const
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
G4ThreeVector GetObjectTranslation() const
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
G4double GetInnerRadiusMinusZ() const
G4double GetTanAlpha1() const
static const G4double alpha
G4double GetDeltaThetaAngle() const
G4double GetOuterRadiusPlusZ() const
G4double GetZHalfLength() const
G4double GetOuterRadius() const
G4double GetDeltaPhiAngle() const
G4double GetYHalfLength() const