CMS 3D CMS Logo

HGCalCLUEAlgo.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_HGCalRecProducers_HGCalCLUEAlgo_h
2 #define RecoLocalCalo_HGCalRecProducers_HGCalCLUEAlgo_h
3 
5 
7 
13 
17 
20 
22 
23 // C/C++ headers
24 #include <set>
25 #include <string>
26 #include <vector>
27 
28 template <typename TILE, typename STRATEGY>
30 public:
33  (HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3),
34  reco::CaloCluster::undefined),
35  vecDeltas_(ps.getParameter<std::vector<double>>("deltac")),
36  kappa_(ps.getParameter<double>("kappa")),
37  ecut_(ps.getParameter<double>("ecut")),
38  dependSensor_(ps.getParameter<bool>("dependSensor")),
39  dEdXweights_(ps.getParameter<std::vector<double>>("dEdXweights")),
40  thicknessCorrection_(ps.getParameter<std::vector<double>>("thicknessCorrection")),
41  sciThicknessCorrection_(ps.getParameter<double>("sciThicknessCorrection")),
42  deltasi_index_regemfac_(ps.getParameter<int>("deltasi_index_regemfac")),
43  maxNumberOfThickIndices_(ps.getParameter<unsigned>("maxNumberOfThickIndices")),
44  fcPerMip_(ps.getParameter<std::vector<double>>("fcPerMip")),
45  fcPerEle_(ps.getParameter<double>("fcPerEle")),
46  nonAgedNoises_(ps.getParameter<std::vector<double>>("noises")),
47  noiseMip_(ps.getParameter<edm::ParameterSet>("noiseMip").getParameter<double>("noise_MIP")),
48  use2x2_(ps.getParameter<bool>("use2x2")),
50 
51  ~HGCalCLUEAlgoT() override {}
52 
53  void getEventSetupPerAlgorithm(const edm::EventSetup& es) override;
54 
55  void populate(const HGCRecHitCollection& hits) override;
56 
57  // this is the method that will start the clusterisation (it is possible to invoke this method
58  // more than once - but make sure it is with different hit collections (or else use reset)
59 
60  void makeClusters() override;
61 
62  // this is the method to get the cluster collection out
63  std::vector<reco::BasicCluster> getClusters(bool) override;
64 
65  void reset() override {
66  clusters_v_.clear();
67  clusters_v_.shrink_to_fit();
68  for (auto& cl : numberOfClustersPerLayer_) {
69  cl = 0;
70  }
71 
72  for (auto& cells : cells_) {
73  cells.clear();
74  cells.shrink_to_fit();
75  }
76  }
77 
78  void computeThreshold();
79 
81  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
82  iDesc.add<double>("positionDeltaRho2", 1.69);
83  iDesc.add<std::vector<double>>("deltac",
84  {
85  1.3,
86  1.3,
87  1.3,
88  0.0315, // for scintillator
89  });
90  iDesc.add<bool>("dependSensor", true);
91  iDesc.add<double>("ecut", 3.0);
92  iDesc.add<double>("kappa", 9.0);
93  iDesc.addUntracked<unsigned int>("verbosity", 3);
94  iDesc.add<std::vector<double>>("dEdXweights", {});
95  iDesc.add<std::vector<double>>("thicknessCorrection", {});
96  iDesc.add<double>("sciThicknessCorrection", 0.9);
97  iDesc.add<int>("deltasi_index_regemfac", 3);
98  iDesc.add<unsigned>("maxNumberOfThickIndices", 6);
99  iDesc.add<std::vector<double>>("fcPerMip", {});
100  iDesc.add<double>("fcPerEle", 0.0);
101  iDesc.add<std::vector<double>>("noises", {});
102  edm::ParameterSetDescription descNestedNoiseMIP;
103  descNestedNoiseMIP.add<bool>("scaleByDose", false);
104  descNestedNoiseMIP.add<unsigned int>("scaleByDoseAlgo", 0);
105  descNestedNoiseMIP.add<double>("scaleByDoseFactor", 1.);
106  descNestedNoiseMIP.add<std::string>("doseMap", "");
107  descNestedNoiseMIP.add<std::string>("sipmMap", "");
108  descNestedNoiseMIP.add<double>("referenceIdark", -1);
109  descNestedNoiseMIP.add<double>("referenceXtalk", -1);
110  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
111  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
112  iDesc.add<bool>("use2x2", true); // use 2x2 or 3x3 scenario for scint density calculation
113  }
114 
117 
118 private:
119  // The two parameters used to identify clusters
120  std::vector<double> vecDeltas_;
121  double kappa_;
122 
123  // The hit energy cutoff
124  double ecut_;
125 
126  // various parameters used for calculating the noise levels for a given sensor (and whether to use
127  // them)
129  std::vector<double> dEdXweights_;
130  std::vector<double> thicknessCorrection_;
134  std::vector<double> fcPerMip_;
135  double fcPerEle_;
136  std::vector<double> nonAgedNoises_;
137  double noiseMip_;
138  std::vector<std::vector<double>> thresholds_;
139  std::vector<std::vector<double>> v_sigmaNoise_;
140 
141  bool use2x2_;
142 
143  // initialization bool
145 
146  float outlierDeltaFactor_ = 2.f;
147 
148  struct CellsOnLayer {
149  std::vector<DetId> detid;
150  std::vector<float> dim1;
151  std::vector<float> dim2;
152 
153  std::vector<float> weight;
154  std::vector<float> rho;
155 
156  std::vector<float> delta;
157  std::vector<int> nearestHigher;
158  std::vector<int> clusterIndex;
159  std::vector<float> sigmaNoise;
160  std::vector<std::vector<int>> followers;
161  std::vector<bool> isSeed;
162 
163  void clear() {
164  detid.clear();
165  dim1.clear();
166  dim2.clear();
167  weight.clear();
168  rho.clear();
169  delta.clear();
170  nearestHigher.clear();
171  clusterIndex.clear();
172  sigmaNoise.clear();
173  followers.clear();
174  isSeed.clear();
175  }
176 
177  void shrink_to_fit() {
178  detid.shrink_to_fit();
179  dim1.shrink_to_fit();
180  dim2.shrink_to_fit();
181  weight.shrink_to_fit();
182  rho.shrink_to_fit();
183  delta.shrink_to_fit();
184  nearestHigher.shrink_to_fit();
185  clusterIndex.shrink_to_fit();
186  sigmaNoise.shrink_to_fit();
187  followers.shrink_to_fit();
188  isSeed.shrink_to_fit();
189  }
190  };
191 
192  std::vector<CellsOnLayer> cells_;
193 
194  std::vector<int> numberOfClustersPerLayer_;
195 
196  inline float distance(const TILE& lt, int cell1, int cell2, int layerId) const { // 2-d distance on the layer (x-y)
197  return std::sqrt(lt.distance2(cells_[layerId].dim1[cell1],
198  cells_[layerId].dim2[cell1],
199  cells_[layerId].dim1[cell2],
200  cells_[layerId].dim2[cell2]));
201  }
202 
203  void prepareDataStructures(const unsigned int layerId);
204  void calculateLocalDensity(const TILE& lt, const unsigned int layerId,
205  float delta); // return max density
206  void calculateLocalDensity(const TILE& lt, const unsigned int layerId, float delta, HGCalSiliconStrategy strategy);
207  void calculateLocalDensity(const TILE& lt,
208  const unsigned int layerId,
209  float delta,
210  HGCalScintillatorStrategy strategy);
211  void calculateDistanceToHigher(const TILE& lt, const unsigned int layerId, float delta);
212  int findAndAssignClusters(const unsigned int layerId, float delta);
213 };
214 
215 // explicit template instantiation
219 
223 
224 #endif
void computeThreshold()
std::vector< double > nonAgedNoises_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
~HGCalCLUEAlgoT() override
Definition: HGCalCLUEAlgo.h:51
HGCalCLUEAlgoT(const edm::ParameterSet &ps)
Definition: HGCalCLUEAlgo.h:31
std::vector< double > vecDeltas_
std::vector< reco::BasicCluster > clusters_v_
Definition: weight.py:1
std::vector< int > clusterIndex
void calculateDistanceToHigher(const TILE &lt, const unsigned int layerId, float delta)
int deltasi_index_regemfac_
void populate(const HGCRecHitCollection &hits) override
float outlierDeltaFactor_
double sciThicknessCorrection_
std::vector< std::vector< int > > followers
std::vector< int > nearestHigher
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< double > thicknessCorrection_
void makeClusters() override
std::vector< bool > isSeed
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< double > fcPerMip_
std::vector< float > weight
float distance(const TILE &lt, int cell1, int cell2, int layerId) const
std::vector< double > dEdXweights_
std::vector< std::vector< double > > v_sigmaNoise_
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
std::vector< int > numberOfClustersPerLayer_
void calculateLocalDensity(const TILE &lt, const unsigned int layerId, float delta)
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
Definition: HGCalCLUEAlgo.h:80
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
unsigned maxNumberOfThickIndices_
std::vector< CellsOnLayer > cells_
math::XYZPoint Point
point in the space
void prepareDataStructures(const unsigned int layerId)
std::vector< float > delta
std::vector< DetId > detid
fixed size matrix
HLT enums.
std::vector< float > dim1
std::vector< float > rho
std::vector< float > sigmaNoise
std::vector< std::vector< double > > thresholds_
void reset() override
Definition: HGCalCLUEAlgo.h:65
int findAndAssignClusters(const unsigned int layerId, float delta)
std::vector< reco::BasicCluster > getClusters(bool) override
std::vector< float > dim2