CMS 3D CMS Logo

PrintSensitive.cc
Go to the documentation of this file.
2 
4 
5 #include "G4Run.hh"
6 #include "G4VPhysicalVolume.hh"
7 #include "G4LogicalVolume.hh"
8 #include "G4TransportationManager.hh"
9 
10 #include <set>
11 #include <map>
12 
13 using namespace CLHEP;
14 
16  name = p.getUntrackedParameter<std::string>("Name","*");
17  nchar = name.find("*");
18  name.assign(name,0,nchar);
19  std::cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
20  << " for names (0-" << nchar << ") = " << name << "\n";
21 }
22 
24 
26 
27  G4VPhysicalVolume * theTopPV = getTopPV();
28  dumpTouch(theTopPV, 0, false, std::cout);
29 }
30 
31 void PrintSensitive::dumpTouch(G4VPhysicalVolume * pv, unsigned int leafDepth,
32  bool printIt, std::ostream & out) {
33 
34  if (leafDepth == 0) fHistory.SetFirstEntry(pv);
35  else fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
36 
37  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().
38  TransformPoint(G4ThreeVector(0,0,0));
39  G4LogicalVolume * lv = pv->GetLogicalVolume();
40 
41  std::string mother = "World";
42  if (pv->GetMotherLogical()) mother = pv->GetMotherLogical()->GetName();
43  std::string lvname = lv->GetName();
44  lvname.assign(lvname,0,nchar);
45  if (lvname == name) printIt = true;
46 
47  if (lv->GetSensitiveDetector() && printIt) {
48  out << leafDepth << " ### VOLUME = " << lv->GetName()
49  << " Copy No " << pv->GetCopyNo() << " in " << mother
50  << " global position of centre " << globalpoint << " (r="
51  << globalpoint.perp() << ", phi=" << globalpoint.phi()/deg
52  << ")\n";
53  }
54 
55  int NoDaughters = lv->GetNoDaughters();
56  while ((NoDaughters--)>0) {
57  G4VPhysicalVolume * pvD = lv->GetDaughter(NoDaughters);
58  if (!pvD->IsReplicated()) dumpTouch(pvD, leafDepth+1, printIt, out);
59  }
60 
61  if (leafDepth > 0) fHistory.BackLevel();
62 }
63 
64 G4VPhysicalVolume * PrintSensitive::getTopPV() {
65  return G4TransportationManager::GetTransportationManager()
66  ->GetNavigatorForTracking()->GetWorldVolume();
67 }
T getUntrackedParameter(std::string const &, T const &) const
PrintSensitive(edm::ParameterSet const &p)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
G4VPhysicalVolume * getTopPV()
def pv(vc)
Definition: MetAnalyzer.py:6
void dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, std::ostream &out=std::cout)
~PrintSensitive() override