CMS 3D CMS Logo

PatternRecognitionbyCA.cc
Go to the documentation of this file.
1 // Author: Felice Pantaleo, Marco Rovere - felice.pantaleo@cern.ch, marco.rovere@cern.ch
2 // Date: 11/2018
3 #include <algorithm>
4 #include <set>
5 #include <vector>
6 
10 
11 #include "TrackstersPCA.h"
15 
16 using namespace ticl;
17 
18 template <typename TILES>
20  : PatternRecognitionAlgoBaseT<TILES>(conf, iC),
21  caloGeomToken_(iC.esConsumes<CaloGeometry, CaloGeometryRecord>()),
22  theGraph_(std::make_unique<HGCGraphT<TILES>>()),
23  oneTracksterPerTrackSeed_(conf.getParameter<bool>("oneTracksterPerTrackSeed")),
24  promoteEmptyRegionToTrackster_(conf.getParameter<bool>("promoteEmptyRegionToTrackster")),
25  out_in_dfs_(conf.getParameter<bool>("out_in_dfs")),
26  max_out_in_hops_(conf.getParameter<int>("max_out_in_hops")),
27  min_cos_theta_(conf.getParameter<double>("min_cos_theta")),
28  min_cos_pointing_(conf.getParameter<double>("min_cos_pointing")),
29  root_doublet_max_distance_from_seed_squared_(
30  conf.getParameter<double>("root_doublet_max_distance_from_seed_squared")),
31  etaLimitIncreaseWindow_(conf.getParameter<double>("etaLimitIncreaseWindow")),
32  skip_layers_(conf.getParameter<int>("skip_layers")),
33  max_missing_layers_in_trackster_(conf.getParameter<int>("max_missing_layers_in_trackster")),
34  check_missing_layers_(max_missing_layers_in_trackster_ < 100),
35  shower_start_max_layer_(conf.getParameter<int>("shower_start_max_layer")),
36  min_layers_per_trackster_(conf.getParameter<int>("min_layers_per_trackster")),
37  filter_on_categories_(conf.getParameter<std::vector<int>>("filter_on_categories")),
38  pid_threshold_(conf.getParameter<double>("pid_threshold")),
39  energy_em_over_total_threshold_(conf.getParameter<double>("energy_em_over_total_threshold")),
40  max_longitudinal_sigmaPCA_(conf.getParameter<double>("max_longitudinal_sigmaPCA")),
41  min_clusters_per_ntuplet_(min_layers_per_trackster_),
42  max_delta_time_(conf.getParameter<double>("max_delta_time")),
43  eidInputName_(conf.getParameter<std::string>("eid_input_name")),
44  eidOutputNameEnergy_(conf.getParameter<std::string>("eid_output_name_energy")),
45  eidOutputNameId_(conf.getParameter<std::string>("eid_output_name_id")),
46  eidMinClusterEnergy_(conf.getParameter<double>("eid_min_cluster_energy")),
47  eidNLayers_(conf.getParameter<int>("eid_n_layers")),
48  eidNClusters_(conf.getParameter<int>("eid_n_clusters")),
49  siblings_maxRSquared_(conf.getParameter<std::vector<double>>("siblings_maxRSquared")){};
50 
51 template <typename TILES>
53 
54 template <typename TILES>
57  std::vector<Trackster> &result,
58  std::unordered_map<int, std::vector<int>> &seedToTracksterAssociation) {
59  // Protect from events with no seeding regions
60  if (input.regions.empty())
61  return;
62 
63  edm::EventSetup const &es = input.es;
64  const CaloGeometry &geom = es.getData(caloGeomToken_);
65  rhtools_.setGeometry(geom);
66 
68  theGraph_->clear();
70  LogDebug("HGCPatternRecoByCA") << "Making Tracksters with CA" << std::endl;
71  }
72 
76 
77  bool isRegionalIter = (input.regions[0].index != -1);
78  std::vector<HGCDoublet::HGCntuplet> foundNtuplets;
79  std::vector<int> seedIndices;
80  std::vector<uint8_t> layer_cluster_usage(input.layerClusters.size(), 0);
81  theGraph_->makeAndConnectDoublets(input.tiles,
82  input.regions,
83  nEtaBin,
84  nPhiBin,
85  input.layerClusters,
86  input.mask,
87  input.layerClustersTime,
88  1,
89  1,
90  min_cos_theta_,
91  min_cos_pointing_,
92  root_doublet_max_distance_from_seed_squared_,
93  etaLimitIncreaseWindow_,
94  skip_layers_,
95  rhtools_.lastLayer(isHFnose),
96  max_delta_time_,
97  rhtools_.lastLayerEE(isHFnose),
98  rhtools_.lastLayerFH(),
99  siblings_maxRSquared_);
100 
101  theGraph_->findNtuplets(foundNtuplets, seedIndices, min_clusters_per_ntuplet_, out_in_dfs_, max_out_in_hops_);
102  //#ifdef FP_DEBUG
103  const auto &doublets = theGraph_->getAllDoublets();
104  int tracksterId = -1;
105 
106  // container for holding tracksters before selection
107  std::vector<Trackster> tmpTracksters;
108  tmpTracksters.reserve(foundNtuplets.size());
109 
110  for (auto const &ntuplet : foundNtuplets) {
111  tracksterId++;
112 
113  std::set<unsigned int> effective_cluster_idx;
114 
115  for (auto const &doublet : ntuplet) {
116  auto innerCluster = doublets[doublet].innerClusterId();
117  auto outerCluster = doublets[doublet].outerClusterId();
118 
119  effective_cluster_idx.insert(innerCluster);
120  effective_cluster_idx.insert(outerCluster);
121 
123  LogDebug("HGCPatternRecoByCA") << " New doublet " << doublet << " for trackster: " << result.size()
124  << " InnerCl " << innerCluster << " " << input.layerClusters[innerCluster].x()
125  << " " << input.layerClusters[innerCluster].y() << " "
126  << input.layerClusters[innerCluster].z() << " OuterCl " << outerCluster << " "
127  << input.layerClusters[outerCluster].x() << " "
128  << input.layerClusters[outerCluster].y() << " "
129  << input.layerClusters[outerCluster].z() << " " << tracksterId << std::endl;
130  }
131  }
132  unsigned showerMinLayerId = 99999;
133  std::vector<unsigned int> uniqueLayerIds;
134  uniqueLayerIds.reserve(effective_cluster_idx.size());
135  std::vector<std::pair<unsigned int, unsigned int>> lcIdAndLayer;
136  lcIdAndLayer.reserve(effective_cluster_idx.size());
137  for (auto const i : effective_cluster_idx) {
138  auto const &haf = input.layerClusters[i].hitsAndFractions();
139  auto layerId = rhtools_.getLayerWithOffset(haf[0].first);
140  showerMinLayerId = std::min(layerId, showerMinLayerId);
141  uniqueLayerIds.push_back(layerId);
142  lcIdAndLayer.emplace_back(i, layerId);
143  }
144  std::sort(uniqueLayerIds.begin(), uniqueLayerIds.end());
145  uniqueLayerIds.erase(std::unique(uniqueLayerIds.begin(), uniqueLayerIds.end()), uniqueLayerIds.end());
146  unsigned int numberOfLayersInTrackster = uniqueLayerIds.size();
147  if (check_missing_layers_) {
148  int numberOfMissingLayers = 0;
149  unsigned int j = showerMinLayerId;
150  unsigned int indexInVec = 0;
151  for (const auto &layer : uniqueLayerIds) {
152  if (layer != j) {
153  numberOfMissingLayers++;
154  j++;
155  if (numberOfMissingLayers > max_missing_layers_in_trackster_) {
156  numberOfLayersInTrackster = indexInVec;
157  for (auto &llpair : lcIdAndLayer) {
158  if (llpair.second >= layer) {
159  effective_cluster_idx.erase(llpair.first);
160  }
161  }
162  break;
163  }
164  }
165  indexInVec++;
166  j++;
167  }
168  }
169  if ((numberOfLayersInTrackster >= min_layers_per_trackster_) and (showerMinLayerId <= shower_start_max_layer_)) {
170  // Put back indices, in the form of a Trackster, into the results vector
171  Trackster tmp;
172  tmp.vertices().reserve(effective_cluster_idx.size());
173  tmp.vertex_multiplicity().resize(effective_cluster_idx.size(), 1);
174  //regions and seedIndices can have different size
175  //if a seeding region does not lead to any trackster
176  tmp.setSeed(input.regions[0].collectionID, seedIndices[tracksterId]);
177 
178  std::copy(std::begin(effective_cluster_idx), std::end(effective_cluster_idx), std::back_inserter(tmp.vertices()));
179  tmpTracksters.push_back(tmp);
180  }
181  }
182  ticl::assignPCAtoTracksters(tmpTracksters,
183  input.layerClusters,
184  input.layerClustersTime,
185  rhtools_.getPositionLayer(rhtools_.lastLayerEE(isHFnose), isHFnose).z());
186 
187  // run energy regression and ID
188  energyRegressionAndID(input.layerClusters, input.tfSession, tmpTracksters);
189  // Filter results based on PID criteria or EM/Total energy ratio.
190  // We want to **keep** tracksters whose cumulative
191  // probability summed up over the selected categories
192  // is greater than the chosen threshold. Therefore
193  // the filtering function should **discard** all
194  // tracksters **below** the threshold.
195  auto filter_on_pids = [&](Trackster &t) -> bool {
196  auto cumulative_prob = 0.;
197  for (auto index : filter_on_categories_) {
198  cumulative_prob += t.id_probabilities(index);
199  }
200  return (cumulative_prob <= pid_threshold_) &&
201  (t.raw_em_energy() < energy_em_over_total_threshold_ * t.raw_energy());
202  };
203 
204  std::vector<unsigned int> selectedTrackstersIds;
205  for (unsigned i = 0; i < tmpTracksters.size(); ++i) {
206  if (!filter_on_pids(tmpTracksters[i]) and tmpTracksters[i].sigmasPCA()[0] < max_longitudinal_sigmaPCA_) {
207  selectedTrackstersIds.push_back(i);
208  }
209  }
210 
211  result.reserve(selectedTrackstersIds.size());
212 
213  for (unsigned i = 0; i < selectedTrackstersIds.size(); ++i) {
214  const auto &t = tmpTracksters[selectedTrackstersIds[i]];
215  for (auto const lcId : t.vertices()) {
216  layer_cluster_usage[lcId]++;
218  LogDebug("HGCPatternRecoByCA") << "LayerID: " << lcId << " count: " << (int)layer_cluster_usage[lcId]
219  << std::endl;
220  }
221  if (isRegionalIter) {
222  seedToTracksterAssociation[t.seedIndex()].push_back(i);
223  }
224  result.push_back(t);
225  }
226 
227  for (auto &trackster : result) {
228  assert(trackster.vertices().size() <= trackster.vertex_multiplicity().size());
229  for (size_t i = 0; i < trackster.vertices().size(); ++i) {
230  trackster.vertex_multiplicity()[i] = layer_cluster_usage[trackster.vertices(i)];
232  LogDebug("HGCPatternRecoByCA") << "LayerID: " << trackster.vertices(i)
233  << " count: " << (int)trackster.vertex_multiplicity(i) << std::endl;
234  }
235  }
236  // Now decide if the tracksters from the track-based iterations have to be merged
237  if (oneTracksterPerTrackSeed_) {
238  std::vector<Trackster> tmp;
239  mergeTrackstersTRK(result, input.layerClusters, tmp, seedToTracksterAssociation);
240  tmp.swap(result);
241  }
242 
244  input.layerClusters,
245  input.layerClustersTime,
246  rhtools_.getPositionLayer(rhtools_.lastLayerEE(isHFnose), isHFnose).z());
247 
248  // run energy regression and ID
249  energyRegressionAndID(input.layerClusters, input.tfSession, result);
250 
251  // now adding dummy tracksters from seeds not connected to any shower in the result collection
252  // these are marked as charged hadrons with probability 1.
253  if (promoteEmptyRegionToTrackster_) {
254  emptyTrackstersFromSeedsTRK(result, seedToTracksterAssociation, input.regions[0].collectionID);
255  }
256 
258  for (auto &trackster : result) {
259  LogDebug("HGCPatternRecoByCA") << "Trackster characteristics: " << std::endl;
260  LogDebug("HGCPatternRecoByCA") << "Size: " << trackster.vertices().size() << std::endl;
261  auto counter = 0;
262  for (auto const &p : trackster.id_probabilities()) {
263  LogDebug("HGCPatternRecoByCA") << counter++ << ": " << p << std::endl;
264  }
265  }
266  }
267  theGraph_->clear();
268 }
269 
270 template <typename TILES>
272  const std::vector<Trackster> &input,
273  const std::vector<reco::CaloCluster> &layerClusters,
274  std::vector<Trackster> &output,
275  std::unordered_map<int, std::vector<int>> &seedToTracksterAssociation) const {
276  output.reserve(input.size());
277  for (auto &thisSeed : seedToTracksterAssociation) {
278  auto &tracksters = thisSeed.second;
279  if (!tracksters.empty()) {
280  auto numberOfTrackstersInSeed = tracksters.size();
281  output.emplace_back(input[tracksters[0]]);
282  auto &outTrackster = output.back();
283  tracksters[0] = output.size() - 1;
284  auto updated_size = outTrackster.vertices().size();
285  for (unsigned int j = 1; j < numberOfTrackstersInSeed; ++j) {
286  auto &thisTrackster = input[tracksters[j]];
287  updated_size += thisTrackster.vertices().size();
289  LogDebug("HGCPatternRecoByCA") << "Updated size: " << updated_size << std::endl;
290  }
291  outTrackster.vertices().reserve(updated_size);
292  outTrackster.vertex_multiplicity().reserve(updated_size);
293  std::copy(std::begin(thisTrackster.vertices()),
294  std::end(thisTrackster.vertices()),
295  std::back_inserter(outTrackster.vertices()));
296  std::copy(std::begin(thisTrackster.vertex_multiplicity()),
297  std::end(thisTrackster.vertex_multiplicity()),
298  std::back_inserter(outTrackster.vertex_multiplicity()));
299  }
300  tracksters.resize(1);
301 
302  // Find duplicate LCs
303  auto &orig_vtx = outTrackster.vertices();
304  auto vtx_sorted{orig_vtx};
305  std::sort(std::begin(vtx_sorted), std::end(vtx_sorted));
306  for (unsigned int iLC = 1; iLC < vtx_sorted.size(); ++iLC) {
307  if (vtx_sorted[iLC] == vtx_sorted[iLC - 1]) {
308  // Clean up duplicate LCs
309  const auto lcIdx = vtx_sorted[iLC];
310  const auto firstEl = std::find(orig_vtx.begin(), orig_vtx.end(), lcIdx);
311  const auto firstPos = std::distance(std::begin(orig_vtx), firstEl);
312  auto iDup = std::find(std::next(firstEl), orig_vtx.end(), lcIdx);
313  while (iDup != orig_vtx.end()) {
314  orig_vtx.erase(iDup);
315  outTrackster.vertex_multiplicity().erase(outTrackster.vertex_multiplicity().begin() +
316  std::distance(std::begin(orig_vtx), iDup));
317  outTrackster.vertex_multiplicity()[firstPos] -= 1;
318  iDup = std::find(std::next(firstEl), orig_vtx.end(), lcIdx);
319  };
320  }
321  }
322  }
323  }
324  output.shrink_to_fit();
325 }
326 
327 template <typename TILES>
329  std::vector<Trackster> &tracksters,
330  std::unordered_map<int, std::vector<int>> &seedToTracksterAssociation,
331  const edm::ProductID &collectionID) const {
332  for (auto &thisSeed : seedToTracksterAssociation) {
333  if (thisSeed.second.empty()) {
334  Trackster t;
335  t.setRegressedEnergy(0.f);
336  t.zeroProbabilities();
338  t.setSeed(collectionID, thisSeed.first);
339  tracksters.emplace_back(t);
340  thisSeed.second.emplace_back(tracksters.size() - 1);
341  }
342  }
343 }
344 
345 template <typename TILES>
346 void PatternRecognitionbyCA<TILES>::energyRegressionAndID(const std::vector<reco::CaloCluster> &layerClusters,
347  const tensorflow::Session *eidSession,
348  std::vector<Trackster> &tracksters) {
349  // Energy regression and particle identification strategy:
350  //
351  // 1. Set default values for regressed energy and particle id for each trackster.
352  // 2. Store indices of tracksters whose total sum of cluster energies is above the
353  // eidMinClusterEnergy_ (GeV) threshold. Inference is not applied for soft tracksters.
354  // 3. When no trackster passes the selection, return.
355  // 4. Create input and output tensors. The batch dimension is determined by the number of
356  // selected tracksters.
357  // 5. Fill input tensors with layer cluster features. Per layer, clusters are ordered descending
358  // by energy. Given that tensor data is contiguous in memory, we can use pointer arithmetic to
359  // fill values, even with batching.
360  // 6. Zero-fill features for empty clusters in each layer.
361  // 7. Batched inference.
362  // 8. Assign the regressed energy and id probabilities to each trackster.
363  //
364  // Indices used throughout this method:
365  // i -> batch element / trackster
366  // j -> layer
367  // k -> cluster
368  // l -> feature
369 
370  // set default values per trackster, determine if the cluster energy threshold is passed,
371  // and store indices of hard tracksters
372  std::vector<int> tracksterIndices;
373  for (int i = 0; i < (int)tracksters.size(); i++) {
374  // calculate the cluster energy sum (2)
375  // note: after the loop, sumClusterEnergy might be just above the threshold which is enough to
376  // decide whether to run inference for the trackster or not
377  float sumClusterEnergy = 0.;
378  for (const unsigned int &vertex : tracksters[i].vertices()) {
379  sumClusterEnergy += (float)layerClusters[vertex].energy();
380  // there might be many clusters, so try to stop early
381  if (sumClusterEnergy >= eidMinClusterEnergy_) {
382  // set default values (1)
383  tracksters[i].setRegressedEnergy(0.f);
384  tracksters[i].zeroProbabilities();
385  tracksterIndices.push_back(i);
386  break;
387  }
388  }
389  }
390 
391  // do nothing when no trackster passes the selection (3)
392  int batchSize = (int)tracksterIndices.size();
393  if (batchSize == 0) {
394  return;
395  }
396 
397  // create input and output tensors (4)
398  tensorflow::TensorShape shape({batchSize, eidNLayers_, eidNClusters_, eidNFeatures_});
399  tensorflow::Tensor input(tensorflow::DT_FLOAT, shape);
400  tensorflow::NamedTensorList inputList = {{eidInputName_, input}};
401 
402  std::vector<tensorflow::Tensor> outputs;
403  std::vector<std::string> outputNames;
404  if (!eidOutputNameEnergy_.empty()) {
405  outputNames.push_back(eidOutputNameEnergy_);
406  }
407  if (!eidOutputNameId_.empty()) {
408  outputNames.push_back(eidOutputNameId_);
409  }
410 
411  // fill input tensor (5)
412  for (int i = 0; i < batchSize; i++) {
413  const Trackster &trackster = tracksters[tracksterIndices[i]];
414 
415  // per layer, we only consider the first eidNClusters_ clusters in terms of energy, so in order
416  // to avoid creating large / nested structures to do the sorting for an unknown number of total
417  // clusters, create a sorted list of layer cluster indices to keep track of the filled clusters
418  std::vector<int> clusterIndices(trackster.vertices().size());
419  for (int k = 0; k < (int)trackster.vertices().size(); k++) {
420  clusterIndices[k] = k;
421  }
422  sort(clusterIndices.begin(), clusterIndices.end(), [&layerClusters, &trackster](const int &a, const int &b) {
423  return layerClusters[trackster.vertices(a)].energy() > layerClusters[trackster.vertices(b)].energy();
424  });
425 
426  // keep track of the number of seen clusters per layer
427  std::vector<int> seenClusters(eidNLayers_);
428 
429  // loop through clusters by descending energy
430  for (const int &k : clusterIndices) {
431  // get features per layer and cluster and store the values directly in the input tensor
432  const reco::CaloCluster &cluster = layerClusters[trackster.vertices(k)];
433  int j = rhtools_.getLayerWithOffset(cluster.hitsAndFractions()[0].first) - 1;
434  if (j < eidNLayers_ && seenClusters[j] < eidNClusters_) {
435  // get the pointer to the first feature value for the current batch, layer and cluster
436  float *features = &input.tensor<float, 4>()(i, j, seenClusters[j], 0);
437 
438  // fill features
439  *(features++) = float(cluster.energy() / float(trackster.vertex_multiplicity(k)));
440  *(features++) = float(std::abs(cluster.eta()));
441  *(features) = float(cluster.phi());
442 
443  // increment seen clusters
444  seenClusters[j]++;
445  }
446  }
447 
448  // zero-fill features of empty clusters in each layer (6)
449  for (int j = 0; j < eidNLayers_; j++) {
450  for (int k = seenClusters[j]; k < eidNClusters_; k++) {
451  float *features = &input.tensor<float, 4>()(i, j, k, 0);
452  for (int l = 0; l < eidNFeatures_; l++) {
453  *(features++) = 0.f;
454  }
455  }
456  }
457  }
458 
459  // run the inference (7)
461 
462  // store regressed energy per trackster (8)
463  if (!eidOutputNameEnergy_.empty()) {
464  // get the pointer to the energy tensor, dimension is batch x 1
465  float *energy = outputs[0].flat<float>().data();
466 
467  for (const int &i : tracksterIndices) {
468  tracksters[i].setRegressedEnergy(*(energy++));
469  }
470  }
471 
472  // store id probabilities per trackster (8)
473  if (!eidOutputNameId_.empty()) {
474  // get the pointer to the id probability tensor, dimension is batch x id_probabilities.size()
475  int probsIdx = eidOutputNameEnergy_.empty() ? 0 : 1;
476  float *probs = outputs[probsIdx].flat<float>().data();
477 
478  for (const int &i : tracksterIndices) {
479  tracksters[i].setProbabilities(probs);
480  probs += tracksters[i].id_probabilities().size();
481  }
482  }
483 }
484 
485 template <typename TILES>
487  iDesc.add<int>("algo_verbosity", 0);
488  iDesc.add<bool>("oneTracksterPerTrackSeed", false);
489  iDesc.add<bool>("promoteEmptyRegionToTrackster", false);
490  iDesc.add<bool>("out_in_dfs", true);
491  iDesc.add<int>("max_out_in_hops", 10);
492  iDesc.add<double>("min_cos_theta", 0.915);
493  iDesc.add<double>("min_cos_pointing", -1.);
494  iDesc.add<double>("root_doublet_max_distance_from_seed_squared", 9999);
495  iDesc.add<double>("etaLimitIncreaseWindow", 2.1);
496  iDesc.add<int>("skip_layers", 0);
497  iDesc.add<int>("max_missing_layers_in_trackster", 9999);
498  iDesc.add<int>("shower_start_max_layer", 9999)->setComment("make default such that no filtering is applied");
499  iDesc.add<int>("min_layers_per_trackster", 10);
500  iDesc.add<std::vector<int>>("filter_on_categories", {0});
501  iDesc.add<double>("pid_threshold", 0.)->setComment("make default such that no filtering is applied");
502  iDesc.add<double>("energy_em_over_total_threshold", -1.)
503  ->setComment("make default such that no filtering is applied");
504  iDesc.add<double>("max_longitudinal_sigmaPCA", 9999);
505  iDesc.add<double>("max_delta_time", 3.)->setComment("nsigma");
506  iDesc.add<std::string>("eid_input_name", "input");
507  iDesc.add<std::string>("eid_output_name_energy", "output/regressed_energy");
508  iDesc.add<std::string>("eid_output_name_id", "output/id_probabilities");
509  iDesc.add<double>("eid_min_cluster_energy", 1.);
510  iDesc.add<int>("eid_n_layers", 50);
511  iDesc.add<int>("eid_n_clusters", 10);
512  iDesc.add<std::vector<double>>("siblings_maxRSquared", {6e-4, 6e-4, 6e-4});
513 }
514 
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:31
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void mergeTrackstersTRK(const std::vector< Trackster > &, const std::vector< reco::CaloCluster > &, std::vector< Trackster > &, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) const
void emptyTrackstersFromSeedsTRK(std::vector< Trackster > &tracksters, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation, const edm::ProductID &collectionID) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:183
void makeTracksters(const typename PatternRecognitionAlgoBaseT< TILES >::Inputs &input, std::vector< Trackster > &result, std::unordered_map< int, std::vector< int >> &seedToTracksterAssociation) override
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)
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
def unique(seq, keepstr=True)
Definition: tier0.py:24
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:281
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
string inputList
Definition: crabTemplate.py:6
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double energy() const
cluster energy
Definition: CaloCluster.h:148
std::vector< unsigned int > & vertices()
Definition: Trackster.h:56
void energyRegressionAndID(const std::vector< reco::CaloCluster > &layerClusters, const tensorflow::Session *, std::vector< Trackster > &result)
PatternRecognitionbyCA(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
std::vector< float > & vertex_multiplicity()
Definition: Trackster.h:57
double b
Definition: hdecay.h:120
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
double a
Definition: hdecay.h:121
Definition: Common.h:8
Definition: output.py:1
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:180
tmp
align.sh
Definition: createJobs.py:716
#define LogDebug(id)