CMS 3D CMS Logo

PrintSensitive.cc
Go to the documentation of this file.
5 
6 #include "G4Run.hh"
7 #include "G4VPhysicalVolume.hh"
8 #include "G4LogicalVolume.hh"
9 #include "G4NavigationHistory.hh"
10 #include "G4TransportationManager.hh"
11 
12 #include <set>
13 #include <map>
14 #include <iostream>
15 #include <string>
16 
17 class PrintSensitive : public SimWatcher, public Observer<const BeginOfRun *> {
18 public:
20  ~PrintSensitive() override;
21 
22 private:
23  void update(const BeginOfRun *run) override;
24  int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out = G4cout);
25  G4VPhysicalVolume *getTopPV();
26 
27 private:
29  int nchar_;
30  G4NavigationHistory fHistory;
31 };
32 
34  name_ = p.getUntrackedParameter<std::string>("Name", "*");
35  nchar_ = name_.find('*');
36  name_.assign(name_, 0, nchar_);
37  G4cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
38  << " for names (0-" << nchar_ << ") = " << name_ << G4endl;
39 }
40 
42 
44  G4VPhysicalVolume *theTopPV = getTopPV();
45  int nsens = dumpTouch(theTopPV, 0, false, 0, G4cout);
46  G4cout << "\nTotal number of sensitive detector volumes for " << name_ << " is " << nsens << G4endl;
47 }
48 
49 int PrintSensitive::dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out) {
50  if (leafDepth == 0)
51  fHistory.SetFirstEntry(pv);
52  else
53  fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
54 
55  int nsens(ns);
56  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0, 0, 0));
57  G4LogicalVolume *lv = pv->GetLogicalVolume();
58 
59  std::string mother = "World";
60  if (pv->GetMotherLogical())
61  mother = pv->GetMotherLogical()->GetName();
62  std::string lvname = lv->GetName();
63  lvname.assign(lvname, 0, nchar_);
64  if (lvname == name_)
65  printIt = true;
66 
67  if (lv->GetSensitiveDetector() && printIt) {
68  ++nsens;
69  out << nsens << ":" << leafDepth << " ### VOLUME = " << lv->GetName() << " Copy No " << pv->GetCopyNo() << " in "
70  << mother << " global position of centre " << globalpoint << " (r=" << globalpoint.perp()
71  << ", phi=" << globalpoint.phi() / CLHEP::deg << ")\n";
72  }
73 
74  int NoDaughters = lv->GetNoDaughters();
75  while ((NoDaughters--) > 0) {
76  G4VPhysicalVolume *pvD = lv->GetDaughter(NoDaughters);
77  if (!pvD->IsReplicated())
78  nsens = dumpTouch(pvD, leafDepth + 1, printIt, nsens, out);
79  }
80 
81  if (leafDepth > 0)
82  fHistory.BackLevel();
83  return nsens;
84 }
85 
86 G4VPhysicalVolume *PrintSensitive::getTopPV() {
87  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
88 }
89 
92 
#define DEFINE_SIMWATCHER(type)
G4NavigationHistory fHistory
PrintSensitive(edm::ParameterSet const &p)
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
G4VPhysicalVolume * getTopPV()
std::string name_
int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out=G4cout)
~PrintSensitive() override
__host__ __device__ void printIt(C *m, const char *prefix="")
Definition: FitUtils.h:53