CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HGCalCLUEAlgoT< TILE > Class Template Reference

#include <HGCalCLUEAlgo.h>

Inheritance diagram for HGCalCLUEAlgoT< TILE >:
HGCalClusteringAlgoBase

Classes

struct  CellsOnLayer
 

Public Types

typedef math::XYZPoint Point
 point in the space More...
 
- Public Types inherited from HGCalClusteringAlgoBase
enum  VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 }
 

Public Member Functions

void computeThreshold ()
 
std::vector< reco::BasicClustergetClusters (bool) override
 
Density getDensity () override
 
void getEventSetupPerAlgorithm (const edm::EventSetup &es) override
 
 HGCalCLUEAlgoT (const edm::ParameterSet &ps, edm::ConsumesCollector iC)
 
void makeClusters () override
 
void populate (const HGCRecHitCollection &hits) override
 
void reset () override
 
 ~HGCalCLUEAlgoT () override
 
- Public Member Functions inherited from HGCalClusteringAlgoBase
void getEventSetup (const edm::EventSetup &es)
 
 HGCalClusteringAlgoBase (VerbosityLevel v, reco::CaloCluster::AlgoId algo, edm::ConsumesCollector iC)
 
void setAlgoId (reco::CaloCluster::AlgoId algo, bool isNose=false)
 
void setVerbosity (VerbosityLevel the_verbosity)
 
virtual ~HGCalClusteringAlgoBase ()
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &iDesc)
 

Private Member Functions

void calculateDistanceToHigher (const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
 
void calculateLocalDensity (const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
 
math::XYZPoint calculatePosition (const std::vector< int > &v, const unsigned int layerId) const
 
float distance (int cell1, int cell2, int layerId, bool isEtaPhi) const
 
float distance2 (int cell1, int cell2, int layerId, bool isEtaPhi) const
 
int findAndAssignClusters (const unsigned int layerId, float delta_c, float delta_r)
 
void prepareDataStructures (const unsigned int layerId)
 
void setDensity (const unsigned int layerId)
 

Private Attributes

std::vector< CellsOnLayercells_
 
std::vector< double > dEdXweights_
 
int deltasi_index_regemfac_
 
Density density_
 
bool dependSensor_
 
double ecut_
 
double fcPerEle_
 
std::vector< double > fcPerMip_
 
bool initialized_
 
double kappa_
 
unsigned maxNumberOfThickIndices_
 
double noiseMip_
 
std::vector< double > nonAgedNoises_
 
std::vector< int > numberOfClustersPerLayer_
 
float outlierDeltaFactor_ = 2.f
 
const double positionDeltaRho2_
 
double sciThicknessCorrection_
 
std::vector< double > thicknessCorrection_
 
std::vector< std::vector< double > > thresholds_
 
std::vector< double > thresholdW0_
 
bool use2x2_
 
std::vector< std::vector< double > > v_sigmaNoise_
 
std::vector< double > vecDeltas_
 

Additional Inherited Members

- Public Attributes inherited from HGCalClusteringAlgoBase
unsigned int firstLayerBH_
 
bool isNose_
 
unsigned int lastLayerEE_
 
unsigned int lastLayerFH_
 
unsigned int maxlayer_
 
int scintMaxIphi_
 
- Protected Attributes inherited from HGCalClusteringAlgoBase
reco::CaloCluster::AlgoId algoId_
 
edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeomToken_
 
std::vector< reco::BasicClusterclusters_v_
 
hgcal::RecHitTools rhtools_
 
VerbosityLevel verbosity_
 

Detailed Description

template<typename TILE>
class HGCalCLUEAlgoT< TILE >

Definition at line 30 of file HGCalCLUEAlgo.h.

Member Typedef Documentation

◆ Point

template<typename TILE >
typedef math::XYZPoint HGCalCLUEAlgoT< TILE >::Point

point in the space

Definition at line 123 of file HGCalCLUEAlgo.h.

Constructor & Destructor Documentation

◆ HGCalCLUEAlgoT()

template<typename TILE >
HGCalCLUEAlgoT< TILE >::HGCalCLUEAlgoT ( const edm::ParameterSet ps,
edm::ConsumesCollector  iC 
)
inline

Definition at line 32 of file HGCalCLUEAlgo.h.

34  (HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3),
36  iC),
37  thresholdW0_(ps.getParameter<std::vector<double>>("thresholdW0")),
38  positionDeltaRho2_(ps.getParameter<double>("positionDeltaRho2")),
39  vecDeltas_(ps.getParameter<std::vector<double>>("deltac")),
40  kappa_(ps.getParameter<double>("kappa")),
41  ecut_(ps.getParameter<double>("ecut")),
42  dependSensor_(ps.getParameter<bool>("dependSensor")),
43  dEdXweights_(ps.getParameter<std::vector<double>>("dEdXweights")),
44  thicknessCorrection_(ps.getParameter<std::vector<double>>("thicknessCorrection")),
45  sciThicknessCorrection_(ps.getParameter<double>("sciThicknessCorrection")),
46  deltasi_index_regemfac_(ps.getParameter<int>("deltasi_index_regemfac")),
47  maxNumberOfThickIndices_(ps.getParameter<unsigned>("maxNumberOfThickIndices")),
48  fcPerMip_(ps.getParameter<std::vector<double>>("fcPerMip")),
49  fcPerEle_(ps.getParameter<double>("fcPerEle")),
50  nonAgedNoises_(ps.getParameter<std::vector<double>>("noises")),
51  noiseMip_(ps.getParameter<edm::ParameterSet>("noiseMip").getParameter<double>("noise_MIP")),
52  use2x2_(ps.getParameter<bool>("use2x2")),
53  initialized_(false) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const double positionDeltaRho2_
HGCalClusteringAlgoBase(VerbosityLevel v, reco::CaloCluster::AlgoId algo, edm::ConsumesCollector iC)
std::vector< double > thresholdW0_
unsigned maxNumberOfThickIndices_
int deltasi_index_regemfac_
std::vector< double > dEdXweights_
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > vecDeltas_
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_
std::vector< double > thicknessCorrection_
double sciThicknessCorrection_

◆ ~HGCalCLUEAlgoT()

template<typename TILE >
HGCalCLUEAlgoT< TILE >::~HGCalCLUEAlgoT ( )
inlineoverride

Definition at line 55 of file HGCalCLUEAlgo.h.

55 {}

Member Function Documentation

◆ calculateDistanceToHigher()

template<typename TILE >
void HGCalCLUEAlgoT< T >::calculateDistanceToHigher ( const TILE &  lt,
const unsigned int  layerId,
float  delta_c,
float  delta_r 
)
private

Definition at line 346 of file HGCalCLUEAlgo.cc.

References dumpMFGeometry_cfg::delta, hitfit::delta_r(), HLT_2022v15_cff::distance, mps_fire::i, dqmiolumiharvest::j, LogDebug, SiStripPI::max, allConversions_cfi::maxDelta, FastTimerService_cff::range, photonAnalyzer_cfi::xBin, and photonAnalyzer_cfi::yBin.

349  {
350  auto& cellsOnLayer = cells_[layerId];
351  unsigned int numberOfCells = cellsOnLayer.detid.size();
352  bool isOnlySi(false);
353  if (rhtools_.isOnlySilicon(layerId))
354  isOnlySi = true;
355 
356  for (unsigned int i = 0; i < numberOfCells; i++) {
357  bool isSi = isOnlySi || cellsOnLayer.isSi[i];
358  // initialize delta and nearest higher for i
360  float i_delta = maxDelta;
361  int i_nearestHigher = -1;
362  if (isSi) {
363  float delta = delta_c;
364  // get search box for ith hit
365  // guarantee to cover a range "outlierDeltaFactor_*delta_c"
367  std::array<int, 4> search_box = lt.searchBox(
368  cellsOnLayer.x[i] - range, cellsOnLayer.x[i] + range, cellsOnLayer.y[i] - range, cellsOnLayer.y[i] + range);
369  // loop over all bins in the search box
370  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
371  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
372  // get the id of this bin
373  size_t binId = lt.getGlobalBinByBin(xBin, yBin);
374  // get the size of this bin
375  size_t binSize = lt[binId].size();
376 
377  // loop over all hits in this bin
378  for (unsigned int j = 0; j < binSize; j++) {
379  unsigned int otherId = lt[binId][j];
380  bool otherSi = isOnlySi || cellsOnLayer.isSi[otherId];
381  if (otherSi) { //silicon cells cannot talk to scintillator cells
382  float dist = distance(i, otherId, layerId, false);
383  bool foundHigher = (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[i]) ||
384  (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[i] &&
385  cellsOnLayer.detid[otherId] > cellsOnLayer.detid[i]);
386  // if dist == i_delta, then last comer being the nearest higher
387  if (foundHigher && dist <= i_delta) {
388  // update i_delta
389  i_delta = dist;
390  // update i_nearestHigher
391  i_nearestHigher = otherId;
392  }
393  }
394  }
395  }
396  }
397 
398  bool foundNearestHigherInSearchBox = (i_delta != maxDelta);
399  if (foundNearestHigherInSearchBox) {
400  cellsOnLayer.delta[i] = i_delta;
401  cellsOnLayer.nearestHigher[i] = i_nearestHigher;
402  } else {
403  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
404  // we can safely maximize delta to be maxDelta
405  cellsOnLayer.delta[i] = maxDelta;
406  cellsOnLayer.nearestHigher[i] = -1;
407  }
408  } else {
409  //similar to silicon
410  float delta = delta_r;
412  std::array<int, 4> search_box = lt.searchBoxEtaPhi(cellsOnLayer.eta[i] - range,
413  cellsOnLayer.eta[i] + range,
414  cellsOnLayer.phi[i] - range,
415  cellsOnLayer.phi[i] + range);
416  // loop over all bins in the search box
417  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
418  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
419  // get the id of this bin
420  int phi = (yBin % T::type::nRowsPhi);
421  size_t binId = lt.getGlobalBinByBinEtaPhi(xBin, phi);
422  // get the size of this bin
423  size_t binSize = lt[binId].size();
424 
425  // loop over all hits in this bin
426  for (unsigned int j = 0; j < binSize; j++) {
427  unsigned int otherId = lt[binId][j];
428  if (!cellsOnLayer.isSi[otherId]) { //scintillator cells cannot talk to silicon cells
429  float dist = distance(i, otherId, layerId, true);
430  bool foundHigher = (cellsOnLayer.rho[otherId] > cellsOnLayer.rho[i]) ||
431  (cellsOnLayer.rho[otherId] == cellsOnLayer.rho[i] &&
432  cellsOnLayer.detid[otherId] > cellsOnLayer.detid[i]);
433  // if dist == i_delta, then last comer being the nearest higher
434  if (foundHigher && dist <= i_delta) {
435  // update i_delta
436  i_delta = dist;
437  // update i_nearestHigher
438  i_nearestHigher = otherId;
439  }
440  }
441  }
442  }
443  }
444 
445  bool foundNearestHigherInSearchBox = (i_delta != maxDelta);
446  if (foundNearestHigherInSearchBox) {
447  cellsOnLayer.delta[i] = i_delta;
448  cellsOnLayer.nearestHigher[i] = i_nearestHigher;
449  } else {
450  // otherwise delta is guaranteed to be larger outlierDeltaFactor_*delta_c
451  // we can safely maximize delta to be maxDelta
452  cellsOnLayer.delta[i] = maxDelta;
453  cellsOnLayer.nearestHigher[i] = -1;
454  }
455  }
456  LogDebug("HGCalCLUEAlgo") << "Debugging calculateDistanceToHigher: \n"
457  << " cell: " << i << " isSilicon: " << cellsOnLayer.isSi[i]
458  << " eta: " << cellsOnLayer.eta[i] << " phi: " << cellsOnLayer.phi[i]
459  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i]
460  << " nearest higher: " << cellsOnLayer.nearestHigher[i]
461  << " distance: " << cellsOnLayer.delta[i] << "\n";
462  }
463 }
std::vector< CellsOnLayer > cells_
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:238
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:435
float distance(int cell1, int cell2, int layerId, bool isEtaPhi) const
float outlierDeltaFactor_
#define LogDebug(id)

◆ calculateLocalDensity()

template<typename TILE >
void HGCalCLUEAlgoT< T >::calculateLocalDensity ( const TILE &  lt,
const unsigned int  layerId,
float  delta_c,
float  delta_r 
)
private

Definition at line 247 of file HGCalCLUEAlgo.cc.

References funct::abs(), python.cmstools::all(), dumpMFGeometry_cfg::delta, hitfit::delta_r(), HLT_2022v15_cff::distance, muonRecoAnalyzer_cfi::etaBin, mps_fire::i, HGCScintillatorDetId::ieta(), l1tTowerCalibrationProducer_cfi::iEta, HGCScintillatorDetId::iphi(), dqmiolumiharvest::j, LogDebug, SiStripPI::max, BeamMonitor_cff::phiBin, photonAnalyzer_cfi::xBin, and photonAnalyzer_cfi::yBin.

247  {
248  auto& cellsOnLayer = cells_[layerId];
249  unsigned int numberOfCells = cellsOnLayer.detid.size();
250  bool isOnlySi(false);
251  if (rhtools_.isOnlySilicon(layerId))
252  isOnlySi = true;
253 
254  for (unsigned int i = 0; i < numberOfCells; i++) {
255  bool isSi = isOnlySi || cellsOnLayer.isSi[i];
256  if (isSi) {
257  float delta = delta_c;
258  std::array<int, 4> search_box = lt.searchBox(
259  cellsOnLayer.x[i] - delta, cellsOnLayer.x[i] + delta, cellsOnLayer.y[i] - delta, cellsOnLayer.y[i] + delta);
260 
261  for (int xBin = search_box[0]; xBin < search_box[1] + 1; ++xBin) {
262  for (int yBin = search_box[2]; yBin < search_box[3] + 1; ++yBin) {
263  int binId = lt.getGlobalBinByBin(xBin, yBin);
264  size_t binSize = lt[binId].size();
265 
266  for (unsigned int j = 0; j < binSize; j++) {
267  unsigned int otherId = lt[binId][j];
268  bool otherSi = isOnlySi || cellsOnLayer.isSi[otherId];
269  if (otherSi) { //silicon cells cannot talk to scintillator cells
270  if (distance(i, otherId, layerId, false) < delta) {
271  cellsOnLayer.rho[i] += (i == otherId ? 1.f : 0.5f) * cellsOnLayer.weight[otherId];
272  }
273  }
274  }
275  }
276  }
277  } else {
278  float delta = delta_r;
279  std::array<int, 4> search_box = lt.searchBoxEtaPhi(cellsOnLayer.eta[i] - delta,
280  cellsOnLayer.eta[i] + delta,
281  cellsOnLayer.phi[i] - delta,
282  cellsOnLayer.phi[i] + delta);
283  cellsOnLayer.rho[i] += cellsOnLayer.weight[i];
284  float northeast(0), northwest(0), southeast(0), southwest(0), all(0);
285 
286  for (int etaBin = search_box[0]; etaBin < search_box[1] + 1; ++etaBin) {
287  for (int phiBin = search_box[2]; phiBin < search_box[3] + 1; ++phiBin) {
288  int phi = (phiBin % T::type::nRowsPhi);
289  int binId = lt.getGlobalBinByBinEtaPhi(etaBin, phi);
290  size_t binSize = lt[binId].size();
291 
292  for (unsigned int j = 0; j < binSize; j++) {
293  unsigned int otherId = lt[binId][j];
294  if (!cellsOnLayer.isSi[otherId]) { //scintillator cells cannot talk to silicon cells
295  if (distance(i, otherId, layerId, true) < delta) {
296  int iPhi = HGCScintillatorDetId(cellsOnLayer.detid[i]).iphi();
297  int otherIPhi = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).iphi();
298  int iEta = HGCScintillatorDetId(cellsOnLayer.detid[i]).ieta();
299  int otherIEta = HGCScintillatorDetId(cellsOnLayer.detid[otherId]).ieta();
300  int dIPhi = otherIPhi - iPhi;
301  dIPhi += abs(dIPhi) < 2 ? 0
302  : dIPhi < 0 ? scintMaxIphi_
303  : -scintMaxIphi_; // cells with iPhi=288 and iPhi=1 should be neiboring cells
304  int dIEta = otherIEta - iEta;
305  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
306  << " cell: " << otherId << " energy: " << cellsOnLayer.weight[otherId]
307  << " otherIPhi: " << otherIPhi << " iPhi: " << iPhi
308  << " otherIEta: " << otherIEta << " iEta: " << iEta << "\n";
309 
310  if (otherId != i) {
311  auto neighborCellContribution = 0.5f * cellsOnLayer.weight[otherId];
312  all += neighborCellContribution;
313  if (dIPhi >= 0 && dIEta >= 0)
314  northeast += neighborCellContribution;
315  if (dIPhi <= 0 && dIEta >= 0)
316  southeast += neighborCellContribution;
317  if (dIPhi >= 0 && dIEta <= 0)
318  northwest += neighborCellContribution;
319  if (dIPhi <= 0 && dIEta <= 0)
320  southwest += neighborCellContribution;
321  }
322  LogDebug("HGCalCLUEAlgo") << " Debugging calculateLocalDensity for Scintillator: \n"
323  << " northeast: " << northeast << " southeast: " << southeast
324  << " northwest: " << northwest << " southwest: " << southwest << "\n";
325  }
326  }
327  }
328  }
329  }
330  float neighborsval = (std::max(northeast, northwest) > std::max(southeast, southwest))
331  ? std::max(northeast, northwest)
332  : std::max(southeast, southwest);
333  if (use2x2_)
334  cellsOnLayer.rho[i] += neighborsval;
335  else
336  cellsOnLayer.rho[i] += all;
337  }
338  LogDebug("HGCalCLUEAlgo") << "Debugging calculateLocalDensity: \n"
339  << " cell: " << i << " isSilicon: " << cellsOnLayer.isSi[i]
340  << " eta: " << cellsOnLayer.eta[i] << " phi: " << cellsOnLayer.phi[i]
341  << " energy: " << cellsOnLayer.weight[i] << " density: " << cellsOnLayer.rho[i] << "\n";
342  }
343 }
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
int iphi() const
get the phi index
std::vector< CellsOnLayer > cells_
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:238
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:435
float distance(int cell1, int cell2, int layerId, bool isEtaPhi) const
#define LogDebug(id)

◆ calculatePosition()

template<typename T >
math::XYZPoint HGCalCLUEAlgoT< T >::calculatePosition ( const std::vector< int > &  v,
const unsigned int  layerId 
) const
private

Definition at line 187 of file HGCalCLUEAlgo.cc.

References f, mps_fire::i, dqm-mbProfile::log, SiStripPI::max, findQualityFiles::v, and x.

187  {
188  float total_weight = 0.f;
189  float x = 0.f;
190  float y = 0.f;
191 
192  unsigned int maxEnergyIndex = 0;
193  float maxEnergyValue = 0.f;
194 
195  auto& cellsOnLayer = cells_[layerId];
196 
197  // loop over hits in cluster candidate
198  // determining the maximum energy hit
199  for (auto i : v) {
200  total_weight += cellsOnLayer.weight[i];
201  if (cellsOnLayer.weight[i] > maxEnergyValue) {
202  maxEnergyValue = cellsOnLayer.weight[i];
203  maxEnergyIndex = i;
204  }
205  }
206 
207  // Si cell or Scintillator. Used to set approach and parameters
208  auto thick = rhtools_.getSiThickIndex(cellsOnLayer.detid[maxEnergyIndex]);
209  bool isSiliconCell = (thick != -1);
210 
211  // TODO: this is recomputing everything twice and overwriting the position with log weighting position
212  if (isSiliconCell) {
213  float total_weight_log = 0.f;
214  float x_log = 0.f;
215  float y_log = 0.f;
216  for (auto i : v) {
217  //for silicon only just use 1+6 cells = 1.3cm for all thicknesses
218  if (distance2(i, maxEnergyIndex, layerId, false) > positionDeltaRho2_)
219  continue;
220  float rhEnergy = cellsOnLayer.weight[i];
221  float Wi = std::max(thresholdW0_[thick] + std::log(rhEnergy / total_weight), 0.);
222  x_log += cellsOnLayer.x[i] * Wi;
223  y_log += cellsOnLayer.y[i] * Wi;
224  total_weight_log += Wi;
225  }
226 
227  total_weight = total_weight_log;
228  x = x_log;
229  y = y_log;
230  } else {
231  for (auto i : v) {
232  float rhEnergy = cellsOnLayer.weight[i];
233 
234  x += cellsOnLayer.x[i] * rhEnergy;
235  y += cellsOnLayer.y[i] * rhEnergy;
236  }
237  }
238  if (total_weight != 0.) {
239  float inv_tot_weight = 1.f / total_weight;
240  return math::XYZPoint(
241  x * inv_tot_weight, y * inv_tot_weight, rhtools_.getPosition(cellsOnLayer.detid[maxEnergyIndex]).z());
242  } else
243  return math::XYZPoint(0.f, 0.f, 0.f);
244 }
const double positionDeltaRho2_
std::vector< double > thresholdW0_
T z() const
Definition: PV3DBase.h:61
std::vector< CellsOnLayer > cells_
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
double f[11][100]
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:205
float distance2(int cell1, int cell2, int layerId, bool isEtaPhi) const

◆ computeThreshold()

template<typename T >
void HGCalCLUEAlgoT< T >::computeThreshold ( )

Definition at line 514 of file HGCalCLUEAlgo.cc.

References LogDebug.

514  {
515  // To support the TDR geometry and also the post-TDR one (v9 onwards), we
516  // need to change the logic of the vectors containing signal to noise and
517  // thresholds. The first 3 indices will keep on addressing the different
518  // thicknesses of the Silicon detectors in CE_E , the next 3 indices will address
519  // the thicknesses of the Silicon detectors in CE_H, while the last one, number 6 (the
520  // seventh) will address the Scintillators. This change will support both
521  // geometries at the same time.
522 
523  if (initialized_)
524  return; // only need to calculate thresholds once
525 
526  initialized_ = true;
527 
528  std::vector<double> dummy;
529 
530  dummy.resize(maxNumberOfThickIndices_ + !isNose_, 0); // +1 to accomodate for the Scintillators
531  thresholds_.resize(maxlayer_, dummy);
532  v_sigmaNoise_.resize(maxlayer_, dummy);
533 
534  for (unsigned ilayer = 1; ilayer <= maxlayer_; ++ilayer) {
535  for (unsigned ithick = 0; ithick < maxNumberOfThickIndices_; ++ithick) {
536  float sigmaNoise = 0.001f * fcPerEle_ * nonAgedNoises_[ithick] * dEdXweights_[ilayer] /
537  (fcPerMip_[ithick] * thicknessCorrection_[ithick]);
538  thresholds_[ilayer - 1][ithick] = sigmaNoise * ecut_;
539  v_sigmaNoise_[ilayer - 1][ithick] = sigmaNoise;
540  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " nonAgedNoises: " << nonAgedNoises_[ithick]
541  << " fcPerEle: " << fcPerEle_ << " fcPerMip: " << fcPerMip_[ithick]
542  << " noiseMip: " << fcPerEle_ * nonAgedNoises_[ithick] / fcPerMip_[ithick]
543  << " sigmaNoise: " << sigmaNoise << "\n";
544  }
545 
546  if (!isNose_) {
547  float scintillators_sigmaNoise = 0.001f * noiseMip_ * dEdXweights_[ilayer] / sciThicknessCorrection_;
548  thresholds_[ilayer - 1][maxNumberOfThickIndices_] = ecut_ * scintillators_sigmaNoise;
549  v_sigmaNoise_[ilayer - 1][maxNumberOfThickIndices_] = scintillators_sigmaNoise;
550  LogDebug("HGCalCLUEAlgo") << "ilayer: " << ilayer << " noiseMip: " << noiseMip_
551  << " scintillators_sigmaNoise: " << scintillators_sigmaNoise << "\n";
552  }
553  }
554 }
unsigned maxNumberOfThickIndices_
std::vector< std::vector< double > > thresholds_
std::vector< double > dEdXweights_
std::vector< double > fcPerMip_
std::vector< std::vector< double > > v_sigmaNoise_
std::vector< double > nonAgedNoises_
std::vector< double > thicknessCorrection_
double sciThicknessCorrection_
#define LogDebug(id)

◆ distance()

template<typename TILE >
float HGCalCLUEAlgoT< TILE >::distance ( int  cell1,
int  cell2,
int  layerId,
bool  isEtaPhi 
) const
inlineprivate

Definition at line 231 of file HGCalCLUEAlgo.h.

References hgcalTopologyTester_cfi::cell1, hgcalTopologyTester_cfi::cell2, HGCalCLUEAlgoT< TILE >::distance2(), and mathSSE::sqrt().

231  { // 2-d distance on the layer (x-y)
232  return std::sqrt(distance2(cell1, cell2, layerId, isEtaPhi));
233  }
T sqrt(T t)
Definition: SSEVec.h:19
float distance2(int cell1, int cell2, int layerId, bool isEtaPhi) const

◆ distance2()

template<typename TILE >
float HGCalCLUEAlgoT< TILE >::distance2 ( int  cell1,
int  cell2,
int  layerId,
bool  isEtaPhi 
) const
inlineprivate

Definition at line 219 of file HGCalCLUEAlgo.h.

References hgcalTopologyTester_cfi::cell1, hgcalTopologyTester_cfi::cell2, HGCalCLUEAlgoT< TILE >::cells_, reco::deltaPhi(), PVValHelper::dx, PVValHelper::dy, and phi.

Referenced by HGCalCLUEAlgoT< TILE >::distance().

219  { // distance squared
220  if (isEtaPhi) {
221  const float dphi = reco::deltaPhi(cells_[layerId].phi[cell1], cells_[layerId].phi[cell2]);
222  const float deta = cells_[layerId].eta[cell1] - cells_[layerId].eta[cell2];
223  return (deta * deta + dphi * dphi);
224  } else {
225  const float dx = cells_[layerId].x[cell1] - cells_[layerId].x[cell2];
226  const float dy = cells_[layerId].y[cell1] - cells_[layerId].y[cell2];
227  return (dx * dx + dy * dy);
228  }
229  }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
std::vector< CellsOnLayer > cells_

◆ fillPSetDescription()

template<typename TILE >
static void HGCalCLUEAlgoT< TILE >::fillPSetDescription ( edm::ParameterSetDescription iDesc)
inlinestatic

Definition at line 87 of file HGCalCLUEAlgo.h.

References edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), and AlCaHLTBitMon_QueryRunRegistry::string.

87  {
88  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
89  iDesc.add<double>("positionDeltaRho2", 1.69);
90  iDesc.add<std::vector<double>>("deltac",
91  {
92  1.3,
93  1.3,
94  1.3,
95  0.0315, // for scintillator
96  });
97  iDesc.add<bool>("dependSensor", true);
98  iDesc.add<double>("ecut", 3.0);
99  iDesc.add<double>("kappa", 9.0);
100  iDesc.addUntracked<unsigned int>("verbosity", 3);
101  iDesc.add<std::vector<double>>("dEdXweights", {});
102  iDesc.add<std::vector<double>>("thicknessCorrection", {});
103  iDesc.add<double>("sciThicknessCorrection", 0.9);
104  iDesc.add<int>("deltasi_index_regemfac", 3);
105  iDesc.add<unsigned>("maxNumberOfThickIndices", 6);
106  iDesc.add<std::vector<double>>("fcPerMip", {});
107  iDesc.add<double>("fcPerEle", 0.0);
108  iDesc.add<std::vector<double>>("noises", {});
109  edm::ParameterSetDescription descNestedNoiseMIP;
110  descNestedNoiseMIP.add<bool>("scaleByDose", false);
111  descNestedNoiseMIP.add<unsigned int>("scaleByDoseAlgo", 0);
112  descNestedNoiseMIP.add<double>("scaleByDoseFactor", 1.);
113  descNestedNoiseMIP.add<std::string>("doseMap", "");
114  descNestedNoiseMIP.add<std::string>("sipmMap", "");
115  descNestedNoiseMIP.add<double>("referenceIdark", -1);
116  descNestedNoiseMIP.add<double>("referenceXtalk", -1);
117  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
118  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
119  iDesc.add<bool>("use2x2", true); // use 2x2 or 3x3 scenario for scint density calculation
120  }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ findAndAssignClusters()

template<typename T >
int HGCalCLUEAlgoT< T >::findAndAssignClusters ( const unsigned int  layerId,
float  delta_c,
float  delta_r 
)
private

Definition at line 466 of file HGCalCLUEAlgo.cc.

References dumpMFGeometry_cfg::delta, hitfit::delta_r(), mps_fire::i, and dqmiolumiharvest::j.

466  {
467  // this is called once per layer and endcap...
468  // so when filling the cluster temporary vector of Hexels we resize each time
469  // by the number of clusters found. This is always equal to the number of
470  // cluster centers...
471  unsigned int nClustersOnLayer = 0;
472  auto& cellsOnLayer = cells_[layerId];
473  unsigned int numberOfCells = cellsOnLayer.detid.size();
474  std::vector<int> localStack;
475  // find cluster seeds and outlier
476  for (unsigned int i = 0; i < numberOfCells; i++) {
477  float rho_c = kappa_ * cellsOnLayer.sigmaNoise[i];
478  bool isSi = rhtools_.isOnlySilicon(layerId) || cellsOnLayer.isSi[i];
479  float delta = isSi ? delta_c : delta_r;
480 
481  // initialize clusterIndex
482  cellsOnLayer.clusterIndex[i] = -1;
483  bool isSeed = (cellsOnLayer.delta[i] > delta) && (cellsOnLayer.rho[i] >= rho_c);
484  bool isOutlier = (cellsOnLayer.delta[i] > outlierDeltaFactor_ * delta) && (cellsOnLayer.rho[i] < rho_c);
485  if (isSeed) {
486  cellsOnLayer.clusterIndex[i] = nClustersOnLayer;
487  cellsOnLayer.isSeed[i] = true;
488  nClustersOnLayer++;
489  localStack.push_back(i);
490 
491  } else if (!isOutlier) {
492  cellsOnLayer.followers[cellsOnLayer.nearestHigher[i]].push_back(i);
493  }
494  }
495 
496  // need to pass clusterIndex to their followers
497  while (!localStack.empty()) {
498  int endStack = localStack.back();
499  auto& thisSeed = cellsOnLayer.followers[endStack];
500  localStack.pop_back();
501 
502  // loop over followers
503  for (int j : thisSeed) {
504  // pass id to a follower
505  cellsOnLayer.clusterIndex[j] = cellsOnLayer.clusterIndex[endStack];
506  // push this follower to localStack
507  localStack.push_back(j);
508  }
509  }
510  return nClustersOnLayer;
511 }
std::vector< CellsOnLayer > cells_
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:238
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:435
float outlierDeltaFactor_

◆ getClusters()

template<typename T >
std::vector< reco::BasicCluster > HGCalCLUEAlgoT< T >::getClusters ( bool  )
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 131 of file HGCalCLUEAlgo.cc.

References AlignmentPI::calculatePosition(), haddnano::cl, reco::CaloID::DET_HGCAL_ENDCAP, HCALHighEnergyHPDFilter_cfi::energy, mps_fire::i, SiStripPI::max, unpackBuffers-CaloStage1::offsets, and position.

131  {
132  std::vector<int> offsets(numberOfClustersPerLayer_.size(), 0);
133 
134  int maxClustersOnLayer = numberOfClustersPerLayer_[0];
135 
136  for (unsigned layerId = 1; layerId < offsets.size(); ++layerId) {
137  offsets[layerId] = offsets[layerId - 1] + numberOfClustersPerLayer_[layerId - 1];
138 
139  maxClustersOnLayer = std::max(maxClustersOnLayer, numberOfClustersPerLayer_[layerId]);
140  }
141 
142  auto totalNumberOfClusters = offsets.back() + numberOfClustersPerLayer_.back();
143  clusters_v_.resize(totalNumberOfClusters);
144  std::vector<std::vector<int>> cellsIdInCluster;
145  cellsIdInCluster.reserve(maxClustersOnLayer);
146 
147  for (unsigned int layerId = 0; layerId < 2 * maxlayer_ + 2; ++layerId) {
148  cellsIdInCluster.resize(numberOfClustersPerLayer_[layerId]);
149  auto& cellsOnLayer = cells_[layerId];
150  unsigned int numberOfCells = cellsOnLayer.detid.size();
151  auto firstClusterIdx = offsets[layerId];
152 
153  for (unsigned int i = 0; i < numberOfCells; ++i) {
154  auto clusterIndex = cellsOnLayer.clusterIndex[i];
155  if (clusterIndex != -1)
156  cellsIdInCluster[clusterIndex].push_back(i);
157  }
158 
159  std::vector<std::pair<DetId, float>> thisCluster;
160 
161  for (auto& cl : cellsIdInCluster) {
162  auto position = calculatePosition(cl, layerId);
163  float energy = 0.f;
164  int seedDetId = -1;
165 
166  for (auto cellIdx : cl) {
167  energy += cellsOnLayer.weight[cellIdx];
168  thisCluster.emplace_back(cellsOnLayer.detid[cellIdx], 1.f);
169  if (cellsOnLayer.isSeed[cellIdx]) {
170  seedDetId = cellsOnLayer.detid[cellIdx];
171  }
172  }
173  auto globalClusterIndex = cellsOnLayer.clusterIndex[cl[0]] + firstClusterIdx;
174 
175  clusters_v_[globalClusterIndex] =
177  clusters_v_[globalClusterIndex].setSeed(seedDetId);
178  thisCluster.clear();
179  }
180 
181  cellsIdInCluster.clear();
182  }
183  return clusters_v_;
184 }
std::vector< int > numberOfClustersPerLayer_
CaloCluster BasicCluster
std::vector< reco::BasicCluster > clusters_v_
std::vector< CellsOnLayer > cells_
reco::CaloCluster::AlgoId algoId_
static int position[264][3]
Definition: ReadPGInfo.cc:289
math::XYZPoint calculatePosition(const std::vector< int > &v, const unsigned int layerId) const

◆ getDensity()

template<typename T >
Density HGCalCLUEAlgoT< T >::getDensity ( )
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 565 of file HGCalCLUEAlgo.cc.

565  {
566  return density_;
567 }

◆ getEventSetupPerAlgorithm()

template<typename T >
void HGCalCLUEAlgoT< T >::getEventSetupPerAlgorithm ( const edm::EventSetup es)
overridevirtual

Reimplemented from HGCalClusteringAlgoBase.

Definition at line 19 of file HGCalCLUEAlgo.cc.

19  {
20  cells_.clear();
22  cells_.resize(2 * (maxlayer_ + 1));
23  numberOfClustersPerLayer_.resize(2 * (maxlayer_ + 1), 0);
24 }
std::vector< int > numberOfClustersPerLayer_
std::vector< CellsOnLayer > cells_

◆ makeClusters()

template<typename T >
void HGCalCLUEAlgoT< T >::makeClusters ( )
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 99 of file HGCalCLUEAlgo.cc.

References hitfit::delta_r(), PVValHelper::eta, mps_fire::i, LogDebug, and x.

99  {
100  // assign all hits in each layer to a cluster core
101  tbb::this_task_arena::isolate([&] {
102  tbb::parallel_for(size_t(0), size_t(2 * maxlayer_ + 2), [&](size_t i) {
104  T lt;
105  lt.clear();
106  lt.fill(cells_[i].x, cells_[i].y, cells_[i].eta, cells_[i].phi, cells_[i].isSi);
107  float delta_c; // maximum search distance (critical distance) for local
108  // density calculation
109  if (i % maxlayer_ < lastLayerEE_)
110  delta_c = vecDeltas_[0];
111  else if (i % maxlayer_ < (firstLayerBH_ - 1))
112  delta_c = vecDeltas_[1];
113  else
114  delta_c = vecDeltas_[2];
115  float delta_r = vecDeltas_[3];
116  LogDebug("HGCalCLUEAlgo") << "maxlayer: " << maxlayer_ << " lastLayerEE: " << lastLayerEE_
117  << " firstLayerBH: " << firstLayerBH_ << "\n";
118 
119  calculateLocalDensity(lt, i, delta_c, delta_r);
120  calculateDistanceToHigher(lt, i, delta_c, delta_r);
122  });
123  });
124  //Now that we have the density per point we can store it
125  for (unsigned int i = 0; i < 2 * maxlayer_ + 2; ++i) {
126  setDensity(i);
127  }
128 }
int findAndAssignClusters(const unsigned int layerId, float delta_c, float delta_r)
std::vector< int > numberOfClustersPerLayer_
void setDensity(const unsigned int layerId)
void prepareDataStructures(const unsigned int layerId)
void calculateDistanceToHigher(const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
void calculateLocalDensity(const TILE &lt, const unsigned int layerId, float delta_c, float delta_r)
std::vector< CellsOnLayer > cells_
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:238
std::vector< double > vecDeltas_
long double T
#define LogDebug(id)

◆ populate()

template<typename T >
void HGCalCLUEAlgoT< T >::populate ( const HGCRecHitCollection hits)
overridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 27 of file HGCalCLUEAlgo.cc.

References DetId::det(), CaloRecHit::detid(), CaloRecHit::energy(), DetId::HGCalHSi, HGCHEF, hfClusterShapes_cfi::hits, mps_fire::i, hltrates_dqm_sourceclient-live_cfg::offset, position, and DetId::subdetId().

27  {
28  // loop over all hits and create the Hexel structure, skip energies below ecut
29 
30  if (dependSensor_) {
31  // for each layer and wafer calculate the thresholds (sigmaNoise and energy)
32  // once
34  }
35 
36  for (unsigned int i = 0; i < hits.size(); ++i) {
37  const HGCRecHit& hgrh = hits[i];
38  DetId detid = hgrh.detid();
39  unsigned int layerOnSide = (rhtools_.getLayerWithOffset(detid) - 1);
40 
41  // set sigmaNoise default value 1 to use kappa value directly in case of
42  // sensor-independent thresholds
43  float sigmaNoise = 1.f;
44  if (dependSensor_) {
45  int thickness_index = rhtools_.getSiThickIndex(detid);
46  if (thickness_index == -1)
47  thickness_index = maxNumberOfThickIndices_;
48 
49  double storedThreshold = thresholds_[layerOnSide][thickness_index];
50  if (detid.det() == DetId::HGCalHSi || detid.subdetId() == HGCHEF) {
51  storedThreshold = thresholds_[layerOnSide][thickness_index + deltasi_index_regemfac_];
52  }
53  sigmaNoise = v_sigmaNoise_[layerOnSide][thickness_index];
54 
55  if (hgrh.energy() < storedThreshold)
56  continue; // this sets the ZS threshold at ecut times the sigma noise
57  // for the sensor
58  }
59  if (!dependSensor_ && hgrh.energy() < ecut_)
60  continue;
62  int offset = ((rhtools_.zside(detid) + 1) >> 1) * maxlayer_;
63  int layer = layerOnSide + offset;
64 
65  cells_[layer].detid.emplace_back(detid);
66  cells_[layer].x.emplace_back(position.x());
67  cells_[layer].y.emplace_back(position.y());
69  cells_[layer].isSi.emplace_back(rhtools_.isSilicon(detid));
70  cells_[layer].eta.emplace_back(position.eta());
71  cells_[layer].phi.emplace_back(position.phi());
72  } // else, isSilicon == true and eta phi values will not be used
73  cells_[layer].weight.emplace_back(hgrh.energy());
74  cells_[layer].sigmaNoise.emplace_back(sigmaNoise);
75  }
76 }
unsigned maxNumberOfThickIndices_
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
int deltasi_index_regemfac_
std::vector< std::vector< double > > thresholds_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
int zside(const DetId &id) const
Definition: RecHitTools.cc:163
constexpr float energy() const
Definition: CaloRecHit.h:29
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
std::vector< CellsOnLayer > cells_
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
bool isSilicon(const DetId &) const
Definition: RecHitTools.cc:428
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:435
Definition: DetId.h:17
std::vector< std::vector< double > > v_sigmaNoise_
static int position[264][3]
Definition: ReadPGInfo.cc:289
void computeThreshold()
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:205
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365

◆ prepareDataStructures()

template<typename T >
void HGCalCLUEAlgoT< T >::prepareDataStructures ( const unsigned int  layerId)
private

Definition at line 79 of file HGCalCLUEAlgo.cc.

References f, and cmsLHEtoEOSManager::l.

79  {
80  auto cellsSize = cells_[l].detid.size();
81  cells_[l].rho.resize(cellsSize, 0.f);
82  cells_[l].delta.resize(cellsSize, 9999999);
83  cells_[l].nearestHigher.resize(cellsSize, -1);
84  cells_[l].clusterIndex.resize(cellsSize, -1);
85  cells_[l].followers.resize(cellsSize);
86  cells_[l].isSeed.resize(cellsSize, false);
87  if (rhtools_.isOnlySilicon(l)) {
88  cells_[l].isSi.resize(cellsSize, true);
89  cells_[l].eta.resize(cellsSize, 0.f);
90  cells_[l].phi.resize(cellsSize, 0.f);
91  }
92 }
std::vector< CellsOnLayer > cells_
double f[11][100]
bool isOnlySilicon(const unsigned int layer) const
Definition: RecHitTools.cc:435

◆ reset()

template<typename TILE >
void HGCalCLUEAlgoT< TILE >::reset ( void  )
inlineoverridevirtual

Implements HGCalClusteringAlgoBase.

Definition at line 69 of file HGCalCLUEAlgo.h.

References hgcalTBTopologyTester_cfi::cells, HGCalCLUEAlgoT< TILE >::cells_, haddnano::cl, HGCalClusteringAlgoBase::clusters_v_, HGCalCLUEAlgoT< TILE >::density_, and HGCalCLUEAlgoT< TILE >::numberOfClustersPerLayer_.

69  {
70  clusters_v_.clear();
71  clusters_v_.shrink_to_fit();
72  for (auto& cl : numberOfClustersPerLayer_) {
73  cl = 0;
74  }
75 
76  for (auto& cells : cells_) {
77  cells.clear();
78  cells.shrink_to_fit();
79  }
80  density_.clear();
81  }
std::vector< int > numberOfClustersPerLayer_
std::vector< reco::BasicCluster > clusters_v_
std::vector< CellsOnLayer > cells_

◆ setDensity()

template<typename T >
void HGCalCLUEAlgoT< T >::setDensity ( const unsigned int  layerId)
private

Definition at line 557 of file HGCalCLUEAlgo.cc.

References mps_fire::i.

557  {
558  auto& cellsOnLayer = cells_[layerId];
559  unsigned int numberOfCells = cellsOnLayer.detid.size();
560  for (unsigned int i = 0; i < numberOfCells; ++i)
561  density_[cellsOnLayer.detid[i]] = cellsOnLayer.rho[i];
562 }
std::vector< CellsOnLayer > cells_

Member Data Documentation

◆ cells_

template<typename TILE >
std::vector<CellsOnLayer> HGCalCLUEAlgoT< TILE >::cells_
private

◆ dEdXweights_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::dEdXweights_
private

Definition at line 143 of file HGCalCLUEAlgo.h.

◆ deltasi_index_regemfac_

template<typename TILE >
int HGCalCLUEAlgoT< TILE >::deltasi_index_regemfac_
private

Definition at line 146 of file HGCalCLUEAlgo.h.

◆ density_

template<typename TILE >
Density HGCalCLUEAlgoT< TILE >::density_
private

Definition at line 138 of file HGCalCLUEAlgo.h.

Referenced by HGCalCLUEAlgoT< TILE >::reset().

◆ dependSensor_

template<typename TILE >
bool HGCalCLUEAlgoT< TILE >::dependSensor_
private

Definition at line 142 of file HGCalCLUEAlgo.h.

◆ ecut_

template<typename TILE >
double HGCalCLUEAlgoT< TILE >::ecut_
private

Definition at line 135 of file HGCalCLUEAlgo.h.

◆ fcPerEle_

template<typename TILE >
double HGCalCLUEAlgoT< TILE >::fcPerEle_
private

Definition at line 149 of file HGCalCLUEAlgo.h.

◆ fcPerMip_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::fcPerMip_
private

Definition at line 148 of file HGCalCLUEAlgo.h.

◆ initialized_

template<typename TILE >
bool HGCalCLUEAlgoT< TILE >::initialized_
private

Definition at line 158 of file HGCalCLUEAlgo.h.

◆ kappa_

template<typename TILE >
double HGCalCLUEAlgoT< TILE >::kappa_
private

Definition at line 132 of file HGCalCLUEAlgo.h.

◆ maxNumberOfThickIndices_

template<typename TILE >
unsigned HGCalCLUEAlgoT< TILE >::maxNumberOfThickIndices_
private

Definition at line 147 of file HGCalCLUEAlgo.h.

◆ noiseMip_

template<typename TILE >
double HGCalCLUEAlgoT< TILE >::noiseMip_
private

Definition at line 151 of file HGCalCLUEAlgo.h.

◆ nonAgedNoises_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::nonAgedNoises_
private

Definition at line 150 of file HGCalCLUEAlgo.h.

◆ numberOfClustersPerLayer_

template<typename TILE >
std::vector<int> HGCalCLUEAlgoT< TILE >::numberOfClustersPerLayer_
private

Definition at line 217 of file HGCalCLUEAlgo.h.

Referenced by HGCalCLUEAlgoT< TILE >::reset().

◆ outlierDeltaFactor_

template<typename TILE >
float HGCalCLUEAlgoT< TILE >::outlierDeltaFactor_ = 2.f
private

Definition at line 160 of file HGCalCLUEAlgo.h.

◆ positionDeltaRho2_

template<typename TILE >
const double HGCalCLUEAlgoT< TILE >::positionDeltaRho2_
private

Definition at line 128 of file HGCalCLUEAlgo.h.

◆ sciThicknessCorrection_

template<typename TILE >
double HGCalCLUEAlgoT< TILE >::sciThicknessCorrection_
private

Definition at line 145 of file HGCalCLUEAlgo.h.

◆ thicknessCorrection_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::thicknessCorrection_
private

Definition at line 144 of file HGCalCLUEAlgo.h.

◆ thresholds_

template<typename TILE >
std::vector<std::vector<double> > HGCalCLUEAlgoT< TILE >::thresholds_
private

Definition at line 152 of file HGCalCLUEAlgo.h.

◆ thresholdW0_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::thresholdW0_
private

Definition at line 127 of file HGCalCLUEAlgo.h.

◆ use2x2_

template<typename TILE >
bool HGCalCLUEAlgoT< TILE >::use2x2_
private

Definition at line 155 of file HGCalCLUEAlgo.h.

◆ v_sigmaNoise_

template<typename TILE >
std::vector<std::vector<double> > HGCalCLUEAlgoT< TILE >::v_sigmaNoise_
private

Definition at line 153 of file HGCalCLUEAlgo.h.

◆ vecDeltas_

template<typename TILE >
std::vector<double> HGCalCLUEAlgoT< TILE >::vecDeltas_
private

Definition at line 131 of file HGCalCLUEAlgo.h.