Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteStructure.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$
28//
29// class G4GDMLWriteStructure Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
36
37#include "G4Material.hh"
38#include "G4ReflectedSolid.hh"
39#include "G4DisplacedSolid.hh"
42#include "G4PVDivision.hh"
43#include "G4PVReplica.hh"
44#include "G4OpticalSurface.hh"
47
50{
51}
52
54{
55}
56
57void
58G4GDMLWriteStructure::DivisionvolWrite(xercesc::DOMElement* volumeElement,
59 const G4PVDivision* const divisionvol)
60{
61 EAxis axis = kUndefined;
62 G4int number = 0;
63 G4double width = 0.0;
64 G4double offset = 0.0;
65 G4bool consuming = false;
66
67 divisionvol->GetReplicationData(axis,number,width,offset,consuming);
68 axis = divisionvol->GetDivisionAxis();
69
70 G4String unitString("mm");
71 G4String axisString("kUndefined");
72 if (axis==kXAxis) { axisString = "kXAxis"; }
73 else if (axis==kYAxis) { axisString = "kYAxis"; }
74 else if (axis==kZAxis) { axisString = "kZAxis"; }
75 else if (axis==kRho) { axisString = "kRho"; }
76 else if (axis==kPhi) { axisString = "kPhi"; unitString = "rad"; }
77
78 const G4String name
79 = GenerateName(divisionvol->GetName(),divisionvol);
80 const G4String volumeref
81 = GenerateName(divisionvol->GetLogicalVolume()->GetName(),
82 divisionvol->GetLogicalVolume());
83
84 xercesc::DOMElement* divisionvolElement = NewElement("divisionvol");
85 divisionvolElement->setAttributeNode(NewAttribute("axis",axisString));
86 divisionvolElement->setAttributeNode(NewAttribute("number",number));
87 divisionvolElement->setAttributeNode(NewAttribute("width",width));
88 divisionvolElement->setAttributeNode(NewAttribute("offset",offset));
89 divisionvolElement->setAttributeNode(NewAttribute("unit",unitString));
90 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
91 volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
92 divisionvolElement->appendChild(volumerefElement);
93 volumeElement->appendChild(divisionvolElement);
94}
95
96void G4GDMLWriteStructure::PhysvolWrite(xercesc::DOMElement* volumeElement,
97 const G4VPhysicalVolume* const physvol,
98 const G4Transform3D& T,
99 const G4String& ModuleName)
100{
101 HepGeom::Scale3D scale;
102 HepGeom::Rotate3D rotate;
103 HepGeom::Translate3D translate;
104
105 T.getDecomposition(scale,rotate,translate);
106
107 const G4ThreeVector scl(scale(0,0),scale(1,1),scale(2,2));
108 const G4ThreeVector rot = GetAngles(rotate.getRotation());
109 const G4ThreeVector pos = T.getTranslation();
110
111 const G4String name = GenerateName(physvol->GetName(),physvol);
112
113 xercesc::DOMElement* physvolElement = NewElement("physvol");
114 physvolElement->setAttributeNode(NewAttribute("name",name));
115 volumeElement->appendChild(physvolElement);
116
117 const G4String volumeref
118 = GenerateName(physvol->GetLogicalVolume()->GetName(),
119 physvol->GetLogicalVolume());
120
121 if (ModuleName.empty())
122 {
123 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
124 volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
125 physvolElement->appendChild(volumerefElement);
126 }
127 else
128 {
129 xercesc::DOMElement* fileElement = NewElement("file");
130 fileElement->setAttributeNode(NewAttribute("name",ModuleName));
131 fileElement->setAttributeNode(NewAttribute("volname",volumeref));
132 physvolElement->appendChild(fileElement);
133 }
134
135 if (std::fabs(pos.x()) > kLinearPrecision
136 || std::fabs(pos.y()) > kLinearPrecision
137 || std::fabs(pos.z()) > kLinearPrecision)
138 {
139 PositionWrite(physvolElement,name+"_pos",pos);
140 }
141 if (std::fabs(rot.x()) > kAngularPrecision
142 || std::fabs(rot.y()) > kAngularPrecision
143 || std::fabs(rot.z()) > kAngularPrecision)
144 {
145 RotationWrite(physvolElement,name+"_rot",rot);
146 }
147 if (std::fabs(scl.x()-1.0) > kRelativePrecision
148 || std::fabs(scl.y()-1.0) > kRelativePrecision
149 || std::fabs(scl.z()-1.0) > kRelativePrecision)
150 {
151 ScaleWrite(physvolElement,name+"_scl",scl);
152 }
153}
154
155void G4GDMLWriteStructure::ReplicavolWrite(xercesc::DOMElement* volumeElement,
156 const G4VPhysicalVolume* const replicavol)
157{
158 EAxis axis = kUndefined;
159 G4int number = 0;
160 G4double width = 0.0;
161 G4double offset = 0.0;
162 G4bool consuming = false;
163 G4String unitString("mm");
164
165 replicavol->GetReplicationData(axis,number,width,offset,consuming);
166
167 const G4String volumeref
168 = GenerateName(replicavol->GetLogicalVolume()->GetName(),
169 replicavol->GetLogicalVolume());
170
171 xercesc::DOMElement* replicavolElement = NewElement("replicavol");
172 replicavolElement->setAttributeNode(NewAttribute("number",number));
173 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
174 volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
175 replicavolElement->appendChild(volumerefElement);
176 xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
177 replicavolElement->appendChild(replicateElement);
178
179 xercesc::DOMElement* dirElement = NewElement("direction");
180 if(axis==kXAxis)
181 { dirElement->setAttributeNode(NewAttribute("x","1")); }
182 else if(axis==kYAxis)
183 { dirElement->setAttributeNode(NewAttribute("y","1")); }
184 else if(axis==kZAxis)
185 { dirElement->setAttributeNode(NewAttribute("z","1")); }
186 else if(axis==kRho)
187 { dirElement->setAttributeNode(NewAttribute("rho","1")); }
188 else if(axis==kPhi)
189 { dirElement->setAttributeNode(NewAttribute("phi","1"));
190 unitString="rad"; }
191 replicateElement->appendChild(dirElement);
192
193 xercesc::DOMElement* widthElement = NewElement("width");
194 widthElement->setAttributeNode(NewAttribute("value",width));
195 widthElement->setAttributeNode(NewAttribute("unit",unitString));
196 replicateElement->appendChild(widthElement);
197
198 xercesc::DOMElement* offsetElement = NewElement("offset");
199 offsetElement->setAttributeNode(NewAttribute("value",offset));
200 offsetElement->setAttributeNode(NewAttribute("unit",unitString));
201 replicateElement->appendChild(offsetElement);
202
203 volumeElement->appendChild(replicavolElement);
204}
205
208{
209 if (!bsurf) { return; }
210
211 const G4SurfaceProperty* psurf = bsurf->GetSurfaceProperty();
212
213 // Generate the new element for border-surface
214 //
215 xercesc::DOMElement* borderElement = NewElement("bordersurface");
216 borderElement->setAttributeNode(NewAttribute("name", bsurf->GetName()));
217 borderElement->setAttributeNode(NewAttribute("surfaceproperty",
218 psurf->GetName()));
219
220 const G4String volumeref1 = GenerateName(bsurf->GetVolume1()->GetName(),
221 bsurf->GetVolume1());
222 const G4String volumeref2 = GenerateName(bsurf->GetVolume2()->GetName(),
223 bsurf->GetVolume2());
224 xercesc::DOMElement* volumerefElement1 = NewElement("physvolref");
225 xercesc::DOMElement* volumerefElement2 = NewElement("physvolref");
226 volumerefElement1->setAttributeNode(NewAttribute("ref",volumeref1));
227 volumerefElement2->setAttributeNode(NewAttribute("ref",volumeref2));
228 borderElement->appendChild(volumerefElement1);
229 borderElement->appendChild(volumerefElement2);
230
231 if (FindOpticalSurface(psurf))
232 {
233 const G4OpticalSurface* opsurf =
234 dynamic_cast<const G4OpticalSurface*>(psurf);
235 if (!opsurf)
236 {
237 G4Exception("G4GDMLWriteStructure::BorderSurfaceCache()",
238 "InvalidSetup", FatalException, "No optical surface found!");
239 return;
240 }
242 }
243
244 borderElementVec.push_back(borderElement);
245}
246
248SkinSurfaceCache(const G4LogicalSkinSurface* const ssurf)
249{
250 if (!ssurf) { return; }
251
252 const G4SurfaceProperty* psurf = ssurf->GetSurfaceProperty();
253
254 // Generate the new element for border-surface
255 //
256 xercesc::DOMElement* skinElement = NewElement("skinsurface");
257 skinElement->setAttributeNode(NewAttribute("name", ssurf->GetName()));
258 skinElement->setAttributeNode(NewAttribute("surfaceproperty",
259 psurf->GetName()));
260
261 const G4String volumeref = GenerateName(ssurf->GetLogicalVolume()->GetName(),
262 ssurf->GetLogicalVolume());
263 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
264 volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
265 skinElement->appendChild(volumerefElement);
266
267 if (FindOpticalSurface(psurf))
268 {
269 const G4OpticalSurface* opsurf =
270 dynamic_cast<const G4OpticalSurface*>(psurf);
271 if (!opsurf)
272 {
273 G4Exception("G4GDMLWriteStructure::SkinSurfaceCache()",
274 "InvalidSetup", FatalException, "No optical surface found!");
275 return;
276 }
278 }
279
280 skinElementVec.push_back(skinElement);
281}
282
284{
285 const G4OpticalSurface* osurf = dynamic_cast<const G4OpticalSurface*>(psurf);
286 std::vector<const G4OpticalSurface*>::const_iterator pos;
287 pos = std::find(opt_vec.begin(), opt_vec.end(), osurf);
288 if (pos != opt_vec.end()) { return false; } // item already created!
289
290 opt_vec.push_back(osurf); // cache it for future reference
291 return true;
292}
293
296{
297 G4LogicalSkinSurface* surf = 0;
299 if (nsurf)
300 {
301 const G4LogicalSkinSurfaceTable* stable =
303 std::vector<G4LogicalSkinSurface*>::const_iterator pos;
304 for (pos = stable->begin(); pos != stable->end(); pos++)
305 {
306 if (lvol == (*pos)->GetLogicalVolume())
307 {
308 surf = *pos; break;
309 }
310 }
311 }
312 return surf;
313}
314
317{
318 G4LogicalBorderSurface* surf = 0;
320 if (nsurf)
321 {
322 const G4LogicalBorderSurfaceTable* btable =
324 std::vector<G4LogicalBorderSurface*>::const_iterator pos;
325 for (pos = btable->begin(); pos != btable->end(); pos++)
326 {
327 if (pvol == (*pos)->GetVolume1()) // just the first in the couple
328 { // is enough
329 surf = *pos; break;
330 }
331 }
332 }
333 return surf;
334}
335
337{
338 G4cout << "G4GDML: Writing surfaces..." << G4endl;
339
340 std::vector<xercesc::DOMElement*>::const_iterator pos;
341 for (pos = skinElementVec.begin(); pos != skinElementVec.end(); pos++)
342 {
343 structureElement->appendChild(*pos);
344 }
345 for (pos = borderElementVec.begin(); pos != borderElementVec.end(); pos++)
346 {
347 structureElement->appendChild(*pos);
348 }
349}
350
351void G4GDMLWriteStructure::StructureWrite(xercesc::DOMElement* gdmlElement)
352{
353 G4cout << "G4GDML: Writing structure..." << G4endl;
354
355 structureElement = NewElement("structure");
356 gdmlElement->appendChild(structureElement);
357}
358
360TraverseVolumeTree(const G4LogicalVolume* const volumePtr, const G4int depth)
361{
362 if (VolumeMap().find(volumePtr) != VolumeMap().end())
363 {
364 return VolumeMap()[volumePtr]; // Volume is already processed
365 }
366
367 G4VSolid* solidPtr = volumePtr->GetSolid();
368 G4Transform3D R,invR;
369 G4int trans=0;
370
371 while (true) // Solve possible displacement/reflection
372 { // of the referenced solid!
373 if (trans>maxTransforms)
374 {
375 G4String ErrorMessage = "Referenced solid in volume '"
376 + volumePtr->GetName()
377 + "' was displaced/reflected too many times!";
378 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
379 "InvalidSetup", FatalException, ErrorMessage);
380 }
381
382 if (G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
383 {
384 R = R*refl->GetTransform3D();
385 solidPtr = refl->GetConstituentMovedSolid();
386 trans++;
387 continue;
388 }
389
390 if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
391 {
392 R = R*G4Transform3D(disp->GetObjectRotation(),
393 disp->GetObjectTranslation());
394 solidPtr = disp->GetConstituentMovedSolid();
395 trans++;
396 continue;
397 }
398
399 break;
400 }
401
402 // Only compute the inverse when necessary!
403 //
404 if (trans>0) { invR = R.inverse(); }
405
406 const G4String name
407 = GenerateName(volumePtr->GetName(),volumePtr);
408 const G4String materialref
409 = GenerateName(volumePtr->GetMaterial()->GetName(),
410 volumePtr->GetMaterial());
411 const G4String solidref
412 = GenerateName(solidPtr->GetName(),solidPtr);
413
414 xercesc::DOMElement* volumeElement = NewElement("volume");
415 volumeElement->setAttributeNode(NewAttribute("name",name));
416 xercesc::DOMElement* materialrefElement = NewElement("materialref");
417 materialrefElement->setAttributeNode(NewAttribute("ref",materialref));
418 volumeElement->appendChild(materialrefElement);
419 xercesc::DOMElement* solidrefElement = NewElement("solidref");
420 solidrefElement->setAttributeNode(NewAttribute("ref",solidref));
421 volumeElement->appendChild(solidrefElement);
422
423 const G4int daughterCount = volumePtr->GetNoDaughters();
424
425 for (G4int i=0;i<daughterCount;i++) // Traverse all the children!
426 {
427 const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
428 const G4String ModuleName = Modularize(physvol,depth);
429
430 G4Transform3D daughterR;
431
432 if (ModuleName.empty()) // Check if subtree requested to be
433 { // a separate module!
434 daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(),depth+1);
435 }
436 else
437 {
439 daughterR = writer.Write(ModuleName,physvol->GetLogicalVolume(),
440 SchemaLocation,depth+1);
441 }
442
443 if (const G4PVDivision* const divisionvol
444 = dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
445 {
446 if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
447 {
448 G4String ErrorMessage = "Division volume in '"
449 + name
450 + "' can not be related to reflected solid!";
451 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
452 "InvalidSetup", FatalException, ErrorMessage);
453 }
454 DivisionvolWrite(volumeElement,divisionvol);
455 } else
456 if (physvol->IsParameterised()) // Is it a paramvol?
457 {
458 if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
459 {
460 G4String ErrorMessage = "Parameterised volume in '"
461 + name
462 + "' can not be related to reflected solid!";
463 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
464 "InvalidSetup", FatalException, ErrorMessage);
465 }
466 ParamvolWrite(volumeElement,physvol);
467 } else
468 if (physvol->IsReplicated()) // Is it a replicavol?
469 {
470 if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
471 {
472 G4String ErrorMessage = "Replica volume in '"
473 + name
474 + "' can not be related to reflected solid!";
475 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
476 "InvalidSetup", FatalException, ErrorMessage);
477 }
478 ReplicavolWrite(volumeElement,physvol);
479 }
480 else // Is it a physvol?
481 {
483
484 if (physvol->GetFrameRotation() != 0)
485 {
486 rot = *(physvol->GetFrameRotation());
487 }
488 G4Transform3D P(rot,physvol->GetObjectTranslation());
489 PhysvolWrite(volumeElement,physvol,invR*P*daughterR,ModuleName);
490 }
492 }
493
494 structureElement->appendChild(volumeElement);
495 // Append the volume AFTER traversing the children so that
496 // the order of volumes will be correct!
497
498 VolumeMap()[volumePtr] = R;
499
500 AddExtension(volumeElement, volumePtr);
501 // Add any possible user defined extension attached to a volume
502
503 AddMaterial(volumePtr->GetMaterial());
504 // Add the involved materials and solids!
505
506 AddSolid(solidPtr);
507
509
510 return R;
511}
@ FatalException
std::vector< G4LogicalBorderSurface * > G4LogicalBorderSurfaceTable
std::vector< G4LogicalSkinSurface * > G4LogicalSkinSurfaceTable
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double z() const
double x() const
double y() const
static const G4double kLinearPrecision
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
static const G4double kAngularPrecision
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
G4ThreeVector GetAngles(const G4RotationMatrix &)
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
static const G4double kRelativePrecision
void AddMaterial(const G4Material *const)
virtual void ParamvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
virtual void AddSolid(const G4VSolid *const)
xercesc::DOMElement * solidsElement
static const G4int maxTransforms
void OpticalSurfaceWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
const G4LogicalBorderSurface * GetBorderSurface(const G4VPhysicalVolume *const)
void BorderSurfaceCache(const G4LogicalBorderSurface *const)
std::vector< xercesc::DOMElement * > borderElementVec
std::vector< xercesc::DOMElement * > skinElementVec
void SkinSurfaceCache(const G4LogicalSkinSurface *const)
const G4LogicalSkinSurface * GetSkinSurface(const G4LogicalVolume *const)
void PhysvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const topVol, const G4Transform3D &transform, const G4String &moduleName)
virtual void StructureWrite(xercesc::DOMElement *)
G4Transform3D TraverseVolumeTree(const G4LogicalVolume *const topVol, const G4int depth)
xercesc::DOMElement * structureElement
void DivisionvolWrite(xercesc::DOMElement *, const G4PVDivision *const)
G4bool FindOpticalSurface(const G4SurfaceProperty *)
void ReplicavolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:90
G4String SchemaLocation
Definition: G4GDMLWrite.hh:121
G4String Modularize(const G4VPhysicalVolume *const topvol, const G4int depth)
Definition: G4GDMLWrite.cc:295
virtual void AddExtension(xercesc::DOMElement *, const G4LogicalVolume *const)
Definition: G4GDMLWrite.cc:78
G4Transform3D Write(const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
Definition: G4GDMLWrite.cc:133
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:103
VolumeMapType & VolumeMap()
Definition: G4GDMLWrite.cc:60
static const G4LogicalBorderSurfaceTable * GetSurfaceTable()
const G4VPhysicalVolume * GetVolume2() const
const G4VPhysicalVolume * GetVolume1() const
static size_t GetNumberOfBorderSurfaces()
const G4LogicalVolume * GetLogicalVolume() const
static const G4LogicalSkinSurfaceTable * GetSurfaceTable()
static size_t GetNumberOfSkinSurfaces()
const G4String & GetName() const
G4SurfaceProperty * GetSurfaceProperty() const
G4VSolid * GetSolid() const
G4int GetNoDaughters() const
G4String GetName() const
G4Material * GetMaterial() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
const G4String & GetName() const
Definition: G4Material.hh:177
EAxis GetDivisionAxis() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
const G4String & GetName() const
virtual G4bool IsReplicated() const =0
G4LogicalVolume * GetLogicalVolume() const
const G4RotationMatrix * GetFrameRotation() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
const G4String & GetName() const
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
G4String GetName() const
static const Transform3D Identity
Definition: Transform3D.h:197
CLHEP::HepRotation getRotation() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:174
CLHEP::Hep3Vector getTranslation() const
Transform3D inverse() const
Definition: Transform3D.cc:142
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
@ kUndefined
Definition: geomdefs.hh:54
@ kRho
Definition: geomdefs.hh:54
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41