CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  iC),
38  thresholdW0_(ps.getParameter<std::vector<double>>("thresholdW0")),
39  positionDeltaRho_c_(ps.getParameter<std::vector<double>>("positionDeltaRho_c")),
40  vecDeltas_(ps.getParameter<std::vector<double>>("deltac")),
41  kappa_(ps.getParameter<double>("kappa")),
42  ecut_(ps.getParameter<double>("ecut")),
43  sigma2_(1.0),
44  dependSensor_(ps.getParameter<bool>("dependSensor")),
45  dEdXweights_(ps.getParameter<std::vector<double>>("dEdXweights")),
46  thicknessCorrection_(ps.getParameter<std::vector<double>>("thicknessCorrection")),
47  fcPerMip_(ps.getParameter<std::vector<double>>("fcPerMip")),
48  fcPerEle_(ps.getParameter<double>("fcPerEle")),
49  nonAgedNoises_(ps.getParameter<edm::ParameterSet>("noises").getParameter<std::vector<double>>("values")),
50  noiseMip_(ps.getParameter<edm::ParameterSet>("noiseMip").getParameter<double>("noise_MIP")),
52 
53  ~HGCalImagingAlgo() override {}
54 
55  void getEventSetupPerAlgorithm(const edm::EventSetup &es) override;
56 
57  void populate(const HGCRecHitCollection &hits) override;
58  // 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
59  // different hit collections (or else use reset)
60 
61  void makeClusters() override;
62 
63  // this is the method to get the cluster collection out
64  std::vector<reco::BasicCluster> getClusters(bool) override;
65 
66  // use this if you want to reuse the same cluster object but don't want to accumulate clusters (hardly useful?)
67  void reset() override {
68  clusters_v_.clear();
69  clusters_v_.shrink_to_fit();
70  layerClustersPerLayer_.clear();
71  layerClustersPerLayer_.shrink_to_fit();
72  for (auto &it : points_) {
73  it.clear();
74  it.shrink_to_fit();
75  std::vector<KDNode>().swap(it);
76  }
77  for (unsigned int i = 0; i < minpos_.size(); i++) {
78  minpos_[i][0] = 0.;
79  minpos_[i][1] = 0.;
80  maxpos_[i][0] = 0.;
81  maxpos_[i][1] = 0.;
82  }
83  }
84  void computeThreshold();
85 
86  //getDensity
87  Density getDensity() override;
88 
90  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
91  iDesc.add<std::vector<double>>("positionDeltaRho_c", {1.3, 1.3, 1.3});
92  iDesc.add<std::vector<double>>("deltac",
93  {
94  2.0,
95  2.0,
96  5.0,
97  });
98  iDesc.add<bool>("dependSensor", true);
99  iDesc.add<double>("ecut", 3.0);
100  iDesc.add<double>("kappa", 9.0);
101  iDesc.addUntracked<unsigned int>("verbosity", 3);
102  iDesc.add<std::vector<double>>("dEdXweights", {});
103  iDesc.add<std::vector<double>>("thicknessCorrection", {});
104  iDesc.add<std::vector<double>>("fcPerMip", {});
105  iDesc.add<double>("fcPerEle", 0.0);
106  edm::ParameterSetDescription descNestedNoises;
107  descNestedNoises.add<std::vector<double>>("values", {});
108  iDesc.add<edm::ParameterSetDescription>("noises", descNestedNoises);
109  edm::ParameterSetDescription descNestedNoiseMIP;
110  descNestedNoiseMIP.add<bool>("scaleByDose", false);
111  descNestedNoiseMIP.add<double>("scaleByDoseFactor", 1.);
112  iDesc.add<edm::ParameterSetDescription>("scaleByDose", descNestedNoiseMIP);
113  descNestedNoiseMIP.add<std::string>("doseMap", "");
114  iDesc.add<edm::ParameterSetDescription>("doseMap", descNestedNoiseMIP);
115  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
116  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
117  }
118 
121 
122 private:
123  // To compute the cluster position
124  std::vector<double> thresholdW0_;
125  std::vector<double> positionDeltaRho_c_;
126 
127  // The two parameters used to identify clusters
128  std::vector<double> vecDeltas_;
129  double kappa_;
130 
131  // The hit energy cutoff
132  double ecut_;
133 
134  // for energy sharing
135  double sigma2_; // transverse shower size
136 
137  // The vector of clusters
138  std::vector<reco::BasicCluster> clusters_v_;
139 
140  // For keeping the density per hit
142 
143  // various parameters used for calculating the noise levels for a given sensor (and whether to use them)
145  std::vector<double> dEdXweights_;
146  std::vector<double> thicknessCorrection_;
147  std::vector<double> fcPerMip_;
148  double fcPerEle_;
149  std::vector<double> nonAgedNoises_;
150  double noiseMip_;
151  std::vector<std::vector<double>> thresholds_;
152  std::vector<std::vector<double>> sigmaNoise_;
153 
154  // initialization bool
156 
157  struct Hexel {
158  double x;
159  double y;
160  double z;
162  double weight;
163  double fraction;
165  double rho;
166  double delta;
168  bool isBorder;
169  bool isHalo;
171  float sigmaNoise;
172  float thickness;
174 
176  DetId id_in,
177  bool isHalf,
178  float sigmaNoise_in,
179  float thickness_in,
180  const hgcal::RecHitTools *tools_in)
181  : isHalfCell(isHalf),
182  weight(0.),
183  fraction(1.0),
184  detid(id_in),
185  rho(0.),
186  delta(0.),
187  nearestHigher(-1),
188  isBorder(false),
189  isHalo(false),
190  clusterIndex(-1),
191  sigmaNoise(sigmaNoise_in),
192  thickness(thickness_in),
193  tools(tools_in) {
195  weight = hit.energy();
196  x = position.x();
197  y = position.y();
198  z = position.z();
199  }
201  : x(0.),
202  y(0.),
203  z(0.),
204  isHalfCell(false),
205  weight(0.),
206  fraction(1.0),
207  detid(),
208  rho(0.),
209  delta(0.),
210  nearestHigher(-1),
211  isBorder(false),
212  isHalo(false),
213  clusterIndex(-1),
214  sigmaNoise(0.),
215  thickness(0.),
216  tools(nullptr) {}
217  bool operator>(const Hexel &rhs) const { return (rho > rhs.rho); }
218  };
219 
222 
223  std::vector<std::vector<std::vector<KDNode>>> layerClustersPerLayer_;
224 
225  std::vector<size_t> sort_by_delta(const std::vector<KDNode> &v) const {
226  std::vector<size_t> idx(v.size());
227  std::iota(std::begin(idx), std::end(idx), 0);
228  sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1].data.delta > v[i2].data.delta; });
229  return idx;
230  }
231 
232  std::vector<std::vector<KDNode>> points_; //a vector of vectors of hexels, one for each layer
233  //@@EM todo: the number of layers should be obtained programmatically - the range is 1-n instead of 0-n-1...
234 
235  std::vector<std::array<float, 2>> minpos_;
236  std::vector<std::array<float, 2>> maxpos_;
237 
238  //these functions should be in a helper class.
239  inline double distance2(const Hexel &pt1, const Hexel &pt2) const { //distance squared
240  const double dx = pt1.x - pt2.x;
241  const double dy = pt1.y - pt2.y;
242  return (dx * dx + dy * dy);
243  } //distance squaredq
244  inline double distance(const Hexel &pt1, const Hexel &pt2) const { //2-d distance on the layer (x-y)
245  return std::sqrt(distance2(pt1, pt2));
246  }
247  double calculateLocalDensity(std::vector<KDNode> &, KDTree &, const unsigned int) const; //return max density
248  double calculateDistanceToHigher(std::vector<KDNode> &) const;
249  int findAndAssignClusters(std::vector<KDNode> &,
250  KDTree &,
251  double,
252  KDTreeBox<2> &,
253  const unsigned int,
254  std::vector<std::vector<KDNode>> &) const;
255  math::XYZPoint calculatePosition(std::vector<KDNode> &) const;
256 
257  //For keeping the density information
258  void setDensity(const std::vector<KDNode> &nd);
259 
260  // attempt to find subclusters within a given set of hexels
261  std::vector<unsigned> findLocalMaximaInCluster(const std::vector<KDNode> &);
262  math::XYZPoint calculatePositionWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
263  double calculateEnergyWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
264  // outputs
265  void shareEnergy(const std::vector<KDNode> &, const std::vector<unsigned> &, std::vector<std::vector<double>> &);
266 };
267 
268 #endif
constexpr float energy() const
Definition: CaloRecHit.h:29
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox< 2 > &, const unsigned int, std::vector< std::vector< KDNode >> &) const
void getEventSetupPerAlgorithm(const edm::EventSetup &es) override
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
KDTreeNodeInfo< Hexel > KDNode
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int) const
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 >> &)
KDTreeLinkerAlgo< Hexel > KDTree
math::XYZPoint Point
point in the space
std::vector< reco::BasicCluster > getClusters(bool) override
std::vector< std::array< float, 2 > > minpos_
std::vector< double > positionDeltaRho_c_
double distance2(const Hexel &pt1, const Hexel &pt2) const
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< std::array< float, 2 > > maxpos_
double distance(const Hexel &pt1, const Hexel &pt2) const
T sqrt(T t)
Definition: SSEVec.h:19
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)
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
~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
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_
double calculateDistanceToHigher(std::vector< KDNode > &) const
static int position[264][3]
Definition: ReadPGInfo.cc:289
string end
Definition: dataset.py:937
math::XYZPoint calculatePosition(std::vector< KDNode > &) const
bool operator>(const Hexel &rhs) const
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
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const
HGCalImagingAlgo(const edm::ParameterSet &ps, edm::ConsumesCollector iC)