File indexing completed on 2024-04-06 12:30:30
0001 #include "SimG4Core/Notification/interface/BeginOfRun.h"
0002 #include "SimG4Core/Notification/interface/Observer.h"
0003 #include "SimG4Core/Watcher/interface/SimWatcher.h"
0004 #include "SimG4Core/Geometry/interface/DD4hep2DDDName.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006
0007 #include "G4Run.hh"
0008 #include "G4VPhysicalVolume.hh"
0009 #include "G4LogicalVolume.hh"
0010 #include "G4NavigationHistory.hh"
0011 #include "G4TransportationManager.hh"
0012
0013 #include <iostream>
0014 #include <string>
0015
0016 class PrintSensitive : public SimWatcher, public Observer<const BeginOfRun *> {
0017 public:
0018 PrintSensitive(edm::ParameterSet const &p);
0019 ~PrintSensitive() override;
0020
0021 private:
0022 void update(const BeginOfRun *run) override;
0023 int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out = G4cout);
0024 G4VPhysicalVolume *getTopPV();
0025
0026 private:
0027 std::string name_;
0028 bool dd4hep_;
0029 int nchar_;
0030 G4NavigationHistory fHistory;
0031 };
0032
0033 PrintSensitive::PrintSensitive(const edm::ParameterSet &p) {
0034 name_ = p.getParameter<std::string>("Name");
0035 dd4hep_ = p.getParameter<bool>("DD4hep");
0036 nchar_ = name_.find('*');
0037 name_.assign(name_, 0, nchar_);
0038 G4cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
0039 << " for names (0-" << nchar_ << ") = " << name_ << " dd4hep " << dd4hep_ << G4endl;
0040 }
0041
0042 PrintSensitive::~PrintSensitive() {}
0043
0044 void PrintSensitive::update(const BeginOfRun *run) {
0045 G4VPhysicalVolume *theTopPV = getTopPV();
0046 int nsens = dumpTouch(theTopPV, 0, false, 0, G4cout);
0047 G4cout << "\nTotal number of sensitive detector volumes for " << name_ << " is " << nsens << G4endl;
0048 }
0049
0050 int PrintSensitive::dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out) {
0051 if (leafDepth == 0)
0052 fHistory.SetFirstEntry(pv);
0053 else
0054 fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
0055
0056 int nsens(ns);
0057 G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0, 0, 0));
0058 G4LogicalVolume *lv = pv->GetLogicalVolume();
0059
0060 std::string mother = (pv->GetMotherLogical())
0061 ? (DD4hep2DDDName::nameSolid(
0062 static_cast<std::string>(pv->GetMotherLogical()->GetSolid()->GetName()), dd4hep_))
0063 : "World";
0064 std::string lvname = DD4hep2DDDName::nameSolid(static_cast<std::string>(lv->GetSolid()->GetName()), dd4hep_);
0065 if (nchar_ > 0) {
0066 lvname.assign(lvname, 0, nchar_);
0067 if (lvname == name_)
0068 printIt = true;
0069 } else {
0070 printIt = true;
0071 }
0072
0073 if (lv->GetSensitiveDetector() && printIt) {
0074 ++nsens;
0075 lvname = DD4hep2DDDName::nameSolid(static_cast<std::string>(lv->GetName()), dd4hep_);
0076 out << nsens << ":" << leafDepth << " ### VOLUME = " << lvname << " Copy No " << pv->GetCopyNo() << " in " << mother
0077 << " global position of centre " << globalpoint << " (r=" << globalpoint.perp()
0078 << ", phi=" << globalpoint.phi() / CLHEP::deg << ")\n";
0079 }
0080
0081 int NoDaughters = lv->GetNoDaughters();
0082 while ((NoDaughters--) > 0) {
0083 G4VPhysicalVolume *pvD = lv->GetDaughter(NoDaughters);
0084 if (!pvD->IsReplicated())
0085 nsens = dumpTouch(pvD, leafDepth + 1, printIt, nsens, out);
0086 }
0087
0088 if (leafDepth > 0)
0089 fHistory.BackLevel();
0090 return nsens;
0091 }
0092
0093 G4VPhysicalVolume *PrintSensitive::getTopPV() {
0094 return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
0095 }
0096
0097 #include "SimG4Core/Watcher/interface/SimWatcherFactory.h"
0098 #include "FWCore/PluginManager/interface/ModuleDef.h"
0099
0100 DEFINE_SIMWATCHER(PrintSensitive);