ATLAS Offline Software
Loading...
Searching...
No Matches
SimKernelMT.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
10
11#include "SimKernelMT.h"
12
13// ISF includes
19
20// STL
21#include <memory>
22#include <queue>
23#include <utility>
24
25#undef ISFDEBUG
26
28
29 ATH_CHECK( m_simulationTools.retrieve() );
30 for ( auto& curSimTool : m_simulationTools )
31 {
32 if ( curSimTool )
33 {
34 const auto flavor = curSimTool->simFlavor();
35 auto itr = m_simToolMap.find(flavor);
36 if (itr != m_simToolMap.end() )
37 {
38 ATH_MSG_FATAL("Two ISimulatorTool instances (" << itr->second->name() << "," << curSimTool->name() << ") with the same flavor in this job!\n Check your configuration!");
39 return StatusCode::FAILURE;
40 }
41 // New flavour add it to the map.
42 m_simToolMap[flavor] = &*curSimTool;
43 if (flavor == ISF::ParticleKiller) {
44 m_particleKillerTool = &*curSimTool;
45 }
46 }
47 }
48 // Check that a particle killer simulator tool was provided
50 ATH_MSG_FATAL("No fallback ParticleKiller Simulator Tool provided in SimulationTools, the job will bail out now.");
51 return StatusCode::FAILURE;
52 }
53
54 ATH_MSG_INFO("The following Simulators will be used in this job: \t" << m_simulationTools);
55 // retrieve simulation selectors (i.e. the "routing chain")
56 for ( auto& selectorsToolHandleArray: m_simSelectors )
57 {
58 ATH_CHECK( selectorsToolHandleArray.retrieve() );
59 }
60
61 // info screen output
62 ATH_MSG_INFO( "The following routing chains are defined:" );
64 {
65 auto& localSelectors = m_simSelectors[geoID];
66 ATH_MSG_INFO( AtlasDetDescr::AtlasRegionHelper::getName(geoID) << " (GeoID=" << geoID << "): \t" << localSelectors );
67
68 for ( auto& selector : localSelectors )
69 {
70 const auto flavor = selector->simFlavor();
71 auto itr = m_simToolMap.find(flavor);
72 if (itr == m_simToolMap.end() )
73 {
74 ATH_MSG_WARNING( "Map from SimulationFlavor to SimulatorTool:" );
75 for ( auto& entry : m_simToolMap )
76 {
77 ATH_MSG_WARNING( "SimulationFlavor: " << entry.first << ", SimulatorTool Name: " << entry.second->name() );
78 }
79 ATH_MSG_FATAL( "No SimulationTool with flavor " << flavor << " expected by " << selector->name() << " found in this job!\n Check your configuration!" );
80 return StatusCode::FAILURE;
81 }
82 }
83 }
84
85 ATH_CHECK( m_entryLayerTool.retrieve() );
86
87 if ( not m_orderingTool.empty() ) {
88 ATH_CHECK( m_orderingTool.retrieve() );
89 }
90
91 ATH_CHECK( m_inputEvgenKey.initialize() );
92 ATH_CHECK( m_outputTruthKey.initialize() );
93
94 ATH_CHECK( m_caloEntryLayerKey.initialize() );
95 ATH_CHECK( m_muonEntryLayerKey.initialize() );
96 ATH_CHECK( m_muonExitLayerKey.initialize() );
97
98 ATH_CHECK( m_inputConverter.retrieve() );
99 if ( not m_truthPreselectionTool.empty() ) {
101 }
102
103 ATH_CHECK ( m_truthRecordSvc.retrieve() );
104
105 if ( not m_qspatcher.empty() ) {
106 ATH_CHECK( m_qspatcher.retrieve() );
107 }
108
109 ATH_CHECK( m_geoIDSvc.retrieve() );
110
111 return StatusCode::SUCCESS;
112}
113
114
116
117 const EventContext& ctx = Gaudi::Hive::currentContext();
118 auto hitCollections = std::make_shared<HitCollectionMap>();
119 // Call setupEvent for all simulators (TODO: make the tools do this)
120 for (auto& curSimTool: m_simulationTools) {
121 if ( curSimTool ) {
122 if (auto* g4sim =
123 dynamic_cast<ISF::BaseSimulatorG4Tool*>(curSimTool.get())) {
124 ATH_CHECK(g4sim->setupEvent(ctx, *hitCollections));
125 } else {
126 ATH_CHECK(curSimTool->setupEvent(ctx));
127 }
128 ATH_MSG_DEBUG( "Event setup done for " << curSimTool->name() );
129 }
130 }
131
133 if (!inputEvgen.isValid()) {
134 ATH_MSG_FATAL("Unable to read input GenEvent collection '" << inputEvgen.key() << "'");
135 return StatusCode::FAILURE;
136 }
137
138 // create output Truth collection
140 std::unique_ptr<McEventCollection> shadowTruth{};
141 if (m_useShadowEvent) {
142 outputTruth = std::make_unique<McEventCollection>();
143 // copy input Evgen collection to shadow Truth collection
144 shadowTruth = std::make_unique<McEventCollection>(*inputEvgen);
145 for (HepMC::GenEvent* currentGenEvent : *shadowTruth ) {
146 // Apply QS patch if required
147 if ( not m_qspatcher.empty() ) {
148 ATH_CHECK(m_qspatcher->applyWorkaround(*currentGenEvent));
149 }
150 // Copy GenEvent and remove daughters of quasi-stable particles to be simulated
151 std::unique_ptr<HepMC::GenEvent> outputEvent = m_truthPreselectionTool->filterGenEvent(*currentGenEvent);
152 outputTruth->push_back(outputEvent.release());
153 }
154 }
155 else {
156 // copy input Evgen collection to output Truth collection
157 outputTruth = std::make_unique<McEventCollection>(*inputEvgen);
158 // Apply QS patch if required
159 if ( not m_qspatcher.empty() ) {
160 for (HepMC::GenEvent* currentGenEvent : *outputTruth ) {
161 ATH_CHECK(m_qspatcher->applyWorkaround(*currentGenEvent));
162 }
163 }
164 }
165
166 const int largestGeneratedParticleBC = (outputTruth->empty()) ? HepMC::UNDEFINED_ID
167 : HepMC::maxGeneratedParticleBarcode(outputTruth->at(0)); // TODO make this more robust
168 const int largestGeneratedVertexBC = (outputTruth->empty()) ? HepMC::UNDEFINED_ID
169 : HepMC::maxGeneratedVertexBarcode(outputTruth->at(0)); // TODO make this more robust
170 // tell TruthService we're starting a new event
171 ATH_CHECK( m_truthRecordSvc->initializeTruthCollection(largestGeneratedParticleBC, largestGeneratedVertexBC) );
172
173 // Create TrackRecordCollections and pass them to the entryLayerTool
175 caloEntryLayer = std::make_unique<TrackRecordCollection>(caloEntryLayer.name());
176 ATH_CHECK(m_entryLayerTool->registerTrackRecordCollection(caloEntryLayer.ptr(), fAtlasCaloEntry));
178 muonEntryLayer = std::make_unique<TrackRecordCollection>(muonEntryLayer.name());
179 ATH_CHECK(m_entryLayerTool->registerTrackRecordCollection(muonEntryLayer.ptr(), fAtlasMuonEntry));
181 muonExitLayer = std::make_unique<TrackRecordCollection>(muonExitLayer.name());
182 ATH_CHECK(m_entryLayerTool->registerTrackRecordCollection(muonExitLayer.ptr(), fAtlasMuonExit));
183
184 // read and convert input
185 ISFParticleContainer simParticles;
186 ATH_CHECK( m_inputConverter->convert(*outputTruth, simParticles) );
187
188 // create an ordered queue of particles
189 ISFParticleOrderedQueue particleQueue;
190 for ( auto& particle : simParticles ) {
191
192 if ( m_orderingTool.empty() ) {
193 // Without a defined ordering, preserve old FIFO behaviour
194 particle->setOrder( -particleQueue.size() );
195 } else {
196 m_orderingTool->setOrder( *particle );
197 }
198
199 particleQueue.push( particle );
200 }
201
202 unsigned int loopCounter{0};
203 // loop until there are no more particles to simulate
204 ISF::ISFParticleVector particles{};
205 ISimulatorTool* lastSimulator{};
206 ISFParticleContainer newSecondaries{};
207 while ( particleQueue.size() ) {
208 ++loopCounter;
209 ATH_MSG_VERBOSE("Main Loop pass no. " << loopCounter);
210 ATH_MSG_VERBOSE("Queue starts with " << particleQueue.size() << " particles.");
211 // Create a vector of particles with the same simulator
212 ISFParticleOrderedQueue tempQueue;
213 while ( particleQueue.size() ) {
214 auto particlePtr = particleQueue.top();
215 ISFParticle& curParticle( *particlePtr );
216 particleQueue.pop();
217
218 // Get the geo ID for the particle
219 if ( m_forceGeoIDSvc || !validAtlasRegion( curParticle.nextGeoID() ) ) {
220 m_geoIDSvc->identifyAndRegNextGeoID( curParticle );
221 }
222
223 // Get the simulator using the GeoID
224 auto& simTool = identifySimulator(curParticle);
225
226 // Fill the vector
227 if ( particles.empty() ) {
228 // First particle in the vector defines the simulator
229 particles.push_back(particlePtr);
230 lastSimulator=&simTool;
231 }
232 else if (&simTool!=lastSimulator || particles.size() >= m_maxParticleVectorSize.value() ) {
233 // Change of simulator, end the current vector
234 tempQueue.push(particlePtr);
235 }
236 else {
237 particles.push_back(particlePtr);
238 }
239 }
240 particleQueue = std::move(tempQueue);
241 #ifdef ISFDEBUG
242 if (loopCounter>100 && particles.size()<3) {
243 ATH_MSG_INFO("Main Loop pass no. " << loopCounter);
244 ATH_MSG_INFO("Selected " << particles.size() << " particles to be processed by " << lastSimulator->name());
245 for ( const ISFParticle *particle : particles ) {
246 ATH_MSG_INFO(*particle);
247 }
248 }
249 #endif // ISFDEBUG
250
251 ATH_MSG_VERBOSE("Selected " << particles.size() << " particles to be processed by " << lastSimulator->name());
252 // Run the simulation
253 if (auto* g4sim = dynamic_cast<ISF::BaseSimulatorG4Tool*>(lastSimulator)) {
254 ATH_CHECK(g4sim->simulateVector(ctx, particles, newSecondaries,
255 outputTruth.ptr(), hitCollections,
256 shadowTruth.get()));
257 } else {
258 ATH_CHECK(lastSimulator->simulateVector(ctx, particles, newSecondaries,
259 outputTruth.ptr(),
260 shadowTruth.get()));
261 }
262
263 ATH_MSG_VERBOSE(lastSimulator->name() << " returned " << newSecondaries.size() << " new particles to be added to the queue." );
264 // Register returned particles with the entry layer tool, set their order and enqueue them
265 for ( auto* secondary : newSecondaries ) {
266
267 // Set up particle in ISF
268 m_entryLayerTool->registerParticle( *secondary );
269 if ( m_forceGeoIDSvc || !validAtlasRegion( secondary->nextGeoID() ) ) {
270 m_geoIDSvc->identifyAndRegNextGeoID( *secondary );
271 }
272
273 if ( m_orderingTool.empty() ) {
274 // Without a defined ordering, preserve old FIFO behaviour
275 secondary->setOrder( -particleQueue.size() );
276 } else {
277 m_orderingTool->setOrder( *secondary );
278 }
279
280 particleQueue.push( secondary );
281 }
282 newSecondaries.clear();
283
284 // Delete simulated particles
285 for ( auto usedParticle : particles ) {
286 delete usedParticle;
287 }
288 particles.clear();
289 }
290 ATH_MSG_VERBOSE("Final status: queue contains " << particleQueue.size() << " particles.");
291
292 // Release the event from all simulators (TODO: make the tools do this)
293 for (auto& curSimTool: m_simulationTools) {
294 if ( curSimTool ) {
295 if (auto* g4sim =
296 dynamic_cast<ISF::BaseSimulatorG4Tool*>(curSimTool.get())) {
297 ATH_CHECK(g4sim->releaseEvent(ctx, *hitCollections));
298 } else {
299 ATH_CHECK(curSimTool->releaseEvent(ctx));
300 }
301 ATH_MSG_DEBUG( "releaseEvent() completed for " << curSimTool->name() );
302 }
303 }
304
305 // Remove QS patch if required
306 if(!m_qspatcher.empty()) {
307 for (HepMC::GenEvent* currentGenEvent : *outputTruth ) {
308 ATH_CHECK(m_qspatcher->removeWorkaround(*currentGenEvent));
309 }
310 }
311
312 return StatusCode::SUCCESS;
313}
314
315
317 return StatusCode::SUCCESS;
318}
319
320
323 AtlasDetDescr::AtlasRegion geoID = particle.nextGeoID();
324
325 auto& localSelectors = m_simSelectors[geoID];
326 for (auto& selector: localSelectors) {
327 bool selected = selector->selfSelect(particle);
328 if (selected) {
329 return *m_simToolMap.at(selector->simFlavor());
330 }
331 }
332
333 ATH_MSG_WARNING("No simulator found for particle (" << particle << ")."
334 << " Will send it to " << m_particleKillerTool->name());
335 return *m_particleKillerTool;
336}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
bool validAtlasRegion(AtlasDetDescr::AtlasRegion region)
Check a given AtlasRegion for its validity.
Definition AtlasRegion.h:40
static const char * getName(int region)
The generic ISF particle definition,.
Definition ISFParticle.h:42
AtlasDetDescr::AtlasRegion nextGeoID() const
next geoID the particle will be simulated in
virtual StatusCode simulateVector(const EventContext &ctx, const ISFParticleVector &particles, ISFParticleContainer &secondaries, McEventCollection *mcEventCollection, McEventCollection *shadowTruth=nullptr)=0
Simulation call for vectors of particles.
SG::WriteHandleKey< McEventCollection > m_outputTruthKey
Output Simulation Truth collection.
Definition SimKernelMT.h:91
ToolHandleArray< ISimulatorTool > m_simulationTools
Simulation Tools.
SG::ReadHandleKey< McEventCollection > m_inputEvgenKey
Input Generator Truth collection.
Definition SimKernelMT.h:89
BooleanProperty m_useShadowEvent
Definition SimKernelMT.h:97
SG::WriteHandleKey< TrackRecordCollection > m_muonEntryLayerKey
Definition SimKernelMT.h:94
std::map< ISF::SimulationFlavor, ISimulatorTool * > m_simToolMap
Map of the simulation flavours used in this job to the corresponding Simulation Services.
ToolHandle< IGenEventFilter > m_truthPreselectionTool
ToolHandle< IParticleOrderingTool > m_orderingTool
Tool to set particle ordering.
Gaudi::Property< size_t > m_maxParticleVectorSize
Number of particles simultaneously sent to simulator.
ISimulatorTool & identifySimulator(const ISF::ISFParticle &particle)
Returns the simulator to use for the given particle.
ServiceHandle< IGeoIDSvc > m_geoIDSvc
Service to set particle GeoIDs.
SG::WriteHandleKey< TrackRecordCollection > m_caloEntryLayerKey
Output TrackRecordCollections.
Definition SimKernelMT.h:93
ServiceHandle< Simulation::IZeroLifetimePatcher > m_qspatcher
Quasi-Stable Particle Simulation Patcher.
StatusCode execute() override
Check the validity of the MC truth event in the current Athena event.
ServiceHandle< IInputConverter > m_inputConverter
Input converter service (from Generator->ISF particle types)
ISimulatorTool * m_particleKillerTool
When no appropriate simulator can be found for a given particle, the particle is sent to this "partic...
std::array< PublicToolHandleArray< ISimulationSelector >, AtlasDetDescr::fNumAtlasRegions > m_simSelectors
The simulation selectors defining the "routing chain".
SG::WriteHandleKey< TrackRecordCollection > m_muonExitLayerKey
Definition SimKernelMT.h:95
Gaudi::Property< bool > m_forceGeoIDSvc
Force geoID recalculation for each particle.
Definition SimKernelMT.h:99
PublicToolHandle< IEntryLayerTool > m_entryLayerTool
AthenaTool responsible for writing Calo/Muon Entry/Exit Layer collection.
StatusCode initialize() override
Check user configuration and configure the algorithm state accordingly.
ServiceHandle< ITruthSvc > m_truthRecordSvc
Central truth service.
StatusCode finalize() override
Implements empty finalize (implementation required by AthAlgorithm)
virtual bool isValid() override final
Can the handle be successfully dereferenced?
virtual const std::string & key() const override final
Return the StoreGate ID for the referenced object.
const std::string & name() const
Return the StoreGate ID for the referenced object.
pointer_type ptr()
Dereference the pointer.
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
int maxGeneratedVertexBarcode(const HepMC::GenEvent *genEvent)
Get the maximal absolute value of barcode of vertex present in the event. Returns a negative number.
constexpr int UNDEFINED_ID
int maxGeneratedParticleBarcode(const HepMC::GenEvent *genEvent)
Get the maximal value of barcode of particle present in the event.
std::priority_queue< ISF::ISFParticle *, ISF::ISFParticleVector, ISF::ISFParticleOrdering > ISFParticleOrderedQueue
the actual particle priority_queue
@ fAtlasMuonExit
Definition EntryLayer.h:39
@ fAtlasCaloEntry
Definition EntryLayer.h:37
@ fAtlasMuonEntry
Definition EntryLayer.h:38
@ ParticleKiller
std::list< ISF::ISFParticle * > ISFParticleContainer
generic ISFParticle container (not necessarily a std::list!)
std::vector< ISF::ISFParticle * > ISFParticleVector
ISFParticle vector.