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

Last change on this file was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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