Belle II Software  release-05-02-19
OverlapCheckerModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2011 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Thomas Kuhr *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <geometry/modules/overlapchecker/OverlapCheckerModule.h>
12 #include <geometry/GeometryManager.h>
13 #include <framework/gearbox/Unit.h>
14 
15 #include "G4VPhysicalVolume.hh"
16 #include "G4LogicalVolume.hh"
17 #include "G4UnitsTable.hh"
18 #include "G4VExceptionHandler.hh"
19 #include "G4StateManager.hh"
20 #include "G4AffineTransform.hh"
21 #include "G4VSolid.hh"
22 
23 #include <regex>
24 #include <utility>
25 
26 using namespace Belle2;
27 
28 REG_MODULE(OverlapChecker);
29 
30 namespace {
32  class OverlapHandler: public G4VExceptionHandler {
33  public:
35  explicit OverlapHandler(std::function<void(const std::string&)> callback): m_callback(std::move(callback)) {}
37  bool Notify(const char*, const char*, G4ExceptionSeverity, const char* description) override
38  {
39  m_callback(description);
40  return false;
41  }
42  private:
44  std::function<void(const std::string&)> m_callback;
45  };
46 }
47 
49 {
50  //Set module properties and the description
51  setDescription("Checks the geometry for overlaps.");
52 
53  //Parameter definition
54  addParam("points", m_points, R"DOC(
55 Number of test points we will generate on randomly on the surface of each geometry volume and then check for all of them that
56 
57 1. The points are inside the mother volume
58 2. The points are not inside any neighbor volume
59 
60 The higher the number the more precise the check for overlaps becomes, and the slower it gets.
61 See also https://questions.belle2.org/question/7264/ )DOC",
62  m_points);
63  addParam("tolerance", m_tolerance, "Tolerance of overlap check.", m_tolerance);
64  addParam("maxErrors", m_maxErrors, "Number of overlap errors per volume before continuing with next volume", m_maxErrors);
65  addParam("maxDepth", m_maxDepth, "Maximum depth to go into the geometry tree, 0 means no maximum", m_maxDepth);
66  addParam("prefix", m_prefix, "Prefix path, only volumes starting with the given path are checked", m_prefix);
67 }
68 
70 {
71  m_displayData.registerInDataStore();
72  if (m_maxErrors <= 0) m_maxErrors = m_points + 1;
73 }
74 
76 {
77  if (!m_displayData.isValid()) m_displayData.create();
78  //Check the geometry tree for overlaps
79  G4UnitDefinition::BuildUnitsTable();
80  G4VPhysicalVolume* volume = geometry::GeometryManager::getInstance().getTopVolume();
81  if (!volume) {
82  B2ERROR("No geometry found. => Add the Geometry module to the path before the OverlapChecker module.");
83  return;
84  }
85  m_seen.clear();
86  // reset the navigation history to be in the top level volume
87  m_nav.Reset();
88  // remember the existing G4Exception handler and set our own one to find the position of the intersection
89  G4VExceptionHandler* old = G4StateManager::GetStateManager()->GetExceptionHandler();
90  // set our own exception handler
91  OverlapHandler handler([&](const std::string & message) { handleOverlap(message); });
92  G4StateManager::GetStateManager()->SetExceptionHandler(&handler);
93  // now check for overlaps
94  checkVolume(volume, "");
95  // and reset the exception handler
96  G4StateManager::GetStateManager()->SetExceptionHandler(old);
97 
98  //Print the list of found overlaps
99  for (const auto& m_overlap : m_overlaps) {
100  B2ERROR("Overlaps detected for " << m_overlap);
101  }
102 }
103 
104 void OverlapCheckerModule::handleOverlap(const std::string& geant4Message)
105 {
106  // ok, let's handle the Overlap message
107  G4VPhysicalVolume* volume = m_nav.GetTopVolume();
108  m_nav.BackLevel();
109  B2ERROR(geant4Message);
110  std::regex r(R"((mother)?\s*local point \(([-+0-9eE.]+),([-+0-9eE.]+),([-+0-9eE.]+)\))");
111  std::smatch m;
112  if (std::regex_search(geant4Message, m, r)) {
113  G4ThreeVector posLocal(std::atof(m[2].str().c_str()), std::atof(m[3].str().c_str()), std::atof(m[4].str().c_str()));
114  if (m[1].length() == 0) {
115  // Damn you Geant4, giving me the coordinates in the coordinate system of
116  // the sister volume which I don't know exactly is almost useless. Aaah
117  // darn it, let's check all sisters, see if they have the correct name
118  // and if so check if the point transformed to the mother system is
119  // inside the one were checking. So let's remember the transformation
120  // relative to the top volume of our current volume.
121  G4AffineTransform trans_volume(volume->GetRotation(), volume->GetTranslation());
122  trans_volume.Invert();
123  // And the solid
124  G4VSolid* solid = volume->GetLogicalVolume()->GetSolid();
125  // Now look for the name of the intersecting volume in the exception message (jaaaay)
126  std::regex nameRegex("with (.*) volume's");
127  std::smatch nameMatch;
128  // And if we can find the name we can draw the intersection
129  if (std::regex_search(geant4Message, nameMatch, nameRegex)) {
130  const std::string& name = nameMatch[1].str();
131  // By looping over all sisters
132  for (int i = 0; i < volume->GetMotherLogical()->GetNoDaughters(); ++i) {
133  G4VPhysicalVolume* sister = volume->GetMotherLogical()->GetDaughter(i);
134  // ignoring the ones which don't match the name
135  if (name != sister->GetName()) continue;
136  // now transform the point into the coordinate system of the volume we actually look at
137  G4AffineTransform trans_sister(sister->GetRotation(), sister->GetTranslation());
138  G4ThreeVector posMother = trans_sister.TransformPoint(posLocal);
139  G4ThreeVector posSister = trans_volume.TransformPoint(posMother);
140  // and check if the point is inside our volume
141  if (solid->Inside(posSister) != kOutside) {
142  // if so this is the right volume which is intersecting.
143  B2INFO("Found intersecting volume " << sister->GetName() << "." << sister->GetCopyNo());
144  // so add the point to the display data
145  G4AffineTransform t = m_nav.GetTopTransform().Inverse();
146  G4ThreeVector global = t.TransformPoint(posMother);
147  m_displayData->addPoint(geant4Message, TVector3(global[0], global[1], global[2]) * Unit::mm);
148  }
149  }
150  } else {
151  B2ERROR("Could not find name of intersecting volume");
152  }
153  } else {
154  // in case of overlap with mother the life is so much simpler: just convert the point to global and save it.
155  G4AffineTransform t = m_nav.GetTopTransform().Inverse();
156  G4ThreeVector global = t.TransformPoint(posLocal);
157  m_displayData->addPoint(geant4Message, TVector3(global[0], global[1], global[2]) * Unit::mm);
158  }
159  }
160  m_nav.NewLevel(volume);
161 }
162 
163 bool OverlapCheckerModule::checkVolume(G4VPhysicalVolume* volume, const std::string& path, int depth)
164 {
165  // check if we exceeded maximum recursion depth
166  if (m_maxDepth > 0 && depth >= m_maxDepth) return false;
167  // add the volume to the navigation history
168  m_nav.NewLevel(volume);
169  // remember the path to the volume
170  std::string volumePath = path + "/" + volume->GetName();
171  bool result{false};
172  // check if we have a prefix we have to match. if not prefix or prefix matches, check for overlaps
173  if (m_prefix.empty() || volumePath.substr(0, m_prefix.size()) == m_prefix) {
174  result = volume->CheckOverlaps(m_points, m_tolerance, true, m_maxErrors);
175  if (result) {
176  m_overlaps.push_back(volumePath);
177  }
178  }
179 
180  //Check the daughter volumes for overlaps
181  G4LogicalVolume* logicalVolume = volume->GetLogicalVolume();
182  for (int iDaughter = 0; iDaughter < logicalVolume->GetNoDaughters(); iDaughter++) {
183  G4VPhysicalVolume* daughter = logicalVolume->GetDaughter(iDaughter);
184  // check if we already checked this particular volume, if so skip it
185  auto it = m_seen.insert(daughter);
186  if (!it.second) continue;
187  result |= checkVolume(daughter, volumePath, depth + 1);
188  }
189  m_nav.BackLevel();
190  return result;
191 }
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::OverlapCheckerModule::m_maxErrors
int m_maxErrors
maximum number of errors before skipping current volume
Definition: OverlapCheckerModule.h:73
Belle2::OverlapCheckerModule::initialize
void initialize() override
Initialize the module.
Definition: OverlapCheckerModule.cc:69
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::OverlapCheckerModule::checkVolume
bool checkVolume(G4VPhysicalVolume *volume, const std::string &path, int depth=0)
Check a volume for overlaps.
Definition: OverlapCheckerModule.cc:163
Belle2::OverlapCheckerModule::m_nav
G4NavigationHistory m_nav
navigation history to remember coordinate transformations
Definition: OverlapCheckerModule.h:79
Belle2::OverlapCheckerModule::event
void event() override
event function: this runs the overlap checker for each event
Definition: OverlapCheckerModule.cc:75
Belle2::geometry::GeometryManager::getTopVolume
G4VPhysicalVolume * getTopVolume()
Return a pointer to the top volume.
Definition: GeometryManager.h:59
Belle2::geometry::GeometryManager::getInstance
static GeometryManager & getInstance()
Return a reference to the instance.
Definition: GeometryManager.cc:98
Belle2::OverlapCheckerModule::m_points
int m_points
number of test points
Definition: OverlapCheckerModule.h:72
Belle2::OverlapCheckerModule::m_overlaps
std::vector< std::string > m_overlaps
list of overlapping volumes
Definition: OverlapCheckerModule.h:77
Belle2::OverlapCheckerModule::m_seen
std::set< G4VPhysicalVolume * > m_seen
set of logical volumes we already checked
Definition: OverlapCheckerModule.h:78
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::OverlapCheckerModule::m_prefix
std::string m_prefix
check only volumes beginning with prefix
Definition: OverlapCheckerModule.h:76
Belle2::OverlapCheckerModule::OverlapCheckerModule
OverlapCheckerModule()
Constructor of the module.
Definition: OverlapCheckerModule.cc:48
Belle2::OverlapCheckerModule::m_maxDepth
int m_maxDepth
maximum depth to check
Definition: OverlapCheckerModule.h:74
Belle2::Module::addParam
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:562
Belle2::Unit::mm
static const double mm
[millimeters]
Definition: Unit.h:80
Belle2::OverlapCheckerModule::m_displayData
StoreObjPtr< DisplayData > m_displayData
Pointer to the DisplayData where we add the overlap points for rendering.
Definition: OverlapCheckerModule.h:80
Belle2::OverlapCheckerModule::m_tolerance
double m_tolerance
tolerance of overlap check
Definition: OverlapCheckerModule.h:75
Belle2::OverlapCheckerModule::handleOverlap
void handleOverlap(const std::string &geant4Message)
Handle a G4Exception with the overlap message issued by Geant4.
Definition: OverlapCheckerModule.cc:104