CMS 3D CMS Logo

TracksterInferenceByDNN.cc
Go to the documentation of this file.
8 #include "TrackstersPCA.h"
9 
10 namespace ticl {
11  using namespace cms::Ort; // Use ONNXRuntime namespace
12 
13  // Constructor for TracksterInferenceByDNN
16  id_modelPath_(
17  conf.getParameter<edm::FileInPath>("onnxPIDModelPath").fullPath()), // Path to the PID model CLU3D
18  en_modelPath_(
19  conf.getParameter<edm::FileInPath>("onnxEnergyModelPath").fullPath()), // Path to the Energy model CLU3D
20  inputNames_(conf.getParameter<std::vector<std::string>>("inputNames")), // Define input names for inference
21  output_en_(conf.getParameter<std::vector<std::string>>("output_en")), // Define output energy for inference
22  output_id_(conf.getParameter<std::vector<std::string>>("output_id")), // Define output PID for inference
23  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")), // Minimum cluster energy
24  eidNLayers_(conf.getParameter<int>("eid_n_layers")), // Number of layers
25  eidNClusters_(conf.getParameter<int>("eid_n_clusters")), // Number of clusters
26  doPID_(conf.getParameter<int>("doPID")), // Number of clusters
27  doRegression_(conf.getParameter<int>("doRegression")) // Number of clusters
28  {
29  // Initialize ONNX Runtime sessions for PID and Energy models
30  static std::unique_ptr<cms::Ort::ONNXRuntime> onnxPIDRuntimeInstance =
31  std::make_unique<cms::Ort::ONNXRuntime>(id_modelPath_.c_str());
32  onnxPIDSession_ = onnxPIDRuntimeInstance.get();
33  static std::unique_ptr<cms::Ort::ONNXRuntime> onnxEnergyRuntimeInstance =
34  std::make_unique<cms::Ort::ONNXRuntime>(en_modelPath_.c_str());
35  onnxEnergySession_ = onnxEnergyRuntimeInstance.get();
36  }
37 
38  // Method to process input data and prepare it for inference
39  void TracksterInferenceByDNN::inputData(const std::vector<reco::CaloCluster>& layerClusters,
40  std::vector<Trackster>& tracksters) {
41  tracksterIndices_.clear(); // Clear previous indices
42  for (int i = 0; i < static_cast<int>(tracksters.size()); i++) {
43  float sumClusterEnergy = 0.;
44  for (const unsigned int& vertex : tracksters[i].vertices()) {
45  sumClusterEnergy += static_cast<float>(layerClusters[vertex].energy());
46  if (sumClusterEnergy >= eidMinClusterEnergy_) {
47  tracksters[i].setRegressedEnergy(0.f); // Set regressed energy to 0
48  tracksters[i].zeroProbabilities(); // Zero out probabilities
49  tracksterIndices_.push_back(i); // Add index to the list
50  break;
51  }
52  }
53  }
54 
55  // Prepare input shapes and data for inference
56  batchSize_ = static_cast<int>(tracksterIndices_.size());
57  if (batchSize_ == 0)
58  return; // Exit if no tracksters
59 
60  std::vector<int64_t> inputShape = {batchSize_, eidNLayers_, eidNClusters_, eidNFeatures_};
61  input_shapes_ = {inputShape};
62 
63  input_Data_.clear();
65 
66  for (int i = 0; i < batchSize_; i++) {
67  const Trackster& trackster = tracksters[tracksterIndices_[i]];
68 
69  // Prepare indices and sort clusters based on energy
70  std::vector<int> clusterIndices(trackster.vertices().size());
71  for (int k = 0; k < static_cast<int>(trackster.vertices().size()); k++) {
72  clusterIndices[k] = k;
73  }
74 
75  std::sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int& a, const int& b) {
76  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
77  });
78 
79  std::vector<int> seenClusters(eidNLayers_, 0);
80 
81  // Fill input data with cluster information
82  for (const int& k : clusterIndices) {
83  const reco::CaloCluster& cluster = layerClusters[trackster.vertices(k)];
84  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
85  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
86  auto index = (i * eidNLayers_ + j) * eidNFeatures_ * eidNClusters_ + seenClusters[j] * eidNFeatures_;
87  input_Data_[0][index] =
88  static_cast<float>(cluster.energy() / static_cast<float>(trackster.vertex_multiplicity(k)));
89  input_Data_[0][index + 1] = static_cast<float>(std::abs(cluster.eta()));
90  input_Data_[0][index + 2] = static_cast<float>(cluster.phi());
91  seenClusters[j]++;
92  }
93  }
94  }
95  }
96 
97  // Method to run inference and update tracksters
98  void TracksterInferenceByDNN::runInference(std::vector<Trackster>& tracksters) {
99  if (batchSize_ == 0)
100  return; // Exit if no batch
101 
102  if (doPID_ and doRegression_) {
103  // Run energy model inference
105  auto& energyOutputTensor = result[0];
106  if (!output_en_.empty()) {
107  for (int i = 0; i < static_cast<int>(batchSize_); i++) {
108  const float energy = energyOutputTensor[i];
109  tracksters[tracksterIndices_[i]].setRegressedEnergy(energy); // Update energy
110  }
111  }
112  }
113 
114  if (doPID_) {
115  // Run PID model inference
117  auto pidOutputTensor = pidOutput[0];
118  float* probs = pidOutputTensor.data();
119  if (!output_id_.empty()) {
120  for (int i = 0; i < batchSize_; i++) {
121  tracksters[tracksterIndices_[i]].setProbabilities(probs); // Update probabilities
122  probs += tracksters[tracksterIndices_[i]].id_probabilities().size(); // Move to next set of probabilities
123  }
124  }
125  }
126  }
127  // Method to fill parameter set description for configuration
129  iDesc.add<int>("algo_verbosity", 0);
130  iDesc
131  .add<edm::FileInPath>("onnxPIDModelPath",
132  edm::FileInPath("RecoHGCal/TICL/data/ticlv5/onnx_models/patternrecognition/id_v0.onnx"))
133  ->setComment("Path to ONNX PID model CLU3D");
134  iDesc
136  "onnxEnergyModelPath",
137  edm::FileInPath("RecoHGCal/TICL/data/ticlv5/onnx_models/patternrecognition/energy_v0.onnx"))
138  ->setComment("Path to ONNX Energy model CLU3D");
139  iDesc.add<std::vector<std::string>>("inputNames", {"input"});
140  iDesc.add<std::vector<std::string>>("output_en", {"enreg_output"});
141  iDesc.add<std::vector<std::string>>("output_id", {"pid_output"});
142  iDesc.add<double>("eid_min_cluster_energy", 1.0);
143  iDesc.add<int>("eid_n_layers", 50);
144  iDesc.add<int>("eid_n_clusters", 10);
145  iDesc.add<int>("doPID", 1);
146  iDesc.add<int>("doRegression", 1);
147  }
148 } // namespace ticl
const std::vector< std::string > output_en_
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
std::vector< std::vector< int64_t > > input_shapes_
const std::vector< std::string > inputNames_
void inputData(const std::vector< reco::CaloCluster > &layerClusters, std::vector< Trackster > &tracksters) override
TracksterInferenceByDNN(const edm::ParameterSet &conf)
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
const cms::Ort::ONNXRuntime * onnxEnergySession_
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
const cms::Ort::ONNXRuntime * onnxPIDSession_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void runInference(std::vector< Trackster > &tracksters) override
const std::vector< std::string > output_id_
double f[11][100]
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double energy() const
cluster energy
Definition: CaloCluster.h:149
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
std::vector< float > & vertex_multiplicity()
Definition: Trackster.h:58
double b
Definition: hdecay.h:120
HLT enums.
double a
Definition: hdecay.h:121
Definition: Common.h:10
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
std::vector< std::vector< float > > input_Data_
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381
FloatArrays run(const std::vector< std::string > &input_names, FloatArrays &input_values, const std::vector< std::vector< int64_t >> &input_shapes={}, const std::vector< std::string > &output_names={}, int64_t batch_size=1) const
Definition: ONNXRuntime.cc:87