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 
14  name_ = p.getUntrackedParameter<std::string>("Name", "*");
15  nchar_ = name_.find('*');
16  name_.assign(name_, 0, nchar_);
17  std::cout << "PrintSensitive:: Print position of all Sensitive Touchables: "
18  << " for names (0-" << nchar_ << ") = " << name_ << "\n";
19 }
20 
22 
24  G4VPhysicalVolume *theTopPV = getTopPV();
25  int nsens = dumpTouch(theTopPV, 0, false, 0, std::cout);
26  std::cout << "\nTotal number of sensitive detector volumes for " << name_ << " is " << nsens << std::endl;
27 }
28 
29 int PrintSensitive::dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out) {
30  if (leafDepth == 0)
31  fHistory.SetFirstEntry(pv);
32  else
33  fHistory.NewLevel(pv, kNormal, pv->GetCopyNo());
34 
35  int nsens(ns);
36  G4ThreeVector globalpoint = fHistory.GetTopTransform().Inverse().TransformPoint(G4ThreeVector(0, 0, 0));
37  G4LogicalVolume *lv = pv->GetLogicalVolume();
38 
39  std::string mother = "World";
40  if (pv->GetMotherLogical())
41  mother = pv->GetMotherLogical()->GetName();
42  std::string lvname = lv->GetName();
43  lvname.assign(lvname, 0, nchar_);
44  if (lvname == name_)
45  printIt = true;
46 
47  if (lv->GetSensitiveDetector() && printIt) {
48  ++nsens;
49  out << nsens << ":" << leafDepth << " ### VOLUME = " << lv->GetName() << " Copy No " << pv->GetCopyNo() << " in "
50  << mother << " global position of centre " << globalpoint << " (r=" << globalpoint.perp()
51  << ", phi=" << globalpoint.phi() / CLHEP::deg << ")\n";
52  }
53 
54  int NoDaughters = lv->GetNoDaughters();
55  while ((NoDaughters--) > 0) {
56  G4VPhysicalVolume *pvD = lv->GetDaughter(NoDaughters);
57  if (!pvD->IsReplicated())
58  nsens = dumpTouch(pvD, leafDepth + 1, printIt, nsens, out);
59  }
60 
61  if (leafDepth > 0)
62  fHistory.BackLevel();
63  return nsens;
64 }
65 
66 G4VPhysicalVolume *PrintSensitive::getTopPV() {
67  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
68 }
PrintSensitive::name_
std::string name_
Definition: PrintSensitive.h:27
PrintSensitive::dumpTouch
int dumpTouch(G4VPhysicalVolume *pv, unsigned int leafDepth, bool printIt, int ns, std::ostream &out=std::cout)
Definition: PrintSensitive.cc:29
PrintSensitive.h
PrintSensitive::nchar_
int nchar_
Definition: PrintSensitive.h:28
PrintSensitive::getTopPV
G4VPhysicalVolume * getTopPV()
Definition: PrintSensitive.cc:66
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
gather_cfg.cout
cout
Definition: gather_cfg.py:144
PrintSensitive::PrintSensitive
PrintSensitive(edm::ParameterSet const &p)
Definition: PrintSensitive.cc:13
BeginOfRun.h
PrintSensitive::fHistory
G4NavigationHistory fHistory
Definition: PrintSensitive.h:29
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PrintSensitive::update
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
Definition: PrintSensitive.cc:23
edm::ParameterSet
Definition: ParameterSet.h:47
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
BeginOfRun
Definition: BeginOfRun.h:6
writedatasetfile.run
run
Definition: writedatasetfile.py:27
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
PrintSensitive::~PrintSensitive
~PrintSensitive() override
Definition: PrintSensitive.cc:21