8#include <framework/logging/Logger.h>
9#include <tracking/modules/genfitUtilities/Geant4MaterialInterface.h>
10#include <geometry/GeometryManager.h>
12#include "genfit/Exception.h"
13#include "genfit/RKTrackRep.h"
18#include "G4ThreeVector.hh"
19#include "G4Navigator.hh"
20#include "G4VPhysicalVolume.hh"
21#include "G4LogicalVolume.hh"
22#include "G4Material.hh"
23#include "G4TouchableHistory.hh"
26#include <G4VExceptionHandler.hh>
27#include <G4StateManager.hh>
53 virtual bool Notify(
const char* origin,
const char* code, G4ExceptionSeverity serv,
const char* description)
55 B2WARNING(
"Geant4 exception during material lookup -- origin: "
56 << origin <<
"\n code:"
57 << code <<
"\n description:"
60 if (serv == EventMustBeAborted || serv == RunMustBeAborted) {
64 if (serv == FatalException || serv == FatalErrorInArgument) {
120 nav_.SetWorldVolume(pWorld);
121 worldsolid_ = pWorld->GetLogicalVolume()->GetSolid();
129 const G4ThreeVector* direction = 0,
130 const G4bool pRelativeSearch =
true,
131 const G4bool ignoreDirection =
true);
137 const G4ThreeVector& pDirection,
138 const G4double pCurrentProposedStepLength,
139 G4double& pNewSafety);
146 const G4ThreeVector& direction,
147 const G4TouchableHistory& h);
159 return nav_.CreateTouchableHistory();
170 otherHandler = G4StateManager::GetStateManager()->GetExceptionHandler();
181 G4StateManager::GetStateManager()->SetExceptionHandler(
otherHandler);
185 genfit::Exception exc(
"Geant4MaterialInterface::findNextBoundary ==> Geant4 exception during geometry navigation", __LINE__,
210 nav_.SetGeometricallyLimitedStep();
214 const G4ThreeVector* direction,
215 const G4bool pRelativeSearch,
216 const G4bool ignoreDirection)
224 G4VPhysicalVolume* volume =
nav_.LocateGlobalPointAndSetup(point, direction, pRelativeSearch, ignoreDirection);
228 volume =
nav_.GetWorldVolume();
237 const G4ThreeVector& direction,
238 const G4TouchableHistory& h)
245 G4VPhysicalVolume* volume =
nav_.ResetHierarchyAndLocate(point, direction, h);
249 volume =
nav_.GetWorldVolume();
258 const G4ThreeVector& direction,
259 const G4double pCurrentProposedStepLength,
260 G4double& pNewSafety)
265 return worldsolid_->DistanceToIn(point, direction);
269 const auto distance =
nav_.CheckNextStep(point, direction, pCurrentProposedStepLength, pNewSafety);
276Geant4MaterialInterface::Geant4MaterialInterface()
280 nav_->SetWorldVolume(world);
283Geant4MaterialInterface::~Geant4MaterialInterface()
290 double dirX,
double dirY,
double dirZ)
292 G4ThreeVector pos(posX * CLHEP::cm, posY * CLHEP::cm, posZ * CLHEP::cm);
293 G4ThreeVector dir(dirX, dirY, dirZ);
297 nav_->SetGeometricallyLimitedStep();
300 const G4VPhysicalVolume* newVolume =
nav_->LocateGlobalPointAndSetup(pos, &dir);
313 const G4Material* mat =
currentVolume_->GetLogicalVolume()->GetMaterial();
315 double density, Z, A, radiationLength, mEE;
316 if (mat->GetNumberOfElements() == 1) {
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];
329 density = mat->GetDensity() / CLHEP::g * CLHEP::cm3;
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);
340 const genfit::M1x7& stateOrig,
347 const double delta(1.E-2);
348 const double epsilon(1.E-1);
351 genfit::M1x7 state7, oldState7;
352 oldState7 = stateOrig;
354 int stepSign(sMax < 0 ? -1 : 1);
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]);
365 const unsigned maxIt = 300;
372 double slDist =
nav_->CheckNextStep(pointOld, dirOld, fabs(sMax) * CLHEP::cm, safety);
373 if (slDist == kInfinity) {
382 if (safety > fabs(sMax)) {
383 B2DEBUG(20,
" next boundary is farther away than sMax ");
384 return stepSign * safety;
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;
394 double step = slDist;
398 genfit::Exception exc(
"Geant4MaterialInterface::findNextBoundary ==> maximum number of iterations exceeded", __LINE__, __FILE__);
404 return stepSign * (s + step);
412 rep->RKPropagate(state7,
nullptr, SA, stepSign * (s + step), varField);
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]);
419 double dist2 = (pow(state7[0] - oldState7[0], 2)
420 + pow(state7[1] - oldState7[1], 2)
421 + pow(state7[2] - oldState7[2], 2));
425 if (dist2 > safety * safety) {
434 double maxDeviation2 = 0.25 * (step * step - dist2);
435 if (maxDeviation2 > epsilon * epsilon) {
440 step = safety + 0.5 * (step - safety);
451 nav_->SetGeometricallyLimitedStep();
454 std::unique_ptr<G4TouchableHistory> hist(
nav_->CreateTouchableHistory());
455 G4VPhysicalVolume* newVolume =
nav_->LocateGlobalPointAndSetup(pos, &dir);
466 return stepSign * (s + step);
470 G4ThreeVector dirCloser(pos - pointOld);
471 dirCloser.setMag(1.);
472 nav_->ResetHierarchyAndLocate(pointOld, dirCloser, *hist);
477 double secantDist =
nav_->CheckNextStep(pointOld, dirCloser,
478 step * CLHEP::cm, safety) / CLHEP::cm;
480 if (secantDist >= step) {
484 step = std::max(0.9 * step, safety);
486 slDist = step = std::max(secantDist, safety);
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);
510 nav_->SetGeometricallyLimitedStep();
512 nav_->LocateGlobalPointAndSetup(pos, &dir);
513 slDist =
nav_->CheckNextStep(pos, dir,
514 (fabs(sMax) - s) * CLHEP::cm, safety);
515 if (slDist == kInfinity) {
517 slDist = fabs(sMax) - s;
525 if (s + safety > fabs(sMax)) {
526 B2DEBUG(20,
" next boundary is farther away than sMax");
527 return stepSign * (s + safety);
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);
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 occurred.
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 extrapolation is currently 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.