ATLAS Offline Software
Loading...
Searching...
No Matches
ActsFatrasSimTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4#include <algorithm>
5#include <random>
6
7#include "ActsFatrasSimTool.h"
8#include "Acts/ActsVersion.hpp"
9#include <Acts/Utilities/StringHelpers.hpp>
10
11#include "CLHEP/Random/RandFlat.h"
12#include "CLHEP/Random/RandomEngine.h"
13
16
17using namespace Acts::UnitLiterals;
18
20 const std::string& name,
21 const IInterface* parent)
22 : BaseSimulatorTool(type, name, parent) {}
23
25
28 ATH_MSG_INFO("ISF::ActsFatrasSimTool update with ACTS version: v"
29 << Acts::VersionMajor << "." << Acts::VersionMinor << "."
30 << Acts::VersionPatch << " [" << Acts::CommitHash.value_or("unknown hash") << "]");
31 // Retrieve particle filter
32 if (!m_particleFilter.empty()) ATH_CHECK(m_particleFilter.retrieve());
33
34 // setup logger
35 m_logger = makeActsAthenaLogger(this, std::string("ActsFatras"),std::string("ActsFatrasSimTool"));
36
37 // retrieve tracking geo tool
39 m_trackingGeometry = m_trackingGeometryTool->trackingGeometry();
40
41 //retrieve Magnetfield tool
42 ATH_MSG_VERBOSE("Using ATLAS magnetic field service");
44
45 // Random number service
46 if (m_rngSvc.retrieve().isFailure()) {
47 ATH_MSG_FATAL("Could not retrieve " << m_rngSvc);
48 return StatusCode::FAILURE;
49 }
50 // Get own engine with own seeds
51 m_randomEngine = m_rngSvc->getEngine(this, m_randomEngineName.value());
52 if (!m_randomEngine) {
53 ATH_MSG_FATAL("Could not get random engine '" << m_randomEngineName.value() << "'");
54 return StatusCode::FAILURE;
55 }
56
57 // ISF truth service
58 ATH_CHECK (m_truthRecordSvc.retrieve());
59 ATH_MSG_DEBUG( "- Using ISF TruthRecordSvc : " << m_truthRecordSvc.typeAndName() );
60 return StatusCode::SUCCESS;
61}
62
63StatusCode ISF::ActsFatrasSimTool::simulate(const EventContext& ctx,
64 ISFParticle& isp, ISFParticleContainer& secondaries,
65 McEventCollection* mcEventCollection) {
66 ATH_MSG_VERBOSE("Particle " << isp << " received for simulation.");
67 // Check if particle passes filter, if there is one
68 if (!m_particleFilter.empty() && !m_particleFilter->passFilter(isp)) {
69 ATH_MSG_VERBOSE("ISFParticle " << isp << " does not pass selection. Ignoring.");
70 return StatusCode::SUCCESS;
71 }
72 // Process ParticleState from particle stack
73 // Wrap the input ISFParticle in an STL vector with size of 1
74 const ISF::ISFParticleVector ispVector(1, &isp);
75 ATH_CHECK(this->simulateVector(ctx, ispVector, secondaries, mcEventCollection));
76 ATH_MSG_VERBOSE("Simulation done");
77 return StatusCode::SUCCESS;
78}
79
81 const EventContext& ctx,
82 const ISFParticleVector& particles,
83 ISFParticleContainer& secondaries,
84 McEventCollection* /*mcEventCollection*/, McEventCollection *) {
85
86 m_randomEngine->setSeed(m_randomEngineName, ctx);
87 CLHEP::HepRandomEngine* randomEngine = m_randomEngine->getEngine(ctx);
88 Generator generator(CLHEP::RandFlat::shoot(randomEngine->flat()));
89 ATH_MSG_VERBOSE(name() << " RNG seed " << CLHEP::RandFlat::shoot(randomEngine->flat()));
90 ATH_MSG_VERBOSE(name() << " received vector of size "
91 << particles.size() << " particles for simulation.");
92
93 // construct the ACTS simulator
94 Acts::Navigator navigator( Acts::Navigator::Config{ m_trackingGeometry }, m_logger);
95 auto bField = std::make_shared<ATLASMagneticFieldWrapper>();
96 auto chargedStepper = ChargedStepper(std::move(bField));
97 auto neutralStepper = NeutralStepper();
98 auto chargedPropagator = ChargedPropagator(chargedStepper, navigator, m_logger);
99 auto neutralPropagator = NeutralPropagator(neutralStepper, navigator, m_logger);
100 ChargedSimulation simulatorCharged(std::move(chargedPropagator), m_logger);
101 NeutralSimulation simulatorNeutral(std::move(neutralPropagator), m_logger);
102 Simulation simulator=Simulation(std::move(simulatorCharged),std::move(simulatorNeutral));
103 ATH_MSG_VERBOSE(name() << " Min pT for interaction " << m_interact_minPt * Acts::UnitConstants::MeV << " GeV");
104 // Acts propagater options
105 simulator.charged.maxStepSize = m_maxStepSize;
106 simulator.charged.maxStep = m_maxStep;
107 simulator.charged.pathLimit = m_pathLimit;
108 simulator.charged.maxRungeKuttaStepTrials = m_maxRungeKuttaStepTrials;
109 simulator.charged.loopProtection = m_loopProtection;
110 simulator.charged.loopFraction = m_loopFraction;
111 simulator.charged.targetTolerance = m_tolerance;
112 simulator.charged.stepSizeCutOff = m_stepSizeCutOff;
113 // Create interaction list
114 simulator.charged.interactions = ActsFatras::makeStandardChargedElectroMagneticInteractions(m_interact_minPt * Acts::UnitConstants::MeV);
115 // get Geo and Mag map
116 ATH_MSG_VERBOSE(name() << " Getting per event Geo and Mag map");
117 Acts::MagneticFieldContext mctx = getMagneticFieldContext(ctx);
118 const ActsTrk::GeometryContext& gctx = m_trackingGeometryTool->getNominalGeometryContext();
119 auto anygctx = gctx.context();
120 // Loop over ISFParticleVector and process each separately
121 ATH_MSG_VERBOSE(name() << " Processing particles in ISFParticleVector.");
122 for (const auto isfp : particles) {
123 // ====ACTSFatras Simulation====
124 // //
125 // input/output particle and hits containers
126 // Convert to ActsFatras::Particle
127 // ISF: Energy, mass, and momentum are in MeV, position in mm
128 // Acts: Energy, mass, and momentum are in GeV, position in mm
129 ATH_MSG_DEBUG(name() << " Convert ISF::Particle(mass) " << isfp->id()<<"|" << isfp<<"(" << isfp->mass() << ")");
130 std::vector<ActsFatras::Particle> input = std::vector<ActsFatras::Particle>{
131 ActsFatras::Particle(ActsFatras::Barcode().withVertexPrimary(0).withParticle(isfp->id()), static_cast<Acts::PdgParticle>(isfp->pdgCode()),
132 isfp->charge(),isfp->mass() * Acts::UnitConstants::MeV)
133 .setDirection(Acts::makeDirectionFromPhiEta(isfp->momentum().phi(), isfp->momentum().eta()))
134 .setAbsoluteMomentum(isfp->momentum().mag() * Acts::UnitConstants::MeV)
135 .setPosition4(ActsTrk::convertPosToActs(isfp->position(), isfp->timeStamp()))};
136 ATH_MSG_DEBUG(name() << " Propagating ActsFatras::Particle vertex|particle|generation|subparticle, " << input[0]);
137 std::vector<ActsFatras::Particle> simulatedInitial;
138 std::vector<ActsFatras::Particle> simulatedFinal;
139 std::vector<ActsFatras::Hit> hits;
140 // simulate
141 auto result=simulator.simulate(anygctx, mctx, generator, input, simulatedInitial, simulatedFinal, hits);
142 auto simulatedFailure=result.value();
143 if (simulatedFailure.size()>0){
144 for (const auto& simfail : simulatedFailure){
145 auto errCode = Acts::make_error_code(Acts::PropagatorError(simfail.error.value()));
146 ATH_MSG_WARNING(name() << " Particle id " <<simfail.particle.particleId()<< ": fail to be simulated during Propagation: " << errCode.message());
147 ATH_MSG_WARNING(name() << " Particle vertex|particle|generation|subparticle"<<simfail.particle << " starts from position" << Acts::toString(simfail.particle.position()) << " and direction " << Acts::toString(simfail.particle.direction()));
148 return StatusCode::SUCCESS;
149 }
150 }
151
152 ATH_MSG_DEBUG(name() << " initial particle " << simulatedInitial[0]);
153 ATH_MSG_DEBUG(name() << " ActsFatras simulator hits: " << hits.size());
154 int i = 0;
155 for (const auto& hit : hits) {
156 ATH_MSG_DEBUG(name() << " hit pos: " << hit.position() );
157 ++i;
158 if (i>5) break;
159 }
160 ATH_MSG_DEBUG(name() << " No. of particles after ActsFatras simulator: " << simulatedFinal.size());
161 if (!simulatedFinal.empty()){
162 ATH_MSG_DEBUG(name() << " start procesing secondaries");
163 auto itr = simulatedFinal.begin();
164 // Save hits of isfp
165 std::vector<ActsFatras::Hit> particle_hits;
166 if (itr->numberOfHits() > 0) {
167 std::copy(hits.begin(), hits.begin()+itr->numberOfHits(), std::back_inserter(particle_hits));
169 }
170 // Process secondaries
171 auto isKilled = !itr->isAlive();
172 int maxGeneration = simulatedFinal.back().particleId().generation();
173 ATH_MSG_DEBUG(name() << " maxGeneration: "<< maxGeneration);
174 for (int gen = 0; gen <= maxGeneration; ++gen){
175 ATH_MSG_DEBUG(name() << " start with generation "<< gen << "|" << maxGeneration << ": "<< *itr);
176 auto vecsecisfp = std::make_unique<ISF::ISFParticleVector>();
177 while (itr != simulatedFinal.end() && static_cast<int>(itr->particleId().generation()) == gen) {
178 ATH_MSG_DEBUG(name() << " genration "<< gen << "|" << maxGeneration << ": "<< *itr);
179 if(itr->isSecondary()){
180 // convert final particles to ISF::particle
181 const auto pos = ActsTrk::convertPosFromActs(itr->fourPosition()).first;
182 const auto mom = ActsTrk::convertMomFromActs(itr->fourMomentum()).first;
183 double mass = itr->mass() / Acts::UnitConstants::MeV;
184 double charge = itr->charge();
185 int pdgid = itr->pdg();
186 auto properTime = ActsTrk::timeToAthena(itr->time());
187 const int status = 1 + HepMC::SIM_STATUS_THRESHOLD;
188 const int id = HepMC::UNDEFINED_ID;
189 auto secisfp = new ISF::ISFParticle (pos,mom,mass,charge,pdgid,status,properTime,*isfp,id);
190 ATH_MSG_DEBUG(name() <<" secondaries particle (ACTS): "<<*itr<< "("<<itr->momentum()<<")|time "<<itr->time()<<"|process "<< getATLASProcessCode(itr->process()));
191 ATH_MSG_DEBUG(name() <<" secondaries particle (ISF): " << *secisfp << " time "<<secisfp->timeStamp());
192 vecsecisfp->push_back(secisfp);
193 }
194 else{
195 ATH_MSG_DEBUG(name() <<" primary particle found with generation "<< gen);
196 }
197 ++itr;
198 }
199 if (!vecsecisfp->empty()) {
200 ISF::ISFTruthIncident truth(*isfp,
201 *vecsecisfp,
202 getATLASProcessCode((itr-1)->process()),
203 isfp->nextGeoID(),
204 isKilled&&gen==maxGeneration ? ISF::fKillsPrimary : ISF::fPrimarySurvives
205 );
206 ATH_MSG_DEBUG(name() << " Truth incident parentPt2(MinPt2) " << truth.parentPt2() <<" (100 MeV)");
207 ATH_MSG_DEBUG(name() << " Truth incident ChildPt2(MinPt2) " << truth.childrenPt2Pass(300) <<" (300 MeV)");
208 m_truthRecordSvc->registerTruthIncident(truth, true);
211 for (auto *secisfp : *vecsecisfp){
212 if (secisfp->getTruthBinding()) {
213 secondaries.push_back(secisfp);
214 }
215 else {
216 ATH_MSG_WARNING("Secondary particle not written out to truth.\n Parent (" << isfp << ")\n Secondary (" << *secisfp <<")");
217 }
218 } // end of truth binding
219 }// end of store truth bind secondaries
220 }
221 }// end of secondaries
222 ATH_MSG_VERBOSE(name() << " No. of secondaries: " << secondaries.size());
223 ATH_MSG_DEBUG(name() << " End of particle " << isfp->id());
224
225 std::vector<ActsFatras::Particle>().swap(input);
226 std::vector<ActsFatras::Particle>().swap(simulatedInitial);
227 std::vector<ActsFatras::Particle>().swap(simulatedFinal);
228 std::vector<ActsFatras::Hit>().swap(hits);
229 } // end of isfp loop
230 return StatusCode::SUCCESS;
231}
232
233Acts::MagneticFieldContext ISF::ActsFatrasSimTool::getMagneticFieldContext(const EventContext& ctx) const {
235 if (!readHandle.isValid()) {
236 ATH_MSG_ERROR(name() + ": Failed to retrieve magnetic field condition data " + m_fieldCacheCondObjInputKey.key() + ".");
237 }
238 else ATH_MSG_DEBUG(name() << "retrieved magnetic field condition data "<< m_fieldCacheCondObjInputKey.key());
239 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
240
241 return Acts::MagneticFieldContext(fieldCondObj);
242}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#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)
double charge(const T &p)
Definition AtlasPID.h:997
std::unique_ptr< const Acts::Logger > makeActsAthenaLogger(IMessageSvc *svc, const std::string &name, int level, std::optional< std::string > parent_name)
Acts::GeometryContext context() const
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
Gaudi::Property< double > m_maxStepSize
Acts::EigenStepper< Acts::EigenStepperDefaultExtension > ChargedStepper
PublicToolHandle< ActsTrk::ITrackingGeometryTool > m_trackingGeometryTool
Gaudi::Property< double > m_tolerance
Gaudi::Property< double > m_interact_minPt
ServiceHandle< IAthRNGSvc > m_rngSvc
SingleParticleSimulation< NeutralPropagator, NeutralInteractions, ActsFatras::NoSurface, ActsFatras::NoDecay > NeutralSimulation
Acts::Propagator< NeutralStepper, Navigator > NeutralPropagator
virtual StatusCode simulate(const EventContext &ctx, ISFParticle &isp, ISFParticleContainer &, McEventCollection *) override
std::shared_ptr< const Acts::TrackingGeometry > m_trackingGeometry
std::shared_ptr< const Acts::Logger > m_logger
Gaudi::Property< double > m_maxStep
int getATLASProcessCode(ActsFatras::ProcessType actspt)
Gaudi::Property< double > m_pathLimit
Gaudi::Property< double > m_stepSizeCutOff
Gaudi::Property< bool > m_loopProtection
ActsFatrasSimTool(const std::string &type, const std::string &name, const IInterface *parent)
ServiceHandle< ISF::ITruthSvc > m_truthRecordSvc
SiHitCollection m_pixelSiHits
SingleParticleSimulation< ChargedPropagator, ChargedInteractions, HitSurfaceSelector, ActsFatras::NoDecay > ChargedSimulation
virtual StatusCode initialize() override
PublicToolHandle< ISF::IParticleFilter > m_particleFilter
ToolHandle< ActsFatrasWriteHandler > m_ActsFatrasWriteHandler
virtual StatusCode simulateVector(const EventContext &ctx, const ISFParticleVector &particles, ISFParticleContainer &secondaries, McEventCollection *mcEventCollection, McEventCollection *shadowTruth=nullptr) override
Simulation call for vectors of particles.
Gaudi::Property< double > m_maxRungeKuttaStepTrials
virtual Acts::MagneticFieldContext getMagneticFieldContext(const EventContext &) const
Gaudi::Property< double > m_loopFraction
Acts::Propagator< ChargedStepper, Navigator > ChargedPropagator
Gaudi::Property< std::string > m_randomEngineName
Acts::StraightLineStepper NeutralStepper
virtual StatusCode initialize() override
BaseSimulatorTool(const std::string &type, const std::string &name, const IInterface *parent)
The generic ISF particle definition,.
Definition ISFParticle.h:42
Interface class for all truth incidents handled by the ISF.
void updateChildParticleProperties()
Update the id and particleLink properties of the child particles (to be called after registerTruthInc...
double parentPt2() const override final
Return pT^2 of the parent particle.
void updateParentAfterIncidentProperties()
Update the id and particleLink properties of the parentAfterIncident (to be called after registerTrut...
bool childrenPt2Pass(double pt2cut)
Return true if at least one child particle passes the given pT^2 cut (= at least one child with pT^2 ...
This defines the McEventCollection, which is really just an ObjectVector of McEvent objectsFile: Gene...
const std::string process
constexpr double timeToAthena(const double actsT)
Converts a time unit from Acts to Athena units.
std::pair< Amg::Vector3D, double > convertMomFromActs(const Acts::Vector4 &actsMom)
Converts an Acts four-momentum vector into an pair of an Athena three-momentum and the paritcle's ene...
std::pair< Amg::Vector3D, double > convertPosFromActs(const Acts::Vector4 &actsPos)
Converts an Acts 4-vector into a pair of an Athena spatial vector and the passed time.
Acts::Vector4 convertPosToActs(const Amg::Vector3D &athenaPos, const double athenaTime=0.)
Converts a position vector & time from Athena units into Acts units.
constexpr int UNDEFINED_ID
constexpr int SIM_STATUS_THRESHOLD
Constant definiting the status threshold for simulated particles, eg. can be used to separate generat...
@ fKillsPrimary
@ fPrimarySurvives
std::list< ISF::ISFParticle * > ISFParticleContainer
generic ISFParticle container (not necessarily a std::list!)
std::vector< ISF::ISFParticle * > ISFParticleVector
ISFParticle vector.