BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSLinkEventAction.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSAuxiliaryNavigator.hh"
20#include "BDSDebug.hh"
21#include "BDSGlobalConstants.hh"
22#include "BDSHitSamplerLink.hh"
23#include "BDSLinkEventAction.hh"
24#include "BDSLinkEventInfo.hh"
25#include "BDSLinkRunAction.hh"
26#include "BDSOutput.hh"
27#include "BDSSDSampler.hh"
28#include "BDSSDSamplerLink.hh"
29#include "BDSSDManager.hh"
30#include "BDSUtilities.hh"
31
32#include "globals.hh"
33#include "G4Event.hh"
34#include "G4EventManager.hh"
35#include "G4HCofThisEvent.hh"
36#include "G4PropagatorInField.hh"
37#include "G4Run.hh"
38#include "G4SDManager.hh"
39#include "G4THitsMap.hh"
40#include "G4TransportationManager.hh"
41#include "G4VUserEventInformation.hh"
42
43#include <map>
44
45
46BDSLinkEventAction::BDSLinkEventAction(BDSOutput* outputIn,
47 BDSLinkRunAction* runActionIn,
48 G4bool debugIn):
49 output(outputIn),
50 runAction(runActionIn),
51 debug(debugIn),
52 collIDSamplerLink(-1),
53 collIDSampler(-1),
54 currentEventIndex(0),
55 primaryAbsorbedInCollimator(false)
56{
58 verboseEventBDSIM = globals->VerboseEventBDSIM();
59 verboseEventStart = globals->VerboseEventStart();
60 verboseEventStop = BDS::VerboseEventStop(verboseEventStart, globals->VerboseEventContinueFor());
61 printModulo = globals->PrintModuloEvents();
62}
63
64BDSLinkEventAction::~BDSLinkEventAction()
65{;}
66
67void BDSLinkEventAction::BeginOfEventAction(const G4Event* evt)
68{
69 currentEventIndex = evt->GetEventID();
70 primaryAbsorbedInCollimator = false;
71
72 // reset navigators to ensure no mis-navigating and that events are truly independent
73 BDSAuxiliaryNavigator::ResetNavigatorStates();
74 // make further attempts to clear Geant4's tracking history between events to make them
75 // truly independent.
76 G4Navigator* trackingNavigator = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
77 trackingNavigator->ResetStackAndState();
78 G4TransportationManager* tm = G4TransportationManager::GetTransportationManager();
79 int i = 0;
80 for (auto it = tm->GetActiveNavigatorsIterator(); i < (int)tm->GetNoActiveNavigators(); it++)
81 {(*it)->ResetStackAndState(); i++;}
82 tm->GetPropagatorInField()->ClearPropagatorState(); // <- this one really makes a difference
83
84 // number feedback
85 if (debug)
86 {
87 if (currentEventIndex % printModulo == 0)
88 {G4cout << "---> Begin of event: " << currentEventIndex << G4endl;}
89 }
90 if (verboseEventBDSIM) // always print this out
91 {G4cout << __METHOD_NAME__ << "event #" << currentEventIndex << G4endl;}
92
93 // cache hit collection IDs for quicker access
94 if (collIDSamplerLink < 0)
95 {// if one is -1 then all need initialised.
96 G4SDManager* g4SDMan = G4SDManager::GetSDMpointer();
97 BDSSDManager* bdsSDMan = BDSSDManager::Instance();
98 collIDSamplerLink = g4SDMan->GetCollectionID(bdsSDMan->SamplerLink()->GetName());
99 collIDSampler = g4SDMan->GetCollectionID(bdsSDMan->SamplerPlane()->GetName());
100 }
101}
102
103void BDSLinkEventAction::EndOfEventAction(const G4Event* evt)
104{
105 // Get the hits collection of this event - all hits from different SDs.
106 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
107 typedef BDSHitsCollectionSamplerLink slhc;
108 slhc* samplerLink = HCE ? dynamic_cast<slhc*>(HCE->GetHC(collIDSamplerLink)) : nullptr;
109 typedef BDSHitsCollectionSampler shc;
110 shc* sampHC = HCE ? dynamic_cast<shc*>(HCE->GetHC(collIDSampler)) : nullptr;
111 std::vector<shc*> allSamplerHits = {sampHC};
112
113 G4VUserEventInformation* evtInfoG4 = evt->GetUserInformation();
114 BDSLinkEventInfo* evtInfo = dynamic_cast<BDSLinkEventInfo*>(evtInfoG4);
115 G4int primaryExternalParticleID = 0;
116 G4int primaryExternalParentID = 0;
117 if (evtInfo)
118 {
119 primaryExternalParticleID = evtInfo->externalParticleIDofPrimary;
120 primaryExternalParentID = evtInfo->externalParentIDofPrimary;
121 }
122
123 if (!samplerLink)
124 {return;}
125 if (samplerLink->entries() <= 0)
126 {return;}
127 else
128 {runAction->AppendHits(currentEventIndex, primaryExternalParticleID, primaryExternalParentID, samplerLink);}
129
130 output->FillEvent(nullptr,
131 evt->GetPrimaryVertex(),
132 allSamplerHits,
133 std::vector<BDSHitsCollectionSamplerCylinder*>(),
134 std::vector<BDSHitsCollectionSamplerSphere*>(),
135 samplerLink,
136 nullptr,
137 nullptr,
138 nullptr,
139 nullptr,
140 nullptr,
141 nullptr,
142 nullptr,
143 std::vector<const BDSTrajectoryPointHit*>(),
144 std::vector<const BDSTrajectoryPointHit*>(),
145 nullptr,
146 nullptr,
147 nullptr,
148 std::map<G4String, G4THitsMap<G4double>*>(),
149 BDSGlobalConstants::Instance()->TurnsTaken());
150}
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
BDSOutput * output
Cache of output instance. Not owned by this class.
Simple extension to cache extra variables through an event.
Simplified run action to hold link hits.
Output base class that defines interface for all output types.
Definition: BDSOutput.hh:73
void FillEvent(const BDSEventInfo *info, const G4PrimaryVertex *vertex, const std::vector< BDSHitsCollectionSampler * > &samplerHitsPlane, const std::vector< BDSHitsCollectionSamplerCylinder * > &samplerHitsCylinder, const std::vector< BDSHitsCollectionSamplerSphere * > &samplerHitsSphere, const BDSHitsCollectionSamplerLink *samplerHitsLink, const BDSHitsCollectionEnergyDeposition *energyLoss, const BDSHitsCollectionEnergyDeposition *energyLossFull, const BDSHitsCollectionEnergyDeposition *energyLossVacuum, const BDSHitsCollectionEnergyDeposition *energyLossTunnel, const BDSHitsCollectionEnergyDepositionGlobal *energyLossWorld, const BDSHitsCollectionEnergyDepositionGlobal *energyLossWorldContents, const BDSHitsCollectionEnergyDepositionGlobal *worldExitHits, const std::vector< const BDSTrajectoryPointHit * > &primaryHits, const std::vector< const BDSTrajectoryPointHit * > &primaryLosses, const BDSTrajectoriesToStore *trajectories, const BDSHitsCollectionCollimator *collimatorHits, const BDSHitsCollectionApertureImpacts *apertureImpactHits, const std::map< G4String, G4THitsMap< G4double > * > &scorerHitsMap, const G4int turnsTaken)
Copy event information from Geant4 simulation structures to output structures.
Definition: BDSOutput.cc:297
A singleton class that holds all required sensitive detector class instances.
Definition: BDSSDManager.hh:68
BDSSDSamplerLink * SamplerLink() const
SD for link samplers.
Definition: BDSSDManager.hh:95
BDSSDSampler * SamplerPlane() const
SD for samplers (plane type). See also SamplerPlaneWithFilter below.
Definition: BDSSDManager.hh:86
G4int VerboseEventStop(G4int verboseEventStart, G4int verboseEventContinueFor)