ATLAS Offline Software
Loading...
Searching...
No Matches
NnClusterizationFactory.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5 #ifndef SICLUSTERIZATIONTOOL_NnClusterizationFactory_C
6 #define SICLUSTERIZATIONTOOL_NnClusterizationFactory_C
7
8 /******************************************************
9 @class NnClusterizationFactory
10 @author Giacinto Piacquadio (PH-ADE-ID)
11 Package : SiClusterizationTool
12 Created : January 2011
13 DESCRIPTION: Load neural networks used for clustering
14 and deal with:
15 1) estimate of number of particles
16 2) estimate of position of cluster / subclusters
17 3) estimate of error in all cases
18 (full PDF or just RMS)
19 ********************************************************/
20
21
22#include "GaudiKernel/ToolHandle.h"
23#include "GaudiKernel/ServiceHandle.h"
25
26
27//this is a typedef: no forward decl possible
37
38#include <RtypesCore.h> //Double_t
39#include <Eigen/Dense>
40#include <vector>
41#include <array>
42#include <string>
43#include <string_view>
44#include <regex>
45
46
47class TTrainedNetwork;
48class TH1;
49
50namespace lwt {
51 class NanReplacer;
52 class LightweightGraph;
53}
54
55namespace Trk {
56 class NeuralNetworkToHistoTool;
57 class Surface;
58}
59
60namespace InDetDD{
61 class SiLocalPosition;
62}
63
64namespace InDet {
65
66 class PixelCluster;
67
68 struct NNinput{
69 operator bool() const {
70 return !matrixOfToT.empty();
71 }
72 int sizeX = 0;
73 int sizeY = 0;
74 std::vector<std::vector<float> > matrixOfToT;
75 std::vector<float> vectorOfPitchesY;
78 float phi = 0;
79 float theta = 0;
80 float etaModule = 0;
81 bool useTrackInfo = 0;
84 };
85
86 static const InterfaceID IID_NnClusterizationFactory("InDet::NnClusterizationFactory", 1, 0);
87
89
90 public:
91
93 static const InterfaceID& interfaceID() { return IID_NnClusterizationFactory; };
94
95 NnClusterizationFactory(const std::string& name,
96 const std::string& n, const IInterface* p);
98
99 virtual StatusCode initialize() override;
100 virtual StatusCode finalize() override { return StatusCode::SUCCESS; };
101
102 std::vector<double> estimateNumberOfParticles(const InDet::PixelCluster& pCluster,
103 Amg::Vector3D & beamSpotPosition) const;
104
105 std::vector<double> estimateNumberOfParticles(const InDet::PixelCluster& pCluster,
106 const Trk::Surface& pixelSurface,
107 const Trk::TrackParameters& trackParsAtSurface) const;
108
109 /* Public-facing method 1: no track parameters */
110 std::vector<Amg::Vector2D> estimatePositions(const InDet::PixelCluster& pCluster,
111 Amg::Vector3D & beamSpotPosition,
112 std::vector<Amg::MatrixX> & errors,
113 int numberSubClusters) const;
114
115 /* Public-facing method 1: with track parameters */
116 std::vector<Amg::Vector2D> estimatePositions(const InDet::PixelCluster& pCluster,
117 const Trk::Surface& pixelSurface,
118 const Trk::TrackParameters& trackParsAtSurface,
119 std::vector<Amg::MatrixX> & errors,
120 int numberSubClusters) const;
121
122 private:
123
124 // Handling lwtnn inputs
125 typedef std::vector<Eigen::VectorXd> InputVector;
126
127 /* Estimate number of particles for both with and w/o tracks */
128 /* Method 1: using older TTrainedNetworks */
129 std::vector<double> estimateNumberOfParticlesTTN(const TTrainedNetworkCollection &nn_collection,
130 const std::vector<double>& inputData) const;
131
132 /* Estimate number of particles for both with and w/o tracks */
133 /* Method 2: using lwtnn for more flexible interfacing with an ordered vector
134 * Vector order MUST match variable order. */
136
137 /* Estimate position for both with and w/o tracks */
138 /* Method 1: using older TTrainedNetworks */
139 std::vector<Amg::Vector2D> estimatePositionsTTN(
140 const TTrainedNetworkCollection &nn_collection,
141 const std::vector<double>& inputData,
142 const NNinput& input,
143 const InDet::PixelCluster& pCluster,
144 int numberSubClusters,
145 std::vector<Amg::MatrixX> & errors) const;
146
147 /* Estimate position for both with and w/o tracks */
148 /* Method 2: using lwtnn for more flexible interfacing with an ordered vector
149 * Vector order MUST match variable order. */
150 std::vector<Amg::Vector2D> estimatePositionsLWTNN(
152 NNinput& rawInput,
153 const InDet::PixelCluster& pCluster,
154 int numberSubClusters,
155 std::vector<Amg::MatrixX> & errors) const;
156
157 // For error formatting in lwtnn cases
158 static double correctedRMSX(double posPixels) ;
159
160 double correctedRMSY(double posPixels, std::vector<float>& pitches) const;
161
162 /* algorithmic component */
164 Amg::Vector3D & beamSpotPosition,
165 double & tanl) const;
166
167 void addTrackInfoToInput(NNinput& input,
168 const Trk::Surface& pixelSurface,
169 const Trk::TrackParameters& trackParsAtSurface,
170 const double tanl) const;
171
172
173 std::vector<double> assembleInputRunI(NNinput& input) const;
174
175
176 std::vector<double> assembleInputRunII(NNinput& input) const;
177
178 InputVector eigenInput(NNinput & input) const;
179
180 std::vector<Amg::Vector2D> getPositionsFromOutput(std::vector<double> & output,
181 const NNinput & input,
182 const InDet::PixelCluster& pCluster) const;
183
184
185 void getErrorMatrixFromOutput(std::vector<double>& outputX,
186 std::vector<double>& outputY,
187 std::vector<Amg::MatrixX>& errorMatrix,
188 int nParticles) const;
189
190
191 Gaudi::Property< std::vector<std::string> > m_nnOrder
192 {this, "NetworkOrder", {
193 "NumberParticles",
194 "ImpactPoints1P",
195 "ImpactPoints2P",
196 "ImpactPoints3P",
197 "ImpactPointErrorsX1",
198 "ImpactPointErrorsX2",
199 "ImpactPointErrorsX3",
200 "ImpactPointErrorsY1",
201 "ImpactPointErrorsY2",
202 "ImpactPointErrorsY3"},
203 "The order in which the networks will appear in the TTrainedNetworkCollection"};
204
210 static constexpr std::array<std::string_view, kNNetworkTypes> s_nnTypeNames{
211 "NumberParticlesNN",
212 "PositionNN",
213 "ErrorXNN",
214 "ErrorYNN" };
215 static constexpr std::array<unsigned int, kNNetworkTypes> m_nParticleGroup{0U,1U,1U,1U}; // unsigned int
216 static const std::array<std::regex, kNNetworkTypes> m_nnNames;
217
218 unsigned int m_nParticleNNId{};
219 std::vector< std::vector<unsigned int> > m_NNId{};
220
221
222 // Function to be called to assemble the inputs
224
225 //Calculate flat vector dimension, according to input
226 size_t calculateVectorDimension(const bool useTrackInfo) const;
227
228 // Function to be called to compute the output
229 using ReturnType = std::vector<Double_t>;
230 using InputType = std::vector<Double_t>;
231 //the following declares a member variable m_calculateOutput which is a function pointer
232 //to a member function of the TTrainedNetwork. Note to anyone brave enough to update this to C++17using std::function:
233 //TTrainedNetwork::calculateNormalized is overloaded so template resolution does not work trivially.
235
236 ToolHandle<ISiLorentzAngleTool> m_pixelLorentzAngleTool
237 {this, "PixelLorentzAngleTool", "SiLorentzAngleTool/PixelLorentzAngleTool", "Tool to retreive Lorentz angle of Pixel"};
238
240 {this, "PixelChargeCalibCondData", "PixelChargeCalibCondData", "Output key"};
241
243 {this, "NnCollectionReadKey", "PixelClusterNN", "The conditions store key for the pixel cluster NNs"};
244
246 {this, "NnCollectionWithTrackReadKey", "PixelClusterNNWithTrack",
247 "The conditions store key for the pixel cluster NNs which needs tracks as input"};
248
250 {this, "NnCollectionJSONReadKey", "PixelClusterNNJSON",
251 "The conditions key for the pixel cluster NNs configured via JSON file and accessed with lwtnn"};
252
253 // this is written into the JSON config "node_index"
254 // this can be found from the LWTNN GraphConfig object used to initalize the collection objects
255 // option size_t index = graph_config.outputs.at("output_node_name").node_index
256 //
257 Gaudi::Property< std::size_t > m_outputNodesPos1
258 {this, "OutputNodePos1", 7,
259 "Output node for the 1 position networks (LWTNN)"};
260
261 Gaudi::Property< std::vector<std::size_t> > m_outputNodesPos2
262 {this, "OutputNodePos2", { 10, 11 },
263 "List of output nodes for the 2 position network (LWTNN)"};
264
265 Gaudi::Property< std::vector<std::size_t> > m_outputNodesPos3
266 {this, "OutputNodePos3", { 13, 14, 15 },
267 "List of output nodes for the 3 position networks (LWTNN)"};
268
269 Gaudi::Property<unsigned int> m_maxSubClusters
270 {this, "MaxSubClusters", 3, "Maximum number of sub cluster supported by the networks." };
271
273 {this, "correctLorShiftBarrelWithoutTracks",0.,"Lorentz shift correction factor when evaluating NN without track input."};
274
275 Gaudi::Property<double> m_correctLorShiftBarrelWithTracks
276 {this, "correctLorShiftBarrelWithTracks",0.,"Lorentz shift correction factor when evaluating NN with track input."};
277
278 Gaudi::Property<bool> m_useToT
279 {this, "useToT",true,"Use Tot rather than charge." }; // @TODO toggle mode depending on whether a PxielCalibSvc is set ?
280
281 Gaudi::Property<bool> m_addIBL
282 {this, "addIBL", false, "Also apply to clusters in IBL." };
283
284 Gaudi::Property<bool> m_doRunI
285 {this, "doRunI", false, "Use runI style network (outputs are not normalised; add pitches; use charge if not m_useToT)"};
286
287 Gaudi::Property<bool> m_useTTrainedNetworks
288 {this, "useTTrainedNetworks", false, "Use earlier (release-21-like) neural networks stored in ROOT files and accessed via TTrainedNetowrk."};
289
290 Gaudi::Property<bool> m_useRecenteringNNWithouTracks
291 {this, "useRecenteringNNWithoutTracks",false,"Recenter x position when evaluating NN without track input."};
292
293 Gaudi::Property<bool> m_useRecenteringNNWithTracks
294 {this, "useRecenteringNNWithTracks",false,"Recenter x position when evaluating NN with track input."};
295
296 Gaudi::Property<unsigned int> m_sizeX
297 {this, "sizeX",7,"Size of pixel matrix along X"};
298
299 Gaudi::Property<unsigned int> m_sizeY
300 {this, "sizeY",7,"Size of pixel matrix along Y"};
301
302 };
303
304 }//end InDet namespace
305
306 #endif
Store pixel constant parameters in PixelModuleData.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Class to represent a position in the natural frame of a silicon sensor, for Pixel and SCT For Pixel: ...
std::vector< double > assembleInputRunII(NNinput &input) const
void addTrackInfoToInput(NNinput &input, const Trk::Surface &pixelSurface, const Trk::TrackParameters &trackParsAtSurface, const double tanl) const
Gaudi::Property< unsigned int > m_maxSubClusters
SG::ReadCondHandleKey< PixelChargeCalibCondData > m_chargeDataKey
std::vector< double > estimateNumberOfParticlesLWTNN(NnClusterizationFactory::InputVector &input) const
virtual StatusCode finalize() override
Gaudi::Property< unsigned int > m_sizeX
double correctedRMSY(double posPixels, std::vector< float > &pitches) const
ReturnType(::TTrainedNetwork::* m_calculateOutput)(const InputType &input) const
NNinput createInput(const InDet::PixelCluster &pCluster, Amg::Vector3D &beamSpotPosition, double &tanl) const
virtual StatusCode initialize() override
Gaudi::Property< double > m_correctLorShiftBarrelWithoutTracks
static const InterfaceID & interfaceID()
AlgTool interface methods.
Gaudi::Property< std::size_t > m_outputNodesPos1
Gaudi::Property< std::vector< std::size_t > > m_outputNodesPos2
ToolHandle< ISiLorentzAngleTool > m_pixelLorentzAngleTool
Gaudi::Property< std::vector< std::size_t > > m_outputNodesPos3
std::vector< Amg::Vector2D > estimatePositionsLWTNN(NnClusterizationFactory::InputVector &input, NNinput &rawInput, const InDet::PixelCluster &pCluster, int numberSubClusters, std::vector< Amg::MatrixX > &errors) const
std::vector< double > assembleInputRunI(NNinput &input) const
SG::ReadCondHandleKey< LWTNNCollection > m_readKeyJSON
std::vector< Amg::Vector2D > estimatePositions(const InDet::PixelCluster &pCluster, Amg::Vector3D &beamSpotPosition, std::vector< Amg::MatrixX > &errors, int numberSubClusters) const
Gaudi::Property< std::vector< std::string > > m_nnOrder
Gaudi::Property< double > m_correctLorShiftBarrelWithTracks
Gaudi::Property< bool > m_useTTrainedNetworks
std::vector< double > estimateNumberOfParticlesTTN(const TTrainedNetworkCollection &nn_collection, const std::vector< double > &inputData) const
static constexpr std::array< unsigned int, kNNetworkTypes > m_nParticleGroup
SG::ReadCondHandleKey< TTrainedNetworkCollection > m_readKeyWithoutTrack
std::vector< double >(InDet::NnClusterizationFactory::* m_assembleInput)(NNinput &input) const
std::vector< Eigen::VectorXd > InputVector
NnClusterizationFactory(const std::string &name, const std::string &n, const IInterface *p)
std::vector< Amg::Vector2D > estimatePositionsTTN(const TTrainedNetworkCollection &nn_collection, const std::vector< double > &inputData, const NNinput &input, const InDet::PixelCluster &pCluster, int numberSubClusters, std::vector< Amg::MatrixX > &errors) const
static const std::array< std::regex, kNNetworkTypes > m_nnNames
Gaudi::Property< bool > m_useRecenteringNNWithouTracks
static constexpr std::array< std::string_view, kNNetworkTypes > s_nnTypeNames
InputVector eigenInput(NNinput &input) const
void getErrorMatrixFromOutput(std::vector< double > &outputX, std::vector< double > &outputY, std::vector< Amg::MatrixX > &errorMatrix, int nParticles) const
Gaudi::Property< bool > m_useRecenteringNNWithTracks
std::vector< std::vector< unsigned int > > m_NNId
std::vector< Amg::Vector2D > getPositionsFromOutput(std::vector< double > &output, const NNinput &input, const InDet::PixelCluster &pCluster) const
static double correctedRMSX(double posPixels)
size_t calculateVectorDimension(const bool useTrackInfo) const
Gaudi::Property< unsigned int > m_sizeY
SG::ReadCondHandleKey< TTrainedNetworkCollection > m_readKeyWithTrack
std::vector< double > estimateNumberOfParticles(const InDet::PixelCluster &pCluster, Amg::Vector3D &beamSpotPosition) const
Abstract Base Class for tracking surfaces.
Eigen::Matrix< double, 3, 1 > Vector3D
Message Stream Member.
Primary Vertex Finder.
static const InterfaceID IID_NnClusterizationFactory("InDet::NnClusterizationFactory", 1, 0)
Ensure that the ATLAS eigen extensions are properly loaded.
ParametersBase< TrackParametersDim, Charged > TrackParameters
std::vector< std::vector< float > > matrixOfToT
std::vector< float > vectorOfPitchesY