CMS 3D CMS Logo

HGCalImagingAlgo.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_HGCalRecProducers_HGCalImagingAlgo_h
2 #define RecoLocalCalo_HGCalRecProducers_HGCalImagingAlgo_h
3 
5 
7 
17 
20 
23 
24 // C/C++ headers
25 #include <string>
26 #include <vector>
27 #include <set>
28 
30 
32 public:
35  (HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3),
36  reco::CaloCluster::undefined),
37  thresholdW0_(ps.getParameter<std::vector<double>>("thresholdW0")),
38  positionDeltaRho_c_(ps.getParameter<std::vector<double>>("positionDeltaRho_c")),
39  vecDeltas_(ps.getParameter<std::vector<double>>("deltac")),
40  kappa_(ps.getParameter<double>("kappa")),
41  ecut_(ps.getParameter<double>("ecut")),
42  sigma2_(1.0),
43  dependSensor_(ps.getParameter<bool>("dependSensor")),
44  dEdXweights_(ps.getParameter<std::vector<double>>("dEdXweights")),
45  thicknessCorrection_(ps.getParameter<std::vector<double>>("thicknessCorrection")),
46  fcPerMip_(ps.getParameter<std::vector<double>>("fcPerMip")),
47  fcPerEle_(ps.getParameter<double>("fcPerEle")),
48  nonAgedNoises_(ps.getParameter<edm::ParameterSet>("noises").getParameter<std::vector<double>>("values")),
49  noiseMip_(ps.getParameter<edm::ParameterSet>("noiseMip").getParameter<double>("noise_MIP")),
51 
52  ~HGCalImagingAlgo() override {}
53 
54  void getEventSetupPerAlgorithm(const edm::EventSetup &es) override;
55 
56  void populate(const HGCRecHitCollection &hits) override;
57  // this is the method that will start the clusterisation (it is possible to invoke this method more than once - but make sure it is with
58  // 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  // use this if you want to reuse the same cluster object but don't want to accumulate clusters (hardly useful?)
66  void reset() override {
67  clusters_v_.clear();
68  clusters_v_.shrink_to_fit();
69  layerClustersPerLayer_.clear();
70  layerClustersPerLayer_.shrink_to_fit();
71  for (auto &it : points_) {
72  it.clear();
73  it.shrink_to_fit();
74  std::vector<KDNode>().swap(it);
75  }
76  for (unsigned int i = 0; i < minpos_.size(); i++) {
77  minpos_[i][0] = 0.;
78  minpos_[i][1] = 0.;
79  maxpos_[i][0] = 0.;
80  maxpos_[i][1] = 0.;
81  }
82  }
83  void computeThreshold();
84 
85  //getDensity
86  Density getDensity() override;
87 
89  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
90  iDesc.add<std::vector<double>>("positionDeltaRho_c", {1.3, 1.3, 1.3});
91  iDesc.add<std::vector<double>>("deltac",
92  {
93  2.0,
94  2.0,
95  5.0,
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<std::vector<double>>("fcPerMip", {});
104  iDesc.add<double>("fcPerEle", 0.0);
105  edm::ParameterSetDescription descNestedNoises;
106  descNestedNoises.add<std::vector<double>>("values", {});
107  iDesc.add<edm::ParameterSetDescription>("noises", descNestedNoises);
108  edm::ParameterSetDescription descNestedNoiseMIP;
109  descNestedNoiseMIP.add<bool>("scaleByDose", false);
110  descNestedNoiseMIP.add<double>("scaleByDoseFactor", 1.);
111  iDesc.add<edm::ParameterSetDescription>("scaleByDose", descNestedNoiseMIP);
112  descNestedNoiseMIP.add<std::string>("doseMap", "");
113  iDesc.add<edm::ParameterSetDescription>("doseMap", descNestedNoiseMIP);
114  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
115  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
116  }
117 
120 
121 private:
122  // To compute the cluster position
123  std::vector<double> thresholdW0_;
124  std::vector<double> positionDeltaRho_c_;
125 
126  // The two parameters used to identify clusters
127  std::vector<double> vecDeltas_;
128  double kappa_;
129 
130  // The hit energy cutoff
131  double ecut_;
132 
133  // for energy sharing
134  double sigma2_; // transverse shower size
135 
136  // The vector of clusters
137  std::vector<reco::BasicCluster> clusters_v_;
138 
139  // For keeping the density per hit
141 
142  // various parameters used for calculating the noise levels for a given sensor (and whether to use them)
144  std::vector<double> dEdXweights_;
145  std::vector<double> thicknessCorrection_;
146  std::vector<double> fcPerMip_;
147  double fcPerEle_;
148  std::vector<double> nonAgedNoises_;
149  double noiseMip_;
150  std::vector<std::vector<double>> thresholds_;
151  std::vector<std::vector<double>> sigmaNoise_;
152 
153  // initialization bool
155 
156  struct Hexel {
157  double x;
158  double y;
159  double z;
161  double weight;
162  double fraction;
164  double rho;
165  double delta;
167  bool isBorder;
168  bool isHalo;
170  float sigmaNoise;
171  float thickness;
173 
175  DetId id_in,
176  bool isHalf,
177  float sigmaNoise_in,
178  float thickness_in,
179  const hgcal::RecHitTools *tools_in)
180  : isHalfCell(isHalf),
181  weight(0.),
182  fraction(1.0),
183  detid(id_in),
184  rho(0.),
185  delta(0.),
186  nearestHigher(-1),
187  isBorder(false),
188  isHalo(false),
189  clusterIndex(-1),
190  sigmaNoise(sigmaNoise_in),
191  thickness(thickness_in),
192  tools(tools_in) {
193  const GlobalPoint position(tools->getPosition(detid));
194  weight = hit.energy();
195  x = position.x();
196  y = position.y();
197  z = position.z();
198  }
200  : x(0.),
201  y(0.),
202  z(0.),
203  isHalfCell(false),
204  weight(0.),
205  fraction(1.0),
206  detid(),
207  rho(0.),
208  delta(0.),
209  nearestHigher(-1),
210  isBorder(false),
211  isHalo(false),
212  clusterIndex(-1),
213  sigmaNoise(0.),
214  thickness(0.),
215  tools(nullptr) {}
216  bool operator>(const Hexel &rhs) const { return (rho > rhs.rho); }
217  };
218 
221 
222  std::vector<std::vector<std::vector<KDNode>>> layerClustersPerLayer_;
223 
224  std::vector<size_t> sort_by_delta(const std::vector<KDNode> &v) const {
225  std::vector<size_t> idx(v.size());
226  std::iota(std::begin(idx), std::end(idx), 0);
227  sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1].data.delta > v[i2].data.delta; });
228  return idx;
229  }
230 
231  std::vector<std::vector<KDNode>> points_; //a vector of vectors of hexels, one for each layer
232  //@@EM todo: the number of layers should be obtained programmatically - the range is 1-n instead of 0-n-1...
233 
234  std::vector<std::array<float, 2>> minpos_;
235  std::vector<std::array<float, 2>> maxpos_;
236 
237  //these functions should be in a helper class.
238  inline double distance2(const Hexel &pt1, const Hexel &pt2) const { //distance squared
239  const double dx = pt1.x - pt2.x;
240  const double dy = pt1.y - pt2.y;
241  return (dx * dx + dy * dy);
242  } //distance squaredq
243  inline double distance(const Hexel &pt1, const Hexel &pt2) const { //2-d distance on the layer (x-y)
244  return std::sqrt(distance2(pt1, pt2));
245  }
246  double calculateLocalDensity(std::vector<KDNode> &, KDTree &, const unsigned int) const; //return max density
247  double calculateDistanceToHigher(std::vector<KDNode> &) const;
248  int findAndAssignClusters(std::vector<KDNode> &,
249  KDTree &,
250  double,
251  KDTreeBox<2> &,
252  const unsigned int,
253  std::vector<std::vector<KDNode>> &) const;
254  math::XYZPoint calculatePosition(std::vector<KDNode> &) const;
255 
256  //For keeping the density information
257  void setDensity(const std::vector<KDNode> &nd);
258 
259  // attempt to find subclusters within a given set of hexels
260  std::vector<unsigned> findLocalMaximaInCluster(const std::vector<KDNode> &);
261  math::XYZPoint calculatePositionWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
262  double calculateEnergyWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
263  // outputs
264  void shareEnergy(const std::vector<KDNode> &, const std::vector<unsigned> &, std::vector<std::vector<double>> &);
265 };
266 
267 #endif
math::XYZPoint calculatePosition(std::vector< KDNode > &) const
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
bool operator>(const Hexel &rhs) const
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
KDTreeNodeInfo< Hexel > KDNode
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
std::vector< std::vector< double > > thresholds_
void shareEnergy(const std::vector< KDNode > &, const std::vector< unsigned > &, std::vector< std::vector< double >> &)
double distance(const Hexel &pt1, const Hexel &pt2) const
KDTreeLinkerAlgo< Hexel > KDTree
math::XYZPoint Point
point in the space
Definition: weight.py:1
double distance2(const Hexel &pt1, const Hexel &pt2) const
std::vector< reco::BasicCluster > getClusters(bool) override
std::vector< std::array< float, 2 > > minpos_
std::vector< double > positionDeltaRho_c_
std::vector< double > thresholdW0_
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
void makeClusters() override
std::map< DetId, float > Density
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const
std::vector< std::array< float, 2 > > maxpos_
T sqrt(T t)
Definition: SSEVec.h:19
Definition: tools.py:1
std::vector< reco::BasicCluster > clusters_v_
std::vector< std::vector< KDNode > > points_
std::vector< double > vecDeltas_
std::vector< double > dEdXweights_
std::vector< std::vector< double > > sigmaNoise_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int) const
void setDensity(const std::vector< KDNode > &nd)
std::vector< unsigned > findLocalMaximaInCluster(const std::vector< KDNode > &)
Density getDensity() override
void populate(const HGCRecHitCollection &hits) override
const hgcal::RecHitTools * tools
Definition: DetId.h:17
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox< 2 > &, const unsigned int, std::vector< std::vector< KDNode >> &) const
~HGCalImagingAlgo() override
std::vector< double > thicknessCorrection_
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_
HGCalImagingAlgo(const edm::ParameterSet &ps)
fixed size matrix
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:289
void reset() override
Hexel(const HGCRecHit &hit, DetId id_in, bool isHalf, float sigmaNoise_in, float thickness_in, const hgcal::RecHitTools *tools_in)
hgcal_clustering::Density Density
double calculateDistanceToHigher(std::vector< KDNode > &) const