ATLAS Offline Software
Loading...
Searching...
No Matches
LooperKiller.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "LooperKiller.h"
6
7#include "G4RunManagerKernel.hh"
8#include "G4TransportationManager.hh"
9#include "G4Navigator.hh"
10#include "G4PropagatorInField.hh"
11#include "G4TrackingManager.hh"
12#include "G4SteppingManager.hh"
13#include "G4StackManager.hh"
14#include "G4EventManager.hh"
15#include "G4Event.hh"
16#include "MCTruth/TrackHelper.h"
22
23#include "GaudiKernel/Bootstrap.h"
24#include "GaudiKernel/ISvcLocator.h"
25#include "GaudiKernel/IMessageSvc.h"
26
27namespace G4UA
28{
29
30 //---------------------------------------------------------------------------
32 : AthMessaging(Gaudi::svcLocator()->service<IMessageSvc>("MessageSvc"),
33 "LooperKiller"),
34 m_evtStore("StoreGateSvc/StoreGateSvc", "LooperKiller"),
35 m_detStore("StoreGateSvc/DetectorStore", "LooperKiller"),
37 {
38 }
39
40 //---------------------------------------------------------------------------
41 void LooperKiller::UserSteppingAction(const G4Step* aStep)
42 {
43
44 if (aStep->GetTrack()->GetCurrentStepNumber() < m_config.MaxSteps) {
45 if (m_count_steps==0) return;
46 // Track recovered...
47 ATH_MSG_WARNING("Track finished on its own. Congrats. Moving on with the event.");
48 m_count_steps = 0;
49 G4TransportationManager *tm = G4TransportationManager::GetTransportationManager();
50 tm->GetNavigatorForTracking()->SetVerboseLevel(0);
51 tm->GetPropagatorInField()->SetVerboseLevel(0);
52 G4RunManagerKernel *rmk = G4RunManagerKernel::GetRunManagerKernel();
53 rmk->GetTrackingManager()->SetVerboseLevel(0);
54 rmk->GetTrackingManager()->GetSteppingManager()->SetVerboseLevel(0);
55 rmk->GetStackManager()->SetVerboseLevel(0);
56 return;
57 } else if (aStep->GetTrack()->GetCurrentStepNumber() == m_config.MaxSteps) {
58 ATH_MSG_WARNING("LooperKiller triggered!! Hold on to your hats!!!!!!!!" );
59 }
60
61 G4TransportationManager *tm = G4TransportationManager::GetTransportationManager();
62 tm->GetNavigatorForTracking()->SetVerboseLevel(m_config.VerboseLevel);
63 tm->GetPropagatorInField()->SetVerboseLevel(m_config.VerboseLevel);
64
65 G4RunManagerKernel *rmk = G4RunManagerKernel::GetRunManagerKernel();
66 rmk->GetTrackingManager()->SetVerboseLevel(m_config.VerboseLevel);
67 rmk->GetTrackingManager()->GetSteppingManager()->SetVerboseLevel(m_config.VerboseLevel);
68 rmk->GetStackManager()->SetVerboseLevel(m_config.VerboseLevel);
69
71
72 if (m_count_steps>m_config.PrintSteps) {
73 m_count_steps = 0;
74 m_report.killed_tracks++;
75 aStep->GetTrack()->SetTrackStatus(fStopAndKill);
76 tm->GetNavigatorForTracking()->SetVerboseLevel(0);
77 tm->GetPropagatorInField()->SetVerboseLevel(0);
78 rmk->GetTrackingManager()->SetVerboseLevel(0);
79 rmk->GetTrackingManager()->GetSteppingManager()->SetVerboseLevel(0);
80 rmk->GetStackManager()->SetVerboseLevel(0);
81 int pdg_id{0};
82 TrackHelper trackHelper(aStep->GetTrack());
83 if ( m_config.BSM_Only && (trackHelper.IsPrimary() || trackHelper.IsRegisteredSecondary()) ) {
85 if (part) { pdg_id = part->pdg_id(); }
86 }
87 if ( !m_config.BSM_Only || MC::isBSM(pdg_id)) { // Sometimes we may ony want to bail out for BSM particles.
88 // Bail out...
89 if (m_config.AbortEvent){
90 rmk->GetEventManager()->AbortCurrentEvent();
91 rmk->GetEventManager()->GetNonconstCurrentEvent()->SetEventAborted();
92 }
93 if (m_config.SetError){
94
95 // Set error state in eventInfo
96 SG::ReadHandle<xAOD::EventInfo> eic("McEventInfo");
97 if (! eic.isValid()){
98 ATH_MSG_WARNING( "Failed to retrieve EventInfo" );
99 } else {
100 eic->updateErrorState(xAOD::EventInfo::Core,xAOD::EventInfo::Error);
101 ATH_MSG_WARNING( "Set error state in event info!" );
102 }
103 } // End of set error
104 } // End of BSM-only check
105 const std::string name = aStep->GetTrack()->GetDefinition()->GetParticleName();
106 if ( ( m_config.BSM_Only && !MC::isBSM(pdg_id) ) || !m_config.AbortEvent ) {
107 ATH_MSG_INFO ("Quietly stopped tracking looping " << name
108 << " (trackID " << aStep->GetTrack()->GetTrackID()
109 << ", track pos: "<<aStep->GetTrack()->GetPosition()
110 << ", mom: "<<aStep->GetTrack()->GetMomentum()
111 << ", parentID " << aStep->GetTrack()->GetParentID() << ")");
112 } else {
113 ATH_MSG_WARNING ("Stopped tracking looping " << name
114 << " (trackID " << aStep->GetTrack()->GetTrackID()
115 << ", track pos: "<<aStep->GetTrack()->GetPosition()
116 << ", mom: "<<aStep->GetTrack()->GetMomentum()
117 << ", parentID " << aStep->GetTrack()->GetParentID() << "). The event will abort now.");
118 }
119 } // End of handling end of error time
120 }
121
122} // namespace G4UA
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
ATLAS-specific HepMC functions.
Handle class for reading from StoreGate.
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
Config m_config
My configuration options.
StoreGateSvc_t m_evtStore
Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
Pointer to StoreGate (detector store by default)
LooperKiller(const Config &config)
virtual void UserSteppingAction(const G4Step *) override
virtual bool isValid() override final
Can the handle be successfully dereferenced?
TrackInformation * GetTrackInformation()
Definition TrackHelper.h:28
bool IsPrimary() const
bool IsRegisteredSecondary() const
virtual HepMC::ConstGenParticlePtr GetCurrentGenParticle() const override
return a pointer to the GenParticle corresponding to the current G4Track (if there is one).
@ Core
Core flags describing the event.
@ Error
The sub-detector issued an error.
=============================================================================
GenParticle * GenParticlePtr
Definition GenParticle.h:37
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.