ATLAS Offline Software
Loading...
Searching...
No Matches
PunchThroughTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef ISF_PUNCHTHROUGHTOOLS_SRC_PUNCHTHROUGHTOOL_H
6#define ISF_PUNCHTHROUGHTOOLS_SRC_PUNCHTHROUGHTOOL_H 1
7
9#include <string>
10
12
13// Athena Base
15
16//Geometry
18
20
21// Gaudi & StoreGate
22#include "GaudiKernel/IPartPropSvc.h"
23
25
27
29
30class XMLCoreNode;
31
32/*-------------------------------------------------------------------------
33 * Forward declarations
34 *-------------------------------------------------------------------------*/
35
36class TFile;
37
38namespace HepPDT {
39 class ParticleDataTable;
40}
41
42namespace ISF {
44 class PDFcreator;
45}
46
47namespace ISF {
48
49 class PunchThroughTool : public extends<AthAlgTool, IPunchThroughTool>
50 {
51 public:
53 PunchThroughTool(const std::string&,const std::string&,const IInterface*);
54
56 virtual ~PunchThroughTool () = default;
57
59 virtual StatusCode initialize();
61 virtual StatusCode finalize ();
63 const ISF::ISFParticleVector* computePunchThroughParticles(const ISF::ISFParticle &isfp, const TFCSSimulationState& simulstate, CLHEP::HepRandomEngine* rndmEngine) const;
64
65 private:
66 struct InfoMap
67 {
68 InfoMap (const XMLCoreNode& node);
69
70 std::string name;
71 std::vector<double> etaMins;
72 std::vector<double> etaMaxs;
73 std::vector<int> pidStr;
74 };
75
76 /*---------------------------------------------------------------------
77 * Private member functions
78 *---------------------------------------------------------------------*/
80 StatusCode registerParticle(int pdgID, bool doAntiparticle = false,
81 double minEnergy = 0., int maxNumParticles = -1,
82 double numParticlesFactor = 1., double energyFactor = 1.,
83 double posAngleFactor = 1.,
84 double momAngleFactor = 1.);
87 StatusCode registerCorrelation(int pdgID1, int pdgID2,double minCorrEnergy = 0., double fullCorrEnergy = 0.);
88
90 std::unique_ptr<ISF::PDFcreator> readLookuptablePDF(int pdgID, const std::string& folderName);
91
97 int getAllParticles(const ISF::ISFParticle &isfp, ISFParticleVector& isfpCont, CLHEP::HepRandomEngine* rndmEngine, int pdg, double interpEnergy, double interpEta, int numParticles = -1) const;
98
102 int getCorrelatedParticles(const ISF::ISFParticle &isfp, ISFParticleVector& isfpCont, int doPdg, int corrParticles, CLHEP::HepRandomEngine* rndmEngine, double interpEnergy, double interpEta) const;
103
105 ISF::ISFParticle *getOneParticle(const ISF::ISFParticle &isfp, int pdg, CLHEP::HepRandomEngine* rndmEngine, double interpEnergy, double interpEta) const;
106
108 ISF::ISFParticle *createExitPs(const ISF::ISFParticle &isfp, int PDGcode, double energy, double theta, double phi, double momTheta, double momPhi) const;
109
111 double getFloatAfterPatternInStr(const char *str, const char *pattern);
113 Amg::Vector3D propagator(double theta, double phi) const;
114
115 //apply the inverse PCA transform
116 std::vector<double> inversePCA(int pcaCdfIterator, std::vector<double> &variables) const;
117
118 //apply the inverse CDF trainsform
119 static double inverseCdfTransform(double variable, const std::map<double, double>& inverse_cdf_map) ;
120
121 //dot product between matrix and vector, used to inverse PCA
122 static std::vector<double> dotProduct(const std::vector<std::vector<double>> &m, const std::vector<double> &v) ;
123
124 //returns normal cdf based on normal gaussian value
125 static double normal_cdf(double x) ;
126
127 //apply energy interpolation
128 double interpolateEnergy(const double &energy, CLHEP::HepRandomEngine* rndmEngine) const;
129
130 //apply eta interpolation
131 double interpolateEta(const double &eta, CLHEP::HepRandomEngine* rndmEngine) const;
132
133 //get the infoMap from xml file based on the xmlpathname and also name of mainNode
134 std::vector<InfoMap> getInfoMap(const std::string& mainNode, const XMLCoreNode& doc);
135
136 //decide the pca / cdf part to read based on pdgId and eta
137 int passedParamIterator(int pid, double eta, const std::vector<InfoMap> &mapvect) const;
138
139 //load inverse quantile transformer from XML
140 StatusCode initializeInverseCDF(const std::string & quantileTransformerConfigFile);
141
142 //get CDF mapping for individual XML node
143 static std::map<double, double> getVariableCDFmappings(const XMLCoreNode* node);
144
145 //load inverse PCA from XML
146 StatusCode initializeInversePCA(const std::string & inversePCAConfigFile);
147
148
149 /*---------------------------------------------------------------------
150 * Private members
151 *---------------------------------------------------------------------*/
152
154 std::vector<double> m_energyPoints;
155 std::vector<double> m_etaPoints;
156
158 double m_R1{0.};
159 double m_R2{0.};
160 double m_z1{0.};
161 double m_z2{0.};
162
164 const HepPDT::ParticleDataTable* m_particleDataTable{nullptr};
165
167 TFile* m_fileLookupTable{nullptr};
168
170 std::map<int, PunchThroughParticle*> m_particles;
171
172 /*---------------------------------------------------------------------
173 * Properties
174 *---------------------------------------------------------------------*/
175
177 StringProperty m_filenameLookupTable{this, "FilenameLookupTable", "CaloPunchThroughParametrisation.root", "holds the filename of the lookup table"};
178 StringProperty m_filenameInverseCDF{this, "FilenameInverseCdf", "", "holds the filename of inverse quantile transformer config"};
179 StringProperty m_filenameInversePCA{this, "FilenameInversePca", "", "holds the filename of inverse PCA config"};
180
181 PublicToolHandle<IPunchThroughClassifier> m_punchThroughClassifier{this, "PunchThroughClassifier", "ISF_PunchThroughClassifier", ""};
182 IntegerArrayProperty m_pdgInitiators{this, "PunchThroughInitiators", {}, "vector of punch-through initiator pgds"};
183 IntegerArrayProperty m_initiatorsMinEnergy{this, "InitiatorsMinEnergy", {}, "vector of punch-through initiator min energies to create punch through"};
184 DoubleArrayProperty m_initiatorsEtaRange{this, "InitiatorsEtaRange", {}, "vector of min and max abs eta range to allow punch through initiators"};
185 IntegerArrayProperty m_punchThroughParticles{this, "PunchThroughParticles", {}, "vector of pdgs of the particles produced in punch-throughs"};
186 BooleanArrayProperty m_doAntiParticles{this, "DoAntiParticles", {}, "vector of bools to determine if anti-particles are created for each punch-through particle type"};
187 IntegerArrayProperty m_correlatedParticle{this, "CorrelatedParticle", {}, "holds the pdg of the correlated particle for each given pdg"};
188 DoubleArrayProperty m_minCorrEnergy{this, "MinCorrelationEnergy", {}, "holds the energy threshold below which no particle correlation is computed"};
189 DoubleArrayProperty m_fullCorrEnergy{this, "FullCorrelationEnergy", {}, "holds the energy threshold above which a particle correlation is fully developed"};
190 DoubleArrayProperty m_posAngleFactor{this, "ScalePosDeflectionAngles", {}, "tuning parameter to scale the position deflection angles"};
191 DoubleArrayProperty m_momAngleFactor{this, "ScaleMomDeflectionAngles", {}, "tuning parameter to scale the momentum deflection angles"};
192 DoubleArrayProperty m_minEnergy{this, "MinEnergy", {}, "punch-through particles minimum energies"};
193 IntegerArrayProperty m_maxNumParticles{this, "MaxNumParticles", {}, "maximum number of punch-through particles for each particle type"};
194 DoubleArrayProperty m_numParticlesFactor{this, "NumParticlesFactor", {}, "scale the number of punch-through particles"};
195 DoubleArrayProperty m_energyFactor{this, "EnergyFactor", {}, "scale the energy of the punch-through particles"};
196
197
198
199 /*---------------------------------------------------------------------
200 * ServiceHandles
201 *---------------------------------------------------------------------*/
202 ServiceHandle<IPartPropSvc> m_particlePropSvc{this, "PartPropSvc", "PartPropSvc", "particle properties svc"};
203 ServiceHandle<IGeoIDSvc> m_geoIDSvc{this, "GeoIDSvc", "ISF::GeoIDSvc"};
204 ServiceHandle<IEnvelopeDefSvc> m_envDefSvc{this, "EnvelopeDefSvc", "AtlasGeometry_EnvelopeDefSvc"};
205
207 DoubleProperty m_beamPipe{this, "BeamPipeRadius", 500.};
208
210 std::vector<std::vector<std::vector<double>>> m_inverse_PCA_matrix;
211 std::vector<std::vector<double>> m_PCA_means;
212
214 std::vector<InfoMap> m_xml_info_pca;
215 std::vector<InfoMap> m_xml_info_cdf;
216
218 std::vector<std::map<double, double>> m_variable0_inverse_cdf;
219 std::vector<std::map<double, double>> m_variable1_inverse_cdf;
220 std::vector<std::map<double, double>> m_variable2_inverse_cdf;
221 std::vector<std::map<double, double>> m_variable3_inverse_cdf;
222 std::vector<std::map<double, double>> m_variable4_inverse_cdf;
223 };
224}
225
226#endif
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
#define x
The generic ISF particle definition,.
Definition ISFParticle.h:42
Creates random numbers with a distribution given as ROOT TF1.
Definition PDFcreator.h:37
This class holds information for different properties of a punch-through particle (energy,...
std::vector< double > inversePCA(int pcaCdfIterator, std::vector< double > &variables) const
std::vector< std::vector< std::vector< double > > > m_inverse_PCA_matrix
pca vectors
const ISF::ISFParticleVector * computePunchThroughParticles(const ISF::ISFParticle &isfp, const TFCSSimulationState &simulstate, CLHEP::HepRandomEngine *rndmEngine) const
interface function: fill a vector with the punch-through particles
static std::vector< double > dotProduct(const std::vector< std::vector< double > > &m, const std::vector< double > &v)
StringProperty m_filenameInverseCDF
IntegerArrayProperty m_correlatedParticle
DoubleArrayProperty m_minCorrEnergy
static std::map< double, double > getVariableCDFmappings(const XMLCoreNode *node)
std::vector< std::map< double, double > > m_variable1_inverse_cdf
std::unique_ptr< ISF::PDFcreator > readLookuptablePDF(int pdgID, const std::string &folderName)
reads out the lookuptable for the given type of particle
static double inverseCdfTransform(double variable, const std::map< double, double > &inverse_cdf_map)
ServiceHandle< IPartPropSvc > m_particlePropSvc
StringProperty m_filenameLookupTable
Properties.
StringProperty m_filenameInversePCA
std::vector< std::map< double, double > > m_variable3_inverse_cdf
PublicToolHandle< IPunchThroughClassifier > m_punchThroughClassifier
BooleanArrayProperty m_doAntiParticles
double getFloatAfterPatternInStr(const char *str, const char *pattern)
get the floating point number in a string, after the given pattern
int getAllParticles(const ISF::ISFParticle &isfp, ISFParticleVector &isfpCont, CLHEP::HepRandomEngine *rndmEngine, int pdg, double interpEnergy, double interpEta, int numParticles=-1) const
create the right number of punch-through particles for the given pdg and return the number of particl...
DoubleArrayProperty m_energyFactor
StatusCode registerCorrelation(int pdgID1, int pdgID2, double minCorrEnergy=0., double fullCorrEnergy=0.)
register a correlation for the two given types of punch-through particles with a given energy thresho...
Amg::Vector3D propagator(double theta, double phi) const
get particle through the calorimeter
DoubleProperty m_beamPipe
beam pipe radius
DoubleArrayProperty m_momAngleFactor
static double normal_cdf(double x)
DoubleArrayProperty m_initiatorsEtaRange
DoubleArrayProperty m_fullCorrEnergy
std::map< int, PunchThroughParticle * > m_particles
needed to create punch-through particles with the right distributions
const HepPDT::ParticleDataTable * m_particleDataTable
ParticleDataTable needed to get connection pdg_code <-> charge.
std::vector< InfoMap > m_xml_info_cdf
StatusCode registerParticle(int pdgID, bool doAntiparticle=false, double minEnergy=0., int maxNumParticles=-1, double numParticlesFactor=1., double energyFactor=1., double posAngleFactor=1., double momAngleFactor=1.)
registers a type of punch-through particles which will be simulated
std::vector< std::vector< double > > m_PCA_means
DoubleArrayProperty m_numParticlesFactor
PunchThroughTool(const std::string &, const std::string &, const IInterface *)
Constructor.
ISF::ISFParticle * getOneParticle(const ISF::ISFParticle &isfp, int pdg, CLHEP::HepRandomEngine *rndmEngine, double interpEnergy, double interpEta) const
create exactly one punch-through particle with the given pdg and the given max energy
std::vector< std::map< double, double > > m_variable2_inverse_cdf
IntegerArrayProperty m_initiatorsMinEnergy
StatusCode initializeInverseCDF(const std::string &quantileTransformerConfigFile)
std::vector< std::map< double, double > > m_variable0_inverse_cdf
(vector of map) for CDF mappings
virtual StatusCode initialize()
AlgTool initialize method.
int passedParamIterator(int pid, double eta, const std::vector< InfoMap > &mapvect) const
std::vector< double > m_energyPoints
energy and eta points in param
ServiceHandle< IGeoIDSvc > m_geoIDSvc
std::vector< InfoMap > m_xml_info_pca
infoMaps
int getCorrelatedParticles(const ISF::ISFParticle &isfp, ISFParticleVector &isfpCont, int doPdg, int corrParticles, CLHEP::HepRandomEngine *rndmEngine, double interpEnergy, double interpEta) const
get the right number of particles for the given pdg while considering the correlation to an other par...
IntegerArrayProperty m_pdgInitiators
DoubleArrayProperty m_minEnergy
IntegerArrayProperty m_maxNumParticles
double m_R1
calo-MS borders
std::vector< std::map< double, double > > m_variable4_inverse_cdf
StatusCode initializeInversePCA(const std::string &inversePCAConfigFile)
double interpolateEta(const double &eta, CLHEP::HepRandomEngine *rndmEngine) const
TFile * m_fileLookupTable
ROOT objects.
virtual StatusCode finalize()
AlgTool finalize method.
ISF::ISFParticle * createExitPs(const ISF::ISFParticle &isfp, int PDGcode, double energy, double theta, double phi, double momTheta, double momPhi) const
create a ISF Particle state at the MS entrace containing a particle with the given properties
ServiceHandle< IEnvelopeDefSvc > m_envDefSvc
virtual ~PunchThroughTool()=default
Destructor.
DoubleArrayProperty m_posAngleFactor
double interpolateEnergy(const double &energy, CLHEP::HepRandomEngine *rndmEngine) const
std::vector< InfoMap > getInfoMap(const std::string &mainNode, const XMLCoreNode &doc)
IntegerArrayProperty m_punchThroughParticles
std::vector< double > m_etaPoints
Simple DOM-like node structure to hold the result of XML parsing.
Definition XMLCoreNode.h:45
Definition node.h:24
Eigen::Matrix< double, 3, 1 > Vector3D
ISFParticleOrderedQueue.
std::vector< ISF::ISFParticle * > ISFParticleVector
ISFParticle vector.
std::vector< double > etaMins
InfoMap(const XMLCoreNode &node)
std::vector< double > etaMaxs