CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  DDValue ddv(attribute,sd,0);
117  DDFilteredView fv(*pDD);
118  std::cout << "PrintGeomMatInfo:: Get Filtered view for "
119  << attribute << " = " << sd << std::endl;
120  fv.addFilter(filter);
121  bool dodet = fv.firstChild();
122 
123  std::string spaces = spacesFromLeafDepth(1);
124 
125  while (dodet) {
126  const DDLogicalPart & log = fv.logicalPart();
127  std::string lvname = DDSplit(log.name()).first;
128  DDTranslation tran = fv.translation();
129  std::vector<int> copy = fv.copyNumbers();
130 
131  unsigned int leafDepth = copy.size();
132  std::cout << leafDepth << spaces << "### VOLUME = " << lvname
133  << " Copy No";
134  for (int k=leafDepth-1; k>=0; k--) std::cout << " " << copy[k];
135  std::cout << " Centre at " << tran << " (r = " << tran.Rho()
136  << ", phi = " << tran.phi()/deg << ")" << std::endl;
137  dodet = fv.next();
138  }
139  }
140  }
141 }
142 
144 {
145  theTopPV = getTopPV();
146 
147  if (_dumpSummary) dumpSummary(std::cout);
148  if (_dumpLVTree) dumpG4LVTree(std::cout);
149  if (_dumpLVMatBudget) dumpG4LVMatBudget(std::cout);
150 
151  //---------- Dump list of objects of each class with detail of parameters
152  if (_dumpMaterial) dumpMaterialList(std::cout);
153  if (_dumpLVList) dumpG4LVList(std::cout);
154 
155  //---------- Dump LV and PV information
156  if (_dumpLV || _dumpPV || _dumpTouch) dumpHierarchyTreePVLV(std::cout);
157 }
158 
159 void PrintGeomMatInfo::dumpSummary(std::ostream & out)
160 {
161  //---------- Dump number of objects of each class
162  out << " @@@@@@@@@@@@@@@@@@ Dumping G4 geometry objects Summary " << std::endl;
163  if(theTopPV == 0)
164  {
165  out << " No volume created " << std::endl;
166  return;
167  }
168  out << " @@@ Geometry built inside world volume: " << theTopPV->GetName() << std::endl;
169  // Get number of solids (< # LV if several LV share a solid)
170  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
171  std::vector<G4LogicalVolume *>::const_iterator lvcite;
172  std::set<G4VSolid *> theSolids;
173  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
174  theSolids.insert((*lvcite)->GetSolid());
175  out << " Number of G4VSolid's: " << theSolids.size() << std::endl;
176  out << " Number of G4LogicalVolume's: " << lvs->size() << std::endl;
177  const G4PhysicalVolumeStore * pvs = G4PhysicalVolumeStore::GetInstance();
178  out << " Number of G4VPhysicalVolume's: " << pvs->size() << std::endl;
179  out << " Number of Touchable's: " << countNoTouchables() << std::endl;
180  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
181  out << " Number of G4Material's: " << matTab->size() << std::endl;
182 }
183 
185 {
186  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's List " << std::endl;
187  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
188  std::vector<G4LogicalVolume*>::const_iterator lvcite;
189  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
190  out << "LV:" << (*lvcite)->GetName() << "\tMaterial: " << (*lvcite)->GetMaterial()->GetName() << std::endl;
191 }
192 
194 {
195  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's Tree " << std::endl;
196  G4LogicalVolume * lv = getTopLV();
197  dumpG4LVLeaf(lv,0,1,out);
198 }
199 
201 {
202  out << " @@@@@@@@@@@@@@@@ DUMPING G4Material List ";
203  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
204  out << " with " << matTab->size() << " materials " << std::endl;
205  std::vector<G4Material*>::const_iterator matite;
206  for (matite = matTab->begin(); matite != matTab->end(); matite++)
207  out << "Material: " << (*matite) << std::endl;
208 }
209 
210 void PrintGeomMatInfo::dumpG4LVLeaf(G4LogicalVolume * lv, unsigned int leafDepth, unsigned int count, std::ostream & out)
211 {
212  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
213  out << " LV:(" << leafDepth << ") " << lv->GetName() << " (" << count
214  << ")" << std::endl;
215  //--- If a volume is placed n types as daughter of this LV, it should only be counted once
216  std::map<G4LogicalVolume*, unsigned int> lvCount;
217  std::map<G4LogicalVolume*, unsigned int>::const_iterator cite;
218  for (int ii = 0; ii < lv->GetNoDaughters(); ii++) {
219  cite = lvCount.find(lv->GetDaughter(ii)->GetLogicalVolume());
220  if (cite != lvCount.end()) lvCount[cite->first] = (cite->second) + 1;
221  else lvCount.insert(std::pair< G4LogicalVolume*,unsigned int>(lv->GetDaughter(ii)->GetLogicalVolume(),1));
222  }
223  for (cite = lvCount.begin(); cite != lvCount.end(); cite++)
224  dumpG4LVLeaf((cite->first), leafDepth+1, (cite->second), out);
225 }
226 
228 {
229  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's Material Budget Tree " << std::endl;
230  G4LogicalVolume * lv = getTopLV();
231  dumpG4LVLeafWithMat(lv,0,1,out);
232 }
233 
234 void PrintGeomMatInfo::dumpG4LVLeafWithMat(G4LogicalVolume * lv, unsigned int leafDepth, unsigned int count, std::ostream & out)
235 {
236  // switch off dumping at the next same level as the dump
237  if(_dumpIt && _level2Dump == leafDepth) {
238  _dumpIt = false;
239 // std::cout << " stopped dumping " << std::endl;
240  }
241  // loop over volumes to dump: volumes cannot not overlap in hierarchy!!
242  for (unsigned int i=0; i<_lvNames2Dump.size(); i++) {
243  if(_lvNames2Dump[i].compare(lv->GetName()) == 0) {
244  _dumpIt = true;
245  _dumpIndex = i;
246  _level2Dump = leafDepth;
247 // std::cout << " start dumping " << _lvNames2Dump[i] << " at level " << _level2Dump << std::endl;
248  break;
249  }
250  }
251 
252  if(_dumpIt) {
253  if(leafDepth < _maxLevelsCounted) _countsPerLevel[leafDepth] = count;
254  unsigned int total_multipler = 1;
255  for (unsigned int ii=_level2Dump; ii <= leafDepth; ii++) total_multipler *= _countsPerLevel[ii];
256  double thick = (lv->GetSolid()->GetCubicVolume() * total_multipler)/_areaLayer[_dumpIndex];
257  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
258  // print out Level, Logical volume name, volume in mm**3, material name, rad len of mat in mm
259  // total number of volumes from dump level start, equivalent thickness when spread over
260  // a cylinder of _radiusLayer and length _zLayer; equivalent thick in rad len
261  out << " LV::" << leafDepth << ": " << lv->GetName() << " :" << count
262  << ": " << lv->GetSolid()->GetName() << " :" << lv->GetSolid()->GetCubicVolume()
263  << ": "<< lv->GetMaterial()->GetName() << " :" << lv->GetMaterial()->GetRadlen()
264  << ":" << " :" << total_multipler << ":"
265  << " thk :" << thick << ": x/X0 :" << thick/lv->GetMaterial()->GetRadlen()
266  << ":" << " Kg : " << lv->GetMass()/kg << std::endl;
267  } else {
268  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
269  out << " LV:(" << leafDepth << ") " << lv->GetName() << " (" << count << ")" << std::endl;
270  }
271 
272  //--- If a volume is placed n types as daughter of this LV, it should only be counted once
273  std::map<G4LogicalVolume*, unsigned int> lvCount;
274  std::map<G4LogicalVolume*, unsigned int>::const_iterator cite;
275  for (int ii = 0; ii < lv->GetNoDaughters(); ii++) {
276  cite = lvCount.find(lv->GetDaughter(ii)->GetLogicalVolume());
277  if (cite != lvCount.end()) lvCount[cite->first] = (cite->second) + 1;
278  else lvCount.insert(std::pair< G4LogicalVolume*,unsigned int>(lv->GetDaughter(ii)->GetLogicalVolume(),1));
279  }
280  for (cite = lvCount.begin(); cite != lvCount.end(); cite++)
281  dumpG4LVLeafWithMat((cite->first), leafDepth+1, (cite->second), out);
282 }
283 
284 
286 {
287  int nTouch = 0;
288  G4LogicalVolume * lv = getTopLV();
289  add1touchable(lv, nTouch);
290  return nTouch;
291 }
292 
293 
294 void PrintGeomMatInfo::add1touchable(G4LogicalVolume * lv, int & nTouch)
295 {
296  int siz = lv->GetNoDaughters();
297  for(int ii = 0; ii < siz; ii++)
298  add1touchable(lv->GetDaughter(ii)->GetLogicalVolume(), ++nTouch);
299 }
300 
302 {
303  //dumps in the following order:
304  // 1) a LV with details
305  // 2) list of PVs daughters of this LV with details
306  // 3) list of LVs daughters of this LV and for each go to 1)
307 
308  //----- Get top PV
309  G4LogicalVolume* topLV = getTopLV();
310 
311  //----- Dump this leaf (it will recursively dump all the tree)
312  dumpHierarchyLeafPVLV(topLV, 0, out);
313  dumpPV(theTopPV, 0, out);
314 
315  //----- Dump the touchables (it will recursively dump all the tree)
316  if (_dumpTouch) dumpTouch(theTopPV, 0, out);
317 }
318 
319 void PrintGeomMatInfo::dumpHierarchyLeafPVLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out)
320 {
321  //----- Dump this LV
322  dumpLV(lv, leafDepth, out);
323 
324  //----- Get LV daughters from list of PV daughters
325  mmlvpv lvpvDaughters;
326  std::set< G4LogicalVolume * > lvDaughters;
327  int NoDaughters = lv->GetNoDaughters();
328  while ((NoDaughters--)>0)
329  {
330  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
331  lvpvDaughters.insert(mmlvpv::value_type(pvD->GetLogicalVolume(), pvD));
332  lvDaughters.insert(pvD->GetLogicalVolume());
333  }
334 
335  std::set< G4LogicalVolume * >::const_iterator scite;
336  mmlvpv::const_iterator mmcite;
337 
338  //----- Dump daughters PV and LV
339  for (scite = lvDaughters.begin(); scite != lvDaughters.end(); scite++)
340  {
341  std::pair< mmlvpv::iterator, mmlvpv::iterator > mmER = lvpvDaughters.equal_range(*scite);
342  //----- Dump daughters PV of this LV
343  for (mmcite = mmER.first ; mmcite != mmER.second; mmcite++)
344  dumpPV((*mmcite).second, leafDepth+1, out);
345  //----- Dump daughters LV
346  dumpHierarchyLeafPVLV(*scite, leafDepth+1, out );
347  }
348 }
349 
350 void PrintGeomMatInfo::dumpLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out)
351 {
352  std::string spaces = spacesFromLeafDepth(leafDepth);
353 
354  //----- dump name
355  if (_dumpLV) {
356  out << leafDepth << spaces << "$$$ VOLUME = " << lv->GetName()
357  << " Solid: " << lv->GetSolid()->GetName() << " MATERIAL: "
358  << lv->GetMaterial()->GetName() << std::endl;
359  if (_dumpSolid)
360  dumpSolid(lv->GetSolid(), leafDepth, out); //----- dump solid
361 
362  //----- dump LV info
363  //--- material
364  if (_dumpAtts) {
365  //--- Visualisation attributes
366  const G4VisAttributes * fVA = lv->GetVisAttributes();
367  if (fVA!=0) {
368  out << spaces << " VISUALISATION ATTRIBUTES: " << std::endl;
369  out << spaces << " IsVisible " << fVA->IsVisible() << std::endl;
370  out << spaces << " IsDaughtersInvisible " << fVA->IsDaughtersInvisible() << std::endl;
371  out << spaces << " Colour " << fVA->GetColour() << std::endl;
372  out << spaces << " LineStyle " << fVA->GetLineStyle() << std::endl;
373  out << spaces << " LineWidth " << fVA->GetLineWidth() << std::endl;
374  out << spaces << " IsForceDrawingStyle " << fVA->IsForceDrawingStyle() << std::endl;
375  out << spaces << " ForcedDrawingStyle " << fVA->GetForcedDrawingStyle() << std::endl;
376  }
377 
378  //--- User Limits
379  G4UserLimits * fUL = lv->GetUserLimits();
380  G4Track dummy;
381  if (fUL!=0) {
382  out << spaces << " MaxAllowedStep " << fUL->GetMaxAllowedStep(dummy) << std::endl;
383  out << spaces << " UserMaxTrackLength " << fUL->GetUserMaxTrackLength(dummy) << std::endl;
384  out << spaces << " UserMaxTime " << fUL->GetUserMaxTime(dummy) << std::endl;
385  out << spaces << " UserMinEkine " << fUL->GetUserMinEkine(dummy) << std::endl;
386  out << spaces << " UserMinRange " << fUL->GetUserMinRange(dummy) << std::endl;
387  }
388 
389  //--- other LV info
390  if (lv->GetSensitiveDetector())
391  out << spaces << " IS SENSITIVE DETECTOR " << std::endl;
392  if (lv->GetFieldManager())
393  out << spaces << " FIELD ON " << std::endl;
394 
395  // Pointer (possibly NULL) to optimisation info objects.
396  out << spaces
397  << " Quality for optimisation, average number of voxels to be spent per content "
398  << lv->GetSmartless() << std::endl;
399 
400  // Pointer (possibly NULL) to G4FastSimulationManager object.
401  if (lv->GetFastSimulationManager())
402  out << spaces << " Logical Volume is an envelope for a FastSimulationManager "
403  << std::endl;
404  out << spaces << " Weight used in the event biasing technique = "
405  << lv->GetBiasWeight() << std::endl;
406  }
407  }
408 }
409 
410 void PrintGeomMatInfo::dumpPV(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out)
411 {
412  std::string spaces = spacesFromLeafDepth(leafDepth);
413 
414  //----- PV info
415  if (_dumpPV)
416  {
417  std::string mother = "World";
418  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
419  out << leafDepth << spaces << "### VOLUME = " << pv->GetName()
420  << " Copy No " << pv->GetCopyNo() << " in " << mother
421  << " at " << pv->GetTranslation();
422  }
423  if (!pv->IsReplicated())
424  {
425  if (_dumpPV)
426  {
427  if(pv->GetRotation() == 0) out << " with no rotation" << std::endl;
428  else if(!_dumpRotation) out << " with rotation" << std::endl; //just rotation name
429  else out << " with rotation " << *(pv->GetRotation()) << std::endl;
430  }
431  }
432  else
433  {
434  if (_dumpReplica )
435  {
436  out << spaces << " It is replica: " << std::endl;
437  EAxis axis;
438  int nReplicas;
439  double width;
440  double offset;
441  bool consuming;
442  pv->GetReplicationData(axis, nReplicas, width, offset, consuming);
443  out << spaces << " axis " << axis << std::endl
444  << spaces << " nReplicas " << nReplicas << std::endl;
445  if (pv->GetParameterisation() != 0)
446  out << spaces << " It is parameterisation " << std::endl;
447  else
448  out << spaces << " width " << width << std::endl
449  << spaces << " offset " << offset << std::endl
450  << spaces << " consuming" << consuming << std::endl;
451  if (pv->GetParameterisation() != 0)
452  out << spaces << " It is parameterisation " << std::endl;
453  }
454  }
455 }
456 
457 void PrintGeomMatInfo::dumpTouch(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out)
458 {
459  std::string spaces = spacesFromLeafDepth(leafDepth);
460  if (leafDepth == 0) fHistory.SetFirstEntry(pv);
461  else fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
462 
463  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().
464  TransformPoint(G4ThreeVector(0,0,0));
465  G4LogicalVolume * lv = pv->GetLogicalVolume();
466 
467  std::string mother = "World";
468  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
469  std::string lvname = lv->GetName();
470  lvname.assign(lvname,0,nchar);
471  if (lvname == name)
472  out << leafDepth << spaces << "### VOLUME = " << lv->GetName()
473  << " Copy No " << pv->GetCopyNo() << " in " << mother
474  << " global position of centre " << globalpoint << " (r = "
475  << globalpoint.perp() << ", phi = " << globalpoint.phi()/deg
476  << ")" << std::endl;
477 
478  int NoDaughters = lv->GetNoDaughters();
479  while ((NoDaughters--)>0)
480  {
481  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
482  if (!pvD->IsReplicated()) dumpTouch(pvD, leafDepth+1, out);
483  }
484 
485  if (leafDepth > 0) fHistory.BackLevel();
486 }
487 
489 {
490  std::string spaces;
491  unsigned int ii;
492  for(ii = 0; ii < leafDepth; ii++) { spaces += " "; }
493  return spaces;
494 }
495 
496 void PrintGeomMatInfo::dumpSolid(G4VSolid * sol, unsigned int leafDepth, std::ostream & out)
497 {
498  std::string spaces = spacesFromLeafDepth(leafDepth);
499  out << spaces << *(sol) << std::endl;
500 }
501 
502 G4VPhysicalVolume * PrintGeomMatInfo::getTopPV()
503 {
504  return G4TransportationManager::GetTransportationManager()
505  ->GetNavigatorForTracking()->GetWorldVolume();
506 }
507 
508 G4LogicalVolume * PrintGeomMatInfo::getTopLV()
509 { return theTopPV->GetLogicalVolume(); }
510 
511 
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
G4VPhysicalVolume * getTopPV()
void dumpG4LVTree(std::ostream &out=std::cout)
void addFilter(const DDFilter &, log_op op=AND)
const N & name() const
Definition: DDBase.h:82
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
const reco::GenParticle * mother(const reco::GenParticle &p, unsigned int imoth=0)
void dumpMaterialList(std::ostream &out=std::cout)
int ii
Definition: cuy.py:588
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)
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:88
void dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, std::ostream &out=std::cout)
Container::value_type value_type
void dumpSummary(std::ostream &out=std::cout)
tuple out
Definition: dbtoconf.py:99
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)
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
tuple cout
Definition: gather_cfg.py:121
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.
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37