CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PrintGeomInfoAction.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 
33 
34  _dumpSummary = p.getUntrackedParameter<bool>("DumpSummary", true);
35  _dumpLVTree = p.getUntrackedParameter<bool>("DumpLVTree", true);
36  _dumpMaterial= p.getUntrackedParameter<bool>("DumpMaterial",false);
37  _dumpLVList = p.getUntrackedParameter<bool>("DumpLVList", false);
38  _dumpLV = p.getUntrackedParameter<bool>("DumpLV", false);
39  _dumpSolid = p.getUntrackedParameter<bool>("DumpSolid", false);
40  _dumpAtts = p.getUntrackedParameter<bool>("DumpAttributes", false);
41  _dumpPV = p.getUntrackedParameter<bool>("DumpPV", false);
42  _dumpRotation= p.getUntrackedParameter<bool>("DumpRotation",false);
43  _dumpReplica = p.getUntrackedParameter<bool>("DumpReplica", false);
44  _dumpTouch = p.getUntrackedParameter<bool>("DumpTouch", false);
45  _dumpSense = p.getUntrackedParameter<bool>("DumpSense", false);
46  name = p.getUntrackedParameter<std::string>("Name","*");
47  nchar = name.find("*");
48  name.assign(name,0,nchar);
49  names = p.getUntrackedParameter<std::vector<std::string> >("Names");
50  std::cout << "PrintGeomInfoAction:: initialised with verbosity levels:"
51  << " Summary " << _dumpSummary << " LVTree " << _dumpLVTree
52  << " LVList " << _dumpLVList << " Material " << _dumpMaterial
53  << "\n "
54  << " LV " << _dumpLV << " Solid " << _dumpSolid
55  << " Attribs " << _dumpAtts
56  << "\n "
57  << " PV " << _dumpPV << " Rotation " << _dumpRotation
58  << " Replica " << _dumpReplica
59  << "\n "
60  << " Touchable " << _dumpTouch << " for names (0-" << nchar
61  << ") = " << name
62  << "\n "
63  << " Sensitive " << _dumpSense << " for " << names.size()
64  << " namess";
65  for (unsigned int i=0; i<names.size(); i++) std::cout << " " << names[i];
66  std::cout << std::endl;
67 }
68 
70 
72 
73  if (_dumpSense) {
75  (*job)()->get<IdealGeometryRecord>().get(pDD);
76 
77  std::cout << "PrintGeomInfoAction::Get Printout of Sensitive Volumes "
78  << "for " << names.size() << " Readout Units" << std::endl;
79  for (unsigned int i=0; i<names.size(); i++) {
80  std::string attribute = "ReadOutName";
81  std::string sd = names[i];
83  DDValue ddv(attribute,sd,0);
85  DDFilteredView fv(*pDD);
86  std::cout << "PrintGeomInfoAction:: Get Filtered view for "
87  << attribute << " = " << sd << std::endl;
88  fv.addFilter(filter);
89  bool dodet = fv.firstChild();
90 
91  std::string spaces = spacesFromLeafDepth(1);
92 
93  while (dodet) {
94  const DDLogicalPart & log = fv.logicalPart();
95  std::string lvname = log.name().name();
96  DDTranslation tran = fv.translation();
97  std::vector<int> copy = fv.copyNumbers();
98 
99  unsigned int leafDepth = copy.size();
100  std::cout << leafDepth << spaces << "### VOLUME = " << lvname
101  << " Copy No";
102  for (int k=leafDepth-1; k>=0; k--) std::cout << " " << copy[k];
103  std::cout << " Centre at " << tran << " (r = " << tran.Rho()
104  << ", phi = " << tran.phi()/deg << ")" << std::endl;
105  dodet = fv.next();
106  }
107  }
108  }
109 }
110 
112 
113  theTopPV = getTopPV();
114 
117 
118  //---------- Dump list of objects of each class with detail of parameters
121 
122  //---------- Dump LV and PV information
124 }
125 
126 void PrintGeomInfoAction::dumpSummary(std::ostream & out) {
127 
128  //---------- Dump number of objects of each class
129  out << " @@@@@@@@@@@@@@@@@@ Dumping G4 geometry objects Summary " << std::endl;
130  if (theTopPV == 0) {
131  out << " No volume created " << std::endl;
132  return;
133  }
134  out << " @@@ Geometry built inside world volume: " << theTopPV->GetName() << std::endl;
135  // Get number of solids (< # LV if several LV share a solid)
136  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
137  std::vector<G4LogicalVolume *>::const_iterator lvcite;
138  std::set<G4VSolid *> theSolids;
139  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
140  theSolids.insert((*lvcite)->GetSolid());
141  out << " Number of G4VSolid's: " << theSolids.size() << std::endl;
142  out << " Number of G4LogicalVolume's: " << lvs->size() << std::endl;
143  const G4PhysicalVolumeStore * pvs = G4PhysicalVolumeStore::GetInstance();
144  out << " Number of G4VPhysicalVolume's: " << pvs->size() << std::endl;
145  out << " Number of Touchable's: " << countNoTouchables() << std::endl;
146  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
147  out << " Number of G4Material's: " << matTab->size() << std::endl;
148 }
149 
151 
152  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's List " << std::endl;
153  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
154  std::vector<G4LogicalVolume*>::const_iterator lvcite;
155  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
156  out << "LV:" << (*lvcite)->GetName() << "\tMaterial: " << (*lvcite)->GetMaterial()->GetName() << std::endl;
157 }
158 
160 
161  out << " @@@@@@@@@@@@@@@@ DUMPING G4LogicalVolume's Tree " << std::endl;
162  G4LogicalVolume * lv = getTopLV();
163  dumpG4LVLeaf(lv,0,1,out);
164 }
165 
167 
168  out << " @@@@@@@@@@@@@@@@ DUMPING G4Material List ";
169  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
170  out << " with " << matTab->size() << " materials " << std::endl;
171  std::vector<G4Material*>::const_iterator matite;
172  for (matite = matTab->begin(); matite != matTab->end(); matite++)
173  out << "Material: " << (*matite) << std::endl;
174 }
175 
176 void PrintGeomInfoAction::dumpG4LVLeaf(G4LogicalVolume * lv, unsigned int leafDepth, unsigned int count, std::ostream & out) {
177 
178  for (unsigned int ii=0; ii < leafDepth; ii++) out << " ";
179  out << " LV:(" << leafDepth << ") " << lv->GetName() << " (" << count
180  << ")" << std::endl;
181  //--- If a volume is placed n types as daughter of this LV, it should only be counted once
182  std::map<G4LogicalVolume*, unsigned int> lvCount;
183  std::map<G4LogicalVolume*, unsigned int>::const_iterator cite;
184  for (int ii = 0; ii < lv->GetNoDaughters(); ii++) {
185  cite = lvCount.find(lv->GetDaughter(ii)->GetLogicalVolume());
186  if (cite != lvCount.end()) lvCount[cite->first] = (cite->second) + 1;
187  else lvCount.insert(std::pair< G4LogicalVolume*,unsigned int>(lv->GetDaughter(ii)->GetLogicalVolume(),1));
188  }
189  for (cite = lvCount.begin(); cite != lvCount.end(); cite++)
190  dumpG4LVLeaf((cite->first), leafDepth+1, (cite->second), out);
191 }
192 
194 
195  int nTouch = 0;
196  G4LogicalVolume * lv = getTopLV();
197  add1touchable(lv, nTouch);
198  return nTouch;
199 }
200 
201 void PrintGeomInfoAction::add1touchable(G4LogicalVolume * lv, int & nTouch) {
202 
203  int siz = lv->GetNoDaughters();
204  for(int ii = 0; ii < siz; ii++)
205  add1touchable(lv->GetDaughter(ii)->GetLogicalVolume(), ++nTouch);
206 }
207 
209 
210  //dumps in the following order:
211  // 1) a LV with details
212  // 2) list of PVs daughters of this LV with details
213  // 3) list of LVs daughters of this LV and for each go to 1)
214 
215  //----- Get top PV
216  G4LogicalVolume* topLV = getTopLV();
217 
218  //----- Dump this leaf (it will recursively dump all the tree)
219  dumpHierarchyLeafPVLV(topLV, 0, out);
220  dumpPV(theTopPV, 0, out);
221 
222  //----- Dump the touchables (it will recursively dump all the tree)
223  if (_dumpTouch) dumpTouch(theTopPV, 0, out);
224 }
225 
226 void PrintGeomInfoAction::dumpHierarchyLeafPVLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out) {
227 
228  //----- Dump this LV
229  dumpLV(lv, leafDepth, out);
230 
231  //----- Get LV daughters from list of PV daughters
232  mmlvpv lvpvDaughters;
233  std::set< G4LogicalVolume * > lvDaughters;
234  int NoDaughters = lv->GetNoDaughters();
235  while ((NoDaughters--)>0) {
236  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
237  lvpvDaughters.insert(mmlvpv::value_type(pvD->GetLogicalVolume(), pvD));
238  lvDaughters.insert(pvD->GetLogicalVolume());
239  }
240 
241  std::set< G4LogicalVolume * >::const_iterator scite;
242  mmlvpv::const_iterator mmcite;
243 
244  //----- Dump daughters PV and LV
245  for (scite = lvDaughters.begin(); scite != lvDaughters.end(); scite++) {
246  std::pair< mmlvpv::iterator, mmlvpv::iterator > mmER = lvpvDaughters.equal_range(*scite);
247  //----- Dump daughters PV of this LV
248  for (mmcite = mmER.first ; mmcite != mmER.second; mmcite++)
249  dumpPV((*mmcite).second, leafDepth+1, out);
250  //----- Dump daughters LV
251  dumpHierarchyLeafPVLV(*scite, leafDepth+1, out );
252  }
253 }
254 
255 void PrintGeomInfoAction::dumpLV(G4LogicalVolume * lv, unsigned int leafDepth, std::ostream & out) {
256 
257  std::string spaces = spacesFromLeafDepth(leafDepth);
258 
259  //----- dump name
260  if (_dumpLV) {
261  out << leafDepth << spaces << "$$$ VOLUME = " << lv->GetName()
262  << " Solid: " << lv->GetSolid()->GetName() << " MATERIAL: "
263  << lv->GetMaterial()->GetName() << std::endl;
264  if (_dumpSolid)
265  dumpSolid(lv->GetSolid(), leafDepth, out); //----- dump solid
266 
267  //----- dump LV info
268  //--- material
269  if (_dumpAtts) {
270  //--- Visualisation attributes
271  const G4VisAttributes * fVA = lv->GetVisAttributes();
272  if (fVA!=0) {
273  out << spaces << " VISUALISATION ATTRIBUTES: " << std::endl;
274  out << spaces << " IsVisible " << fVA->IsVisible() << std::endl;
275  out << spaces << " IsDaughtersInvisible " << fVA->IsDaughtersInvisible() << std::endl;
276  out << spaces << " Colour " << fVA->GetColour() << std::endl;
277  out << spaces << " LineStyle " << fVA->GetLineStyle() << std::endl;
278  out << spaces << " LineWidth " << fVA->GetLineWidth() << std::endl;
279  out << spaces << " IsForceDrawingStyle " << fVA->IsForceDrawingStyle() << std::endl;
280  out << spaces << " ForcedDrawingStyle " << fVA->GetForcedDrawingStyle() << std::endl;
281  }
282 
283  //--- User Limits
284  G4UserLimits * fUL = lv->GetUserLimits();
285  G4Track dummy;
286  if (fUL!=0) {
287  out << spaces << " MaxAllowedStep " << fUL->GetMaxAllowedStep(dummy) << std::endl;
288  out << spaces << " UserMaxTrackLength " << fUL->GetUserMaxTrackLength(dummy) << std::endl;
289  out << spaces << " UserMaxTime " << fUL->GetUserMaxTime(dummy) << std::endl;
290  out << spaces << " UserMinEkine " << fUL->GetUserMinEkine(dummy) << std::endl;
291  out << spaces << " UserMinRange " << fUL->GetUserMinRange(dummy) << std::endl;
292  }
293 
294  //--- other LV info
295  if (lv->GetSensitiveDetector())
296  out << spaces << " IS SENSITIVE DETECTOR " << std::endl;
297  if (lv->GetFieldManager())
298  out << spaces << " FIELD ON " << std::endl;
299 
300  // Pointer (possibly NULL) to optimisation info objects.
301  out << spaces
302  << " Quality for optimisation, average number of voxels to be spent per content "
303  << lv->GetSmartless() << std::endl;
304 
305  // Pointer (possibly NULL) to G4FastSimulationManager object.
306  if (lv->GetFastSimulationManager())
307  out << spaces << " Logical Volume is an envelope for a FastSimulationManager "
308  << std::endl;
309  out << spaces << " Weight used in the event biasing technique = "
310  << lv->GetBiasWeight() << std::endl;
311  }
312  }
313 }
314 
315 void PrintGeomInfoAction::dumpPV(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out) {
316 
317  std::string spaces = spacesFromLeafDepth(leafDepth);
318 
319  //----- PV info
320  if (_dumpPV) {
321  std::string mother = "World";
322  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
323  out << leafDepth << spaces << "### VOLUME = " << pv->GetName()
324  << " Copy No " << pv->GetCopyNo() << " in " << mother
325  << " at " << pv->GetTranslation();
326  }
327  if (!pv->IsReplicated()) {
328  if (_dumpPV) {
329  if(pv->GetRotation() == 0) out << " with no rotation" << std::endl;
330  else if(!_dumpRotation) out << " with rotation" << std::endl; //just rotation name
331  else out << " with rotation " << *(pv->GetRotation()) << std::endl;
332  }
333  } else {
334  if (_dumpReplica ) {
335  out << spaces << " It is replica: " << std::endl;
336  EAxis axis;
337  int nReplicas;
338  double width;
339  double offset;
340  bool consuming;
341  pv->GetReplicationData(axis, nReplicas, width, offset, consuming);
342  out << spaces << " axis " << axis << std::endl
343  << spaces << " nReplicas " << nReplicas << std::endl;
344  if (pv->GetParameterisation() != 0)
345  out << spaces << " It is parameterisation " << std::endl;
346  else
347  out << spaces << " width " << width << std::endl
348  << spaces << " offset " << offset << std::endl
349  << spaces << " consuming" << consuming << std::endl;
350  if (pv->GetParameterisation() != 0)
351  out << spaces << " It is parameterisation " << std::endl;
352  }
353  }
354 }
355 
356 void PrintGeomInfoAction::dumpTouch(G4VPhysicalVolume * pv, unsigned int leafDepth, std::ostream & out) {
357 
358  std::string spaces = spacesFromLeafDepth(leafDepth);
359  if (leafDepth == 0) fHistory.SetFirstEntry(pv);
360  else fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
361 
362  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0,0,0));
363  G4LogicalVolume * lv = pv->GetLogicalVolume();
364 
365  std::string mother = "World";
366  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
367  std::string lvname = lv->GetName();
368  lvname.assign(lvname,0,nchar);
369  if (lvname == name)
370  out << leafDepth << spaces << "### VOLUME = " << lv->GetName()
371  << " Copy No " << pv->GetCopyNo() << " in " << mother
372  << " global position of centre " << globalpoint << " (r = "
373  << globalpoint.perp() << ", phi = " << globalpoint.phi()/deg
374  << ")" << std::endl;
375 
376  int NoDaughters = lv->GetNoDaughters();
377  while ((NoDaughters--)>0) {
378  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
379  if (!pvD->IsReplicated()) dumpTouch(pvD, leafDepth+1, out);
380  }
381 
382  if (leafDepth > 0) fHistory.BackLevel();
383 }
384 
386 
387  std::string spaces;
388  unsigned int ii;
389  for(ii = 0; ii < leafDepth; ii++) { spaces += " "; }
390  return spaces;
391 }
392 
393 void PrintGeomInfoAction::dumpSolid(G4VSolid * sol, unsigned int leafDepth, std::ostream & out) {
394 
395  std::string spaces = spacesFromLeafDepth(leafDepth);
396  out << spaces << *(sol) << std::endl;
397 }
398 
399 G4VPhysicalVolume * PrintGeomInfoAction::getTopPV() {
400 
401  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
402 }
403 
404 G4LogicalVolume * PrintGeomInfoAction::getTopLV() {
405  return theTopPV->GetLogicalVolume();
406 }
407 
408 
void add1touchable(G4LogicalVolume *lv, int &nTouch)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void dumpMaterialList(std::ostream &out=std::cout)
void dumpSummary(std::ostream &out=std::cout)
PrintGeomInfoAction(edm::ParameterSet const &p)
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
void addFilter(const DDFilter &, log_op op=AND)
const N & name() const
Definition: DDBase.h:82
void dumpPV(G4VPhysicalVolume *pv, unsigned int leafDepth, std::ostream &out=std::cout)
std::multimap< G4LogicalVolume *, G4VPhysicalVolume *, std::less< G4LogicalVolume * > > mmlvpv
nav_type copyNumbers() const
return the stack of copy numbers
G4VPhysicalVolume * theTopPV
void dumpHierarchyTreePVLV(std::ostream &out=std::cout)
int ii
Definition: cuy.py:588
std::vector< std::string > names
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DDTranslation
Definition: DDTranslation.h:7
void update(const BeginOfJob *job)
This routine will be called when the appropriate signal arrives.
void dumpG4LVTree(std::ostream &out=std::cout)
bool next()
set current node to the next node in the filtered tree
void dumpSolid(G4VSolid *sol, unsigned int leafDepth, std::ostream &out=std::cout)
void dumpLV(G4LogicalVolume *lv, unsigned int leafDepth, std::ostream &out=std::cout)
G4NavigationHistory fHistory
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:88
std::string spacesFromLeafDepth(unsigned int leafDepth)
unsigned int offset(bool)
Container::value_type value_type
int k[5][pyjets_maxn]
tuple out
Definition: dbtoconf.py:99
G4VPhysicalVolume * getTopPV()
void dumpHierarchyLeafPVLV(G4LogicalVolume *lv, unsigned int leafDepth, std::ostream &out=std::cout)
double sd
G4LogicalVolume * getTopLV()
void dumpG4LVList(std::ostream &out=std::cout)
void dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, std::ostream &out=std::cout)
bool firstChild()
set the current node to the first child ...
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
const DDTranslation & translation() const
The absolute translation of the current node.
const std::string & name() const
Returns the name.
Definition: DDName.cc:87
void dumpG4LVLeaf(G4LogicalVolume *lv, unsigned int leafDepth, unsigned int count, std::ostream &out=std::cout)
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37