CMS 3D CMS Logo

TrackerG4SimHitNumberingScheme.cc
Go to the documentation of this file.
3 
5 #include "DD4hep/Filter.h"
6 #include "G4TransportationManager.hh"
7 #include "G4Navigator.hh"
8 #include "G4VTouchable.hh"
9 #include "G4TouchableHistory.hh"
10 #include "G4VSensitiveDetector.hh"
11 
12 //#define DEBUG
13 
15  : alreadySet_(false), geomDet_(&det) {}
16 
18  if (alreadySet_)
19  return;
20  alreadySet_ = true;
21 
22  G4Navigator* theStdNavigator = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
23  G4Navigator theNavigator;
24  theNavigator.SetWorldVolume(theStdNavigator->GetWorldVolume());
25 
26  std::vector<const GeometricDet*> allSensitiveDets;
27  geomDet_->deepComponents(allSensitiveDets);
28  edm::LogVerbatim("TrackerSimInfoNumbering")
29  << " TouchableTo History: got " << allSensitiveDets.size() << " sensitive detectors from GeometricDet.";
30 
31  for (auto& theSD : allSensitiveDets) {
32  auto const& t = theSD->translation();
33  theNavigator.LocateGlobalPointAndSetup(G4ThreeVector(t.x(), t.y(), t.z()));
34  G4TouchableHistory* hist = theNavigator.CreateTouchableHistory();
36  touchToNavStory(hist, st);
37 
38  directMap_[st] = theSD->geographicalId();
39 
40  LogDebug("TrackerSimDebugNumbering") << " INSERTING LV " << hist->GetVolume()->GetLogicalVolume()->GetName()
41  << " SD: "
42  << hist->GetVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName()
43  << " Now size is " << directMap_.size();
44  delete hist;
45  }
46  edm::LogVerbatim("TrackerSimInfoNumbering")
47  << " TrackerG4SimHitNumberingScheme: mapped " << directMap_.size() << " detectors to Geant4.";
48 
49  if (directMap_.size() != allSensitiveDets.size()) {
50  edm::LogError("TrackerSimInfoNumbering") << " ERROR: GeomDet sensitive detectors do not match Geant4 ones.";
51  throw cms::Exception("TrackerG4SimHitNumberingScheme::buildAll")
52  << " cannot resolve structure of tracking sensitive detectors";
53  }
54 }
55 
58 #ifdef DEBUG
59  std::vector<int> debugint;
60  std::vector<std::string> debugstring;
61 #endif
62  int levels = v->GetHistoryDepth();
63 
64  for (int k = 0; k <= levels; ++k) {
65  if (dd4hep::dd::noNamespace(v->GetVolume(k)->GetLogicalVolume()->GetName()) != "TOBInactive") {
66  st.emplace_back(
67  std::pair<int, std::string>(v->GetVolume(k)->GetCopyNo(), v->GetVolume(k)->GetLogicalVolume()->GetName()));
68 #ifdef DEBUG
69  debugint.emplace_back(v->GetVolume(k)->GetCopyNo());
70  debugstring.emplace_back(v->GetVolume(k)->GetLogicalVolume()->GetName());
71 #endif
72  }
73  }
74 #ifdef DEBUG
75  LogDebug("TrackerSimDebugNumbering") << " G4 TrackerG4SimHitNumberingScheme " << debugint;
76  for (u_int32_t jj = 0; jj < debugstring.size(); jj++)
77  LogDebug("TrackerSimDebugNumbering") << " " << debugstring[jj];
78 #endif
79 }
80 
81 unsigned int TrackerG4SimHitNumberingScheme::g4ToNumberingScheme(const G4VTouchable* v) {
82  if (alreadySet_ == false) {
83  buildAll();
84  }
86  touchToNavStory(v, st);
87 
88 #ifdef DEBUG
89  dumpG4VPV(v);
90  LogDebug("TrackerSimDebugNumbering") << " Returning: " << directMap_[st];
91 #endif
92 
93  return directMap_[st];
94 }
95 
96 void TrackerG4SimHitNumberingScheme::dumpG4VPV(const G4VTouchable* v) {
97  int levels = v->GetHistoryDepth();
98 
99  LogDebug("TrackerSimDebugNumbering") << " NAME : " << v->GetVolume()->GetLogicalVolume()->GetName();
100  for (int k = 0; k <= levels; k++) {
101  LogDebug("TrackerSimInfoNumbering") << " Hist: " << v->GetVolume(k)->GetLogicalVolume()->GetName() << " Copy "
102  << v->GetVolume(k)->GetCopyNo();
103  }
104 }
TrackerG4SimHitNumberingScheme::buildAll
void buildAll()
Definition: TrackerG4SimHitNumberingScheme.cc:17
MessageLogger.h
funct::false
false
Definition: Factorize.h:29
TrackerG4SimHitNumberingScheme::Nav_Story
std::vector< std::pair< int, std::string > > Nav_Story
Definition: TrackerG4SimHitNumberingScheme.h:15
TrackerG4SimHitNumberingScheme::directMap_
DirectMapType directMap_
Definition: TrackerG4SimHitNumberingScheme.h:28
TrackerG4SimHitNumberingScheme::alreadySet_
bool alreadySet_
Definition: TrackerG4SimHitNumberingScheme.h:29
findQualityFiles.v
v
Definition: findQualityFiles.py:179
TrackerG4SimHitNumberingScheme::TrackerG4SimHitNumberingScheme
TrackerG4SimHitNumberingScheme(const GeometricDet &)
Definition: TrackerG4SimHitNumberingScheme.cc:14
GeometricDet
Definition: GeometricDet.h:31
dqmdumpme.k
k
Definition: dqmdumpme.py:60
TrackerG4SimHitNumberingScheme::touchToNavStory
void touchToNavStory(const G4VTouchable *, Nav_Story &)
Definition: TrackerG4SimHitNumberingScheme.cc:56
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
gpuVertexFinder::hist
__shared__ Hist hist
Definition: gpuClusterTracksDBSCAN.h:48
TrackerG4SimHitNumberingScheme::dumpG4VPV
void dumpG4VPV(const G4VTouchable *)
Definition: TrackerG4SimHitNumberingScheme.cc:96
TrackerG4SimHitNumberingScheme::geomDet_
const GeometricDet * geomDet_
Definition: TrackerG4SimHitNumberingScheme.h:30
TrackerG4SimHitNumberingScheme::g4ToNumberingScheme
unsigned int g4ToNumberingScheme(const G4VTouchable *)
Definition: TrackerG4SimHitNumberingScheme.cc:81
TrackerG4SimHitNumberingScheme.h
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
GeometricDet::deepComponents
ConstGeometricDetContainer deepComponents() const
Definition: GeometricDet.cc:252
GeometricDet.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
Exception
Definition: hltDiff.cc:245
findQualityFiles.jj
string jj
Definition: findQualityFiles.py:188
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
jets_cff.levels
levels
Definition: jets_cff.py:21