Belle II Software development
Geant4MaterialInterface.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8#include <framework/logging/Logger.h>
9#include <tracking/modules/genfitUtilities/Geant4MaterialInterface.h>
10#include <geometry/GeometryManager.h>
11
12#include "genfit/Exception.h"
13#include "genfit/RKTrackRep.h"
14
15#include <assert.h>
16#include <math.h>
17
18#include "G4ThreeVector.hh"
19#include "G4Navigator.hh"
20#include "G4VPhysicalVolume.hh"
21#include "G4LogicalVolume.hh"
22#include "G4Material.hh"
23#include "G4TouchableHistory.hh"
24
25
26#include <G4VExceptionHandler.hh>
27#include <G4StateManager.hh>
28
29
30using namespace Belle2;
31
32namespace Belle2 {
43 class Geant4MaterialInterfaceExceptioHandler : public G4VExceptionHandler {
44 public:
49
53 virtual bool Notify(const char* origin, const char* code, G4ExceptionSeverity serv, const char* description)
54 {
55 B2WARNING("Geant4 exception during material lookup -- origin: "
56 << origin << "\n code:"
57 << code << "\n description:"
58 << description);
59
60 if (serv == EventMustBeAborted || serv == RunMustBeAborted) {
61 m_inFailureState = true;
62 }
63
64 if (serv == FatalException || serv == FatalErrorInArgument) {
65 // on fatal exceptions, instruct Geant4 to create core dump
66 // and crash
67 return true;
68 }
69
70 // returning false will continue the program execution
71 // G4ExceptionSeverity = JustWarning will be covered by this, too
72 return false;
73 }
74
78 bool isInFailureState() const
79 {
80 return m_inFailureState;
81 }
82
88 {
89 m_inFailureState = false;
90 }
91
92 private:
96 bool m_inFailureState = false;
97 };
98
106 public:
107 G4SafeNavigator() = default;
108 ~G4SafeNavigator() = default;
109
113 G4VPhysicalVolume* GetWorldVolume() const { return nav_.GetWorldVolume(); }
114
118 void SetWorldVolume(G4VPhysicalVolume* pWorld)
119 {
120 nav_.SetWorldVolume(pWorld);
121 worldsolid_ = pWorld->GetLogicalVolume()->GetSolid();
122 }
123
128 G4VPhysicalVolume* LocateGlobalPointAndSetup(const G4ThreeVector& point,
129 const G4ThreeVector* direction = 0,
130 const G4bool pRelativeSearch = true,
131 const G4bool ignoreDirection = true);
132
136 G4double CheckNextStep(const G4ThreeVector& pGlobalPoint,
137 const G4ThreeVector& pDirection,
138 const G4double pCurrentProposedStepLength,
139 G4double& pNewSafety);
140
145 G4VPhysicalVolume* ResetHierarchyAndLocate(const G4ThreeVector& point,
146 const G4ThreeVector& direction,
147 const G4TouchableHistory& h);
148
153
157 G4TouchableHistory* CreateTouchableHistory() const
158 {
159 return nav_.CreateTouchableHistory();
160 }
161
162 private:
163
169 {
170 otherHandler = G4StateManager::GetStateManager()->GetExceptionHandler();
172 G4StateManager::GetStateManager()->SetExceptionHandler(&exceptionHandler);
173 }
174
180 {
181 G4StateManager::GetStateManager()->SetExceptionHandler(otherHandler);
182
183 // was the a problem in the last usage ?
185 genfit::Exception exc("Geant4MaterialInterface::findNextBoundary ==> Geant4 exception during geometry navigation", __LINE__,
186 __FILE__);
187 exc.setFatal();
188 throw exc;
189 }
190 }
191
193 G4ThreeVector lastpoint_;
195 G4VPhysicalVolume* lastvolume_{0};
197 G4Navigator nav_;
199 const G4VSolid* worldsolid_ {0};
201 G4VExceptionHandler* otherHandler = nullptr;
204 };
206}
207
209{
210 nav_.SetGeometricallyLimitedStep();
211}
212
213G4VPhysicalVolume* G4SafeNavigator::LocateGlobalPointAndSetup(const G4ThreeVector& point,
214 const G4ThreeVector* direction,
215 const G4bool pRelativeSearch,
216 const G4bool ignoreDirection)
217{
218 if (point == lastpoint_ && lastvolume_) {
219 return lastvolume_;
220 }
221 //B2INFO("### init: " << point);
222
224 G4VPhysicalVolume* volume = nav_.LocateGlobalPointAndSetup(point, direction, pRelativeSearch, ignoreDirection);
226
227 if (!volume) {
228 volume = nav_.GetWorldVolume();
229 }
230 // remember last point to speed up setup if possible
231 lastpoint_ = point;
232 lastvolume_ = volume;
233 return volume;
234}
235
236G4VPhysicalVolume* G4SafeNavigator::ResetHierarchyAndLocate(const G4ThreeVector& point,
237 const G4ThreeVector& direction,
238 const G4TouchableHistory& h)
239{
240 if (point == lastpoint_ && lastvolume_) {
241 return lastvolume_;
242 }
243 //B2INFO("### reset: " << point);
245 G4VPhysicalVolume* volume = nav_.ResetHierarchyAndLocate(point, direction, h);
247
248 if (!volume) {
249 volume = nav_.GetWorldVolume();
250 }
251 // remember last point to speed up setup if possible
252 lastpoint_ = point;
253 lastvolume_ = volume;
254 return volume;
255}
256
257G4double G4SafeNavigator::CheckNextStep(const G4ThreeVector& point,
258 const G4ThreeVector& direction,
259 const G4double pCurrentProposedStepLength,
260 G4double& pNewSafety)
261{
262 //make sure we're inside the world volume
263 if (worldsolid_->Inside(point) == kOutside) {
264 pNewSafety = worldsolid_->DistanceToIn(point);
265 return worldsolid_->DistanceToIn(point, direction);
266 }
267
269 const auto distance = nav_.CheckNextStep(point, direction, pCurrentProposedStepLength, pNewSafety);
271
272 return distance;
273}
274
275
276Geant4MaterialInterface::Geant4MaterialInterface()
277 : nav_(new G4SafeNavigator()), currentVolume_(0)
278{
279 G4VPhysicalVolume* world = geometry::GeometryManager::getInstance().getTopVolume();
280 nav_->SetWorldVolume(world);
281}
282
283Geant4MaterialInterface::~Geant4MaterialInterface()
284{
285}
286
287
288bool
289Geant4MaterialInterface::initTrack(double posX, double posY, double posZ,
290 double dirX, double dirY, double dirZ)
291{
292 G4ThreeVector pos(posX * CLHEP::cm, posY * CLHEP::cm, posZ * CLHEP::cm);
293 G4ThreeVector dir(dirX, dirY, dirZ);
294
295 if (m_takingFullStep) {
296 m_takingFullStep = false;
297 nav_->SetGeometricallyLimitedStep();
298 }
299
300 const G4VPhysicalVolume* newVolume = nav_->LocateGlobalPointAndSetup(pos, &dir);
301 bool volChanged = newVolume != currentVolume_;
302 currentVolume_ = newVolume;
303
304 return volChanged;
305}
306
307
308genfit::Material
310{
311 assert(currentVolume_);
312
313 const G4Material* mat = currentVolume_->GetLogicalVolume()->GetMaterial();
314
315 double density, Z, A, radiationLength, mEE;
316 if (mat->GetNumberOfElements() == 1) {
317 Z = mat->GetZ();
318 A = mat->GetA();
319 } else {
320 // Calculate weight-averaged A, Z
321 A = Z = 0;
322 for (unsigned i = 0; i < mat->GetNumberOfElements(); ++i) {
323 const G4Element* element = (*mat->GetElementVector())[i];
324 Z += element->GetZ() * mat->GetFractionVector()[i];
325 A += element->GetA() * mat->GetFractionVector()[i];
326 }
327 }
328
329 density = mat->GetDensity() / CLHEP::g * CLHEP::cm3;
330 // Z has correct units
331 A *= CLHEP::mole / CLHEP::g;
332 radiationLength = mat->GetRadlen() / CLHEP::cm;
333 mEE = mat->GetIonisation()->GetMeanExcitationEnergy() / CLHEP::eV;
334 return genfit::Material(density, Z, A, radiationLength, mEE);
335}
336
337
338double
339Geant4MaterialInterface::findNextBoundary(const genfit::RKTrackRep* rep,
340 const genfit::M1x7& stateOrig,
341 double sMax, // signed
342 bool varField)
343{
344 // This assumes that sMax is small enough to take only a single RK
345 // step. This restriction comes about because RKPropagate only
346 // takes one step.
347 const double delta(1.E-2); // cm, distance limit beneath which straight-line steps are taken.
348 const double epsilon(1.E-1); // cm, allowed upper bound on arch deviation from straight line
349
350 genfit::M1x3 SA;
351 genfit::M1x7 state7, oldState7;
352 oldState7 = stateOrig;
353
354 int stepSign(sMax < 0 ? -1 : 1);
355
356 G4ThreeVector pointOld(stateOrig[0] * CLHEP::cm,
357 stateOrig[1] * CLHEP::cm,
358 stateOrig[2] * CLHEP::cm);
359 G4ThreeVector dirOld(stepSign * stateOrig[3],
360 stepSign * stateOrig[4],
361 stepSign * stateOrig[5]);
362
363 double s = 0; // trajectory length to boundary
364
365 const unsigned maxIt = 300;
366 unsigned it = 0;
367
368 // Initialize the geometry to the current location (set by caller).
369 double safety;
370
371 m_takingFullStep = false;
372 double slDist = nav_->CheckNextStep(pointOld, dirOld, fabs(sMax) * CLHEP::cm, safety);
373 if (slDist == kInfinity) {
374 slDist = fabs(sMax);
375 m_takingFullStep = true;
376 } else {
377 slDist /= CLHEP::cm;
378 }
379 safety /= CLHEP::cm;
380
381 // No boundary in sight?
382 if (safety > fabs(sMax)) {
383 B2DEBUG(20, " next boundary is farther away than sMax ");
384 return stepSign * safety; // sMax
385 }
386
387 // Are we at the boundary?
388 if (slDist < delta) {
389 B2DEBUG(20, " very close to the boundary -> return @ it " << it
390 << " stepSign*slDist = "
391 << stepSign << "*" << slDist);
392 return stepSign * slDist;
393 }
394 double step = slDist;
395
396 while (1) {
397 if (++it > maxIt) {
398 genfit::Exception exc("Geant4MaterialInterface::findNextBoundary ==> maximum number of iterations exceeded", __LINE__, __FILE__);
399 exc.setFatal();
400 throw exc;
401 }
402
403 if (step < delta)
404 return stepSign * (s + step);
405
406 // We have to find whether there's any boundary on our path.
407
408 // Follow curved arch, then see if we may have missed a boundary.
409 // Always propagate complete way from original start to avoid
410 // inconsistent extrapolations. This is always a single RK step.
411 state7 = stateOrig;
412 rep->RKPropagate(state7, nullptr, SA, stepSign * (s + step), varField);
413
414 G4ThreeVector pos(state7[0] * CLHEP::cm, state7[1] * CLHEP::cm, state7[2] * CLHEP::cm);
415 G4ThreeVector dir(stepSign * state7[3], stepSign * state7[4], stepSign * state7[5]);
416
417 // Straight line distance between extrapolation finish and
418 // the end of the previously determined safe segment.
419 double dist2 = (pow(state7[0] - oldState7[0], 2)
420 + pow(state7[1] - oldState7[1], 2)
421 + pow(state7[2] - oldState7[2], 2));
422
423 // If we moved less than safety, the volume cannot possibly have
424 // changed, so we skip further checks.
425 if (dist2 > safety * safety) {
426
427 // Do we need to try again with a shorter step? There are two
428 // possible reasons:
429 //
430 // 1. Too much curvature?
431
432 // Maximal lateral deviation² of the curved path from the
433 // straight line connecting beginning and end.
434 double maxDeviation2 = 0.25 * (step * step - dist2);
435 if (maxDeviation2 > epsilon * epsilon) {
436 // Need to take a shorter step to reliably estimate material,
437 // but only if we didn't move by safety.
438
439 // Take a shorter step, but never shorter than safety.
440 step = safety + 0.5 * (step - safety);
441
442 continue;
443 }
444
445 // 2. Volume changed?
446 //
447 // Where are we after the step?
448
449 if (m_takingFullStep) {
450 m_takingFullStep = false;
451 nav_->SetGeometricallyLimitedStep();
452 }
453
454 std::unique_ptr<G4TouchableHistory> hist(nav_->CreateTouchableHistory());
455 G4VPhysicalVolume* newVolume = nav_->LocateGlobalPointAndSetup(pos, &dir);
456
457 if (newVolume != currentVolume_) {
458
459 // Volume changed during the extrapolation.
460
461 // Extrapolation may not take the exact step length we asked
462 // for, so it can happen that a requested step < safety takes
463 // us across the boundary. This is then the best estimate we
464 // can get of the distance to the boundary with the stepper.
465 if (step <= safety)
466 return stepSign * (s + step);
467
468 // Move back to last good point, but looking in the actual
469 // direction of the step.
470 G4ThreeVector dirCloser(pos - pointOld);
471 dirCloser.setMag(1.);
472 nav_->ResetHierarchyAndLocate(pointOld, dirCloser, *hist);
473
474 // Look along the secant of the actual trajectory instead of
475 // the original direction. There should be a crossing within
476 // distance step.
477 double secantDist = nav_->CheckNextStep(pointOld, dirCloser,
478 step * CLHEP::cm, safety) / CLHEP::cm;
479 safety /= CLHEP::cm;
480 if (secantDist >= step) {
481 // Cannot be. Just take a shorter step, and hope that this
482 // works.
483 slDist = secantDist;
484 step = std::max(0.9 * step, safety);
485 } else {
486 slDist = step = std::max(secantDist, safety);
487 }
488
489 // Are we at the boundary?
490 if (slDist < delta) {
491 B2DEBUG(20, " very close to the boundary -> return @ it " << it
492 << " stepSign*(s + slDist) = "
493 << stepSign << "*(" << s + slDist << ")");
494 return stepSign * (s + slDist);
495 }
496
497 continue;
498 }
499 }
500
501 // We're in the new place, the step was safe, advance.
502 s += step;
503
504 oldState7 = state7;
505 pointOld = pos;
506 if (m_takingFullStep) {
507 m_takingFullStep = false;
508 // inform the navigator that the full geometrical step was taken. This is required for
509 // some Geant4 volume enter/exit optimizations to work.
510 nav_->SetGeometricallyLimitedStep();
511 }
512 nav_->LocateGlobalPointAndSetup(pos, &dir);
513 slDist = nav_->CheckNextStep(pos, dir,
514 (fabs(sMax) - s) * CLHEP::cm, safety);
515 if (slDist == kInfinity) {
516 m_takingFullStep = true;
517 slDist = fabs(sMax) - s;
518 } else {
519 slDist /= CLHEP::cm;
520 }
521 safety /= CLHEP::cm;
522 step = slDist;
523
524 // No boundary in sight?
525 if (s + safety > fabs(sMax)) {
526 B2DEBUG(20, " next boundary is farther away than sMax");
527 return stepSign * (s + safety); // sMax
528 }
529
530 // Are we at the boundary?
531 if (slDist < delta) {
532 B2DEBUG(20, " very close to the boundary -> return @ it " << it
533 << " stepSign*(s + slDist) = "
534 << stepSign << "*(" << s + slDist << ")");
535 return stepSign * (s + slDist);
536 }
537 }
538}
Guards against leaving the physical volume.
void SetGeometricallyLimitedStep()
Calls Geant4's SetGeometricallyLimitedStep.
G4Navigator nav_
Geant4's navigator which calls are forwarded to.
G4VExceptionHandler * otherHandler
Stores the pointer to exception handler which this class temporarily replaces for most calls.
G4TouchableHistory * CreateTouchableHistory() const
Call Geant4's CreateTouchableHistory.
G4VPhysicalVolume * GetWorldVolume() const
Returns the Geant4 world volume.
void uninstallAndCheckOurExceptionHandler()
Reinstate the previously used exception handler and check whether a critical exception was recorded b...
G4ThreeVector lastpoint_
the last point which has been queried with G4
G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Use Geant4's LocateGlobalPointAndSetup to get a G4VPhysicalVolume or use cached values from previous ...
const G4VSolid * worldsolid_
The topmost solid of the G4 world.
G4double CheckNextStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Check if within world volume and call Geant4's CheckNextStep.
Geant4MaterialInterfaceExceptioHandler exceptionHandler
Custom exception handler to handle stuck tracks properly (and abort)
G4VPhysicalVolume * lastvolume_
the last volume which has been queried
void installOurExceptionHandler()
Install our specific exception handler and cache the one which was set before.
G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Call Geant4's ResetHierarchyAndLocate to get a G4VPhysicalVolume or use cached values from previous c...
void SetWorldVolume(G4VPhysicalVolume *pWorld)
Sets the Geant4 world volume.
This class implements a custom exception handler for Geant4 which is used to record whether a critica...
void resetFailureState()
Reset the recorded failure state to be ready for the next calls to Geant4.
bool isInFailureState() const
Returns true if a problematic exception was encountered.
virtual ~Geant4MaterialInterfaceExceptioHandler()=default
virtual destructor for proper resource de-allocation
bool m_inFailureState
Stores whether a problematic exception occured.
virtual bool Notify(const char *origin, const char *code, G4ExceptionSeverity serv, const char *description)
G4VExceptionHandler method called when an exception is raised.
bool initTrack(double posX, double posY, double posZ, double dirX, double dirY, double dirZ) override
Initialize the navigator at given position and with given direction.
bool m_takingFullStep
stores whether to call SetGeometricallyLimitedStep() because the full step length was taken.
genfit::Material getMaterialParameters() override
Get material parameters in current material.
const class G4VPhysicalVolume * currentVolume_
the volume the extraoplation is currenly located in
double findNextBoundary(const genfit::RKTrackRep *rep, const genfit::M1x7 &state7, double sMax, bool varField=true) override
Make a step (following the curvature) until step length sMax or the next boundary is reached.
std::unique_ptr< class G4SafeNavigator > nav_
holds a object of G4SafeNavigator, which is located in Geant4MaterialInterface.cc
G4VPhysicalVolume * getTopVolume()
Return a pointer to the top volume.
static GeometryManager & getInstance()
Return a reference to the instance.
Abstract base class for different kinds of events.