CMS 3D CMS Logo

TracksterInferenceByCNNv4.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 TracksterInferenceByCNNv4
16  modelPath_(conf.getParameter<edm::FileInPath>("onnxModelPath").fullPath()), // Path to the PID model CLU3D
17  inputNames_(conf.getParameter<std::vector<std::string>>("inputNames")), // Define input names for inference
18  outputNames_(conf.getParameter<std::vector<std::string>>("outputNames")), // Define output names for inference
19  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")), // Minimum cluster energy
20  eidNLayers_(conf.getParameter<int>("eid_n_layers")), // Number of layers
21  eidNClusters_(conf.getParameter<int>("eid_n_clusters")), // Number of clusters
22  doPID_(conf.getParameter<int>("doPID")), // Number of clusters
23  doRegression_(conf.getParameter<int>("doRegression")) // Number of clusters
24  {
25  // Initialize ONNX Runtime sessions for PID and Energy models
26  static std::unique_ptr<cms::Ort::ONNXRuntime> onnxRuntimeInstance =
27  std::make_unique<cms::Ort::ONNXRuntime>(modelPath_.c_str());
28  onnxSession_ = onnxRuntimeInstance.get();
29  }
30 
31  // Method to process input data and prepare it for inference
32  void TracksterInferenceByCNNv4::inputData(const std::vector<reco::CaloCluster>& layerClusters,
33  std::vector<Trackster>& tracksters) {
34  tracksterIndices_.clear(); // Clear previous indices
35  for (int i = 0; i < static_cast<int>(tracksters.size()); i++) {
36  float sumClusterEnergy = 0.;
37  for (const unsigned int& vertex : tracksters[i].vertices()) {
38  sumClusterEnergy += static_cast<float>(layerClusters[vertex].energy());
39  if (sumClusterEnergy >= eidMinClusterEnergy_) {
40  tracksters[i].setRegressedEnergy(0.f); // Set regressed energy to 0
41  tracksters[i].zeroProbabilities(); // Zero out probabilities
42  tracksterIndices_.push_back(i); // Add index to the list
43  break;
44  }
45  }
46  }
47 
48  // Prepare input shapes and data for inference
49  batchSize_ = static_cast<int>(tracksterIndices_.size());
50  if (batchSize_ == 0)
51  return; // Exit if no tracksters
52 
53  std::vector<int64_t> inputShape = {batchSize_, eidNLayers_, eidNClusters_, eidNFeatures_};
54  input_shapes_ = {inputShape};
55 
56  input_Data_.clear();
58 
59  for (int i = 0; i < batchSize_; i++) {
60  const Trackster& trackster = tracksters[tracksterIndices_[i]];
61 
62  // Prepare indices and sort clusters based on energy
63  std::vector<int> clusterIndices(trackster.vertices().size());
64  for (int k = 0; k < static_cast<int>(trackster.vertices().size()); k++) {
65  clusterIndices[k] = k;
66  }
67 
68  std::sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int& a, const int& b) {
69  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
70  });
71 
72  std::vector<int> seenClusters(eidNLayers_, 0);
73 
74  // Fill input data with cluster information
75  for (const int& k : clusterIndices) {
76  const reco::CaloCluster& cluster = layerClusters[trackster.vertices(k)];
77  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
78  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
79  auto index = (i * eidNLayers_ + j) * eidNFeatures_ * eidNClusters_ + seenClusters[j] * eidNFeatures_;
80  input_Data_[0][index] =
81  static_cast<float>(cluster.energy() / static_cast<float>(trackster.vertex_multiplicity(k)));
82  input_Data_[0][index + 1] = static_cast<float>(std::abs(cluster.eta()));
83  input_Data_[0][index + 2] = static_cast<float>(cluster.phi());
84  seenClusters[j]++;
85  }
86  }
87  }
88  }
89 
90  // Method to run inference and update tracksters
91  void TracksterInferenceByCNNv4::runInference(std::vector<Trackster>& tracksters) {
92  if (batchSize_ == 0)
93  return; // Exit if no batch
94 
95  std::vector<std::vector<float>> outputTensors;
97  if (doPID_ and doRegression_) {
98  // Run energy model inference
99  if (!outputNames_.empty()) {
100  for (int i = 0; i < static_cast<int>(batchSize_); i++) {
101  const float energy = outputTensors[0][i];
102  tracksters[tracksterIndices_[i]].setRegressedEnergy(energy); // Update energy
103  }
104  }
105  }
106 
107  if (doPID_) {
108  // Run PID model inference
109  if (!outputNames_.empty()) {
110  int probsIdx = outputNames_.empty() ? 0 : 1;
111  std::vector<float> vec = outputTensors[probsIdx];
112  float* probs = vec.data();
113  for (int i = 0; i < batchSize_; i++) {
114  tracksters[tracksterIndices_[i]].setProbabilities(probs); // Update probabilities
115  probs += tracksters[tracksterIndices_[i]].id_probabilities().size(); // Move to next set of probabilities
116  }
117  }
118  }
119  }
120 
121  // Method to fill parameter set description for configuration
123  iDesc.add<int>("algo_verbosity", 0);
124  iDesc
125  .add<edm::FileInPath>("onnxModelPath",
126  edm::FileInPath("RecoHGCal/TICL/data/ticlv4/onnx_models/energy_id_v0.onnx"))
127  ->setComment("Path to ONNX PID model CLU3D");
128  iDesc.add<std::vector<std::string>>("inputNames", {"input:0"});
129  iDesc.add<std::vector<std::string>>("outputNames", {"output/regressed_energy:0", "output/id_probabilities:0"});
130  iDesc.add<double>("eid_min_cluster_energy", 1.0);
131  iDesc.add<int>("eid_n_layers", 50);
132  iDesc.add<int>("eid_n_clusters", 10);
133  iDesc.add<int>("doPID", 1);
134  iDesc.add<int>("doRegression", 0);
135  }
136 } // namespace ticl
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
const std::vector< std::string > outputNames_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const cms::Ort::ONNXRuntime * onnxSession_
void inputData(const std::vector< reco::CaloCluster > &layerClusters, std::vector< Trackster > &tracksters) override
std::vector< std::vector< int64_t > > input_shapes_
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
std::vector< std::vector< float > > input_Data_
double b
Definition: hdecay.h:120
HLT enums.
double a
Definition: hdecay.h:121
TracksterInferenceByCNNv4(const edm::ParameterSet &conf)
void runInference(std::vector< Trackster > &tracksters) override
Definition: Common.h:10
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
const std::vector< std::string > inputNames_
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