source: trunk/source/persistency/gdml/src/G4GDMLWriteStructure.cc @ 1102

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

fichiers manquants

File size: 13.1 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: G4GDMLWriteStructure.cc,v 1.74 2008/11/13 16:48:19 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// class G4GDMLWriteStructure Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLWriteStructure.hh"
37
38void
39G4GDMLWriteStructure::DivisionvolWrite(xercesc::DOMElement* volumeElement,
40                                       const G4PVDivision* const divisionvol)
41{
42   EAxis axis = kUndefined;
43   G4int number = 0;
44   G4double width = 0.0;
45   G4double offset = 0.0;
46   G4bool consuming = false;
47
48   divisionvol->GetReplicationData(axis,number,width,offset,consuming);
49
50   G4String unitString("mm");
51   G4String axisString("kUndefined");
52   if (axis==kXAxis) { axisString = "kXAxis"; } else
53   if (axis==kYAxis) { axisString = "kYAxis"; } else
54   if (axis==kZAxis) { axisString = "kZAxis"; } else
55   if (axis==kRho) { axisString = "kRho";     } else
56   if (axis==kPhi) { axisString = "kPhi"; unitString = "degree"; }
57
58   const G4String name
59         = GenerateName(divisionvol->GetName(),divisionvol);
60   const G4String volumeref
61         = GenerateName(divisionvol->GetLogicalVolume()->GetName(),
62                        divisionvol->GetLogicalVolume());
63
64   xercesc::DOMElement* divisionvolElement = NewElement("divisionvol");
65   divisionvolElement->setAttributeNode(NewAttribute("axis",axisString));
66   divisionvolElement->setAttributeNode(NewAttribute("number",number));
67   divisionvolElement->setAttributeNode(NewAttribute("width",width));
68   divisionvolElement->setAttributeNode(NewAttribute("offset",offset));
69   divisionvolElement->setAttributeNode(NewAttribute("unit",unitString));
70   xercesc::DOMElement* volumerefElement = NewElement("volumeref");
71   volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
72   divisionvolElement->appendChild(volumerefElement);
73   volumeElement->appendChild(divisionvolElement);
74}
75
76void G4GDMLWriteStructure::PhysvolWrite(xercesc::DOMElement* volumeElement,
77                                        const G4VPhysicalVolume* const physvol,
78                                        const G4Transform3D& T,
79                                        const G4String& ModuleName)
80{
81   HepGeom::Scale3D scale;
82   HepGeom::Rotate3D rotate;
83   HepGeom::Translate3D translate;
84
85   T.getDecomposition(scale,rotate,translate);
86
87   const G4ThreeVector scl(scale(0,0),scale(1,1),scale(2,2));
88   const G4ThreeVector rot = GetAngles(rotate.getRotation());
89   const G4ThreeVector pos = T.getTranslation();
90
91   const G4String name = GenerateName(physvol->GetName(),physvol);
92
93   xercesc::DOMElement* physvolElement = NewElement("physvol");
94   physvolElement->setAttributeNode(NewAttribute("name",name));
95   volumeElement->appendChild(physvolElement);
96
97   const G4String volumeref
98         = GenerateName(physvol->GetLogicalVolume()->GetName(),
99                        physvol->GetLogicalVolume());
100
101   if (ModuleName.empty())
102   {
103      xercesc::DOMElement* volumerefElement = NewElement("volumeref");
104      volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
105      physvolElement->appendChild(volumerefElement);
106   }
107   else
108   {
109      xercesc::DOMElement* fileElement = NewElement("file");
110      fileElement->setAttributeNode(NewAttribute("name",ModuleName));
111      fileElement->setAttributeNode(NewAttribute("volname",volumeref));
112      physvolElement->appendChild(fileElement);
113   }
114
115   if (std::fabs(pos.x()) > kLinearPrecision
116    || std::fabs(pos.y()) > kLinearPrecision
117    || std::fabs(pos.z()) > kLinearPrecision)
118   {
119     PositionWrite(physvolElement,name+"_pos",pos);
120   }
121   if (std::fabs(rot.x()) > kAngularPrecision
122    || std::fabs(rot.y()) > kAngularPrecision
123    || std::fabs(rot.z()) > kAngularPrecision)
124   {
125     RotationWrite(physvolElement,name+"_rot",rot);
126   }
127   if (std::fabs(scl.x()-1.0) > kRelativePrecision
128    || std::fabs(scl.y()-1.0) > kRelativePrecision
129    || std::fabs(scl.z()-1.0) > kRelativePrecision)
130   {
131     ScaleWrite(physvolElement,name+"_scl",scl);
132   }
133}
134
135void G4GDMLWriteStructure::ReplicavolWrite(xercesc::DOMElement* volumeElement,
136                                     const G4VPhysicalVolume* const replicavol)
137{
138   EAxis axis = kUndefined;
139   G4int number = 0;
140   G4double width = 0.0;
141   G4double offset = 0.0;
142   G4bool consuming = false;
143   G4String unitString("mm");
144
145   replicavol->GetReplicationData(axis,number,width,offset,consuming);
146
147   const G4String volumeref
148         = GenerateName(replicavol->GetLogicalVolume()->GetName(),
149                        replicavol->GetLogicalVolume());
150
151   xercesc::DOMElement* replicavolElement = NewElement("replicavol");
152   replicavolElement->setAttributeNode(NewAttribute("number",number));
153   xercesc::DOMElement* volumerefElement = NewElement("volumeref");
154   volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
155   replicavolElement->appendChild(volumerefElement);
156
157   xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
158   replicavolElement->appendChild(replicateElement);
159
160   xercesc::DOMElement* dirElement = NewElement("direction");
161   if(axis==kXAxis)dirElement->setAttributeNode(NewAttribute("x","1"));
162   if(axis==kYAxis)dirElement->setAttributeNode(NewAttribute("y","1"));
163   if(axis==kZAxis)dirElement->setAttributeNode(NewAttribute("z","1"));
164   if(axis==kRho)dirElement->setAttributeNode(NewAttribute("rho","1"));
165   if(axis==kPhi)dirElement->setAttributeNode(NewAttribute("phi","1"));
166   replicateElement->appendChild(dirElement);
167
168   xercesc::DOMElement* widthElement = NewElement("width");
169   widthElement->setAttributeNode(NewAttribute("value",width));
170   widthElement->setAttributeNode(NewAttribute("unit",unitString));
171   replicateElement->appendChild(widthElement);
172
173   xercesc::DOMElement* offsetElement = NewElement("offset");
174   offsetElement->setAttributeNode(NewAttribute("value",offset));
175   offsetElement->setAttributeNode(NewAttribute("unit",unitString));
176   replicateElement->appendChild(offsetElement);
177
178   volumeElement->appendChild(replicavolElement);
179}
180
181void G4GDMLWriteStructure::StructureWrite(xercesc::DOMElement* gdmlElement)
182{
183   G4cout << "G4GDML: Writing structure..." << G4endl;
184
185   structureElement = NewElement("structure");
186   gdmlElement->appendChild(structureElement);
187}
188
189G4Transform3D G4GDMLWriteStructure::
190TraverseVolumeTree(const G4LogicalVolume* const volumePtr, const G4int depth)
191{
192   if (VolumeMap().find(volumePtr) != VolumeMap().end())
193   {
194     return VolumeMap()[volumePtr]; // Volume is already processed
195   }
196
197   G4VSolid* solidPtr = volumePtr->GetSolid();
198   G4Transform3D R,invR;
199   G4int reflected = 0;
200
201   while (true) // Solve possible displacement/reflection
202   {            // of the referenced solid!
203      if (reflected>maxReflections)
204      {
205        G4String ErrorMessage = "Referenced solid in volume '"
206                              + volumePtr->GetName()
207                              + "' was displaced/reflected too many times!";
208        G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
209                    "InvalidSetup", FatalException, ErrorMessage);
210      }
211
212      if (G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
213      {
214         R = R*refl->GetTransform3D();
215         solidPtr = refl->GetConstituentMovedSolid();
216         reflected++;
217         continue;
218      }
219
220      if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
221      {
222         R = R*G4Transform3D(disp->GetObjectRotation(),
223                             disp->GetObjectTranslation());
224         solidPtr = disp->GetConstituentMovedSolid();
225         reflected++;
226         continue;
227      }
228
229      break;
230   }
231
232   if (reflected>0) { invR = R.inverse(); }
233     // Only compute the inverse when necessary!
234
235   const G4String name
236         = GenerateName(volumePtr->GetName(),volumePtr);
237   const G4String materialref
238         = GenerateName(volumePtr->GetMaterial()->GetName(),
239                        volumePtr->GetMaterial());
240   const G4String solidref
241         = GenerateName(solidPtr->GetName(),solidPtr);
242
243   xercesc::DOMElement* volumeElement = NewElement("volume");
244   volumeElement->setAttributeNode(NewAttribute("name",name));
245   xercesc::DOMElement* materialrefElement = NewElement("materialref");
246   materialrefElement->setAttributeNode(NewAttribute("ref",materialref));
247   volumeElement->appendChild(materialrefElement);
248   xercesc::DOMElement* solidrefElement = NewElement("solidref");
249   solidrefElement->setAttributeNode(NewAttribute("ref",solidref));
250   volumeElement->appendChild(solidrefElement);
251
252   const G4int daughterCount = volumePtr->GetNoDaughters();
253
254   for (G4int i=0;i<daughterCount;i++)   // Traverse all the children!
255   {
256      const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
257      const G4String ModuleName = Modularize(physvol,depth);
258
259      G4Transform3D daughterR;
260
261      if (ModuleName.empty())   // Check if subtree requested to be
262      {                         // a separate module!
263         daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(),depth+1);
264      }
265      else
266      {   
267         G4GDMLWriteStructure writer;
268         daughterR = writer.Write(ModuleName,physvol->GetLogicalVolume(),
269                                  SchemaLocation,depth+1);
270      }
271
272      if (const G4PVDivision* const divisionvol
273         = dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
274      {
275         if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
276         {
277            G4String ErrorMessage = "Division volume in '"
278                                  + name
279                                  + "' can not be related to reflected solid!";
280            G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
281                        "InvalidSetup", FatalException, ErrorMessage);
282         }
283         DivisionvolWrite(volumeElement,divisionvol); 
284      } else 
285      if (physvol->IsParameterised())   // Is it a paramvol?
286      {
287         if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
288         {
289            G4String ErrorMessage = "Parameterised volume in '"
290                                  + name
291                                  + "' can not be related to reflected solid!";
292            G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
293                        "InvalidSetup", FatalException, ErrorMessage);
294         }
295         ParamvolWrite(volumeElement,physvol);
296      } else
297      if (physvol->IsReplicated())   // Is it a replicavol?
298      {
299         if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
300         {
301            G4String ErrorMessage = "Replica volume in '"
302                                  + name
303                                  + "' can not be related to reflected solid!";
304            G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
305                        "InvalidSetup", FatalException, ErrorMessage);
306         }
307         ReplicavolWrite(volumeElement,physvol); 
308      }
309      else   // Is it a physvol?
310      {
311         G4RotationMatrix rot;
312
313         if (physvol->GetFrameRotation() != 0)
314         {
315           rot = *(physvol->GetFrameRotation());
316         }
317         G4Transform3D P(rot,physvol->GetObjectTranslation());
318         PhysvolWrite(volumeElement,physvol,invR*P*daughterR,ModuleName);
319      }
320   }
321
322   structureElement->appendChild(volumeElement);
323     // Append the volume AFTER traversing the children so that
324     // the order of volumes will be correct!
325
326   VolumeMap()[volumePtr] = R;
327
328   G4GDMLWriteMaterials::AddMaterial(volumePtr->GetMaterial());
329     // Add the involved materials and solids!
330
331   G4GDMLWriteSolids::AddSolid(solidPtr);
332
333   return R;
334}
Note: See TracBrowser for help on using the repository browser.