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