Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4GDMLWriteStructure.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4GDMLWriteStructure.cc
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: PHG4GDMLWriteStructure.cc 68053 2013-03-13 14:39:51Z gcosmo $
28 //
29 // class PHG4GDMLWriteStructure Implementation
30 //
31 // Original author: Zoltan Torzsok, November 2007
32 //
33 // --------------------------------------------------------------------
34 
36 #include "PHG4GDMLConfig.hh"
37 
38 #include <Geant4/G4DisplacedSolid.hh>
39 #include <Geant4/G4LogicalBorderSurface.hh>
40 #include <Geant4/G4LogicalSkinSurface.hh>
41 #include <Geant4/G4LogicalVolumeStore.hh>
42 #include <Geant4/G4Material.hh>
43 #include <Geant4/G4OpticalSurface.hh>
44 #include <Geant4/G4PVDivision.hh>
45 #include <Geant4/G4PVReplica.hh>
46 #include <Geant4/G4PhysicalVolumeStore.hh>
47 #include <Geant4/G4ReflectedSolid.hh>
48 #include <Geant4/G4ReflectionFactory.hh>
49 #include <Geant4/G4Region.hh>
50 
51 #include <Geant4/G4Electron.hh>
52 #include <Geant4/G4Gamma.hh>
53 #include <Geant4/G4Positron.hh>
54 #include <Geant4/G4ProductionCuts.hh>
55 #include <Geant4/G4ProductionCutsTable.hh>
56 #include <Geant4/G4Proton.hh>
57 #include <Geant4/G4Version.hh>
58 
59 #include <cassert>
60 
63  , structureElement(nullptr)
64  , cexport(false)
65  , config(config_input)
66 {
67  assert(config);
68  reflFactory = G4ReflectionFactory::Instance();
69 }
70 
72 {
73 }
74 
75 void PHG4GDMLWriteStructure::DivisionvolWrite(xercesc::DOMElement* volumeElement,
76  const G4PVDivision* const divisionvol)
77 {
78  EAxis axis = kUndefined;
79  G4int number = 0;
80  G4double width = 0.0;
81  G4double offset = 0.0;
82  G4bool consuming = false;
83 
84  divisionvol->GetReplicationData(axis, number, width, offset, consuming);
85  axis = divisionvol->GetDivisionAxis();
86 
87  G4String unitString("mm");
88  G4String axisString("kUndefined");
89  if (axis == kXAxis)
90  {
91  axisString = "kXAxis";
92  }
93  else if (axis == kYAxis)
94  {
95  axisString = "kYAxis";
96  }
97  else if (axis == kZAxis)
98  {
99  axisString = "kZAxis";
100  }
101  else if (axis == kRho)
102  {
103  axisString = "kRho";
104  }
105  else if (axis == kPhi)
106  {
107  axisString = "kPhi";
108  unitString = "rad";
109  }
110 
111  const G4String name = GenerateName(divisionvol->GetName(), divisionvol);
112  const G4String volumeref = 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 void PHG4GDMLWriteStructure::PhysvolWrite(xercesc::DOMElement* volumeElement,
128  const G4VPhysicalVolume* const physvol,
129  const G4Transform3D& T,
130  const G4String& ModuleName)
131 {
132  HepGeom::Scale3D scale;
133  HepGeom::Rotate3D rotate;
134  HepGeom::Translate3D translate;
135 
136  T.getDecomposition(scale, rotate, translate);
137 
138  const G4ThreeVector scl(scale(0, 0), scale(1, 1), scale(2, 2));
139  const G4ThreeVector rot = GetAngles(rotate.getRotation());
140  const G4ThreeVector pos = T.getTranslation();
141 
142  const G4String name = GenerateName(physvol->GetName(), physvol);
143  const G4int copynumber = physvol->GetCopyNo();
144 
145  xercesc::DOMElement* physvolElement = NewElement("physvol");
146  physvolElement->setAttributeNode(NewAttribute("name", name));
147  if (copynumber) physvolElement->setAttributeNode(NewAttribute("copynumber", copynumber));
148 
149  volumeElement->appendChild(physvolElement);
150 
151  G4LogicalVolume* lv;
152  // Is it reflected?
153  if (reflFactory->IsReflected(physvol->GetLogicalVolume()))
154  {
155  lv = reflFactory->GetConstituentLV(physvol->GetLogicalVolume());
156  }
157  else
158  {
159  lv = physvol->GetLogicalVolume();
160  }
161 
162  const G4String volumeref = GenerateName(lv->GetName(), lv);
163 
164  if (ModuleName.empty())
165  {
166  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
167  volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
168  physvolElement->appendChild(volumerefElement);
169  }
170  else
171  {
172  xercesc::DOMElement* fileElement = NewElement("file");
173  fileElement->setAttributeNode(NewAttribute("name", ModuleName));
174  fileElement->setAttributeNode(NewAttribute("volname", volumeref));
175  physvolElement->appendChild(fileElement);
176  }
177 
178  if (std::fabs(pos.x()) > kLinearPrecision || std::fabs(pos.y()) > kLinearPrecision || std::fabs(pos.z()) > kLinearPrecision)
179  {
180  PositionWrite(physvolElement, name + "_pos", pos);
181  }
182  if (std::fabs(rot.x()) > kAngularPrecision || std::fabs(rot.y()) > kAngularPrecision || std::fabs(rot.z()) > kAngularPrecision)
183  {
184  RotationWrite(physvolElement, name + "_rot", rot);
185  }
186  if (std::fabs(scl.x() - 1.0) > kRelativePrecision || std::fabs(scl.y() - 1.0) > kRelativePrecision || std::fabs(scl.z() - 1.0) > kRelativePrecision)
187  {
188  ScaleWrite(physvolElement, name + "_scl", scl);
189  }
190 }
191 
192 void PHG4GDMLWriteStructure::ReplicavolWrite(xercesc::DOMElement* volumeElement,
193  const G4VPhysicalVolume* const replicavol)
194 {
195  EAxis axis = kUndefined;
196  G4int number = 0;
197  G4double width = 0.0;
198  G4double offset = 0.0;
199  G4bool consuming = false;
200  G4String unitString("mm");
201 
202  replicavol->GetReplicationData(axis, number, width, offset, consuming);
203 
204  const G4String volumeref = GenerateName(replicavol->GetLogicalVolume()->GetName(),
205  replicavol->GetLogicalVolume());
206 
207  xercesc::DOMElement* replicavolElement = NewElement("replicavol");
208  replicavolElement->setAttributeNode(NewAttribute("number", number));
209  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
210  volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
211  replicavolElement->appendChild(volumerefElement);
212  xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
213  replicavolElement->appendChild(replicateElement);
214 
215  xercesc::DOMElement* dirElement = NewElement("direction");
216  if (axis == kXAxis)
217  {
218  dirElement->setAttributeNode(NewAttribute("x", "1"));
219  }
220  else if (axis == kYAxis)
221  {
222  dirElement->setAttributeNode(NewAttribute("y", "1"));
223  }
224  else if (axis == kZAxis)
225  {
226  dirElement->setAttributeNode(NewAttribute("z", "1"));
227  }
228  else if (axis == kRho)
229  {
230  dirElement->setAttributeNode(NewAttribute("rho", "1"));
231  }
232  else if (axis == kPhi)
233  {
234  dirElement->setAttributeNode(NewAttribute("phi", "1"));
235  unitString = "rad";
236  }
237  replicateElement->appendChild(dirElement);
238 
239  xercesc::DOMElement* widthElement = NewElement("width");
240  widthElement->setAttributeNode(NewAttribute("value", width));
241  widthElement->setAttributeNode(NewAttribute("unit", unitString));
242  replicateElement->appendChild(widthElement);
243 
244  xercesc::DOMElement* offsetElement = NewElement("offset");
245  offsetElement->setAttributeNode(NewAttribute("value", offset));
246  offsetElement->setAttributeNode(NewAttribute("unit", unitString));
247  replicateElement->appendChild(offsetElement);
248 
249  volumeElement->appendChild(replicavolElement);
250 }
251 
253  BorderSurfaceCache(const G4LogicalBorderSurface* const bsurf)
254 {
255  if (!bsurf)
256  {
257  return;
258  }
259 
260  const G4SurfaceProperty* psurf = bsurf->GetSurfaceProperty();
261 
262  // Generate the new element for border-surface
263  //
264  xercesc::DOMElement* borderElement = NewElement("bordersurface");
265  borderElement->setAttributeNode(NewAttribute("name", bsurf->GetName()));
266  borderElement->setAttributeNode(NewAttribute("surfaceproperty",
267  psurf->GetName()));
268 
269  const G4String volumeref1 = GenerateName(bsurf->GetVolume1()->GetName(),
270  bsurf->GetVolume1());
271  const G4String volumeref2 = GenerateName(bsurf->GetVolume2()->GetName(),
272  bsurf->GetVolume2());
273  xercesc::DOMElement* volumerefElement1 = NewElement("physvolref");
274  xercesc::DOMElement* volumerefElement2 = NewElement("physvolref");
275  volumerefElement1->setAttributeNode(NewAttribute("ref", volumeref1));
276  volumerefElement2->setAttributeNode(NewAttribute("ref", volumeref2));
277  borderElement->appendChild(volumerefElement1);
278  borderElement->appendChild(volumerefElement2);
279 
280  if (FindOpticalSurface(psurf))
281  {
282  const G4OpticalSurface* opsurf =
283  dynamic_cast<const G4OpticalSurface*>(psurf);
284  if (!opsurf)
285  {
286  G4Exception("PHG4GDMLWriteStructure::BorderSurfaceCache()",
287  "InvalidSetup", FatalException, "No optical surface found!");
288  return;
289  }
291  }
292 
293  borderElementVec.push_back(borderElement);
294 }
295 
297  SkinSurfaceCache(const G4LogicalSkinSurface* const ssurf)
298 {
299  if (!ssurf)
300  {
301  return;
302  }
303 
304  const G4SurfaceProperty* psurf = ssurf->GetSurfaceProperty();
305 
306  // Generate the new element for border-surface
307  //
308  xercesc::DOMElement* skinElement = NewElement("skinsurface");
309  skinElement->setAttributeNode(NewAttribute("name", ssurf->GetName()));
310  skinElement->setAttributeNode(NewAttribute("surfaceproperty",
311  psurf->GetName()));
312 
313  const G4String volumeref = GenerateName(ssurf->GetLogicalVolume()->GetName(),
314  ssurf->GetLogicalVolume());
315  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
316  volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
317  skinElement->appendChild(volumerefElement);
318 
319  if (FindOpticalSurface(psurf))
320  {
321  const G4OpticalSurface* opsurf =
322  dynamic_cast<const G4OpticalSurface*>(psurf);
323  if (!opsurf)
324  {
325  G4Exception("PHG4GDMLWriteStructure::SkinSurfaceCache()",
326  "InvalidSetup", FatalException, "No optical surface found!");
327  return;
328  }
330  }
331 
332  skinElementVec.push_back(skinElement);
333 }
334 
335 G4bool PHG4GDMLWriteStructure::FindOpticalSurface(const G4SurfaceProperty* psurf)
336 {
337  const G4OpticalSurface* osurf = dynamic_cast<const G4OpticalSurface*>(psurf);
338  std::vector<const G4OpticalSurface*>::const_iterator pos;
339  pos = std::find(opt_vec.begin(), opt_vec.end(), osurf);
340  if (pos != opt_vec.end())
341  {
342  return false;
343  } // item already created!
344 
345  opt_vec.push_back(osurf); // cache it for future reference
346  return true;
347 }
348 
349 const G4LogicalSkinSurface*
350 PHG4GDMLWriteStructure::GetSkinSurface(const G4LogicalVolume* const lvol)
351 {
352  G4LogicalSkinSurface* surf = 0;
353  G4int nsurf = G4LogicalSkinSurface::GetNumberOfSkinSurfaces();
354  if (nsurf)
355  {
356  const G4LogicalSkinSurfaceTable* stable =
357  G4LogicalSkinSurface::GetSurfaceTable();
358  std::vector<G4LogicalSkinSurface*>::const_iterator pos;
359  for (pos = stable->begin(); pos != stable->end(); ++pos)
360  {
361  if (lvol == (*pos)->GetLogicalVolume())
362  {
363  surf = *pos;
364  break;
365  }
366  }
367  }
368  return surf;
369 }
370 
371 #if G4VERSION_NUMBER >= 1007
372 
373 const G4LogicalBorderSurface* PHG4GDMLWriteStructure::GetBorderSurface(
374  const G4VPhysicalVolume* const pvol)
375 {
376  G4LogicalBorderSurface* surf = nullptr;
377  G4int nsurf = G4LogicalBorderSurface::GetNumberOfBorderSurfaces();
378  if (nsurf)
379  {
380  const G4LogicalBorderSurfaceTable* btable =
381  G4LogicalBorderSurface::GetSurfaceTable();
382  for (auto pos = btable->cbegin(); pos != btable->cend(); ++pos)
383  {
384  if (pvol == pos->first.first) // just the first in the couple
385  { // could be enough?
386  surf = pos->second; // break;
387  BorderSurfaceCache(surf);
388  }
389  }
390  }
391  return surf;
392 }
393 
394 #else
395 
396 const G4LogicalBorderSurface*
398 {
399  G4LogicalBorderSurface* surf = 0;
400  G4int nsurf = G4LogicalBorderSurface::GetNumberOfBorderSurfaces();
401  if (nsurf)
402  {
403  const G4LogicalBorderSurfaceTable* btable =
404  G4LogicalBorderSurface::GetSurfaceTable();
405  std::vector<G4LogicalBorderSurface*>::const_iterator pos;
406  for (pos = btable->begin(); pos != btable->end(); ++pos)
407  {
408  if (pvol == (*pos)->GetVolume1()) // just the first in the couple
409  { // is enough
410  surf = *pos;
411  break;
412  }
413  }
414  }
415  return surf;
416 }
417 
418 #endif
419 
421 {
422  std::cout << "PHG4GDML: Writing surfaces..." << std::endl;
423 
424  std::vector<xercesc::DOMElement*>::const_iterator pos;
425  for (pos = skinElementVec.begin(); pos != skinElementVec.end(); ++pos)
426  {
427  structureElement->appendChild(*pos);
428  }
429  for (pos = borderElementVec.begin(); pos != borderElementVec.end(); ++pos)
430  {
431  structureElement->appendChild(*pos);
432  }
433 }
434 
435 void PHG4GDMLWriteStructure::StructureWrite(xercesc::DOMElement* gdmlElement)
436 {
437  std::cout << "PHG4GDML: Writing structure..." << std::endl;
438 
439  structureElement = NewElement("structure");
440  gdmlElement->appendChild(structureElement);
441 }
442 
443 G4Transform3D PHG4GDMLWriteStructure::
444  TraverseVolumeTree(const G4LogicalVolume* const volumePtr, const G4int depth)
445 {
446  if (VolumeMap().find(volumePtr) != VolumeMap().end())
447  {
448  return VolumeMap()[volumePtr]; // Volume is already processed
449  }
450 
451  //jump over the exclusions
452  assert(config);
453  if (config->get_excluded_logical_vol().find(volumePtr) != config->get_excluded_logical_vol().end())
454  {
455  return G4Transform3D::Identity;
456  }
457 
458  G4VSolid* solidPtr = volumePtr->GetSolid();
459  G4Transform3D R, invR;
460  G4int trans = 0;
461 
462  std::map<const G4LogicalVolume*, PHG4GDMLAuxListType>::iterator auxiter;
463 
464  while (true) // Solve possible displacement/reflection
465  { // of the referenced solid!
466  if (trans > maxTransforms)
467  {
468  G4String ErrorMessage = "Referenced solid in volume '" + volumePtr->GetName() + "' was displaced/reflected too many times!";
469  G4Exception("PHG4GDMLWriteStructure::TraverseVolumeTree()",
470  "InvalidSetup", FatalException, ErrorMessage);
471  }
472 
473  if (G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
474  {
475  R = R * refl->GetTransform3D();
476  solidPtr = refl->GetConstituentMovedSolid();
477  trans++;
478  continue;
479  }
480 
481  if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
482  {
483  R = R * G4Transform3D(disp->GetObjectRotation(),
484  disp->GetObjectTranslation());
485  solidPtr = disp->GetConstituentMovedSolid();
486  trans++;
487  continue;
488  }
489 
490  break;
491  }
492 
493  //check if it is a reflected volume
494  G4LogicalVolume* tmplv = const_cast<G4LogicalVolume*>(volumePtr);
495 
496  if (reflFactory->IsReflected(tmplv))
497  {
498  tmplv = reflFactory->GetConstituentLV(tmplv);
499  if (VolumeMap().find(tmplv) != VolumeMap().end())
500  {
501  return R; // Volume is already processed
502  }
503  }
504 
505  // Only compute the inverse when necessary!
506  //
507  if (trans > 0)
508  {
509  invR = R.inverse();
510  }
511 
512  const G4String name = GenerateName(tmplv->GetName(), tmplv);
513  const G4String materialref = GenerateName(volumePtr->GetMaterial()->GetName(),
514  volumePtr->GetMaterial());
515  const G4String solidref = GenerateName(solidPtr->GetName(), solidPtr);
516 
517  xercesc::DOMElement* volumeElement = NewElement("volume");
518  volumeElement->setAttributeNode(NewAttribute("name", name));
519  xercesc::DOMElement* materialrefElement = NewElement("materialref");
520  materialrefElement->setAttributeNode(NewAttribute("ref", materialref));
521  volumeElement->appendChild(materialrefElement);
522  xercesc::DOMElement* solidrefElement = NewElement("solidref");
523  solidrefElement->setAttributeNode(NewAttribute("ref", solidref));
524  volumeElement->appendChild(solidrefElement);
525 
526  const G4int daughterCount = volumePtr->GetNoDaughters();
527 
528  for (G4int i = 0; i < daughterCount; i++) // Traverse all the children!
529  {
530  const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
531 
532  //jump over the exclusions
533  assert(config);
534  if (config->get_excluded_physical_vol().find(physvol) != config->get_excluded_physical_vol().end())
535  continue;
536 
537  const G4String ModuleName = Modularize(physvol, depth);
538 
539  G4Transform3D daughterR;
540 
541  if (ModuleName.empty()) // Check if subtree requested to be
542  { // a separate module!
543  daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(), depth + 1);
544  }
545  else
546  {
548  daughterR = writer.Write(ModuleName, physvol->GetLogicalVolume(),
549  SchemaLocation, depth + 1);
550  }
551 
552  if (const G4PVDivision* const divisionvol = dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
553  {
554  if (!G4Transform3D::Identity.isNear(invR * daughterR, kRelativePrecision))
555  {
556  G4String ErrorMessage = "Division volume in '" + name + "' can not be related to reflected solid!";
557  G4Exception("PHG4GDMLWriteStructure::TraverseVolumeTree()",
558  "InvalidSetup", FatalException, ErrorMessage);
559  }
560  DivisionvolWrite(volumeElement, divisionvol);
561  }
562  else if (physvol->IsParameterised()) // Is it a paramvol?
563  {
564  if (!G4Transform3D::Identity.isNear(invR * daughterR, kRelativePrecision))
565  {
566  G4String ErrorMessage = "Parameterised volume in '" + name + "' can not be related to reflected solid!";
567  G4Exception("PHG4GDMLWriteStructure::TraverseVolumeTree()",
568  "InvalidSetup", FatalException, ErrorMessage);
569  }
570  ParamvolWrite(volumeElement, physvol);
571  }
572  else if (physvol->IsReplicated()) // Is it a replicavol?
573  {
574  if (!G4Transform3D::Identity.isNear(invR * daughterR, kRelativePrecision))
575  {
576  G4String ErrorMessage = "Replica volume in '" + name + "' can not be related to reflected solid!";
577  G4Exception("PHG4GDMLWriteStructure::TraverseVolumeTree()",
578  "InvalidSetup", FatalException, ErrorMessage);
579  }
580  ReplicavolWrite(volumeElement, physvol);
581  }
582  else // Is it a physvol?
583  {
584  G4RotationMatrix rot;
585  if (physvol->GetFrameRotation() != 0)
586  {
587  rot = *(physvol->GetFrameRotation());
588  }
589  G4Transform3D P(rot, physvol->GetObjectTranslation());
590 
591  PhysvolWrite(volumeElement, physvol, invR * P * daughterR, ModuleName);
592  }
594  }
595 
596  if (cexport)
597  {
598  ExportEnergyCuts(volumePtr);
599  }
600  // Add optional energy cuts
601 
602  // Here write the auxiliary info
603  //
604  auxiter = auxmap.find(volumePtr);
605  if (auxiter != auxmap.end())
606  {
607  AddAuxInfo(&(auxiter->second), volumeElement);
608  }
609 
610  structureElement->appendChild(volumeElement);
611  // Append the volume AFTER traversing the children so that
612  // the order of volumes will be correct!
613 
614  VolumeMap()[tmplv] = R;
615 
616  AddExtension(volumeElement, volumePtr);
617  // Add any possible user defined extension attached to a volume
618 
619  AddMaterial(volumePtr->GetMaterial());
620  // Add the involved materials and solids!
621 
622  AddSolid(solidPtr);
623 
624  SkinSurfaceCache(GetSkinSurface(volumePtr));
625 
626  return R;
627 }
628 
630  const G4LogicalVolume* const lvol)
631 {
632  std::map<const G4LogicalVolume*,
633  PHG4GDMLAuxListType>::iterator pos = auxmap.find(lvol);
634 
635  if (pos == auxmap.end())
636  {
637  auxmap[lvol] = PHG4GDMLAuxListType();
638  }
639 
640  auxmap[lvol].push_back(myaux);
641 }
642 
644 {
645  cexport = fcuts;
646 }
647 
648 void PHG4GDMLWriteStructure::ExportEnergyCuts(const G4LogicalVolume* const lvol)
649 {
650  // PHG4GDMLEvaluator eval;
651  G4ProductionCuts* pcuts = lvol->GetRegion()->GetProductionCuts();
652  G4ProductionCutsTable* ctab = G4ProductionCutsTable::GetProductionCutsTable();
653  G4Gamma* gamma = G4Gamma::Gamma();
654  G4Electron* eminus = G4Electron::Electron();
655  G4Positron* eplus = G4Positron::Positron();
656  G4Proton* proton = G4Proton::Proton();
657 
658  G4double gamma_cut = ctab->ConvertRangeToEnergy(gamma, lvol->GetMaterial(),
659  pcuts->GetProductionCut("gamma"));
660  G4double eminus_cut = ctab->ConvertRangeToEnergy(eminus, lvol->GetMaterial(),
661  pcuts->GetProductionCut("e-"));
662  G4double eplus_cut = ctab->ConvertRangeToEnergy(eplus, lvol->GetMaterial(),
663  pcuts->GetProductionCut("e+"));
664  G4double proton_cut = ctab->ConvertRangeToEnergy(proton, lvol->GetMaterial(),
665  pcuts->GetProductionCut("proton"));
666 
667  PHG4GDMLAuxStructType gammainfo = {"gammaECut",
668  ConvertToString(gamma_cut), "MeV", 0};
669  PHG4GDMLAuxStructType eminusinfo = {"electronECut",
670  ConvertToString(eminus_cut), "MeV", 0};
671  PHG4GDMLAuxStructType eplusinfo = {"positronECut",
672  ConvertToString(eplus_cut), "MeV", 0};
673  PHG4GDMLAuxStructType protinfo = {"protonECut",
674  ConvertToString(proton_cut), "MeV", 0};
675 
676  AddVolumeAuxiliary(gammainfo, lvol);
677  AddVolumeAuxiliary(eminusinfo, lvol);
678  AddVolumeAuxiliary(eplusinfo, lvol);
679  AddVolumeAuxiliary(protinfo, lvol);
680 }
681 
683 {
684  std::ostringstream os;
685  os << dval;
686  G4String vl = os.str();
687  return vl;
688 }