CMS 3D CMS Logo

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