Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ActsFatras::NuclearInteraction Struct Reference

#include <acts/blob/sPHENIX/Fatras/include/ActsFatras/Physics/NuclearInteraction/NuclearInteraction.hpp>

Public Types

using Scalar = Particle::Scalar
 

Public Member Functions

template<typename generator_t >
std::pair< Scalar, ScalargeneratePathLimits (generator_t &generator, const Particle &particle) const
 
template<typename generator_t >
bool run (generator_t &generator, Particle &particle, std::vector< Particle > &generated) const
 

Public Attributes

detail::MultiParticleNuclearInteractionParametrisation multiParticleParameterisation
 The storage of the parameterisation.
 
unsigned int nMatchingTrials = 100
 The number of trials to match momenta and inveriant masses.
 
unsigned int nMatchingTrialsTotal = 1000
 

Private Member Functions

const
detail::NuclearInteractionParameters
findParameters (double rnd, const detail::NuclearInteractionParametrisation &parametrisation, float particleMomentum) const
 
bool softInteraction (double rnd, float probability) const
 
unsigned int finalStateMultiplicity (double rnd, const detail::NuclearInteractionParameters::CumulativeDistribution &distribution) const
 
template<typename generator_t >
std::vector< int > samplePdgIds (generator_t &generator, const detail::NuclearInteractionParameters::PdgMap &pdgMap, unsigned int multiplicity, int particlePdg, bool soft) const
 
template<typename generator_t >
Acts::ActsDynamicVector sampleInvariantMasses (generator_t &generator, const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity &parametrisation) const
 
template<typename generator_t >
Acts::ActsDynamicVector sampleMomenta (generator_t &generator, const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity &parametrisation, float initialMomentum) const
 
bool match (const Acts::ActsDynamicVector &momenta, const Acts::ActsDynamicVector &invariantMasses, float parametrizedMomentum) const
 
template<typename generator_t >
std::optional< std::pair
< Acts::ActsDynamicVector,
Acts::ActsDynamicVector > > 
sampleKinematics (generator_t &generator, const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity &parameters, float momentum) const
 
std::pair
< ActsFatras::Particle::Scalar,
ActsFatras::Particle::Scalar
globalAngle (ActsFatras::Particle::Scalar phi1, ActsFatras::Particle::Scalar theta1, float phi2, float theta2) const
 
template<typename generator_t >
std::vector< ParticleconvertParametersToParticles (generator_t &generator, const std::vector< int > &pdgId, const Acts::ActsDynamicVector &momenta, const Acts::ActsDynamicVector &invariantMasses, Particle &initialParticle, float parametrizedMomentum, bool soft) const
 
unsigned int sampleDiscreteValues (double rnd, const detail::NuclearInteractionParameters::CumulativeDistribution &distribution) const
 
Scalar sampleContinuousValues (double rnd, const detail::NuclearInteractionParameters::CumulativeDistribution &distribution, bool interpolate=false) const
 

Detailed Description

This class provides a parametrised nuclear interaction. The thereby required parametrisation needs to be set and is not provided by default.

Note
This class differs between two different processes labelled as nuclear interaction. Either the initial particle survives (soft) or it gets destroyed (hard) by this process.

Definition at line 39 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 39 of file NuclearInteraction.hpp

Member Typedef Documentation

Definition at line 40 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 40 of file NuclearInteraction.hpp

Member Function Documentation

template<typename generator_t >
std::vector< Particle > ActsFatras::NuclearInteraction::convertParametersToParticles ( generator_t &  generator,
const std::vector< int > &  pdgId,
const Acts::ActsDynamicVector momenta,
const Acts::ActsDynamicVector invariantMasses,
Particle initialParticle,
float  parametrizedMomentum,
bool  soft 
) const
private

Converter from sampled numbers to a vector of particles

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]pdgIdThe PDG IDs
[in]momentaThe momenta
[in]invariantMassesThe invariant masses
[in]initialParticleThe initial particle
[in]parametrizedMomentumMomentum of the parametrisation
[in]softTreat it as soft or hard nuclear interaction
Returns
Vector containing the final state particles

Definition at line 499 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 499 of file NuclearInteraction.hpp

References ActsFatras::Particle::direction(), ActsFatras::Particle::fourPosition(), globalAngle(), i, ActsFatras::Barcode::makeDescendant(), Acts::makeDirectionFromPhiTheta(), momentum, testing::internal::move(), merge_hashes::p, ActsFatras::Particle::particleId(), ActsTests::PropagationDatasets::phi, Acts::VectorHelpers::phi(), ActsFatras::Particle::referenceSurface(), ActsFatras::Particle::setDirection(), ActsFatras::Particle::setPosition4(), ActsFatras::Particle::setProcess(), ActsFatras::Particle::setReferenceSurface(), size, ActsTests::PropagationDatasets::theta, and Acts::VectorHelpers::theta().

Referenced by run().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int ActsFatras::NuclearInteraction::finalStateMultiplicity ( double  rnd,
const detail::NuclearInteractionParameters::CumulativeDistribution distribution 
) const
private

Evaluates the multiplicity of the final state

Parameters
[in]rndRandom number
[in]distributionThe multiplicity distribution
Returns
The final state multiplicity

Definition at line 102 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 102 of file NuclearInteraction.cpp

References sampleDiscreteValues().

Referenced by run().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const detail::NuclearInteractionParameters & ActsFatras::NuclearInteraction::findParameters ( double  rnd,
const detail::NuclearInteractionParametrisation parametrisation,
float  particleMomentum 
) const
private

Retrieves the parametrisation for the particle

Parameters
[in]rndA random number
[in]parametrisationThe storage of parametrisations
[in]particleMomentumThe particles momentum
Returns
The parametrisation

Definition at line 20 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 20 of file NuclearInteraction.cpp

Referenced by generatePathLimits(), and run().

+ Here is the caller graph for this function:

template<typename generator_t >
std::pair<Scalar, Scalar> ActsFatras::NuclearInteraction::generatePathLimits ( generator_t &  generator,
const Particle particle 
) const
inline

This method evaluates the nuclear interaction length L0.

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]particleThe ingoing particle
Returns
valid X0 limit and no limit on L0

Definition at line 57 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 57 of file NuclearInteraction.hpp

References ActsFatras::Particle::absoluteMomentum(), findParameters(), multiParticleParameterisation, ActsFatras::detail::NuclearInteractionParameters::nuclearInteractionProbability, ActsFatras::Particle::pdg(), and sampleContinuousValues().

+ Here is the call graph for this function:

std::pair< ActsFatras::Particle::Scalar, ActsFatras::Particle::Scalar > ActsFatras::NuclearInteraction::globalAngle ( ActsFatras::Particle::Scalar  phi1,
ActsFatras::Particle::Scalar  theta1,
float  phi2,
float  theta2 
) const
private

Converts relative angles to absolute angles wrt the global coordinate system.

Note
It is assumed that the angles of the first particle are provided in the context of the global coordinate system whereas the angles of the second particle are provided relatively to the first particle.
Parameters
[in]phi1The azimuthal angle of the first particle
[in]theta1The polar angle of the first particle
[in]phi2The azimuthal angle of the second particle
[in]theta2The polar angle of the second particle
Returns
Azimuthal and polar angle of the second particle in the global coordinate system

Definition at line 110 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 110 of file NuclearInteraction.cpp

References ActsTests::PropagationDatasets::phi, and ActsTests::PropagationDatasets::theta.

Referenced by convertParametersToParticles().

+ Here is the caller graph for this function:

bool ActsFatras::NuclearInteraction::match ( const Acts::ActsDynamicVector momenta,
const Acts::ActsDynamicVector invariantMasses,
float  parametrizedMomentum 
) const
private

Tests whether the final state momenta and invariant masses are matching to each other to allow the evaluation of particle directions.

Parameters
[in]momentaThe final state momenta
[in]invariantMassesThe final state invariant masses
[in]parametrizedMomentumThe momentum of the parametrized particle
Returns
Decision whether the parameters can be matched to each other or not.

Definition at line 142 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 142 of file NuclearInteraction.cpp

References i, momentum, and size.

Referenced by sampleKinematics().

+ Here is the caller graph for this function:

template<typename generator_t >
bool ActsFatras::NuclearInteraction::run ( generator_t &  generator,
Particle particle,
std::vector< Particle > &  generated 
) const
inline

This method performs a nuclear interaction.

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in,out]particleThe ingoing particle
[out]generatedAdditional generated particles
Returns
True if the particle was killed, false otherwise

Definition at line 100 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 100 of file NuclearInteraction.hpp

References ActsFatras::Particle::absoluteMomentum(), convertParametersToParticles(), finalStateMultiplicity(), findParameters(), ActsFatras::detail::NuclearInteractionParameters::hardKinematicParameters, ActsFatras::detail::NuclearInteractionParameters::hardMultiplicity, ActsFatras::detail::NuclearInteractionParameters::momentum, multiParticleParameterisation, full_chain_odd::multiplicity, particles, ActsFatras::Particle::pdg(), ActsFatras::detail::NuclearInteractionParameters::pdgMap, sampleKinematics(), samplePdgIds(), ActsFatras::Particle::setAbsoluteMomentum(), softInteraction(), ActsFatras::detail::NuclearInteractionParameters::softInteractionProbability, ActsFatras::detail::NuclearInteractionParameters::softKinematicParameters, and ActsFatras::detail::NuclearInteractionParameters::softMultiplicity.

+ Here is the call graph for this function:

Particle::Scalar ActsFatras::NuclearInteraction::sampleContinuousValues ( double  rnd,
const detail::NuclearInteractionParameters::CumulativeDistribution distribution,
bool  interpolate = false 
) const
private

This function performs an inverse sampling to provide a continuous value from a distribition.

Parameters
[in]rndA random number in [0,1]
[in]distributionThe distribution to sample from
[in]interpolateFlag to steer whether an interpolation between neighbouring bins should be performed instead of a bin lookup
Returns
The sampled value

Definition at line 68 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 68 of file NuclearInteraction.cpp

References distance(), it, Acts::UnitConstants::min, and physmon_simulation::rnd.

Referenced by generatePathLimits(), sampleInvariantMasses(), and sampleMomenta().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int ActsFatras::NuclearInteraction::sampleDiscreteValues ( double  rnd,
const detail::NuclearInteractionParameters::CumulativeDistribution distribution 
) const
private

This function performs an inverse sampling to provide a discrete value from a distribution.

Parameters
[in]rndA random number in [0,1]
[in]distributionThe distribution to sample from
Returns
The sampled value

Definition at line 48 of file NuclearInteraction.cpp.

View newest version in sPHENIX GitHub at line 48 of file NuclearInteraction.cpp

References distance(), it, Acts::UnitConstants::min, and physmon_simulation::rnd.

Referenced by finalStateMultiplicity().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename generator_t >
Acts::ActsDynamicVector ActsFatras::NuclearInteraction::sampleInvariantMasses ( generator_t &  generator,
const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity parametrisation 
) const
private

Evaluates the final state invariant masses

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]parametrisationParametrisation of kinematic properties
Returns
Vector containing the invariant masses

Definition at line 406 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 406 of file NuclearInteraction.hpp

References dist(), ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::eigenvaluesInvariantMass, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::eigenvectorsInvariantMass, i, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::invariantMassDistributions, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::meanInvariantMass, Dataset::parameters, sampleContinuousValues(), and size.

Referenced by sampleKinematics().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename generator_t >
std::optional< std::pair< Acts::ActsDynamicVector, Acts::ActsDynamicVector > > ActsFatras::NuclearInteraction::sampleKinematics ( generator_t &  generator,
const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity parameters,
float  momentum 
) const
private

This method samples the kinematics of the final state particles

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]parametersThe parametrisation
[in]momentumThe momentum of the parametrisation
Returns
The final state momenta and invariant masses

Definition at line 473 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 473 of file NuclearInteraction.hpp

References match(), nMatchingTrials, nMatchingTrialsTotal, sampleInvariantMasses(), and sampleMomenta().

Referenced by run().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename generator_t >
Acts::ActsDynamicVector ActsFatras::NuclearInteraction::sampleMomenta ( generator_t &  generator,
const detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity parametrisation,
float  initialMomentum 
) const
private

Evaluates the final state momenta

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]parametrisationParametrisation of kinematic properties
[in]initialMomentumThe initial momentum
Returns
Vector containing the momenta

Definition at line 435 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 435 of file NuclearInteraction.hpp

References dist(), ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::eigenvaluesMomentum, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::eigenvectorsMomentum, i, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::meanMomentum, ActsFatras::detail::NuclearInteractionParameters::ParametersWithFixedMultiplicity::momentumDistributions, Dataset::parameters, sampleContinuousValues(), size, and sum().

Referenced by sampleKinematics().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

template<typename generator_t >
std::vector< int > ActsFatras::NuclearInteraction::samplePdgIds ( generator_t &  generator,
const detail::NuclearInteractionParameters::PdgMap pdgMap,
unsigned int  multiplicity,
int  particlePdg,
bool  soft 
) const
private

Evaluates the final state PDG IDs

Template Parameters
generator_tThe random number generator type
Parameters
[in,out]generatorThe random number generator
[in]pdgMapThe branching probability map
[in]multiplicityThe final state multiplicity
[in]particlePdgThe PDG ID of the initial particle
[in]softTreat it as soft or hard nuclear interaction
Returns
Vector containing the PDG IDs

Definition at line 347 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 347 of file NuclearInteraction.hpp

References fixGDML::element, i, full_chain_odd::multiplicity, and physmon_simulation::rnd.

Referenced by run().

+ Here is the caller graph for this function:

bool ActsFatras::NuclearInteraction::softInteraction ( double  rnd,
float  probability 
) const
inlineprivate

Estimates the interaction type

Parameters
[in]rndRandom number
[in]probabilityThe probability for a soft interaction
Returns
True if a soft interaction occurs

Definition at line 194 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 194 of file NuclearInteraction.hpp

Referenced by run().

+ Here is the caller graph for this function:

Member Data Documentation

detail::MultiParticleNuclearInteractionParametrisation ActsFatras::NuclearInteraction::multiParticleParameterisation

The storage of the parameterisation.

Definition at line 43 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 43 of file NuclearInteraction.hpp

Referenced by generatePathLimits(), and run().

unsigned int ActsFatras::NuclearInteraction::nMatchingTrials = 100

The number of trials to match momenta and inveriant masses.

Definition at line 46 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 46 of file NuclearInteraction.hpp

Referenced by sampleKinematics().

unsigned int ActsFatras::NuclearInteraction::nMatchingTrialsTotal = 1000

Definition at line 47 of file NuclearInteraction.hpp.

View newest version in sPHENIX GitHub at line 47 of file NuclearInteraction.hpp

Referenced by sampleKinematics().


The documentation for this struct was generated from the following files: