CMS 3D CMS Logo

PrintGeomMatInfo.cc
Go to the documentation of this file.
2 
6 
16 
17 #include "G4Run.hh"
18 #include "G4PhysicalVolumeStore.hh"
19 #include "G4LogicalVolumeStore.hh"
20 #include "G4VPhysicalVolume.hh"
21 #include "G4LogicalVolume.hh"
22 #include "G4VSolid.hh"
23 #include "G4Material.hh"
24 #include "G4Track.hh"
25 #include "G4VisAttributes.hh"
26 #include "G4UserLimits.hh"
27 #include "G4TransportationManager.hh"
28 
29 #include <set>
30 #include <map>
31 
32 using namespace CLHEP;
33 
35 {
36  _dumpSummary = p.getUntrackedParameter<bool>("DumpSummary", true);
37  _dumpLVTree = p.getUntrackedParameter<bool>("DumpLVTree", true);
38  _dumpLVMatBudget = p.getUntrackedParameter<bool>("DumpLVMatBudget", false);
39  _lvNames2Dump= p.getUntrackedParameter<std::vector<std::string> >("LVNames2Dump");
40  _level2Dump = 0;
41  _dumpIndex = 0;
42  _dumpIt = false;
43  _radiusLayer = p.getUntrackedParameter<std::vector<double> >("Radius2Use");
44  _zLayer = p.getUntrackedParameter<std::vector<double> >("Z2Use");
45  _maxLevelsCounted = 50;
46  _countsPerLevel.assign(_maxLevelsCounted,0);
47  _dumpMaterial= p.getUntrackedParameter<bool>("DumpMaterial",false);
48  _dumpLVList = p.getUntrackedParameter<bool>("DumpLVList", false);
49  _dumpLV = p.getUntrackedParameter<bool>("DumpLV", false);
50  _dumpSolid = p.getUntrackedParameter<bool>("DumpSolid", false);
51  _dumpAtts = p.getUntrackedParameter<bool>("DumpAttributes", false);
52  _dumpPV = p.getUntrackedParameter<bool>("DumpPV", false);
53  _dumpRotation= p.getUntrackedParameter<bool>("DumpRotation",false);
54  _dumpReplica = p.getUntrackedParameter<bool>("DumpReplica", false);
55  _dumpTouch = p.getUntrackedParameter<bool>("DumpTouch", false);
56  _dumpSense = p.getUntrackedParameter<bool>("DumpSense", false);
57  name = p.getUntrackedParameter<std::string>("Name","*");
58  nchar = name.find("*");
59  name.assign(name,0,nchar);
60  names = p.getUntrackedParameter<std::vector<std::string> >("Names");
61  std::cout << "size of _lvNames2Dump = " << _lvNames2Dump.size()
62  << " size of _radiusLayer = " << _radiusLayer.size()
63  << " size of _zLayer = " << _zLayer.size() << std::endl;
64  std::cout << "PrintGeomMatInfo:: initialised with verbosity levels:"
65  << " Summary " << _dumpSummary << " LVTree " << _dumpLVTree
66  << " LVList " << _dumpLVList << " Material " << _dumpMaterial
67  << "\n "
68  << " LVMatBudget " << _dumpLVMatBudget << " for";
69  _areaLayer.reserve(_lvNames2Dump.size());
70  if(_lvNames2Dump.size() == _radiusLayer.size() &&
71  _lvNames2Dump.size() == _zLayer.size()) {
72  for (unsigned int i=0; i<_lvNames2Dump.size(); i++) {
73  _areaLayer[i] = 2.0*3.14159*_radiusLayer[i]*_zLayer[i];
74  std::cout << "\n "
75  << " " << _lvNames2Dump[i] << " radius = " << _radiusLayer[i]
76  << " z = " << _zLayer[i] << " area = " << _areaLayer[i];
77  }
78  } else {
79  _areaLayer.assign(3,0.0);
80  std::cout << "\n "
81  << " Problem with unequal sizes!! Fix and rerun";
82  }
83  std::cout << "\n "
84  << " and max levels for count = " << _maxLevelsCounted;
85  std::cout << "\n "
86  << " LV " << _dumpLV << " Solid " << _dumpSolid
87  << " Attribs " << _dumpAtts
88  << "\n "
89  << " PV " << _dumpPV << " Rotation " << _dumpRotation
90  << " Replica " << _dumpReplica
91  << "\n "
92  << " Touchable " << _dumpTouch << " for names (0-" << nchar
93  << ") = " << name
94  << "\n "
95  << " Sensitive " << _dumpSense << " for " << names.size()
96  << " names";
97  for (unsigned int i=0; i<names.size(); i++) std::cout << " " << names[i];
98  std::cout << std::endl;
99 }
100 
102 
104 {
105  if (_dumpSense) {
107  (*job)()->get<IdealGeometryRecord>().get(pDD);
108 
109  std::cout << "PrintGeomMatInfo::Get Printout of Sensitive Volumes "
110  << "for " << names.size() << " Readout Units" << std::endl;
111  for (unsigned int i=0; i<names.size(); i++) {
112  std::string attribute = "ReadOutName";
113  std::string sd = names[i];
115  DDFilteredView fv(*pDD,filter);
116  std::cout << "PrintGeomMatInfo:: Get Filtered view for "
117  << attribute << " = " << sd << std::endl;
118  bool dodet = fv.firstChild();
119 
120  std::string spaces = spacesFromLeafDepth(1);
121 
122  while (dodet) {
123  const DDLogicalPart & log = fv.logicalPart();
124  std::string lvname = DDSplit(log.name()).first;
125  DDTranslation tran = fv.translation();
126  std::vector<int> copy = fv.copyNumbers();
127 
128  unsigned int leafDepth = copy.size();
129  std::cout << leafDepth << spaces << "### VOLUME = " << lvname
130  << " Copy No";
131  for (int k=leafDepth-1; k>=0; k--) std::cout << " " << copy[k];
132  std::cout << " Centre at " << tran << " (r = " << tran.Rho()
133  << ", phi = " << tran.phi()/deg << ")" << std::endl;
134  dodet = fv.next();
135  }
136  }
137  }
138 }
139 
141 {
142  theTopPV = getTopPV();
143 
144  if (_dumpSummary) dumpSummary(std::cout);
145  if (_dumpLVTree) dumpG4LVTree(std::cout);
146  if (_dumpLVMatBudget) dumpG4LVMatBudget(std::cout);
147 
148  //---------- Dump list of objects of each class with detail of parameters
149  if (_dumpMaterial) dumpMaterialList(std::cout);
150  if (_dumpLVList) dumpG4LVList(std::cout);
151 
152  //---------- Dump LV and PV information
153  if (_dumpLV || _dumpPV || _dumpTouch) dumpHierarchyTreePVLV(std::cout);
154 }
155 
156 void PrintGeomMatInfo::dumpSummary(std::ostream & out)
157 {
158  //---------- Dump number of objects of each class
159  out << " @@@@@@@@@@@@@@@@@@ Dumping G4 geometry objects Summary " << std::endl;
160  if(theTopPV == 0)
161  {
162  out << " No volume created " << std::endl;
163  return;
164  }
165  out << " @@@ Geometry built inside world volume: " << theTopPV->GetName() << std::endl;
166  // Get number of solids (< # LV if several LV share a solid)
167  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
168  std::vector<G4LogicalVolume *>::const_iterator lvcite;
169  std::set<G4VSolid *> theSolids;
170  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
171  theSolids.insert((*lvcite)->GetSolid());
172  out << " Number of G4VSolid's: " << theSolids.size() << std::endl;
173  out << " Number of G4LogicalVolume's: " << lvs->size() << std::endl;
174  const G4PhysicalVolumeStore * pvs = G4PhysicalVolumeStore::GetInstance();
175  out << " Number of G4VPhysicalVolume's: " << pvs->size() << std::endl;
176  out << " Number of Touchable's: " << countNoTouchables() << std::endl;
177  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
178  out << " Number of G4Material's: " << matTab->size() << std::endl;
179 }
180 
182 {
183  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's List " << std::endl;
184  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
185  std::vector<G4LogicalVolume*>::const_iterator lvcite;
186  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
187  out << "LV:" << (*lvcite)->GetName() << "\tMaterial: " << (*lvcite)->GetMaterial()->GetName() << std::endl;
188 }
189 
191 {
192  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's Tree " << std::endl;
193  G4LogicalVolume * lv = getTopLV();
194  dumpG4LVLeaf(lv,0,1,out);
195 }
196 
198 {
199  out << " @@@@@@@@@@@@@@@@ DUMPING G4Material List ";
200  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
201  out << " with " << matTab->size() << " materials " << std::endl;
202  std::vector<G4Material*>::const_iterator matite;
203  for (matite = matTab->begin(); matite != matTab->end(); matite++)
204  out << "Material: " << (*matite) << std::endl;
205 }
206 
207 void PrintGeomMatInfo::dumpG4LVLeaf(G4LogicalVolume * lv, unsigned int leafDepth, unsigned int count, std::ostream & out)
208 {
209  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
210  out << " LV:(" << leafDepth << ") " << lv->GetName() << " (" << count
211  << ")" << std::endl;
212  //--- If a volume is placed n types as daughter of this LV, it should only be counted once
213  std::map<G4LogicalVolume*, unsigned int> lvCount;
214  std::map<G4LogicalVolume*, unsigned int>::const_iterator cite;
215  for (int ii = 0; ii < lv->GetNoDaughters(); ii++) {
216  cite = lvCount.find(lv->GetDaughter(ii)->GetLogicalVolume());
217  if (cite != lvCount.end()) lvCount[cite->first] = (cite->second) + 1;
218  else lvCount.insert(std::pair< G4LogicalVolume*,unsigned int>(lv->GetDaughter(ii)->GetLogicalVolume(),1));
219  }
220  for (cite = lvCount.begin(); cite != lvCount.end(); cite++)
221  dumpG4LVLeaf((cite->first), leafDepth+1, (cite->second), out);
222 }
223 
225 {
226  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's Material Budget Tree " << std::endl;
227  G4LogicalVolume * lv = getTopLV();
228  dumpG4LVLeafWithMat(lv,0,1,out);
229 }
230 
231 void PrintGeomMatInfo::dumpG4LVLeafWithMat(G4LogicalVolume * lv, unsigned int leafDepth, unsigned int count, std::ostream & out)
232 {
233  // switch off dumping at the next same level as the dump
234  if(_dumpIt && _level2Dump == leafDepth) {
235  _dumpIt = false;
236 // std::cout << " stopped dumping " << std::endl;
237  }
238  // loop over volumes to dump: volumes cannot not overlap in hierarchy!!
239  for (unsigned int i=0; i<_lvNames2Dump.size(); i++) {
240  if(_lvNames2Dump[i].compare(lv->GetName()) == 0) {
241  _dumpIt = true;
242  _dumpIndex = i;
243  _level2Dump = leafDepth;
244 // std::cout << " start dumping " << _lvNames2Dump[i] << " at level " << _level2Dump << std::endl;
245  break;
246  }
247  }
248 
249  if(_dumpIt) {
250  if(leafDepth < _maxLevelsCounted) _countsPerLevel[leafDepth] = count;
251  unsigned int total_multipler = 1;
252  for (unsigned int ii=_level2Dump; ii <= leafDepth; ii++) total_multipler *= _countsPerLevel[ii];
253  double thick = (lv->GetSolid()->GetCubicVolume() * total_multipler)/_areaLayer[_dumpIndex];
254  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
255  // print out Level, Logical volume name, volume in mm**3, material name, rad len of mat in mm
256  // total number of volumes from dump level start, equivalent thickness when spread over
257  // a cylinder of _radiusLayer and length _zLayer; equivalent thick in rad len
258  out << " LV::" << leafDepth << ": " << lv->GetName() << " :" << count
259  << ": " << lv->GetSolid()->GetName() << " :" << lv->GetSolid()->GetCubicVolume()
260  << ": "<< lv->GetMaterial()->GetName() << " :" << lv->GetMaterial()->GetRadlen()
261  << ":" << " :" << total_multipler << ":"
262  << " thk :" << thick << ": x/X0 :" << thick/lv->GetMaterial()->GetRadlen()
263  << ":" << " Kg : " << lv->GetMass()/kg << std::endl;
264  } else {
265  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
266  out << " LV:(" << leafDepth << ") " << lv->GetName() << " (" << count << ")" << std::endl;
267  }
268 
269  //--- If a volume is placed n types as daughter of this LV, it should only be counted once
270  std::map<G4LogicalVolume*, unsigned int> lvCount;
271  std::map<G4LogicalVolume*, unsigned int>::const_iterator cite;
272  for (int ii = 0; ii < lv->GetNoDaughters(); ii++) {
273  cite = lvCount.find(lv->GetDaughter(ii)->GetLogicalVolume());
274  if (cite != lvCount.end()) lvCount[cite->first] = (cite->second) + 1;
275  else lvCount.insert(std::pair< G4LogicalVolume*,unsigned int>(lv->GetDaughter(ii)->GetLogicalVolume(),1));
276  }
277  for (cite = lvCount.begin(); cite != lvCount.end(); cite++)
278  dumpG4LVLeafWithMat((cite->first), leafDepth+1, (cite->second), out);
279 }
280 
281 
283 {
284  int nTouch = 0;
285  G4LogicalVolume * lv = getTopLV();
286  add1touchable(lv, nTouch);
287  return nTouch;
288 }
289 
290 
291 void PrintGeomMatInfo::add1touchable(G4LogicalVolume * lv, int & nTouch)
292 {
293  int siz = lv->GetNoDaughters();
294  for(int ii = 0; ii < siz; ii++)
295  add1touchable(lv->GetDaughter(ii)->GetLogicalVolume(), ++nTouch);
296 }
297 
299 {
300  //dumps in the following order:
301  // 1) a LV with details
302  // 2) list of PVs daughters of this LV with details
303  // 3) list of LVs daughters of this LV and for each go to 1)
304 
305  //----- Get top PV
306  G4LogicalVolume* topLV = getTopLV();
307 
308  //----- Dump this leaf (it will recursively dump all the tree)
309  dumpHierarchyLeafPVLV(topLV, 0, out);
310  dumpPV(theTopPV, 0, out);
311 
312  //----- Dump the touchables (it will recursively dump all the tree)
313  if (_dumpTouch) dumpTouch(theTopPV, 0, out);
314 }
315 
316 void PrintGeomMatInfo::dumpHierarchyLeafPVLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out)
317 {
318  //----- Dump this LV
319  dumpLV(lv, leafDepth, out);
320 
321  //----- Get LV daughters from list of PV daughters
322  mmlvpv lvpvDaughters;
323  std::set< G4LogicalVolume * > lvDaughters;
324  int NoDaughters = lv->GetNoDaughters();
325  while ((NoDaughters--)>0)
326  {
327  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
328  lvpvDaughters.insert(mmlvpv::value_type(pvD->GetLogicalVolume(), pvD));
329  lvDaughters.insert(pvD->GetLogicalVolume());
330  }
331 
332  std::set< G4LogicalVolume * >::const_iterator scite;
333  mmlvpv::const_iterator mmcite;
334 
335  //----- Dump daughters PV and LV
336  for (scite = lvDaughters.begin(); scite != lvDaughters.end(); scite++)
337  {
338  std::pair< mmlvpv::iterator, mmlvpv::iterator > mmER = lvpvDaughters.equal_range(*scite);
339  //----- Dump daughters PV of this LV
340  for (mmcite = mmER.first ; mmcite != mmER.second; mmcite++)
341  dumpPV((*mmcite).second, leafDepth+1, out);
342  //----- Dump daughters LV
343  dumpHierarchyLeafPVLV(*scite, leafDepth+1, out );
344  }
345 }
346 
347 void PrintGeomMatInfo::dumpLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out)
348 {
349  std::string spaces = spacesFromLeafDepth(leafDepth);
350 
351  //----- dump name
352  if (_dumpLV) {
353  out << leafDepth << spaces << "$$$ VOLUME = " << lv->GetName()
354  << " Solid: " << lv->GetSolid()->GetName() << " MATERIAL: "
355  << lv->GetMaterial()->GetName() << std::endl;
356  if (_dumpSolid)
357  dumpSolid(lv->GetSolid(), leafDepth, out); //----- dump solid
358 
359  //----- dump LV info
360  //--- material
361  if (_dumpAtts) {
362  //--- Visualisation attributes
363  const G4VisAttributes * fVA = lv->GetVisAttributes();
364  if (fVA!=0) {
365  out << spaces << " VISUALISATION ATTRIBUTES: " << std::endl;
366  out << spaces << " IsVisible " << fVA->IsVisible() << std::endl;
367  out << spaces << " IsDaughtersInvisible " << fVA->IsDaughtersInvisible() << std::endl;
368  out << spaces << " Colour " << fVA->GetColour() << std::endl;
369  out << spaces << " LineStyle " << fVA->GetLineStyle() << std::endl;
370  out << spaces << " LineWidth " << fVA->GetLineWidth() << std::endl;
371  out << spaces << " IsForceDrawingStyle " << fVA->IsForceDrawingStyle() << std::endl;
372  out << spaces << " ForcedDrawingStyle " << fVA->GetForcedDrawingStyle() << std::endl;
373  }
374 
375  //--- User Limits
376  G4UserLimits * fUL = lv->GetUserLimits();
377  G4Track dummy;
378  if (fUL!=0) {
379  out << spaces << " MaxAllowedStep " << fUL->GetMaxAllowedStep(dummy) << std::endl;
380  out << spaces << " UserMaxTrackLength " << fUL->GetUserMaxTrackLength(dummy) << std::endl;
381  out << spaces << " UserMaxTime " << fUL->GetUserMaxTime(dummy) << std::endl;
382  out << spaces << " UserMinEkine " << fUL->GetUserMinEkine(dummy) << std::endl;
383  out << spaces << " UserMinRange " << fUL->GetUserMinRange(dummy) << std::endl;
384  }
385 
386  //--- other LV info
387  if (lv->GetSensitiveDetector())
388  out << spaces << " IS SENSITIVE DETECTOR " << std::endl;
389  if (lv->GetFieldManager())
390  out << spaces << " FIELD ON " << std::endl;
391 
392  // Pointer (possibly NULL) to optimisation info objects.
393  out << spaces
394  << " Quality for optimisation, average number of voxels to be spent per content "
395  << lv->GetSmartless() << std::endl;
396 
397  // Pointer (possibly NULL) to G4FastSimulationManager object.
398  if (lv->GetFastSimulationManager())
399  out << spaces << " Logical Volume is an envelope for a FastSimulationManager "
400  << std::endl;
401  out << spaces << " Weight used in the event biasing technique = "
402  << lv->GetBiasWeight() << std::endl;
403  }
404  }
405 }
406 
407 void PrintGeomMatInfo::dumpPV(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out)
408 {
409  std::string spaces = spacesFromLeafDepth(leafDepth);
410 
411  //----- PV info
412  if (_dumpPV)
413  {
414  std::string mother = "World";
415  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
416  out << leafDepth << spaces << "### VOLUME = " << pv->GetName()
417  << " Copy No " << pv->GetCopyNo() << " in " << mother
418  << " at " << pv->GetTranslation();
419  }
420  if (!pv->IsReplicated())
421  {
422  if (_dumpPV)
423  {
424  if(pv->GetRotation() == 0) out << " with no rotation" << std::endl;
425  else if(!_dumpRotation) out << " with rotation" << std::endl; //just rotation name
426  else out << " with rotation " << *(pv->GetRotation()) << std::endl;
427  }
428  }
429  else
430  {
431  if (_dumpReplica )
432  {
433  out << spaces << " It is replica: " << std::endl;
434  EAxis axis;
435  int nReplicas;
436  double width;
437  double offset;
438  bool consuming;
439  pv->GetReplicationData(axis, nReplicas, width, offset, consuming);
440  out << spaces << " axis " << axis << std::endl
441  << spaces << " nReplicas " << nReplicas << std::endl;
442  if (pv->GetParameterisation() != 0)
443  out << spaces << " It is parameterisation " << std::endl;
444  else
445  out << spaces << " width " << width << std::endl
446  << spaces << " offset " << offset << std::endl
447  << spaces << " consuming" << consuming << std::endl;
448  if (pv->GetParameterisation() != 0)
449  out << spaces << " It is parameterisation " << std::endl;
450  }
451  }
452 }
453 
454 void PrintGeomMatInfo::dumpTouch(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out)
455 {
456  std::string spaces = spacesFromLeafDepth(leafDepth);
457  if (leafDepth == 0) fHistory.SetFirstEntry(pv);
458  else fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
459 
460  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().
461  TransformPoint(G4ThreeVector(0,0,0));
462  G4LogicalVolume * lv = pv->GetLogicalVolume();
463 
464  std::string mother = "World";
465  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
466  std::string lvname = lv->GetName();
467  lvname.assign(lvname,0,nchar);
468  if (lvname == name)
469  out << leafDepth << spaces << "### VOLUME = " << lv->GetName()
470  << " Copy No " << pv->GetCopyNo() << " in " << mother
471  << " global position of centre " << globalpoint << " (r = "
472  << globalpoint.perp() << ", phi = " << globalpoint.phi()/deg
473  << ")" << std::endl;
474 
475  int NoDaughters = lv->GetNoDaughters();
476  while ((NoDaughters--)>0)
477  {
478  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
479  if (!pvD->IsReplicated()) dumpTouch(pvD, leafDepth+1, out);
480  }
481 
482  if (leafDepth > 0) fHistory.BackLevel();
483 }
484 
486 {
487  std::string spaces;
488  unsigned int ii;
489  for(ii = 0; ii < leafDepth; ii++) { spaces += " "; }
490  return spaces;
491 }
492 
493 void PrintGeomMatInfo::dumpSolid(G4VSolid * sol, unsigned int leafDepth, std::ostream & out)
494 {
495  std::string spaces = spacesFromLeafDepth(leafDepth);
496  out << spaces << *(sol) << std::endl;
497 }
498 
499 G4VPhysicalVolume * PrintGeomMatInfo::getTopPV()
500 {
501  return G4TransportationManager::GetTransportationManager()
502  ->GetNavigatorForTracking()->GetWorldVolume();
503 }
504 
505 G4LogicalVolume * PrintGeomMatInfo::getTopLV()
506 { return theTopPV->GetLogicalVolume(); }
507 
508 
T getUntrackedParameter(std::string const &, T const &) const
bool compare(const P &i, const P &j)
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
G4VPhysicalVolume * getTopPV()
void dumpG4LVTree(std::ostream &out=std::cout)
const N & name() const
Definition: DDBase.h:78
static const HistoName names[]
std::multimap< G4LogicalVolume *, G4VPhysicalVolume *, std::less< G4LogicalVolume * > > mmlvpv
void dumpG4LVLeafWithMat(G4LogicalVolume *lv, unsigned int leafDepth, unsigned int count, std::ostream &out=std::cout)
nav_type copyNumbers() const
return the stack of copy numbers
void dumpMaterialList(std::ostream &out=std::cout)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DDTranslation
Definition: DDTranslation.h:7
void add1touchable(G4LogicalVolume *lv, int &nTouch)
void dumpHierarchyLeafPVLV(G4LogicalVolume *lv, unsigned int leafDepth, std::ostream &out=std::cout)
void dumpG4LVList(std::ostream &out=std::cout)
void dumpG4LVMatBudget(std::ostream &out=std::cout)
void dumpHierarchyTreePVLV(std::ostream &out=std::cout)
void dumpLV(G4LogicalVolume *lv, unsigned int leafDepth, std::ostream &out=std::cout)
bool next()
set current node to the next node in the filtered tree
std::string spacesFromLeafDepth(unsigned int leafDepth)
void dumpPV(G4VPhysicalVolume *pv, unsigned int leafDepth, std::ostream &out=std::cout)
def pv(vc)
Definition: MetAnalyzer.py:6
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:92
void dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, std::ostream &out=std::cout)
ii
Definition: cuy.py:588
void dumpSummary(std::ostream &out=std::cout)
int k[5][pyjets_maxn]
G4LogicalVolume * getTopLV()
double sd
PrintGeomMatInfo(edm::ParameterSet const &p)
void update(const BeginOfJob *job)
This routine will be called when the appropriate signal arrives.
void dumpSolid(G4VSolid *sol, unsigned int leafDepth, std::ostream &out=std::cout)
bool firstChild()
set the current node to the first child ...
void dumpG4LVLeaf(G4LogicalVolume *lv, unsigned int leafDepth, unsigned int count, std::ostream &out=std::cout)
std::pair< std::string, std::string > DDSplit(const std::string &n)
split into (name,namespace), separator = &#39;:&#39;
Definition: DDSplit.cc:4
const DDTranslation & translation() const
The absolute translation of the current node.