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
10#include "CLHEP/Random/RandFlat.h"
11#include "CLHEP/Random/RandomEngine.h"
12
15
16using namespace Acts::UnitLiterals;
17
19 const std::string& name,
20 const IInterface* parent)
21 : BaseSimulatorTool(type, name, parent) {}
22
24
27 ATH_MSG_INFO("ISF::ActsFatrasSimTool update with ACTS version: v"
28 << Acts::VersionMajor << "." << Acts::VersionMinor << "."
29 << Acts::VersionPatch << " [" << Acts::CommitHash.value_or("unknown hash") << "]");
30 // Retrieve particle filter
31 if (!m_particleFilter.empty()) ATH_CHECK(m_particleFilter.retrieve());
32
33 // setup logger
34 m_logger = makeActsAthenaLogger(this, std::string("ActsFatras"),std::string("ActsFatrasSimTool"));
35
36 // retrieve tracking geo tool
38 m_trackingGeometry = m_trackingGeometryTool->trackingGeometry();
39
40 //retrieve Magnetfield tool
41 ATH_MSG_VERBOSE("Using ATLAS magnetic field service");
43
44 // Random number service
45 if (m_rngSvc.retrieve().isFailure()) {
46 ATH_MSG_FATAL("Could not retrieve " << m_rngSvc);
47 return StatusCode::FAILURE;
48 }
49 // Get own engine with own seeds
50 m_randomEngine = m_rngSvc->getEngine(this, m_randomEngineName.value());
51 if (!m_randomEngine) {
52 ATH_MSG_FATAL("Could not get random engine '" << m_randomEngineName.value() << "'");
53 return StatusCode::FAILURE;
54 }
55
56 // ISF truth service
57 ATH_CHECK (m_truthRecordSvc.retrieve());
58 ATH_MSG_DEBUG( "- Using ISF TruthRecordSvc : " << m_truthRecordSvc.typeAndName() );
59 return StatusCode::SUCCESS;
60}
61
62StatusCode ISF::ActsFatrasSimTool::simulate(const EventContext& ctx,
63 ISFParticle& isp, ISFParticleContainer& secondaries,
64 McEventCollection* mcEventCollection) {
65 ATH_MSG_VERBOSE("Particle " << isp << " received for simulation.");
66 // Check if particle passes filter, if there is one
67 if (!m_particleFilter.empty() && !m_particleFilter->passFilter(isp)) {
68 ATH_MSG_VERBOSE("ISFParticle " << isp << " does not pass selection. Ignoring.");
69 return StatusCode::SUCCESS;
70 }
71 // Process ParticleState from particle stack
72 // Wrap the input ISFParticle in an STL vector with size of 1
73 const ISF::ISFParticleVector ispVector(1, &isp);
74 ATH_CHECK(this->simulateVector(ctx, ispVector, secondaries, mcEventCollection));
75 ATH_MSG_VERBOSE("Simulation done");
76 return StatusCode::SUCCESS;
77}
78
80 const EventContext& ctx,
81 const ISFParticleVector& particles,
82 ISFParticleContainer& secondaries,
83 McEventCollection* /*mcEventCollection*/, McEventCollection *) {
84
85 m_randomEngine->setSeed(m_randomEngineName, ctx);
86 CLHEP::HepRandomEngine* randomEngine = m_randomEngine->getEngine(ctx);
87 Generator generator(CLHEP::RandFlat::shoot(randomEngine->flat()));
88 ATH_MSG_VERBOSE(name() << " RNG seed " << CLHEP::RandFlat::shoot(randomEngine->flat()));
89 ATH_MSG_VERBOSE(name() << " received vector of size "
90 << particles.size() << " particles for simulation.");
91
92 // construct the ACTS simulator
93 Acts::Navigator navigator( Acts::Navigator::Config{ m_trackingGeometry }, m_logger);
94 auto bField = std::make_shared<ATLASMagneticFieldWrapper>();
95 auto chargedStepper = ChargedStepper(std::move(bField));
96 auto neutralStepper = NeutralStepper();
97 auto chargedPropagator = ChargedPropagator(chargedStepper, navigator, m_logger);
98 auto neutralPropagator = NeutralPropagator(neutralStepper, navigator, m_logger);
99 ChargedSimulation simulatorCharged(std::move(chargedPropagator), m_logger);
100 NeutralSimulation simulatorNeutral(std::move(neutralPropagator), m_logger);
101 Simulation simulator=Simulation(std::move(simulatorCharged),std::move(simulatorNeutral));
102 ATH_MSG_VERBOSE(name() << " Min pT for interaction " << m_interact_minPt * Acts::UnitConstants::MeV << " GeV");
103 // Acts propagater options
104 simulator.charged.maxStepSize = m_maxStepSize;
105 simulator.charged.maxStep = m_maxStep;
106 simulator.charged.pathLimit = m_pathLimit;
107 simulator.charged.maxRungeKuttaStepTrials = m_maxRungeKuttaStepTrials;
108 simulator.charged.loopProtection = m_loopProtection;
109 simulator.charged.loopFraction = m_loopFraction;
110 simulator.charged.targetTolerance = m_tolerance;
111 simulator.charged.stepSizeCutOff = m_stepSizeCutOff;
112 // Create interaction list
113 simulator.charged.interactions = ActsFatras::makeStandardChargedElectroMagneticInteractions(m_interact_minPt * Acts::UnitConstants::MeV);
114 // get Geo and Mag map
115 ATH_MSG_VERBOSE(name() << " Getting per event Geo and Mag map");
116 Acts::MagneticFieldContext mctx = getMagneticFieldContext(ctx);
117 const ActsTrk::GeometryContext& gctx = m_trackingGeometryTool->getNominalGeometryContext();
118 auto anygctx = gctx.context();
119 // Loop over ISFParticleVector and process each separately
120 ATH_MSG_VERBOSE(name() << " Processing particles in ISFParticleVector.");
121 for (const auto isfp : particles) {
122 // ====ACTSFatras Simulation====
123 // //
124 // input/output particle and hits containers
125 // Convert to ActsFatras::Particle
126 // ISF: Energy, mass, and momentum are in MeV, position in mm
127 // Acts: Energy, mass, and momentum are in GeV, position in mm
128 ATH_MSG_DEBUG(name() << " Convert ISF::Particle(mass) " << isfp->id()<<"|" << isfp<<"(" << isfp->mass() << ")");
129 std::vector<ActsFatras::Particle> input = std::vector<ActsFatras::Particle>{
130 ActsFatras::Particle(ActsFatras::Barcode().withVertexPrimary(0).withParticle(isfp->id()), static_cast<Acts::PdgParticle>(isfp->pdgCode()),
131 isfp->charge(),isfp->mass() * Acts::UnitConstants::MeV)
132 .setDirection(Acts::makeDirectionFromPhiEta(isfp->momentum().phi(), isfp->momentum().eta()))
133 .setAbsoluteMomentum(isfp->momentum().mag() * Acts::UnitConstants::MeV)
134 .setPosition4(ActsTrk::convertPosToActs(isfp->position(), isfp->timeStamp()))};
135 ATH_MSG_DEBUG(name() << " Propagating ActsFatras::Particle vertex|particle|generation|subparticle, " << input[0]);
136 std::vector<ActsFatras::Particle> simulatedInitial;
137 std::vector<ActsFatras::Particle> simulatedFinal;
138 std::vector<ActsFatras::Hit> hits;
139 // simulate
140 auto result=simulator.simulate(anygctx, mctx, generator, input, simulatedInitial, simulatedFinal, hits);
141 auto simulatedFailure=result.value();
142 if (simulatedFailure.size()>0){
143 for (const auto& simfail : simulatedFailure){
144 auto errCode = Acts::make_error_code(Acts::PropagatorError(simfail.error.value()));
145 ATH_MSG_WARNING(name() << " Particle id " <<simfail.particle.particleId()<< ": fail to be simulated during Propagation: " << errCode.message());
146 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()));
147 return StatusCode::SUCCESS;
148 }
149 }
150
151 ATH_MSG_DEBUG(name() << " initial particle " << simulatedInitial[0]);
152 ATH_MSG_DEBUG(name() << " ActsFatras simulator hits: " << hits.size());
153 int i = 0;
154 for (const auto& hit : hits) {
155 ATH_MSG_DEBUG(name() << " hit pos: " << hit.position() );
156 ++i;
157 if (i>5) break;
158 }
159 ATH_MSG_DEBUG(name() << " No. of particles after ActsFatras simulator: " << simulatedFinal.size());
160 if (!simulatedFinal.empty()){
161 ATH_MSG_DEBUG(name() << " start procesing secondaries");
162 auto itr = simulatedFinal.begin();
163 // Save hits of isfp
164 std::vector<ActsFatras::Hit> particle_hits;
165 if (itr->numberOfHits() > 0) {
166 std::copy(hits.begin(), hits.begin()+itr->numberOfHits(), std::back_inserter(particle_hits));
168 }
169 // Process secondaries
170 auto isKilled = !itr->isAlive();
171 int maxGeneration = simulatedFinal.back().particleId().generation();
172 ATH_MSG_DEBUG(name() << " maxGeneration: "<< maxGeneration);
173 for (int gen = 0; gen <= maxGeneration; ++gen){
174 ATH_MSG_DEBUG(name() << " start with generation "<< gen << "|" << maxGeneration << ": "<< *itr);
175 auto vecsecisfp = std::make_unique<ISF::ISFParticleVector>();
176 while (itr != simulatedFinal.end() && static_cast<int>(itr->particleId().generation()) == gen) {
177 ATH_MSG_DEBUG(name() << " genration "<< gen << "|" << maxGeneration << ": "<< *itr);
178 if(itr->isSecondary()){
179 // convert final particles to ISF::particle
180 const auto pos = ActsTrk::convertPosFromActs(itr->fourPosition()).first;
181 const auto mom = ActsTrk::convertMomFromActs(itr->fourMomentum()).first;
182 double mass = itr->mass() / Acts::UnitConstants::MeV;
183 double charge = itr->charge();
184 int pdgid = itr->pdg();
185 auto properTime = ActsTrk::timeToAthena(itr->time());
186 const int status = 1 + HepMC::SIM_STATUS_THRESHOLD;
187 const int id = HepMC::UNDEFINED_ID;
188 auto secisfp = new ISF::ISFParticle (pos,mom,mass,charge,pdgid,status,properTime,*isfp,id);
189 ATH_MSG_DEBUG(name() <<" secondaries particle (ACTS): "<<*itr<< "("<<itr->momentum()<<")|time "<<itr->time()<<"|process "<< getATLASProcessCode(itr->process()));
190 ATH_MSG_DEBUG(name() <<" secondaries particle (ISF): " << *secisfp << " time "<<secisfp->timeStamp());
191 vecsecisfp->push_back(secisfp);
192 }
193 else{
194 ATH_MSG_DEBUG(name() <<" primary particle found with generation "<< gen);
195 }
196 ++itr;
197 }
198 if (!vecsecisfp->empty()) {
199 ISF::ISFTruthIncident truth(*isfp,
200 *vecsecisfp,
201 getATLASProcessCode((itr-1)->process()),
202 isfp->nextGeoID(),
203 isKilled&&gen==maxGeneration ? ISF::fKillsPrimary : ISF::fPrimarySurvives
204 );
205 ATH_MSG_DEBUG(name() << " Truth incident parentPt2(MinPt2) " << truth.parentPt2() <<" (100 MeV)");
206 ATH_MSG_DEBUG(name() << " Truth incident ChildPt2(MinPt2) " << truth.childrenPt2Pass(300) <<" (300 MeV)");
207 m_truthRecordSvc->registerTruthIncident(truth, true);
210 for (auto *secisfp : *vecsecisfp){
211 if (secisfp->getTruthBinding()) {
212 secondaries.push_back(secisfp);
213 }
214 else {
215 ATH_MSG_WARNING("Secondary particle not written out to truth.\n Parent (" << isfp << ")\n Secondary (" << *secisfp <<")");
216 }
217 } // end of truth binding
218 }// end of store truth bind secondaries
219 }
220 }// end of secondaries
221 ATH_MSG_VERBOSE(name() << " No. of secondaries: " << secondaries.size());
222 ATH_MSG_DEBUG(name() << " End of particle " << isfp->id());
223
224 std::vector<ActsFatras::Particle>().swap(input);
225 std::vector<ActsFatras::Particle>().swap(simulatedInitial);
226 std::vector<ActsFatras::Particle>().swap(simulatedFinal);
227 std::vector<ActsFatras::Hit>().swap(hits);
228 } // end of isfp loop
229 return StatusCode::SUCCESS;
230}
231
232Acts::MagneticFieldContext ISF::ActsFatrasSimTool::getMagneticFieldContext(const EventContext& ctx) const {
234 if (!readHandle.isValid()) {
235 ATH_MSG_ERROR(name() + ": Failed to retrieve magnetic field condition data " + m_fieldCacheCondObjInputKey.key() + ".");
236 }
237 else ATH_MSG_DEBUG(name() << "retrieved magnetic field condition data "<< m_fieldCacheCondObjInputKey.key());
238 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
239
240 return Acts::MagneticFieldContext(fieldCondObj);
241}
#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.