![]() |
|
|||
File indexing completed on 2025-07-12 07:52:21
0001 // This file is part of the ACTS project. 0002 // 0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project 0004 // 0005 // This Source Code Form is subject to the terms of the Mozilla Public 0006 // License, v. 2.0. If a copy of the MPL was not distributed with this 0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/. 0008 0009 #pragma once 0010 0011 #include "Acts/Definitions/PdgParticle.hpp" 0012 #include "Acts/Definitions/Units.hpp" 0013 #include "ActsExamples/Framework/RandomNumbers.hpp" 0014 0015 #include <array> 0016 #include <cstddef> 0017 #include <functional> 0018 #include <limits> 0019 #include <memory> 0020 #include <numbers> 0021 #include <optional> 0022 #include <random> 0023 0024 namespace HepMC3 { 0025 class GenEvent; 0026 } 0027 0028 namespace ActsExamples { 0029 0030 /// @brief Generator interface for vertices and particles 0031 struct ParticlesGenerator { 0032 /// @brief Virtual destructor required 0033 virtual ~ParticlesGenerator() = default; 0034 /// @brief Generate vertices and particles for one interaction 0035 /// @note This method cannot be `const` because the Pythia8 generator 0036 /// uses the Pythia8 interfaces, which is non-const 0037 /// 0038 /// @param rng Shared random number generator instance 0039 /// @return The vertex and particle containers 0040 virtual std::shared_ptr<HepMC3::GenEvent> operator()(RandomEngine& rng) = 0; 0041 }; 0042 0043 /// Generate particles from uniform parameter distributions. 0044 /// 0045 /// Generates a single vertex with the given number of tracks. The track 0046 /// direction is drawn from a uniform distribution on the unit sphere (within 0047 /// the given limits). Its absolute momentum is drawn from a uniform 0048 /// distribution. Position and time are always set to zero. 0049 class ParametricParticleGenerator : public ParticlesGenerator { 0050 public: 0051 struct Config { 0052 /// Low, high (exclusive) for the transverse direction angle. 0053 double phiMin = -std::numbers::pi; 0054 double phiMax = std::numbers::pi; 0055 /// Low, high (inclusive) for the longitudinal direction angle. 0056 /// 0057 /// This intentionally uses theta instead of eta so it can represent the 0058 /// full direction space with finite values. 0059 /// 0060 /// @note This is the standard generation, for detector performance 0061 /// classification, where a flat distribution in eta can be useful, 0062 /// this can be set by the etaUniform flag; 0063 /// 0064 double thetaMin = std::numeric_limits<double>::min(); 0065 double thetaMax = std::numbers::pi - std::numeric_limits<double>::epsilon(); 0066 bool etaUniform = false; 0067 /// Low, high (exclusive) for absolute/transverse momentum. 0068 double pMin = 1 * Acts::UnitConstants::GeV; 0069 double pMax = 10 * Acts::UnitConstants::GeV; 0070 /// Indicate if the momentum refers to transverse momentum. 0071 bool pTransverse = false; 0072 /// Indicate if the momentum should be uniformly distributed in log space. 0073 bool pLogUniform = false; 0074 /// (Absolute) PDG particle number to identify the particle type. 0075 Acts::PdgParticle pdg = Acts::PdgParticle::eMuon; 0076 /// Randomize the charge and flip the PDG particle number sign accordingly. 0077 bool randomizeCharge = false; 0078 /// Number of particles. 0079 std::size_t numParticles = 1; 0080 0081 /// Overrides particle charge. 0082 std::optional<double> charge; 0083 /// Overrides particle mass. 0084 std::optional<double> mass; 0085 }; 0086 0087 explicit ParametricParticleGenerator(const Config& cfg); 0088 0089 /// Generate a single primary vertex with the given number of particles. 0090 std::shared_ptr<HepMC3::GenEvent> operator()(RandomEngine& rng) override; 0091 0092 private: 0093 using UniformIndex = std::uniform_int_distribution<std::uint8_t>; 0094 using UniformReal = std::uniform_real_distribution<double>; 0095 0096 Config m_cfg; 0097 0098 // will be automatically set from PDG data tables 0099 double m_mass{}; 0100 0101 // (anti-)particle choice is one random draw but defines two properties 0102 std::array<Acts::PdgParticle, 2> m_pdgChoices{}; 0103 0104 UniformIndex m_particleTypeChoice; 0105 UniformReal m_phiDist; 0106 std::function<std::pair<double, double>(RandomEngine& rng)> m_sinCosThetaDist; 0107 std::function<double(RandomEngine& rng)> m_somePDist; 0108 }; 0109 0110 } // namespace ActsExamples
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |