CMS 3D CMS Logo

PatternRecognitionbyFastJet.cc
Go to the documentation of this file.
1 // Author: Marco Rovere - marco.rovere@cern.ch
2 // Date: 10/2021
3 #include <algorithm>
4 #include <set>
5 #include <vector>
6 
7 #include "tbb/task_arena.h"
8 #include "tbb/tbb.h"
9 
17 
18 #include "TrackstersPCA.h"
22 
23 #include "fastjet/ClusterSequence.hh"
24 
25 using namespace ticl;
26 using namespace fastjet;
27 
28 template <typename TILES>
31  : PatternRecognitionAlgoBaseT<TILES>(conf, iC),
32  caloGeomToken_(iC.esConsumes<CaloGeometry, CaloGeometryRecord>()),
33  antikt_radius_(conf.getParameter<double>("antikt_radius")),
34  minNumLayerCluster_(conf.getParameter<int>("minNumLayerCluster")),
35  eidInputName_(conf.getParameter<std::string>("eid_input_name")),
36  eidOutputNameEnergy_(conf.getParameter<std::string>("eid_output_name_energy")),
37  eidOutputNameId_(conf.getParameter<std::string>("eid_output_name_id")),
38  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")),
39  eidNLayers_(conf.getParameter<int>("eid_n_layers")),
40  eidNClusters_(conf.getParameter<int>("eid_n_clusters")){};
41 
42 template <typename TILES>
43 void PatternRecognitionbyFastJet<TILES>::buildJetAndTracksters(std::vector<PseudoJet> &fjInputs,
44  std::vector<ticl::Trackster> &result) {
46  edm::LogVerbatim("PatternRecogntionbyFastJet")
47  << "Creating FastJet with " << fjInputs.size() << " LayerClusters in input";
48  }
49  fastjet::ClusterSequence sequence(fjInputs, JetDefinition(antikt_algorithm, antikt_radius_));
50  auto jets = fastjet::sorted_by_pt(sequence.inclusive_jets(0));
52  edm::LogVerbatim("PatternRecogntionbyFastJet") << "FastJet produced " << jets.size() << " jets/trackster";
53  }
54 
55  auto trackster_idx = result.size();
56  auto jetsSize = std::count_if(jets.begin(), jets.end(), [=](fastjet::PseudoJet jet) {
57  return jet.constituents().size() > static_cast<unsigned int>(minNumLayerCluster_);
58  });
59  result.resize(trackster_idx + jetsSize);
60 
61  for (const auto &pj : jets) {
62  if (pj.constituents().size() > static_cast<unsigned int>(minNumLayerCluster_)) {
63  for (const auto &component : pj.constituents()) {
64  result[trackster_idx].vertices().push_back(component.user_index());
65  result[trackster_idx].vertex_multiplicity().push_back(1);
67  edm::LogVerbatim("PatternRecogntionbyFastJet")
68  << "Jet has " << pj.constituents().size() << " components that are stored in trackster " << trackster_idx;
69  }
70  }
71  trackster_idx++;
72  } else {
74  edm::LogVerbatim("PatternRecogntionbyFastJet")
75  << "Jet with " << pj.constituents().size() << " constituents discarded since too small wrt "
76  << minNumLayerCluster_;
77  }
78  }
79  }
80  fjInputs.clear();
81 }
82 
83 template <typename TILES>
86  std::vector<Trackster> &result,
87  std::unordered_map<int, std::vector<int>> &seedToTracksterAssociation) {
88  // Protect from events with no seeding regions
89  if (input.regions.empty())
90  return;
91 
92  edm::EventSetup const &es = input.es;
93  const CaloGeometry &geom = es.getData(caloGeomToken_);
94  rhtools_.setGeometry(geom);
95 
96  constexpr auto isHFnose = std::is_same<TILES, TICLLayerTilesHFNose>::value;
97  constexpr int nEtaBin = TILES::constants_type_t::nEtaBins;
98  constexpr int nPhiBin = TILES::constants_type_t::nPhiBins;
99 
100  // We need to partition the two sides of the HGCAL detector
101  auto lastLayerPerSide = static_cast<unsigned int>(rhtools_.lastLayer(isHFnose)) - 1;
102  unsigned int maxLayer = 2 * lastLayerPerSide - 1;
103  std::vector<fastjet::PseudoJet> fjInputs;
104  fjInputs.clear();
105  for (unsigned int currentLayer = 0; currentLayer <= maxLayer; ++currentLayer) {
106  if (currentLayer == lastLayerPerSide) {
107  buildJetAndTracksters(fjInputs, result);
108  }
109  const auto &tileOnLayer = input.tiles[currentLayer];
110  for (int ieta = 0; ieta <= nEtaBin; ++ieta) {
111  auto offset = ieta * nPhiBin;
113  edm::LogVerbatim("PatternRecogntionbyFastJet") << "offset: " << offset;
114  }
115  for (int iphi = 0; iphi <= nPhiBin; ++iphi) {
117  edm::LogVerbatim("PatternRecogntionbyFastJet") << "iphi: " << iphi;
118  edm::LogVerbatim("PatternRecogntionbyFastJet") << "Entries in tileBin: " << tileOnLayer[offset + iphi].size();
119  }
120  for (auto clusterIdx : tileOnLayer[offset + iphi]) {
121  // Skip masked layer clusters
122  if (input.mask[clusterIdx] == 0.) {
124  edm::LogVerbatim("PatternRecogntionbyFastJet") << "Skipping masked layerIdx " << clusterIdx;
125  }
126  continue;
127  }
128  // Should we correct for the position of the PV?
129  auto const &cl = input.layerClusters[clusterIdx];
130  math::XYZVector direction(cl.x(), cl.y(), cl.z());
131  direction = direction.Unit();
132  direction *= cl.energy();
133  auto fpj = fastjet::PseudoJet(direction.X(), direction.Y(), direction.Z(), cl.energy());
134  fpj.set_user_index(clusterIdx);
135  fjInputs.push_back(fpj);
136  } // End of loop on the clusters on currentLayer
137  } // End of loop over phi-bin region
138  } // End of loop over eta-bin region
139  } // End of loop over layers
140 
141  // Collect the jet from the other side wrt to the one taken care of inside the main loop above.
142  buildJetAndTracksters(fjInputs, result);
143 
145  input.layerClusters,
146  input.layerClustersTime,
147  rhtools_.getPositionLayer(rhtools_.lastLayerEE(isHFnose), isHFnose).z());
148 
149  // run energy regression and ID
150  energyRegressionAndID(input.layerClusters, input.tfSession, result);
152  for (auto const &t : result) {
153  edm::LogVerbatim("PatternRecogntionbyFastJet") << "Barycenter: " << t.barycenter();
154  edm::LogVerbatim("PatternRecogntionbyFastJet") << "LCs: " << t.vertices().size();
155  edm::LogVerbatim("PatternRecogntionbyFastJet") << "Energy: " << t.raw_energy();
156  edm::LogVerbatim("PatternRecogntionbyFastJet") << "Regressed: " << t.regressed_energy();
157  }
158  }
159 }
160 
161 template <typename TILES>
163  const tensorflow::Session *eidSession,
164  std::vector<Trackster> &tracksters) {
165  // Energy regression and particle identification strategy:
166  //
167  // 1. Set default values for regressed energy and particle id for each trackster.
168  // 2. Store indices of tracksters whose total sum of cluster energies is above the
169  // eidMinClusterEnergy_ (GeV) treshold. Inference is not applied for soft tracksters.
170  // 3. When no trackster passes the selection, return.
171  // 4. Create input and output tensors. The batch dimension is determined by the number of
172  // selected tracksters.
173  // 5. Fill input tensors with layer cluster features. Per layer, clusters are ordered descending
174  // by energy. Given that tensor data is contiguous in memory, we can use pointer arithmetic to
175  // fill values, even with batching.
176  // 6. Zero-fill features for empty clusters in each layer.
177  // 7. Batched inference.
178  // 8. Assign the regressed energy and id probabilities to each trackster.
179  //
180  // Indices used throughout this method:
181  // i -> batch element / trackster
182  // j -> layer
183  // k -> cluster
184  // l -> feature
185 
186  // set default values per trackster, determine if the cluster energy threshold is passed,
187  // and store indices of hard tracksters
188  std::vector<int> tracksterIndices;
189  for (int i = 0; i < static_cast<int>(tracksters.size()); i++) {
190  // calculate the cluster energy sum (2)
191  // note: after the loop, sumClusterEnergy might be just above the threshold which is enough to
192  // decide whether to run inference for the trackster or not
193  float sumClusterEnergy = 0.;
194  for (const unsigned int &vertex : tracksters[i].vertices()) {
195  sumClusterEnergy += static_cast<float>(layerClusters[vertex].energy());
196  // there might be many clusters, so try to stop early
197  if (sumClusterEnergy >= eidMinClusterEnergy_) {
198  // set default values (1)
199  tracksters[i].setRegressedEnergy(0.f);
200  tracksters[i].zeroProbabilities();
201  tracksterIndices.push_back(i);
202  break;
203  }
204  }
205  }
206 
207  // do nothing when no trackster passes the selection (3)
208  int batchSize = static_cast<int>(tracksterIndices.size());
209  if (batchSize == 0) {
210  return;
211  }
212 
213  // create input and output tensors (4)
214  tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
215  tensorflow::Tensor input(tensorflow::DT_FLOAT, shape);
216  tensorflow::NamedTensorList inputList = {{eidInputName_, input}};
217 
218  std::vector<tensorflow::Tensor> outputs;
219  std::vector<std::string> outputNames;
220  if (!eidOutputNameEnergy_.empty()) {
221  outputNames.push_back(eidOutputNameEnergy_);
222  }
223  if (!eidOutputNameId_.empty()) {
224  outputNames.push_back(eidOutputNameId_);
225  }
226 
227  // fill input tensor (5)
228  for (int i = 0; i < batchSize; i++) {
229  const Trackster &trackster = tracksters[tracksterIndices[i]];
230 
231  // per layer, we only consider the first eidNClusters_ clusters in terms of energy, so in order
232  // to avoid creating large / nested structures to do the sorting for an unknown number of total
233  // clusters, create a sorted list of layer cluster indices to keep track of the filled clusters
234  std::vector<int> clusterIndices(trackster.vertices().size());
235  for (int k = 0; k < (int)trackster.vertices().size(); k++) {
236  clusterIndices[k] = k;
237  }
238  sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int &a, const int &b) {
239  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
240  });
241 
242  // keep track of the number of seen clusters per layer
243  std::vector<int> seenClusters(eidNLayers_);
244 
245  // loop through clusters by descending energy
246  for (const int &k : clusterIndices) {
247  // get features per layer and cluster and store the values directly in the input tensor
248  const reco::CaloCluster &cluster = layerClusters[trackster.vertices(k)];
249  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
250  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
251  // get the pointer to the first feature value for the current batch, layer and cluster
252  float *features = &input.tensor<float, 4>()(i, j, seenClusters[j], 0);
253 
254  // fill features
255  *(features++) = float(cluster.energy() / float(trackster.vertex_multiplicity(k)));
256  *(features++) = float(std::abs(cluster.eta()));
257  *(features) = float(cluster.phi());
258 
259  // increment seen clusters
260  seenClusters[j]++;
261  }
262  }
263 
264  // zero-fill features of empty clusters in each layer (6)
265  for (int j = 0; j < eidNLayers_; j++) {
266  for (int k = seenClusters[j]; k < eidNClusters_; k++) {
267  float *features = &input.tensor<float, 4>()(i, j, k, 0);
268  for (int l = 0; l < eidNFeatures_; l++) {
269  *(features++) = 0.f;
270  }
271  }
272  }
273  }
274 
275  // run the inference (7)
276  tensorflow::run(const_cast<tensorflow::Session *>(eidSession), inputList, outputNames, &outputs);
277 
278  // store regressed energy per trackster (8)
279  if (!eidOutputNameEnergy_.empty()) {
280  // get the pointer to the energy tensor, dimension is batch x 1
281  float *energy = outputs[0].flat<float>().data();
282 
283  for (const int &i : tracksterIndices) {
284  tracksters[i].setRegressedEnergy(*(energy++));
285  }
286  }
287 
288  // store id probabilities per trackster (8)
289  if (!eidOutputNameId_.empty()) {
290  // get the pointer to the id probability tensor, dimension is batch x id_probabilities.size()
291  int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
292  float *probs = outputs[probsIdx].flat<float>().data();
293 
294  for (const int &i : tracksterIndices) {
295  tracksters[i].setProbabilities(probs);
296  probs += tracksters[i].id_probabilities().size();
297  }
298  }
299 }
300 
301 template <typename TILES>
303  iDesc.add<int>("algo_verbosity", 0);
304  iDesc.add<double>("antikt_radius", 0.09)->setComment("Radius to be used while running the Anti-kt clustering");
305  iDesc.add<int>("minNumLayerCluster", 5)->setComment("Not Inclusive");
306  iDesc.add<std::string>("eid_input_name", "input");
307  iDesc.add<std::string>("eid_output_name_energy", "output/regressed_energy");
308  iDesc.add<std::string>("eid_output_name_id", "output/id_probabilities");
309  iDesc.add<double>("eid_min_cluster_energy", 1.);
310  iDesc.add<int>("eid_n_layers", 50);
311  iDesc.add<int>("eid_n_clusters", 10);
312 }
313 
Log< level::Info, true > LogVerbatim
void setComment(std::string const &value)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:30
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:184
static std::string const input
Definition: EdmProvDump.cc:50
void assignPCAtoTracksters(std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, const edm::ValueMap< std::pair< float, float >> &, double, bool energyWeight=true)
PatternRecognitionbyFastJet(const edm::ParameterSet &conf, edm::ConsumesCollector)
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:213
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
bool getData(T &iHolder) const
Definition: EventSetup.h:122
string inputList
Definition: crabTemplate.py:6
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void makeTracksters(const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
double energy() const
cluster energy
Definition: CaloCluster.h:149
void buildJetAndTracksters(std::vector< fastjet::PseudoJet > &, std::vector< ticl::Trackster > &)
std::vector< unsigned int > & vertices()
Definition: Trackster.h:57
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
std::vector< float > & vertex_multiplicity()
Definition: Trackster.h:58
double b
Definition: hdecay.h:118
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:79
double a
Definition: hdecay.h:119
Definition: Common.h:8
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181