Geant4-11
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// G4GDMLWriteStructure implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
32#include "G4GDMLEvaluator.hh"
33
34#include "G4Material.hh"
35#include "G4ReflectedSolid.hh"
36#include "G4DisplacedSolid.hh"
40#include "G4PVDivision.hh"
41#include "G4PVReplica.hh"
42#include "G4Region.hh"
43#include "G4OpticalSurface.hh"
46
47#include "G4ProductionCuts.hh"
49#include "G4Gamma.hh"
50#include "G4Electron.hh"
51#include "G4Positron.hh"
52#include "G4Proton.hh"
53
55#include "G4AssemblyStore.hh"
56#include "G4AssemblyVolume.hh"
57
58G4int G4GDMLWriteStructure::levelNo = 0; // Counter for level being exported
59
60// --------------------------------------------------------------------
63 , maxLevel(INT_MAX)
64{
66}
67
68// --------------------------------------------------------------------
70{
71}
72
73// --------------------------------------------------------------------
75 xercesc::DOMElement* volumeElement, const G4PVDivision* const divisionvol)
76{
77 EAxis axis = kUndefined;
78 G4int number = 0;
79 G4double width = 0.0;
80 G4double offset = 0.0;
81 G4bool consuming = false;
82
83 divisionvol->GetReplicationData(axis, number, width, offset, consuming);
84 axis = divisionvol->GetDivisionAxis();
85
86 G4String unitString("mm");
87 G4String axisString("kUndefined");
88 if(axis == kXAxis)
89 {
90 axisString = "kXAxis";
91 }
92 else if(axis == kYAxis)
93 {
94 axisString = "kYAxis";
95 }
96 else if(axis == kZAxis)
97 {
98 axisString = "kZAxis";
99 }
100 else if(axis == kRho)
101 {
102 axisString = "kRho";
103 }
104 else if(axis == kPhi)
105 {
106 axisString = "kPhi";
107 unitString = "rad";
108 }
109
110 const G4String name = GenerateName(divisionvol->GetName(), divisionvol);
111 const G4String volumeref =
112 GenerateName(divisionvol->GetLogicalVolume()->GetName(),
113 divisionvol->GetLogicalVolume());
114
115 xercesc::DOMElement* divisionvolElement = NewElement("divisionvol");
116 divisionvolElement->setAttributeNode(NewAttribute("axis", axisString));
117 divisionvolElement->setAttributeNode(NewAttribute("number", number));
118 divisionvolElement->setAttributeNode(NewAttribute("width", width));
119 divisionvolElement->setAttributeNode(NewAttribute("offset", offset));
120 divisionvolElement->setAttributeNode(NewAttribute("unit", unitString));
121 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
122 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
123 divisionvolElement->appendChild(volumerefElement);
124 volumeElement->appendChild(divisionvolElement);
125}
126
127// --------------------------------------------------------------------
128void G4GDMLWriteStructure::PhysvolWrite(xercesc::DOMElement* volumeElement,
129 const G4VPhysicalVolume* const physvol,
130 const G4Transform3D& T,
131 const G4String& ModuleName)
132{
133 HepGeom::Scale3D scale;
134 HepGeom::Rotate3D rotate;
135 HepGeom::Translate3D translate;
136
137 T.getDecomposition(scale, rotate, translate);
138
139 const G4ThreeVector scl(scale(0, 0), scale(1, 1), scale(2, 2));
140 const G4ThreeVector rot = GetAngles(rotate.getRotation());
141 const G4ThreeVector pos = T.getTranslation();
142
143 const G4String name = GenerateName(physvol->GetName(), physvol);
144 const G4int copynumber = physvol->GetCopyNo();
145
146 xercesc::DOMElement* physvolElement = NewElement("physvol");
147 physvolElement->setAttributeNode(NewAttribute("name", name));
148 if(copynumber)
149 {
150 physvolElement->setAttributeNode(NewAttribute("copynumber", copynumber));
151 }
152
153 volumeElement->appendChild(physvolElement);
154
155 G4LogicalVolume* lv;
156 // Is it reflected?
158 {
160 }
161 else
162 {
163 lv = physvol->GetLogicalVolume();
164 }
165
166 const G4String volumeref = GenerateName(lv->GetName(), lv);
167
168 if(ModuleName.empty())
169 {
170 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
171 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
172 physvolElement->appendChild(volumerefElement);
173 }
174 else
175 {
176 xercesc::DOMElement* fileElement = NewElement("file");
177 fileElement->setAttributeNode(NewAttribute("name", ModuleName));
178 fileElement->setAttributeNode(NewAttribute("volname", volumeref));
179 physvolElement->appendChild(fileElement);
180 }
181
182 if(std::fabs(pos.x()) > kLinearPrecision ||
183 std::fabs(pos.y()) > kLinearPrecision ||
184 std::fabs(pos.z()) > kLinearPrecision)
185 {
186 PositionWrite(physvolElement, name + "_pos", pos);
187 }
188 if(std::fabs(rot.x()) > kAngularPrecision ||
189 std::fabs(rot.y()) > kAngularPrecision ||
190 std::fabs(rot.z()) > kAngularPrecision)
191 {
192 RotationWrite(physvolElement, name + "_rot", rot);
193 }
194 if(std::fabs(scl.x() - 1.0) > kRelativePrecision ||
195 std::fabs(scl.y() - 1.0) > kRelativePrecision ||
196 std::fabs(scl.z() - 1.0) > kRelativePrecision)
197 {
198 ScaleWrite(physvolElement, name + "_scl", scl);
199 }
200}
201
202// --------------------------------------------------------------------
204 xercesc::DOMElement* volumeElement, const G4VPhysicalVolume* const replicavol)
205{
206 EAxis axis = kUndefined;
207 G4int number = 0;
208 G4double width = 0.0;
209 G4double offset = 0.0;
210 G4bool consuming = false;
211 G4String unitString("mm");
212
213 replicavol->GetReplicationData(axis, number, width, offset, consuming);
214
215 const G4String volumeref = GenerateName(
216 replicavol->GetLogicalVolume()->GetName(), replicavol->GetLogicalVolume());
217
218 xercesc::DOMElement* replicavolElement = NewElement("replicavol");
219 replicavolElement->setAttributeNode(NewAttribute("number", number));
220 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
221 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
222 replicavolElement->appendChild(volumerefElement);
223 xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
224 replicavolElement->appendChild(replicateElement);
225
226 xercesc::DOMElement* dirElement = NewElement("direction");
227 if(axis == kXAxis)
228 {
229 dirElement->setAttributeNode(NewAttribute("x", "1"));
230 }
231 else if(axis == kYAxis)
232 {
233 dirElement->setAttributeNode(NewAttribute("y", "1"));
234 }
235 else if(axis == kZAxis)
236 {
237 dirElement->setAttributeNode(NewAttribute("z", "1"));
238 }
239 else if(axis == kRho)
240 {
241 dirElement->setAttributeNode(NewAttribute("rho", "1"));
242 }
243 else if(axis == kPhi)
244 {
245 dirElement->setAttributeNode(NewAttribute("phi", "1"));
246 unitString = "rad";
247 }
248 replicateElement->appendChild(dirElement);
249
250 xercesc::DOMElement* widthElement = NewElement("width");
251 widthElement->setAttributeNode(NewAttribute("value", width));
252 widthElement->setAttributeNode(NewAttribute("unit", unitString));
253 replicateElement->appendChild(widthElement);
254
255 xercesc::DOMElement* offsetElement = NewElement("offset");
256 offsetElement->setAttributeNode(NewAttribute("value", offset));
257 offsetElement->setAttributeNode(NewAttribute("unit", unitString));
258 replicateElement->appendChild(offsetElement);
259
260 volumeElement->appendChild(replicavolElement);
261}
262
263// --------------------------------------------------------------------
264void G4GDMLWriteStructure::AssemblyWrite(xercesc::DOMElement* volumeElement,
265 const G4int assemblyID)
266{
268 G4AssemblyVolume* myassembly = assemblies->GetAssembly(assemblyID);
269
270 xercesc::DOMElement* assemblyElement = NewElement("assembly");
271 G4String name = "Assembly_" + std::to_string(assemblyID);
272
273 assemblyElement->setAttributeNode(NewAttribute("name", name));
274
275 auto vit = myassembly->GetTripletsIterator();
276
277 G4int depth = 0;
278 const G4String ModuleName;
279
280 for(std::size_t i5 = 0; i5 < myassembly->TotalTriplets(); ++i5)
281 {
282 G4LogicalVolume* lvol = (*vit).GetVolume();
283 if (lvol == nullptr)
284 {
285 G4String message = "Nested assemblies not yet supported for exporting. Sorry!";
286 G4Exception("G4GDMLWriteStructure::AssemblyWrite()", "InvalidSetup",
287 FatalException, message);
288 return;
289 }
290 TraverseVolumeTree(lvol, depth + 1);
291
292 const G4ThreeVector rot = GetAngles((*vit).GetRotation()->inverse());
293 const G4ThreeVector pos = (*vit).GetTranslation();
294
295 const G4String pname =
296 GenerateName((*vit).GetVolume()->GetName() + "_pv", &(*vit));
297
298 xercesc::DOMElement* physvolElement = NewElement("physvol");
299 physvolElement->setAttributeNode(NewAttribute("name", pname));
300
301 assemblyElement->appendChild(physvolElement);
302
303 const G4String volumeref =
304 GenerateName((*vit).GetVolume()->GetName(), (*vit).GetVolume());
305
306 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
307 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
308 physvolElement->appendChild(volumerefElement);
309
310 if(std::fabs(pos.x()) > kLinearPrecision ||
311 std::fabs(pos.y()) > kLinearPrecision ||
312 std::fabs(pos.z()) > kLinearPrecision)
313 {
314 PositionWrite(physvolElement,name+"_position_" + std::to_string(i5), pos);
315 }
316
317 if(std::fabs(rot.x()) > kAngularPrecision ||
318 std::fabs(rot.y()) > kAngularPrecision ||
319 std::fabs(rot.z()) > kAngularPrecision)
320 {
321 RotationWrite(physvolElement,name+"_rotation_" + std::to_string(i5), rot);
322 }
323 ++vit;
324 }
325
326 volumeElement->appendChild(assemblyElement);
327}
328
329// --------------------------------------------------------------------
331 const G4LogicalBorderSurface* const bsurf)
332{
333 if(bsurf == nullptr)
334 {
335 return;
336 }
337
338 const G4SurfaceProperty* psurf = bsurf->GetSurfaceProperty();
339
340 // Generate the new element for border-surface
341 //
342 const G4String& bsname = GenerateName(bsurf->GetName(), bsurf);
343 const G4String& psname = GenerateName(psurf->GetName(), psurf);
344 xercesc::DOMElement* borderElement = NewElement("bordersurface");
345 borderElement->setAttributeNode(NewAttribute("name", bsname));
346 borderElement->setAttributeNode(NewAttribute("surfaceproperty", psname));
347
348 const G4String volumeref1 =
349 GenerateName(bsurf->GetVolume1()->GetName(), bsurf->GetVolume1());
350 const G4String volumeref2 =
351 GenerateName(bsurf->GetVolume2()->GetName(), bsurf->GetVolume2());
352 xercesc::DOMElement* volumerefElement1 = NewElement("physvolref");
353 xercesc::DOMElement* volumerefElement2 = NewElement("physvolref");
354 volumerefElement1->setAttributeNode(NewAttribute("ref", volumeref1));
355 volumerefElement2->setAttributeNode(NewAttribute("ref", volumeref2));
356 borderElement->appendChild(volumerefElement1);
357 borderElement->appendChild(volumerefElement2);
358
359 if(FindOpticalSurface(psurf))
360 {
361 const G4OpticalSurface* opsurf =
362 dynamic_cast<const G4OpticalSurface*>(psurf);
363 if(opsurf == nullptr)
364 {
365 G4Exception("G4GDMLWriteStructure::BorderSurfaceCache()", "InvalidSetup",
366 FatalException, "No optical surface found!");
367 return;
368 }
370 }
371
372 borderElementVec.push_back(borderElement);
373}
374
375// --------------------------------------------------------------------
377 const G4LogicalSkinSurface* const ssurf)
378{
379 if(ssurf == nullptr)
380 {
381 return;
382 }
383
384 const G4SurfaceProperty* psurf = ssurf->GetSurfaceProperty();
385
386 // Generate the new element for border-surface
387 //
388 const G4String& ssname = GenerateName(ssurf->GetName(), ssurf);
389 const G4String& psname = GenerateName(psurf->GetName(), psurf);
390 xercesc::DOMElement* skinElement = NewElement("skinsurface");
391 skinElement->setAttributeNode(NewAttribute("name", ssname));
392 skinElement->setAttributeNode(NewAttribute("surfaceproperty", psname));
393
394 const G4String volumeref = GenerateName(ssurf->GetLogicalVolume()->GetName(),
395 ssurf->GetLogicalVolume());
396 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
397 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
398 skinElement->appendChild(volumerefElement);
399
400 if(FindOpticalSurface(psurf))
401 {
402 const G4OpticalSurface* opsurf =
403 dynamic_cast<const G4OpticalSurface*>(psurf);
404 if(opsurf == nullptr)
405 {
406 G4Exception("G4GDMLWriteStructure::SkinSurfaceCache()", "InvalidSetup",
407 FatalException, "No optical surface found!");
408 return;
409 }
411 }
412
413 skinElementVec.push_back(skinElement);
414}
415
416// --------------------------------------------------------------------
418{
419 const G4OpticalSurface* osurf = dynamic_cast<const G4OpticalSurface*>(psurf);
420 auto pos = std::find(opt_vec.cbegin(), opt_vec.cend(), osurf);
421 if(pos != opt_vec.cend())
422 {
423 return false;
424 } // item already created!
425
426 opt_vec.push_back(osurf); // cache it for future reference
427 return true;
428}
429
430// --------------------------------------------------------------------
432 const G4LogicalVolume* const lvol)
433{
434 G4LogicalSkinSurface* surf = 0;
436 if(nsurf)
437 {
438 const G4LogicalSkinSurfaceTable* stable =
440 for(auto pos = stable->cbegin(); pos != stable->cend(); ++pos)
441 {
442 if(lvol == (*pos)->GetLogicalVolume())
443 {
444 surf = *pos;
445 break;
446 }
447 }
448 }
449 return surf;
450}
451
452// --------------------------------------------------------------------
454 const G4VPhysicalVolume* const pvol)
455{
456 G4LogicalBorderSurface* surf = nullptr;
458 if(nsurf)
459 {
460 const G4LogicalBorderSurfaceTable* btable =
462 for(auto pos = btable->cbegin(); pos != btable->cend(); ++pos)
463 {
464 if(pvol == pos->first.first) // just the first in the couple
465 { // could be enough?
466 surf = pos->second; // break;
467 BorderSurfaceCache(surf);
468 }
469 }
470 }
471 return surf;
472}
473
474// --------------------------------------------------------------------
476{
477#ifdef G4VERBOSE
478 G4cout << "G4GDML: Writing surfaces..." << G4endl;
479#endif
480 for(auto pos = skinElementVec.cbegin();
481 pos != skinElementVec.cend(); ++pos)
482 {
483 structureElement->appendChild(*pos);
484 }
485 for(auto pos = borderElementVec.cbegin();
486 pos != borderElementVec.cend(); ++pos)
487 {
488 structureElement->appendChild(*pos);
489 }
490}
491
492// --------------------------------------------------------------------
493void G4GDMLWriteStructure::StructureWrite(xercesc::DOMElement* gdmlElement)
494{
495#ifdef G4VERBOSE
496 G4cout << "G4GDML: Writing structure..." << G4endl;
497#endif
498
499 // filling the list of phys volumes that are parts of assemblies
500
502
503 for(auto it = assemblies->cbegin(); it != assemblies->cend(); ++it)
504 {
505 auto vit = (*it)->GetVolumesIterator();
506
507 for(std::size_t i5 = 0; i5 < (*it)->TotalImprintedVolumes(); ++i5)
508 {
509 G4String pvname = (*vit)->GetName();
510 std::size_t pos = pvname.find("_impr_") + 6;
511 G4String impID = pvname.substr(pos);
512
513 pos = impID.find("_");
514 impID = impID.substr(0, pos);
515
516 assemblyVolMap[*vit] = (*it)->GetAssemblyID();
517 imprintsMap[*vit] = std::atoi(impID.c_str());
518 ++vit;
519 }
520 }
521
522 structureElement = NewElement("structure");
523 gdmlElement->appendChild(structureElement);
524}
525
526// --------------------------------------------------------------------
528 const G4LogicalVolume* const volumePtr, const G4int depth)
529{
530 if(VolumeMap().find(volumePtr) != VolumeMap().cend())
531 {
532 return VolumeMap()[volumePtr]; // Volume is already processed
533 }
534
535 G4VSolid* solidPtr = volumePtr->GetSolid();
536 G4Transform3D R, invR;
537 G4int trans = 0;
538
539 std::map<const G4LogicalVolume*, G4GDMLAuxListType>::iterator auxiter;
540
541 ++levelNo;
542
543 while(true) // Solve possible displacement/reflection
544 { // of the referenced solid!
545 if(trans > maxTransforms)
546 {
547 G4String ErrorMessage = "Referenced solid in volume '" +
548 volumePtr->GetName() +
549 "' was displaced/reflected too many times!";
550 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()", "InvalidSetup",
551 FatalException, ErrorMessage);
552 }
553
554 if(G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
555 {
556 R = R * refl->GetTransform3D();
557 solidPtr = refl->GetConstituentMovedSolid();
558 ++trans;
559 continue;
560 }
561
562 if(G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
563 {
564 R = R * G4Transform3D(disp->GetObjectRotation(),
565 disp->GetObjectTranslation());
566 solidPtr = disp->GetConstituentMovedSolid();
567 ++trans;
568 continue;
569 }
570
571 break;
572 }
573
574 // check if it is a reflected volume
575 G4LogicalVolume* tmplv = const_cast<G4LogicalVolume*>(volumePtr);
576
577 if(reflFactory->IsReflected(tmplv))
578 {
579 tmplv = reflFactory->GetConstituentLV(tmplv);
580 if(VolumeMap().find(tmplv) != VolumeMap().cend())
581 {
582 return R; // Volume is already processed
583 }
584 }
585
586 // Only compute the inverse when necessary!
587 //
588 if(trans > 0)
589 {
590 invR = R.inverse();
591 }
592
593 const G4String name = GenerateName(tmplv->GetName(), tmplv);
594
595 G4String materialref = "NULL";
596
597 if(volumePtr->GetMaterial())
598 {
599 materialref = GenerateName(volumePtr->GetMaterial()->GetName(),
600 volumePtr->GetMaterial());
601 }
602
603 const G4String solidref = GenerateName(solidPtr->GetName(), solidPtr);
604
605 xercesc::DOMElement* volumeElement = NewElement("volume");
606 volumeElement->setAttributeNode(NewAttribute("name", name));
607 xercesc::DOMElement* materialrefElement = NewElement("materialref");
608 materialrefElement->setAttributeNode(NewAttribute("ref", materialref));
609 volumeElement->appendChild(materialrefElement);
610 xercesc::DOMElement* solidrefElement = NewElement("solidref");
611 solidrefElement->setAttributeNode(NewAttribute("ref", solidref));
612 volumeElement->appendChild(solidrefElement);
613
614 G4int daughterCount = volumePtr->GetNoDaughters();
615
616 if(levelNo == maxLevel) // Stop exporting if reached levels limit
617 {
618 daughterCount = 0;
619 }
620
621 std::map<G4int, std::vector<G4int> > assemblyIDToAddedImprints;
622
623 for(G4int i = 0; i < daughterCount; ++i) // Traverse all the children!
624 {
625 const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
626 const G4String ModuleName = Modularize(physvol, depth);
627
628 G4Transform3D daughterR;
629
630 if(ModuleName.empty()) // Check if subtree requested to be
631 { // a separate module!
632 daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(), depth + 1);
633 }
634 else
635 {
637 daughterR = writer.Write(ModuleName, physvol->GetLogicalVolume(),
638 SchemaLocation, depth + 1);
639 }
640
641 if(const G4PVDivision* const divisionvol =
642 dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
643 {
644 if(!G4Transform3D::Identity.isNear(invR * daughterR, kRelativePrecision))
645 {
646 G4String ErrorMessage = "Division volume in '" + name +
647 "' can not be related to reflected solid!";
648 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
649 "InvalidSetup", FatalException, ErrorMessage);
650 }
651 DivisionvolWrite(volumeElement, divisionvol);
652 }
653 else
654 {
655 if(physvol->IsParameterised()) // Is it a paramvol?
656 {
657 if(!G4Transform3D::Identity.isNear(invR * daughterR,
659 {
660 G4String ErrorMessage = "Parameterised volume in '" + name +
661 "' can not be related to reflected solid!";
662 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
663 "InvalidSetup", FatalException, ErrorMessage);
664 }
665 ParamvolWrite(volumeElement, physvol);
666 }
667 else
668 {
669 if(physvol->IsReplicated()) // Is it a replicavol?
670 {
671 if(!G4Transform3D::Identity.isNear(invR * daughterR,
673 {
674 G4String ErrorMessage = "Replica volume in '" + name +
675 "' can not be related to reflected solid!";
676 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
677 "InvalidSetup", FatalException, ErrorMessage);
678 }
679 ReplicavolWrite(volumeElement, physvol);
680 }
681 else // Is it a physvol or an assembly?
682 {
683 if(assemblyVolMap.find(physvol) != assemblyVolMap.cend())
684 {
685 G4int assemblyID = assemblyVolMap[physvol];
686
687 G4String assemblyref = "Assembly_" + std::to_string(assemblyID);
688
689 // here I need to retrieve the imprint ID
690
691 G4int imprintID = imprintsMap[physvol];
692
693 // there are 2 steps:
694 //
695 // 1) add assembly to the structure if that has not yet been done
696 // (but after the constituents volumes have been added)
697 //
698
699 if(std::find(addedAssemblies.cbegin(), addedAssemblies.cend(),
700 assemblyID) == addedAssemblies.cend())
701 {
702 AssemblyWrite(structureElement, assemblyID);
703 addedAssemblies.push_back(assemblyID);
704 assemblyIDToAddedImprints[assemblyID] = std::vector<G4int>();
705 }
706
707 // 2) add the assembly (as physical volume) to the mother volume
708 // (but only once), using it's original position and rotation.
709 //
710
711 // here I need a check if assembly has been already added to the
712 // mother volume
713 std::vector<G4int>& addedImprints = assemblyIDToAddedImprints[assemblyID];
714 if(std::find(addedImprints.cbegin(), addedImprints.cend(),
715 imprintID) == addedImprints.cend())
716 {
717 G4String imprintname = "Imprint_" + std::to_string(imprintID) + "_";
718 imprintname = GenerateName(imprintname, physvol);
719
720 // I need to get those two from the container of imprints from
721 // the assembly I have the imprint ID, I need to get pos and rot
722 //
724 ->GetAssembly(assemblyID)
725 ->GetImprintTransformation(imprintID);
726
727 HepGeom::Scale3D scale;
728 HepGeom::Rotate3D rotate;
729 HepGeom::Translate3D translate;
730
731 transf.getDecomposition(scale, rotate, translate);
732
733 const G4ThreeVector scl(scale(0, 0), scale(1, 1), scale(2, 2));
734 const G4ThreeVector rot = GetAngles(rotate.getRotation().inverse());
735 const G4ThreeVector pos = transf.getTranslation();
736
737 // here I need a normal physvol referencing to my assemblyref
738
739 xercesc::DOMElement* physvolElement = NewElement("physvol");
740 physvolElement->setAttributeNode(NewAttribute("name", imprintname));
741
742 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
743 volumerefElement->setAttributeNode(NewAttribute("ref", assemblyref));
744 physvolElement->appendChild(volumerefElement);
745
746 if(std::fabs(pos.x()) > kLinearPrecision ||
747 std::fabs(pos.y()) > kLinearPrecision ||
748 std::fabs(pos.z()) > kLinearPrecision)
749 {
750 PositionWrite(physvolElement, imprintname + "_pos", pos);
751 }
752 if(std::fabs(rot.x()) > kAngularPrecision ||
753 std::fabs(rot.y()) > kAngularPrecision ||
754 std::fabs(rot.z()) > kAngularPrecision)
755 {
756 RotationWrite(physvolElement, imprintname + "_rot", rot);
757 }
758 if(std::fabs(scl.x() - 1.0) > kRelativePrecision ||
759 std::fabs(scl.y() - 1.0) > kRelativePrecision ||
760 std::fabs(scl.z() - 1.0) > kRelativePrecision)
761 {
762 ScaleWrite(physvolElement, name + "_scl", scl);
763 }
764
765 volumeElement->appendChild(physvolElement);
766 //
767 addedImprints.push_back(imprintID);
768 }
769 }
770 else // not part of assembly, so a normal physical volume
771 {
773 if(physvol->GetFrameRotation() != nullptr)
774 {
775 rot = *(physvol->GetFrameRotation());
776 }
777 G4Transform3D P(rot, physvol->GetObjectTranslation());
778
779 PhysvolWrite(volumeElement, physvol, invR * P * daughterR,
780 ModuleName);
781 }
782 }
783 }
784 }
785 // BorderSurfaceCache(GetBorderSurface(physvol));
786 GetBorderSurface(physvol);
787 }
788
789 if(cexport)
790 {
791 ExportEnergyCuts(volumePtr);
792 }
793 // Add optional energy cuts
794
795 if(sdexport)
796 {
797 ExportSD(volumePtr);
798 }
799 // Add optional SDs
800
801 // Here write the auxiliary info
802 //
803 auxiter = auxmap.find(volumePtr);
804 if(auxiter != auxmap.cend())
805 {
806 AddAuxInfo(&(auxiter->second), volumeElement);
807 }
808
809 structureElement->appendChild(volumeElement);
810 // Append the volume AFTER traversing the children so that
811 // the order of volumes will be correct!
812
813 VolumeMap()[tmplv] = R;
814
815 AddExtension(volumeElement, volumePtr);
816 // Add any possible user defined extension attached to a volume
817
818 AddMaterial(volumePtr->GetMaterial());
819 // Add the involved materials and solids!
820
821 AddSolid(solidPtr);
822
824
825 return R;
826}
827
828// --------------------------------------------------------------------
830 const G4LogicalVolume* const lvol)
831{
832 auto pos = auxmap.find(lvol);
833
834 if(pos == auxmap.cend())
835 {
836 auxmap[lvol] = G4GDMLAuxListType();
837 }
838
839 auxmap[lvol].push_back(myaux);
840}
841
842// --------------------------------------------------------------------
844{
845 cexport = fcuts;
846}
847
848// --------------------------------------------------------------------
850{
851 G4GDMLEvaluator eval;
852 G4ProductionCuts* pcuts = lvol->GetRegion()->GetProductionCuts();
854 G4Gamma* gamma = G4Gamma::Gamma();
858
859 G4double gamma_cut = ctab->ConvertRangeToEnergy(
860 gamma, lvol->GetMaterial(), pcuts->GetProductionCut("gamma"));
861 G4double eminus_cut = ctab->ConvertRangeToEnergy(
862 eminus, lvol->GetMaterial(), pcuts->GetProductionCut("e-"));
863 G4double eplus_cut = ctab->ConvertRangeToEnergy(
864 eplus, lvol->GetMaterial(), pcuts->GetProductionCut("e+"));
865 G4double proton_cut = ctab->ConvertRangeToEnergy(
866 proton, lvol->GetMaterial(), pcuts->GetProductionCut("proton"));
867
868 G4GDMLAuxStructType gammainfo = { "gammaECut",
869 eval.ConvertToString(gamma_cut), "MeV", 0 };
870 G4GDMLAuxStructType eminusinfo = { "electronECut",
871 eval.ConvertToString(eminus_cut), "MeV",
872 0 };
873 G4GDMLAuxStructType eplusinfo = { "positronECut",
874 eval.ConvertToString(eplus_cut), "MeV", 0 };
875 G4GDMLAuxStructType protinfo = { "protonECut",
876 eval.ConvertToString(proton_cut), "MeV", 0 };
877
878 AddVolumeAuxiliary(gammainfo, lvol);
879 AddVolumeAuxiliary(eminusinfo, lvol);
880 AddVolumeAuxiliary(eplusinfo, lvol);
881 AddVolumeAuxiliary(protinfo, lvol);
882}
883
884// --------------------------------------------------------------------
886{
887 sdexport = fsd;
888}
889
890// --------------------------------------------------------------------
892{
894
895 if(sd != nullptr)
896 {
897 G4String SDname = sd->GetName();
898
899 G4GDMLAuxStructType SDinfo = { "SensDet", SDname, "", 0 };
900 AddVolumeAuxiliary(SDinfo, lvol);
901 }
902}
903
904// --------------------------------------------------------------------
906{
907 return maxLevel;
908}
909
910// --------------------------------------------------------------------
912{
913 if(level <= 0)
914 {
915 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()", "InvalidSetup",
916 FatalException, "Levels to export must be greater than zero!");
917 return;
918 }
919 maxLevel = level;
920 levelNo = 0;
921}
static const G4double pos
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
std::map< std::pair< const G4VPhysicalVolume *, const G4VPhysicalVolume * >, G4LogicalBorderSurface * > G4LogicalBorderSurfaceTable
std::vector< G4LogicalSkinSurface * > G4LogicalSkinSurfaceTable
static constexpr double eplus
Definition: G4SIunits.hh:184
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
HepRotation inverse() const
static G4AssemblyStore * GetInstance()
G4AssemblyVolume * GetAssembly(unsigned int id, G4bool verbose=true) const
G4Transform3D & GetImprintTransformation(unsigned int imprintID)
std::vector< G4AssemblyTriplet >::iterator GetTripletsIterator()
std::size_t TotalTriplets() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4String ConvertToString(G4int ival)
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)
std::map< const G4VPhysicalVolume *, G4int > assemblyVolMap
const G4LogicalBorderSurface * GetBorderSurface(const G4VPhysicalVolume *const)
std::vector< const G4OpticalSurface * > opt_vec
void BorderSurfaceCache(const G4LogicalBorderSurface *const)
std::vector< xercesc::DOMElement * > borderElementVec
G4ReflectionFactory * reflFactory
std::map< const G4VPhysicalVolume *, G4int > imprintsMap
std::map< const G4LogicalVolume *, G4GDMLAuxListType > auxmap
void AssemblyWrite(xercesc::DOMElement *, const int assemblyID)
std::vector< xercesc::DOMElement * > skinElementVec
std::vector< G4int > addedAssemblies
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)
void ExportEnergyCuts(const G4LogicalVolume *const)
virtual void StructureWrite(xercesc::DOMElement *)
G4Transform3D TraverseVolumeTree(const G4LogicalVolume *const topVol, const G4int depth)
xercesc::DOMElement * structureElement
void ExportSD(const G4LogicalVolume *const)
void DivisionvolWrite(xercesc::DOMElement *, const G4PVDivision *const)
G4bool FindOpticalSurface(const G4SurfaceProperty *)
void AddVolumeAuxiliary(G4GDMLAuxStructType myaux, const G4LogicalVolume *const)
void ReplicavolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:192
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
G4String SchemaLocation
Definition: G4GDMLWrite.hh:130
G4String Modularize(const G4VPhysicalVolume *const topvol, const G4int depth)
Definition: G4GDMLWrite.cc:377
virtual void AddExtension(xercesc::DOMElement *, const G4LogicalVolume *const)
Definition: G4GDMLWrite.cc:80
G4Transform3D Write(const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
Definition: G4GDMLWrite.cc:203
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:155
void AddAuxInfo(G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
Definition: G4GDMLWrite.cc:94
VolumeMapType & VolumeMap()
Definition: G4GDMLWrite.cc:60
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static std::size_t GetNumberOfBorderSurfaces()
static const G4LogicalBorderSurfaceTable * GetSurfaceTable()
const G4VPhysicalVolume * GetVolume2() const
const G4VPhysicalVolume * GetVolume1() const
const G4LogicalVolume * GetLogicalVolume() const
static const G4LogicalSkinSurfaceTable * GetSurfaceTable()
static std::size_t GetNumberOfSkinSurfaces()
const G4String & GetName() const
G4SurfaceProperty * GetSurfaceProperty() const
G4VSolid * GetSolid() const
G4VSensitiveDetector * GetSensitiveDetector() const
size_t GetNoDaughters() const
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
const G4String & GetName() const
const G4String & GetName() const
Definition: G4Material.hh:173
EAxis GetDivisionAxis() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4ProductionCutsTable * GetProductionCutsTable()
G4double ConvertRangeToEnergy(const G4ParticleDefinition *particle, const G4Material *material, G4double range)
G4double GetProductionCut(G4int index) const
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4bool IsReflected(G4LogicalVolume *lv) const
static G4ReflectionFactory * Instance()
G4LogicalVolume * GetConstituentLV(G4LogicalVolume *reflLV) const
G4ProductionCuts * GetProductionCuts() 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
virtual G4int GetCopyNo() const =0
const G4String & GetName() const
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
G4String GetName() const
static DLL_API const Transform3D Identity
Definition: Transform3D.h:196
CLHEP::HepRotation getRotation() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:173
CLHEP::Hep3Vector getTranslation() const
Transform3D inverse() const
Definition: Transform3D.cc:141
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:60
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
@ kUndefined
Definition: geomdefs.hh:61
@ kRho
Definition: geomdefs.hh:58
const char * name(G4int ptype)
string pname
Definition: eplot.py:33
static double P[]
#define INT_MAX
Definition: templates.hh:90