source: trunk/source/persistency/gdml/src/G4GDMLWriteParamvol.cc @ 1202

Last change on this file since 1202 was 987, checked in by garnier, 15 years ago

fichiers manquants

File size: 17.2 KB
Line 
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,v 1.23 2008/08/20 08:56:32 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// class G4GDMLParamVol Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLWriteParamvol.hh"
37#include <sstream>
38
39void G4GDMLWriteParamvol::
40Box_dimensionsWrite(xercesc::DOMElement* parametersElement,
41                    const G4Box* const box)
42{
43   xercesc::DOMElement* box_dimensionsElement = NewElement("box_dimensions");
44   box_dimensionsElement->
45     setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
46   box_dimensionsElement->
47     setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
48   box_dimensionsElement->
49     setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
50   box_dimensionsElement->
51     setAttributeNode(NewAttribute("lunit","mm"));
52   parametersElement->appendChild(box_dimensionsElement);
53}
54
55void G4GDMLWriteParamvol::
56Trd_dimensionsWrite(xercesc::DOMElement* parametersElement,
57                    const G4Trd* const trd)
58{
59   xercesc::DOMElement* trd_dimensionsElement = NewElement("trd_dimensions");
60   trd_dimensionsElement->
61     setAttributeNode(NewAttribute("x1",2.0*trd->GetXHalfLength1()/mm));
62   trd_dimensionsElement->
63     setAttributeNode(NewAttribute("x2",2.0*trd->GetXHalfLength2()/mm));
64   trd_dimensionsElement->
65     setAttributeNode(NewAttribute("y1",2.0*trd->GetYHalfLength1()/mm));
66   trd_dimensionsElement->
67     setAttributeNode(NewAttribute("y2",2.0*trd->GetYHalfLength2()/mm));
68   trd_dimensionsElement->
69     setAttributeNode(NewAttribute("z",2.0*trd->GetZHalfLength()/mm));
70   trd_dimensionsElement->
71     setAttributeNode(NewAttribute("lunit","mm"));
72   parametersElement->appendChild(trd_dimensionsElement);
73}
74
75void G4GDMLWriteParamvol::
76Trap_dimensionsWrite(xercesc::DOMElement* parametersElement,
77                     const G4Trap* const trap)
78{
79   const G4ThreeVector simaxis = trap->GetSymAxis();
80   const G4double phi = (simaxis.z() != 1.0)
81                      ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
82   const G4double theta = std::acos(simaxis.z());
83   const G4double alpha1 = std::atan(trap->GetTanAlpha1());
84   const G4double alpha2 = std::atan(trap->GetTanAlpha2());
85
86   xercesc::DOMElement* trap_dimensionsElement = NewElement("trap");
87   trap_dimensionsElement->
88     setAttributeNode(NewAttribute("z",2.0*trap->GetZHalfLength()/mm));
89   trap_dimensionsElement->
90     setAttributeNode(NewAttribute("theta",theta/degree));
91   trap_dimensionsElement->
92     setAttributeNode(NewAttribute("phi",phi/degree));
93   trap_dimensionsElement->
94     setAttributeNode(NewAttribute("y1",2.0*trap->GetYHalfLength1()/mm));
95   trap_dimensionsElement->
96     setAttributeNode(NewAttribute("x1",2.0*trap->GetXHalfLength1()/mm));
97   trap_dimensionsElement->
98     setAttributeNode(NewAttribute("x2",2.0*trap->GetXHalfLength2()/mm));
99   trap_dimensionsElement->
100     setAttributeNode(NewAttribute("alpha1",alpha1/degree));
101   trap_dimensionsElement->
102     setAttributeNode(NewAttribute("y2",2.0*trap->GetYHalfLength2()/mm));
103   trap_dimensionsElement->
104     setAttributeNode(NewAttribute("x3",2.0*trap->GetXHalfLength3()/mm));
105   trap_dimensionsElement->
106     setAttributeNode(NewAttribute("x4",2.0*trap->GetXHalfLength4()/mm));
107   trap_dimensionsElement->
108     setAttributeNode(NewAttribute("alpha2",alpha2/degree));
109   trap_dimensionsElement->
110     setAttributeNode(NewAttribute("aunit","deg"));
111   trap_dimensionsElement->
112     setAttributeNode(NewAttribute("lunit","mm"));
113   parametersElement->appendChild(trap_dimensionsElement);
114}
115
116void G4GDMLWriteParamvol::
117Tube_dimensionsWrite(xercesc::DOMElement* parametersElement,
118                     const G4Tubs* const tube)
119{
120   xercesc::DOMElement* tube_dimensionsElement = NewElement("tube_dimensions");
121   tube_dimensionsElement->
122     setAttributeNode(NewAttribute("InR",tube->GetInnerRadius()/mm));
123   tube_dimensionsElement->
124     setAttributeNode(NewAttribute("OutR",tube->GetOuterRadius()/mm));
125   tube_dimensionsElement->
126     setAttributeNode(NewAttribute("hz",2.0*tube->GetZHalfLength()/mm));
127   tube_dimensionsElement->
128     setAttributeNode(NewAttribute("StartPhi",tube->GetStartPhiAngle()/degree));
129   tube_dimensionsElement->
130     setAttributeNode(NewAttribute("DeltaPhi",tube->GetDeltaPhiAngle()/degree));
131   tube_dimensionsElement->
132     setAttributeNode(NewAttribute("aunit","deg"));
133   tube_dimensionsElement->
134     setAttributeNode(NewAttribute("lunit","mm"));
135   parametersElement->appendChild(tube_dimensionsElement);
136}
137
138
139void G4GDMLWriteParamvol::
140Cone_dimensionsWrite(xercesc::DOMElement* parametersElement,
141                     const G4Cons* const cone)
142{
143   xercesc::DOMElement* cone_dimensionsElement = NewElement("cone_dimensions");
144   cone_dimensionsElement->
145     setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
146   cone_dimensionsElement->
147     setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
148   cone_dimensionsElement->
149     setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
150   cone_dimensionsElement->
151     setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
152   cone_dimensionsElement->
153     setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
154   cone_dimensionsElement->
155     setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
156   cone_dimensionsElement->
157     setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
158   cone_dimensionsElement->
159     setAttributeNode(NewAttribute("aunit","deg"));
160   cone_dimensionsElement->
161     setAttributeNode(NewAttribute("lunit","mm"));
162   parametersElement->appendChild(cone_dimensionsElement);
163}
164
165void G4GDMLWriteParamvol::
166Sphere_dimensionsWrite(xercesc::DOMElement* parametersElement,
167                       const G4Sphere* const sphere)
168{
169   xercesc::DOMElement* sphere_dimensionsElement =
170                        NewElement("sphere_dimensions");
171   sphere_dimensionsElement->setAttributeNode(NewAttribute("rmin",
172                             sphere->GetInsideRadius()/mm));
173   sphere_dimensionsElement->setAttributeNode(NewAttribute("rmax",
174                             sphere->GetOuterRadius()/mm));
175   sphere_dimensionsElement->setAttributeNode(NewAttribute("startphi",
176                             sphere->GetStartPhiAngle()/degree));
177   sphere_dimensionsElement->setAttributeNode(NewAttribute("deltaphi",
178                             sphere->GetDeltaPhiAngle()/degree));
179   sphere_dimensionsElement->setAttributeNode(NewAttribute("starttheta",
180                             sphere->GetStartThetaAngle()/degree));
181   sphere_dimensionsElement->setAttributeNode(NewAttribute("deltatheta",
182                             sphere->GetDeltaThetaAngle()/degree));
183   sphere_dimensionsElement->setAttributeNode(NewAttribute("aunit","deg"));
184   sphere_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
185   parametersElement->appendChild(sphere_dimensionsElement);
186}
187
188void G4GDMLWriteParamvol::
189Orb_dimensionsWrite(xercesc::DOMElement* parametersElement,
190                    const G4Orb* const orb)
191{
192   xercesc::DOMElement* orb_dimensionsElement = NewElement("orb_dimensions");
193   orb_dimensionsElement->setAttributeNode(NewAttribute("r",
194                          orb->GetRadius()/mm));
195   orb_dimensionsElement->setAttributeNode(NewAttribute("lunit","mm"));
196   parametersElement->appendChild(orb_dimensionsElement);
197}
198
199void G4GDMLWriteParamvol::
200Torus_dimensionsWrite(xercesc::DOMElement* parametersElement,
201                      const G4Torus* const torus)
202{
203   xercesc::DOMElement* torus_dimensionsElement =
204                        NewElement("torus_dimensions");
205   torus_dimensionsElement->
206     setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
207   torus_dimensionsElement->
208     setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
209   torus_dimensionsElement->
210     setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
211   torus_dimensionsElement->
212     setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
213   torus_dimensionsElement->
214     setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
215   torus_dimensionsElement->
216     setAttributeNode(NewAttribute("aunit","deg"));
217   torus_dimensionsElement->
218     setAttributeNode(NewAttribute("lunit","mm"));
219   parametersElement->appendChild(torus_dimensionsElement);
220}
221
222void G4GDMLWriteParamvol::
223Para_dimensionsWrite(xercesc::DOMElement* parametersElement,
224                     const G4Para* const para)
225{
226   const G4ThreeVector simaxis = para->GetSymAxis();
227   const G4double alpha = std::atan(para->GetTanAlpha());
228   const G4double theta = std::acos(simaxis.z());
229   const G4double phi = (simaxis.z() != 1.0)
230                      ? (std::atan(simaxis.y()/simaxis.x())) : (0.0);
231
232   xercesc::DOMElement* para_dimensionsElement = NewElement("para_dimensions");
233   para_dimensionsElement->
234     setAttributeNode(NewAttribute("x",2.0*para->GetXHalfLength()/mm));
235   para_dimensionsElement->
236     setAttributeNode(NewAttribute("y",2.0*para->GetYHalfLength()/mm));
237   para_dimensionsElement->
238     setAttributeNode(NewAttribute("z",2.0*para->GetZHalfLength()/mm));
239   para_dimensionsElement->
240     setAttributeNode(NewAttribute("alpha",alpha/degree));
241   para_dimensionsElement->
242     setAttributeNode(NewAttribute("theta",theta/degree));
243   para_dimensionsElement->
244     setAttributeNode(NewAttribute("phi",phi/degree));
245   para_dimensionsElement->
246     setAttributeNode(NewAttribute("aunit","deg"));
247   para_dimensionsElement->
248     setAttributeNode(NewAttribute("lunit","mm"));
249   parametersElement->appendChild(para_dimensionsElement);
250}
251
252void G4GDMLWriteParamvol::
253Hype_dimensionsWrite(xercesc::DOMElement* parametersElement,
254                     const G4Hype* const hype)
255{
256   xercesc::DOMElement* hype_dimensionsElement = NewElement("hype_dimensions");
257   hype_dimensionsElement->
258     setAttributeNode(NewAttribute("rmin",hype->GetInnerRadius()/mm));
259   hype_dimensionsElement->
260     setAttributeNode(NewAttribute("rmax",hype->GetOuterRadius()/mm));
261   hype_dimensionsElement->
262     setAttributeNode(NewAttribute("inst",hype->GetInnerStereo()/degree));
263   hype_dimensionsElement->
264     setAttributeNode(NewAttribute("outst",hype->GetOuterStereo()/degree));
265   hype_dimensionsElement->
266     setAttributeNode(NewAttribute("z",2.0*hype->GetZHalfLength()/mm));
267   hype_dimensionsElement->
268     setAttributeNode(NewAttribute("aunit","deg"));
269   hype_dimensionsElement->
270     setAttributeNode(NewAttribute("lunit","mm"));
271   parametersElement->appendChild(hype_dimensionsElement);
272}
273
274void G4GDMLWriteParamvol::
275ParametersWrite(xercesc::DOMElement* paramvolElement,
276                const G4VPhysicalVolume* const paramvol,const G4int& index)
277{
278   paramvol->GetParameterisation()
279     ->ComputeTransformation(index, const_cast<G4VPhysicalVolume*>(paramvol));
280   G4ThreeVector Angles;
281   G4String name = GenerateName(paramvol->GetName(),paramvol);
282   std::stringstream os; 
283   os.precision(15);
284   os << index;     
285   G4String sncopie = os.str(); 
286
287   xercesc::DOMElement* parametersElement = NewElement("parameters");
288   parametersElement->setAttributeNode(NewAttribute("number",index+1));
289
290   PositionWrite(parametersElement, name+sncopie+"_pos",
291                 paramvol->GetObjectTranslation());
292   Angles=GetAngles(paramvol->GetObjectRotationValue());
293   if (Angles.mag2()>DBL_EPSILON)
294   {
295     RotationWrite(parametersElement, name+sncopie+"_rot",
296                   GetAngles(paramvol->GetObjectRotationValue()));
297   }
298   paramvolElement->appendChild(parametersElement);
299
300   G4VSolid* solid = paramvol->GetLogicalVolume()->GetSolid();
301
302   if (G4Box* box = dynamic_cast<G4Box*>(solid))
303   {
304      paramvol->GetParameterisation()->ComputeDimensions(*box,index,
305                const_cast<G4VPhysicalVolume*>(paramvol));
306      Box_dimensionsWrite(parametersElement,box);
307   } else
308   if (G4Trd* trd = dynamic_cast<G4Trd*>(solid))
309   {
310      paramvol->GetParameterisation()->ComputeDimensions(*trd,index,
311                const_cast<G4VPhysicalVolume*>(paramvol));
312      Trd_dimensionsWrite(parametersElement,trd);
313   } else
314   if (G4Trap* trap = dynamic_cast<G4Trap*>(solid))
315   {
316      paramvol->GetParameterisation()->ComputeDimensions(*trap,index,
317                const_cast<G4VPhysicalVolume*>(paramvol));
318      Trap_dimensionsWrite(parametersElement,trap);
319   } else
320   if (G4Tubs* tube = dynamic_cast<G4Tubs*>(solid))
321   {
322      paramvol->GetParameterisation()->ComputeDimensions(*tube,index,
323                const_cast<G4VPhysicalVolume*>(paramvol));
324      Tube_dimensionsWrite(parametersElement,tube);
325   } else
326   if (G4Cons* cone = dynamic_cast<G4Cons*>(solid))
327   {
328      paramvol->GetParameterisation()->ComputeDimensions(*cone,index,
329                const_cast<G4VPhysicalVolume*>(paramvol));
330      Cone_dimensionsWrite(parametersElement,cone);
331   } else
332   if (G4Sphere* sphere = dynamic_cast<G4Sphere*>(solid))
333   {
334      paramvol->GetParameterisation()->ComputeDimensions(*sphere,index,
335                const_cast<G4VPhysicalVolume*>(paramvol));
336      Sphere_dimensionsWrite(parametersElement,sphere);
337   } else
338   if (G4Orb* orb = dynamic_cast<G4Orb*>(solid))
339   {
340      paramvol->GetParameterisation()->ComputeDimensions(*orb,index,
341                const_cast<G4VPhysicalVolume*>(paramvol));
342      Orb_dimensionsWrite(parametersElement,orb);
343   } else
344   if (G4Torus* torus = dynamic_cast<G4Torus*>(solid))
345   {
346      paramvol->GetParameterisation()->ComputeDimensions(*torus,index,
347                const_cast<G4VPhysicalVolume*>(paramvol));
348      Torus_dimensionsWrite(parametersElement,torus);
349   } else
350   if (G4Para* para = dynamic_cast<G4Para*>(solid))
351   {
352      paramvol->GetParameterisation()->ComputeDimensions(*para,index,
353                const_cast<G4VPhysicalVolume*>(paramvol));
354      Para_dimensionsWrite(parametersElement,para);
355   } else
356   if (G4Hype* hype = dynamic_cast<G4Hype*>(solid))
357   {
358      paramvol->GetParameterisation()->ComputeDimensions(*hype,index,
359                const_cast<G4VPhysicalVolume*>(paramvol));
360      Hype_dimensionsWrite(parametersElement,hype);
361   }
362   else
363   {
364     G4String error_msg = "Solid '" + solid->GetName()
365                        + "' cannot be used in parameterised volume!";
366     G4Exception("G4GDMLWriteParamvol::ParametersWrite()",
367                 "InvalidSetup", FatalException, error_msg);
368   }
369}
370
371void G4GDMLWriteParamvol::
372ParamvolWrite(xercesc::DOMElement* volumeElement,
373              const G4VPhysicalVolume* const paramvol)
374{
375   const G4String volumeref =
376                  GenerateName(paramvol->GetLogicalVolume()->GetName(),
377                               paramvol->GetLogicalVolume());
378   xercesc::DOMElement* paramvolElement = NewElement("paramvol");
379   paramvolElement->setAttributeNode(NewAttribute("ncopies",
380                                     paramvol->GetMultiplicity()));
381   xercesc::DOMElement* volumerefElement = NewElement("volumeref");
382   volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
383
384   xercesc::DOMElement* algorithmElement =
385     NewElement("parameterised_position_size");
386   paramvolElement->appendChild(volumerefElement);
387   paramvolElement->appendChild(algorithmElement);
388   ParamvolAlgorithmWrite(algorithmElement,paramvol);
389   volumeElement->appendChild(paramvolElement);
390}
391
392void G4GDMLWriteParamvol::
393ParamvolAlgorithmWrite(xercesc::DOMElement* paramvolElement,
394                       const G4VPhysicalVolume* const paramvol)
395{
396   const G4String volumeref =
397                  GenerateName(paramvol->GetLogicalVolume()->GetName(),
398                               paramvol->GetLogicalVolume());
399 
400   const G4int parameterCount = paramvol->GetMultiplicity();
401
402   for (G4int i=0; i<parameterCount; i++)
403   {
404     ParametersWrite(paramvolElement,paramvol,i);
405   }
406}
Note: See TracBrowser for help on using the repository browser.