Geant4  10.02.p01
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: G4GDMLWriteStructure.cc 91789 2015-08-06 08:14:46Z gcosmo $
28 //
29 // class G4GDMLWriteStructure Implementation
30 //
31 // Original author: Zoltan Torzsok, November 2007
32 //
33 // --------------------------------------------------------------------
34 
35 #include "G4GDMLWriteStructure.hh"
36 #include "G4GDMLEvaluator.hh"
37 
38 #include "G4Material.hh"
39 #include "G4ReflectedSolid.hh"
40 #include "G4DisplacedSolid.hh"
41 #include "G4LogicalVolumeStore.hh"
42 #include "G4PhysicalVolumeStore.hh"
43 #include "G4ReflectionFactory.hh"
44 #include "G4PVDivision.hh"
45 #include "G4PVReplica.hh"
46 #include "G4Region.hh"
47 #include "G4OpticalSurface.hh"
48 #include "G4LogicalSkinSurface.hh"
50 
51 #include "G4ProductionCuts.hh"
52 #include "G4ProductionCutsTable.hh"
53 #include "G4Gamma.hh"
54 #include "G4Electron.hh"
55 #include "G4Positron.hh"
56 #include "G4Proton.hh"
57 
59  : G4GDMLWriteParamvol(), cexport(false)
60 {
62 }
63 
65 {
66 }
67 
68 void
69 G4GDMLWriteStructure::DivisionvolWrite(xercesc::DOMElement* volumeElement,
70  const G4PVDivision* const divisionvol)
71 {
72  EAxis axis = kUndefined;
73  G4int number = 0;
74  G4double width = 0.0;
75  G4double offset = 0.0;
76  G4bool consuming = false;
77 
78  divisionvol->GetReplicationData(axis,number,width,offset,consuming);
79  axis = divisionvol->GetDivisionAxis();
80 
81  G4String unitString("mm");
82  G4String axisString("kUndefined");
83  if (axis==kXAxis) { axisString = "kXAxis"; }
84  else if (axis==kYAxis) { axisString = "kYAxis"; }
85  else if (axis==kZAxis) { axisString = "kZAxis"; }
86  else if (axis==kRho) { axisString = "kRho"; }
87  else if (axis==kPhi) { axisString = "kPhi"; unitString = "rad"; }
88 
89  const G4String name
90  = GenerateName(divisionvol->GetName(),divisionvol);
91  const G4String volumeref
92  = GenerateName(divisionvol->GetLogicalVolume()->GetName(),
93  divisionvol->GetLogicalVolume());
94 
95  xercesc::DOMElement* divisionvolElement = NewElement("divisionvol");
96  divisionvolElement->setAttributeNode(NewAttribute("axis",axisString));
97  divisionvolElement->setAttributeNode(NewAttribute("number",number));
98  divisionvolElement->setAttributeNode(NewAttribute("width",width));
99  divisionvolElement->setAttributeNode(NewAttribute("offset",offset));
100  divisionvolElement->setAttributeNode(NewAttribute("unit",unitString));
101  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
102  volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
103  divisionvolElement->appendChild(volumerefElement);
104  volumeElement->appendChild(divisionvolElement);
105 }
106 
107 void G4GDMLWriteStructure::PhysvolWrite(xercesc::DOMElement* volumeElement,
108  const G4VPhysicalVolume* const physvol,
109  const G4Transform3D& T,
110  const G4String& ModuleName)
111 {
112  HepGeom::Scale3D scale;
113  HepGeom::Rotate3D rotate;
114  HepGeom::Translate3D translate;
115 
116  T.getDecomposition(scale,rotate,translate);
117 
118  const G4ThreeVector scl(scale(0,0),scale(1,1),scale(2,2));
119  const G4ThreeVector rot = GetAngles(rotate.getRotation());
120  const G4ThreeVector pos = T.getTranslation();
121 
122  const G4String name = GenerateName(physvol->GetName(),physvol);
123  const G4int copynumber = physvol->GetCopyNo();
124 
125  xercesc::DOMElement* physvolElement = NewElement("physvol");
126  physvolElement->setAttributeNode(NewAttribute("name",name));
127  if (copynumber) physvolElement->setAttributeNode(NewAttribute("copynumber", copynumber));
128 
129  volumeElement->appendChild(physvolElement);
130 
131  G4LogicalVolume* lv;
132  // Is it reflected?
133  if (reflFactory->IsReflected(physvol->GetLogicalVolume()))
134  {
136  }
137  else
138  {
139  lv = physvol->GetLogicalVolume();
140  }
141 
142  const G4String volumeref = GenerateName(lv->GetName(), lv);
143 
144  if (ModuleName.empty())
145  {
146  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
147  volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
148  physvolElement->appendChild(volumerefElement);
149  }
150  else
151  {
152  xercesc::DOMElement* fileElement = NewElement("file");
153  fileElement->setAttributeNode(NewAttribute("name",ModuleName));
154  fileElement->setAttributeNode(NewAttribute("volname",volumeref));
155  physvolElement->appendChild(fileElement);
156  }
157 
158  if (std::fabs(pos.x()) > kLinearPrecision
159  || std::fabs(pos.y()) > kLinearPrecision
160  || std::fabs(pos.z()) > kLinearPrecision)
161  {
162  PositionWrite(physvolElement,name+"_pos",pos);
163  }
164  if (std::fabs(rot.x()) > kAngularPrecision
165  || std::fabs(rot.y()) > kAngularPrecision
166  || std::fabs(rot.z()) > kAngularPrecision)
167  {
168  RotationWrite(physvolElement,name+"_rot",rot);
169  }
170  if (std::fabs(scl.x()-1.0) > kRelativePrecision
171  || std::fabs(scl.y()-1.0) > kRelativePrecision
172  || std::fabs(scl.z()-1.0) > kRelativePrecision)
173  {
174  ScaleWrite(physvolElement,name+"_scl",scl);
175  }
176 }
177 
178 void G4GDMLWriteStructure::ReplicavolWrite(xercesc::DOMElement* volumeElement,
179  const G4VPhysicalVolume* const replicavol)
180 {
181  EAxis axis = kUndefined;
182  G4int number = 0;
183  G4double width = 0.0;
184  G4double offset = 0.0;
185  G4bool consuming = false;
186  G4String unitString("mm");
187 
188  replicavol->GetReplicationData(axis,number,width,offset,consuming);
189 
190  const G4String volumeref
191  = GenerateName(replicavol->GetLogicalVolume()->GetName(),
192  replicavol->GetLogicalVolume());
193 
194  xercesc::DOMElement* replicavolElement = NewElement("replicavol");
195  replicavolElement->setAttributeNode(NewAttribute("number",number));
196  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
197  volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
198  replicavolElement->appendChild(volumerefElement);
199  xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
200  replicavolElement->appendChild(replicateElement);
201 
202  xercesc::DOMElement* dirElement = NewElement("direction");
203  if(axis==kXAxis)
204  { dirElement->setAttributeNode(NewAttribute("x","1")); }
205  else if(axis==kYAxis)
206  { dirElement->setAttributeNode(NewAttribute("y","1")); }
207  else if(axis==kZAxis)
208  { dirElement->setAttributeNode(NewAttribute("z","1")); }
209  else if(axis==kRho)
210  { dirElement->setAttributeNode(NewAttribute("rho","1")); }
211  else if(axis==kPhi)
212  { dirElement->setAttributeNode(NewAttribute("phi","1"));
213  unitString="rad"; }
214  replicateElement->appendChild(dirElement);
215 
216  xercesc::DOMElement* widthElement = NewElement("width");
217  widthElement->setAttributeNode(NewAttribute("value",width));
218  widthElement->setAttributeNode(NewAttribute("unit",unitString));
219  replicateElement->appendChild(widthElement);
220 
221  xercesc::DOMElement* offsetElement = NewElement("offset");
222  offsetElement->setAttributeNode(NewAttribute("value",offset));
223  offsetElement->setAttributeNode(NewAttribute("unit",unitString));
224  replicateElement->appendChild(offsetElement);
225 
226  volumeElement->appendChild(replicavolElement);
227 }
228 
231 {
232  if (!bsurf) { return; }
233 
234  const G4SurfaceProperty* psurf = bsurf->GetSurfaceProperty();
235 
236  // Generate the new element for border-surface
237  //
238  xercesc::DOMElement* borderElement = NewElement("bordersurface");
239  borderElement->setAttributeNode(NewAttribute("name", bsurf->GetName()));
240  borderElement->setAttributeNode(NewAttribute("surfaceproperty",
241  psurf->GetName()));
242 
243  const G4String volumeref1 = GenerateName(bsurf->GetVolume1()->GetName(),
244  bsurf->GetVolume1());
245  const G4String volumeref2 = GenerateName(bsurf->GetVolume2()->GetName(),
246  bsurf->GetVolume2());
247  xercesc::DOMElement* volumerefElement1 = NewElement("physvolref");
248  xercesc::DOMElement* volumerefElement2 = NewElement("physvolref");
249  volumerefElement1->setAttributeNode(NewAttribute("ref",volumeref1));
250  volumerefElement2->setAttributeNode(NewAttribute("ref",volumeref2));
251  borderElement->appendChild(volumerefElement1);
252  borderElement->appendChild(volumerefElement2);
253 
254  if (FindOpticalSurface(psurf))
255  {
256  const G4OpticalSurface* opsurf =
257  dynamic_cast<const G4OpticalSurface*>(psurf);
258  if (!opsurf)
259  {
260  G4Exception("G4GDMLWriteStructure::BorderSurfaceCache()",
261  "InvalidSetup", FatalException, "No optical surface found!");
262  return;
263  }
265  }
266 
267  borderElementVec.push_back(borderElement);
268 }
269 
272 {
273  if (!ssurf) { return; }
274 
275  const G4SurfaceProperty* psurf = ssurf->GetSurfaceProperty();
276 
277  // Generate the new element for border-surface
278  //
279  xercesc::DOMElement* skinElement = NewElement("skinsurface");
280  skinElement->setAttributeNode(NewAttribute("name", ssurf->GetName()));
281  skinElement->setAttributeNode(NewAttribute("surfaceproperty",
282  psurf->GetName()));
283 
284  const G4String volumeref = GenerateName(ssurf->GetLogicalVolume()->GetName(),
285  ssurf->GetLogicalVolume());
286  xercesc::DOMElement* volumerefElement = NewElement("volumeref");
287  volumerefElement->setAttributeNode(NewAttribute("ref",volumeref));
288  skinElement->appendChild(volumerefElement);
289 
290  if (FindOpticalSurface(psurf))
291  {
292  const G4OpticalSurface* opsurf =
293  dynamic_cast<const G4OpticalSurface*>(psurf);
294  if (!opsurf)
295  {
296  G4Exception("G4GDMLWriteStructure::SkinSurfaceCache()",
297  "InvalidSetup", FatalException, "No optical surface found!");
298  return;
299  }
301  }
302 
303  skinElementVec.push_back(skinElement);
304 }
305 
307 {
308  const G4OpticalSurface* osurf = dynamic_cast<const G4OpticalSurface*>(psurf);
309  std::vector<const G4OpticalSurface*>::const_iterator pos;
310  pos = std::find(opt_vec.begin(), opt_vec.end(), osurf);
311  if (pos != opt_vec.end()) { return false; } // item already created!
312 
313  opt_vec.push_back(osurf); // cache it for future reference
314  return true;
315 }
316 
319 {
320  G4LogicalSkinSurface* surf = 0;
322  if (nsurf)
323  {
324  const G4LogicalSkinSurfaceTable* stable =
326  std::vector<G4LogicalSkinSurface*>::const_iterator pos;
327  for (pos = stable->begin(); pos != stable->end(); pos++)
328  {
329  if (lvol == (*pos)->GetLogicalVolume())
330  {
331  surf = *pos; break;
332  }
333  }
334  }
335  return surf;
336 }
337 
340 {
341  G4LogicalBorderSurface* surf = 0;
343  if (nsurf)
344  {
345  const G4LogicalBorderSurfaceTable* btable =
347  std::vector<G4LogicalBorderSurface*>::const_iterator pos;
348  for (pos = btable->begin(); pos != btable->end(); pos++)
349  {
350  if (pvol == (*pos)->GetVolume1()) // just the first in the couple
351  { // is enough
352  surf = *pos; break;
353  }
354  }
355  }
356  return surf;
357 }
358 
360 {
361  G4cout << "G4GDML: Writing surfaces..." << G4endl;
362 
363  std::vector<xercesc::DOMElement*>::const_iterator pos;
364  for (pos = skinElementVec.begin(); pos != skinElementVec.end(); pos++)
365  {
366  structureElement->appendChild(*pos);
367  }
368  for (pos = borderElementVec.begin(); pos != borderElementVec.end(); pos++)
369  {
370  structureElement->appendChild(*pos);
371  }
372 }
373 
374 void G4GDMLWriteStructure::StructureWrite(xercesc::DOMElement* gdmlElement)
375 {
376  G4cout << "G4GDML: Writing structure..." << G4endl;
377 
378  structureElement = NewElement("structure");
379  gdmlElement->appendChild(structureElement);
380 }
381 
383 TraverseVolumeTree(const G4LogicalVolume* const volumePtr, const G4int depth)
384 {
385  if (VolumeMap().find(volumePtr) != VolumeMap().end())
386  {
387  return VolumeMap()[volumePtr]; // Volume is already processed
388  }
389 
390  G4VSolid* solidPtr = volumePtr->GetSolid();
391  G4Transform3D R,invR;
392  G4int trans=0;
393 
394  std::map<const G4LogicalVolume*, G4GDMLAuxListType>::iterator auxiter;
395 
396  while (true) // Solve possible displacement/reflection
397  { // of the referenced solid!
398  if (trans>maxTransforms)
399  {
400  G4String ErrorMessage = "Referenced solid in volume '"
401  + volumePtr->GetName()
402  + "' was displaced/reflected too many times!";
403  G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
404  "InvalidSetup", FatalException, ErrorMessage);
405  }
406 
407  if (G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
408  {
409  R = R*refl->GetTransform3D();
410  solidPtr = refl->GetConstituentMovedSolid();
411  trans++;
412  continue;
413  }
414 
415  if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
416  {
417  R = R*G4Transform3D(disp->GetObjectRotation(),
418  disp->GetObjectTranslation());
419  solidPtr = disp->GetConstituentMovedSolid();
420  trans++;
421  continue;
422  }
423 
424  break;
425  }
426 
427  //check if it is a reflected volume
428  G4LogicalVolume* tmplv = const_cast<G4LogicalVolume*>(volumePtr);
429 
430  if (reflFactory->IsReflected(tmplv))
431  {
432  tmplv = reflFactory->GetConstituentLV(tmplv);
433  if (VolumeMap().find(tmplv) != VolumeMap().end())
434  {
435  return R; // Volume is already processed
436  }
437  }
438 
439  // Only compute the inverse when necessary!
440  //
441  if (trans>0) { invR = R.inverse(); }
442 
443  const G4String name
444  = GenerateName(tmplv->GetName(), tmplv);
445  const G4String materialref
446  = GenerateName(volumePtr->GetMaterial()->GetName(),
447  volumePtr->GetMaterial());
448  const G4String solidref
449  = GenerateName(solidPtr->GetName(),solidPtr);
450 
451  xercesc::DOMElement* volumeElement = NewElement("volume");
452  volumeElement->setAttributeNode(NewAttribute("name",name));
453  xercesc::DOMElement* materialrefElement = NewElement("materialref");
454  materialrefElement->setAttributeNode(NewAttribute("ref",materialref));
455  volumeElement->appendChild(materialrefElement);
456  xercesc::DOMElement* solidrefElement = NewElement("solidref");
457  solidrefElement->setAttributeNode(NewAttribute("ref",solidref));
458  volumeElement->appendChild(solidrefElement);
459 
460  const G4int daughterCount = volumePtr->GetNoDaughters();
461 
462  for (G4int i=0;i<daughterCount;i++) // Traverse all the children!
463  {
464  const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
465  const G4String ModuleName = Modularize(physvol,depth);
466 
467  G4Transform3D daughterR;
468 
469  if (ModuleName.empty()) // Check if subtree requested to be
470  { // a separate module!
471  daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(),depth+1);
472  }
473  else
474  {
475  G4GDMLWriteStructure writer;
476  daughterR = writer.Write(ModuleName,physvol->GetLogicalVolume(),
477  SchemaLocation,depth+1);
478  }
479 
480  if (const G4PVDivision* const divisionvol
481  = dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
482  {
483  if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
484  {
485  G4String ErrorMessage = "Division volume in '" + name
486  + "' can not be related to reflected solid!";
487  G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
488  "InvalidSetup", FatalException, ErrorMessage);
489  }
490  DivisionvolWrite(volumeElement,divisionvol);
491  } else
492  if (physvol->IsParameterised()) // Is it a paramvol?
493  {
494  if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
495  {
496  G4String ErrorMessage = "Parameterised volume in '" + name
497  + "' can not be related to reflected solid!";
498  G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
499  "InvalidSetup", FatalException, ErrorMessage);
500  }
501  ParamvolWrite(volumeElement,physvol);
502  } else
503  if (physvol->IsReplicated()) // Is it a replicavol?
504  {
505  if (!G4Transform3D::Identity.isNear(invR*daughterR,kRelativePrecision))
506  {
507  G4String ErrorMessage = "Replica volume in '" + name
508  + "' can not be related to reflected solid!";
509  G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
510  "InvalidSetup", FatalException, ErrorMessage);
511  }
512  ReplicavolWrite(volumeElement,physvol);
513  }
514  else // Is it a physvol?
515  {
516  G4RotationMatrix rot;
517  if (physvol->GetFrameRotation() != 0)
518  {
519  rot = *(physvol->GetFrameRotation());
520  }
521  G4Transform3D P(rot,physvol->GetObjectTranslation());
522 
523  PhysvolWrite(volumeElement,physvol,invR*P*daughterR,ModuleName);
524  }
526  }
527 
528  if (cexport) { ExportEnergyCuts(volumePtr); }
529  // Add optional energy cuts
530 
531  // Here write the auxiliary info
532  //
533  auxiter = auxmap.find(volumePtr);
534  if (auxiter != auxmap.end())
535  {
536  AddAuxInfo(&(auxiter->second), volumeElement);
537  }
538 
539  structureElement->appendChild(volumeElement);
540  // Append the volume AFTER traversing the children so that
541  // the order of volumes will be correct!
542 
543  VolumeMap()[tmplv] = R;
544 
545  AddExtension(volumeElement, volumePtr);
546  // Add any possible user defined extension attached to a volume
547 
548  AddMaterial(volumePtr->GetMaterial());
549  // Add the involved materials and solids!
550 
551  AddSolid(solidPtr);
552 
553  SkinSurfaceCache(GetSkinSurface(volumePtr));
554 
555  return R;
556 }
557 
558 void
560  const G4LogicalVolume* const lvol)
561 {
562  std::map<const G4LogicalVolume*,
563  G4GDMLAuxListType>::iterator pos = auxmap.find(lvol);
564 
565  if (pos == auxmap.end()) { auxmap[lvol] = G4GDMLAuxListType(); }
566 
567  auxmap[lvol].push_back(myaux);
568 }
569 
570 void
572 {
573  cexport = fcuts;
574 }
575 
576 void
578 {
579  G4GDMLEvaluator eval;
580  G4ProductionCuts* pcuts = lvol->GetRegion()->GetProductionCuts();
582  G4Gamma* gamma = G4Gamma::Gamma();
583  G4Electron* eminus = G4Electron::Electron();
586 
587  G4double gamma_cut = ctab->ConvertRangeToEnergy(gamma, lvol->GetMaterial(),
588  pcuts->GetProductionCut("gamma"));
589  G4double eminus_cut = ctab->ConvertRangeToEnergy(eminus, lvol->GetMaterial(),
590  pcuts->GetProductionCut("e-"));
591  G4double eplus_cut = ctab->ConvertRangeToEnergy(eplus, lvol->GetMaterial(),
592  pcuts->GetProductionCut("e+"));
593  G4double proton_cut = ctab->ConvertRangeToEnergy(proton, lvol->GetMaterial(),
594  pcuts->GetProductionCut("proton"));
595 
596  G4GDMLAuxStructType gammainfo = {"gammaECut",
597  eval.ConvertToString(gamma_cut), "MeV", 0};
598  G4GDMLAuxStructType eminusinfo = {"electronECut",
599  eval.ConvertToString(eminus_cut), "MeV", 0};
600  G4GDMLAuxStructType eplusinfo = {"positronECut",
601  eval.ConvertToString(eplus_cut), "MeV", 0};
602  G4GDMLAuxStructType protinfo = {"protonECut",
603  eval.ConvertToString(proton_cut), "MeV", 0};
604 
605  AddVolumeAuxiliary(gammainfo, lvol);
606  AddVolumeAuxiliary(eminusinfo, lvol);
607  AddVolumeAuxiliary(eplusinfo, lvol);
608  AddVolumeAuxiliary(protinfo, lvol);
609 }
static const G4double kAngularPrecision
Definition: geomdefs.hh:54
G4String GetName() const
G4ReflectionFactory * reflFactory
G4ProductionCuts * GetProductionCuts() const
void AddAuxInfo(G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
Definition: G4GDMLWrite.cc:90
std::vector< G4LogicalBorderSurface * > G4LogicalBorderSurfaceTable
G4Transform3D TraverseVolumeTree(const G4LogicalVolume *const topVol, const G4int depth)
xercesc::DOMElement * structureElement
const G4LogicalVolume * GetLogicalVolume() const
CLHEP::Hep3Vector G4ThreeVector
G4String SchemaLocation
Definition: G4GDMLWrite.hh:131
CLHEP::HepRotation G4RotationMatrix
virtual G4bool IsReplicated() const =0
G4Material * GetMaterial() const
const G4VPhysicalVolume * GetVolume2() const
G4double GetProductionCut(G4int index) const
const G4LogicalSkinSurface * GetSkinSurface(const G4LogicalVolume *const)
G4String name
Definition: TRTMaterials.hh:40
const G4String & GetName() const
Definition: G4Material.hh:178
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4LogicalVolume * GetConstituentLV(G4LogicalVolume *reflLV) const
#define width
const G4String & GetName() const
EAxis GetDivisionAxis() const
void AddMaterial(const G4Material *const)
virtual void ParamvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
G4Region * GetRegion() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
int G4int
Definition: G4Types.hh:78
std::vector< const G4OpticalSurface * > opt_vec
G4bool IsReflected(G4LogicalVolume *lv) const
static size_t GetNumberOfSkinSurfaces()
static double P[]
static G4ReflectionFactory * Instance()
static size_t GetNumberOfBorderSurfaces()
const G4RotationMatrix * GetFrameRotation() const
static const G4int maxTransforms
G4String Modularize(const G4VPhysicalVolume *const topvol, const G4int depth)
Definition: G4GDMLWrite.cc:329
static const G4LogicalSkinSurfaceTable * GetSurfaceTable()
G4Transform3D Write(const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
Definition: G4GDMLWrite.cc:166
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
static const G4LogicalBorderSurfaceTable * GetSurfaceTable()
G4String ConvertToString(G4int ival)
static const G4double kRelativePrecision
bool G4bool
Definition: G4Types.hh:79
virtual void StructureWrite(xercesc::DOMElement *)
void OpticalSurfaceWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
static G4Proton * Proton()
Definition: G4Proton.cc:93
HepGeom::Transform3D G4Transform3D
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void BorderSurfaceCache(const G4LogicalBorderSurface *const)
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsParameterised() const =0
void ReplicavolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
const G4VPhysicalVolume * GetVolume1() const
void PhysvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const topVol, const G4Transform3D &transform, const G4String &moduleName)
static G4ProductionCutsTable * GetProductionCutsTable()
static G4Positron * Positron()
Definition: G4Positron.cc:94
std::map< const G4LogicalVolume *, G4GDMLAuxListType > auxmap
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4LogicalVolume * GetLogicalVolume() const
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
static const G4double kLinearPrecision
EAxis
Definition: geomdefs.hh:54
void AddVolumeAuxiliary(G4GDMLAuxStructType myaux, const G4LogicalVolume *const)
G4ThreeVector GetAngles(const G4RotationMatrix &)
virtual G4int GetCopyNo() const =0
void SkinSurfaceCache(const G4LogicalSkinSurface *const)
const G4LogicalBorderSurface * GetBorderSurface(const G4VPhysicalVolume *const)
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
virtual void AddExtension(xercesc::DOMElement *, const G4LogicalVolume *const)
Definition: G4GDMLWrite.cc:78
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4SurfaceProperty * GetSurfaceProperty() const
G4double ConvertRangeToEnergy(const G4ParticleDefinition *particle, const G4Material *material, G4double range)
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
G4ThreeVector GetObjectTranslation() const
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
static const double eplus
Definition: G4SIunits.hh:196
std::vector< G4LogicalSkinSurface * > G4LogicalSkinSurfaceTable
void DivisionvolWrite(xercesc::DOMElement *, const G4PVDivision *const)
std::vector< xercesc::DOMElement * > skinElementVec
void ExportEnergyCuts(const G4LogicalVolume *const)
xercesc::DOMElement * solidsElement
static const G4double pos
VolumeMapType & VolumeMap()
Definition: G4GDMLWrite.cc:60
G4VSolid * GetSolid() const
virtual void AddSolid(const G4VSolid *const)
std::vector< xercesc::DOMElement * > borderElementVec
G4bool FindOpticalSurface(const G4SurfaceProperty *)