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  layerClustersPerLayer_.clear();
69  for (auto &it : points_) {
70  it.clear();
71  std::vector<KDNode>().swap(it);
72  }
73  for (unsigned int i = 0; i < minpos_.size(); i++) {
74  minpos_[i][0] = 0.;
75  minpos_[i][1] = 0.;
76  maxpos_[i][0] = 0.;
77  maxpos_[i][1] = 0.;
78  }
79  }
80  void computeThreshold();
81 
82  //getDensity
83  Density getDensity() override;
84 
86  iDesc.add<std::vector<double>>("thresholdW0", {2.9, 2.9, 2.9});
87  iDesc.add<std::vector<double>>("positionDeltaRho_c", {1.3, 1.3, 1.3});
88  iDesc.add<std::vector<double>>("deltac",
89  {
90  2.0,
91  2.0,
92  5.0,
93  });
94  iDesc.add<bool>("dependSensor", true);
95  iDesc.add<double>("ecut", 3.0);
96  iDesc.add<double>("kappa", 9.0);
97  iDesc.addUntracked<unsigned int>("verbosity", 3);
98  iDesc.add<std::vector<double>>("dEdXweights", {});
99  iDesc.add<std::vector<double>>("thicknessCorrection", {});
100  iDesc.add<std::vector<double>>("fcPerMip", {});
101  iDesc.add<double>("fcPerEle", 0.0);
102  edm::ParameterSetDescription descNestedNoises;
103  descNestedNoises.add<std::vector<double>>("values", {});
104  iDesc.add<edm::ParameterSetDescription>("noises", descNestedNoises);
105  edm::ParameterSetDescription descNestedNoiseMIP;
106  descNestedNoiseMIP.add<bool>("scaleByDose", false);
107  iDesc.add<edm::ParameterSetDescription>("scaleByDose", descNestedNoiseMIP);
108  descNestedNoiseMIP.add<std::string>("doseMap", "");
109  iDesc.add<edm::ParameterSetDescription>("doseMap", descNestedNoiseMIP);
110  descNestedNoiseMIP.add<double>("noise_MIP", 1. / 100.);
111  iDesc.add<edm::ParameterSetDescription>("noiseMip", descNestedNoiseMIP);
112  }
113 
116 
117 private:
118  // To compute the cluster position
119  std::vector<double> thresholdW0_;
120  std::vector<double> positionDeltaRho_c_;
121 
122  // The two parameters used to identify clusters
123  std::vector<double> vecDeltas_;
124  double kappa_;
125 
126  // The hit energy cutoff
127  double ecut_;
128 
129  // for energy sharing
130  double sigma2_; // transverse shower size
131 
132  // The vector of clusters
133  std::vector<reco::BasicCluster> clusters_v_;
134 
135  // For keeping the density per hit
137 
138  // various parameters used for calculating the noise levels for a given sensor (and whether to use them)
140  std::vector<double> dEdXweights_;
141  std::vector<double> thicknessCorrection_;
142  std::vector<double> fcPerMip_;
143  double fcPerEle_;
144  std::vector<double> nonAgedNoises_;
145  double noiseMip_;
146  std::vector<std::vector<double>> thresholds_;
147  std::vector<std::vector<double>> sigmaNoise_;
148 
149  // initialization bool
151 
152  struct Hexel {
153  double x;
154  double y;
155  double z;
157  double weight;
158  double fraction;
160  double rho;
161  double delta;
163  bool isBorder;
164  bool isHalo;
166  float sigmaNoise;
167  float thickness;
169 
171  DetId id_in,
172  bool isHalf,
173  float sigmaNoise_in,
174  float thickness_in,
175  const hgcal::RecHitTools *tools_in)
176  : isHalfCell(isHalf),
177  weight(0.),
178  fraction(1.0),
179  detid(id_in),
180  rho(0.),
181  delta(0.),
182  nearestHigher(-1),
183  isBorder(false),
184  isHalo(false),
185  clusterIndex(-1),
186  sigmaNoise(sigmaNoise_in),
187  thickness(thickness_in),
188  tools(tools_in) {
189  const GlobalPoint position(tools->getPosition(detid));
190  weight = hit.energy();
191  x = position.x();
192  y = position.y();
193  z = position.z();
194  }
196  : x(0.),
197  y(0.),
198  z(0.),
199  isHalfCell(false),
200  weight(0.),
201  fraction(1.0),
202  detid(),
203  rho(0.),
204  delta(0.),
205  nearestHigher(-1),
206  isBorder(false),
207  isHalo(false),
208  clusterIndex(-1),
209  sigmaNoise(0.),
210  thickness(0.),
211  tools(nullptr) {}
212  bool operator>(const Hexel &rhs) const { return (rho > rhs.rho); }
213  };
214 
217 
218  std::vector<std::vector<std::vector<KDNode>>> layerClustersPerLayer_;
219 
220  std::vector<size_t> sort_by_delta(const std::vector<KDNode> &v) const {
221  std::vector<size_t> idx(v.size());
222  std::iota(std::begin(idx), std::end(idx), 0);
223  sort(idx.begin(), idx.end(), [&v](size_t i1, size_t i2) { return v[i1].data.delta > v[i2].data.delta; });
224  return idx;
225  }
226 
227  std::vector<std::vector<KDNode>> points_; //a vector of vectors of hexels, one for each layer
228  //@@EM todo: the number of layers should be obtained programmatically - the range is 1-n instead of 0-n-1...
229 
230  std::vector<std::array<float, 2>> minpos_;
231  std::vector<std::array<float, 2>> maxpos_;
232 
233  //these functions should be in a helper class.
234  inline double distance2(const Hexel &pt1, const Hexel &pt2) const { //distance squared
235  const double dx = pt1.x - pt2.x;
236  const double dy = pt1.y - pt2.y;
237  return (dx * dx + dy * dy);
238  } //distance squaredq
239  inline double distance(const Hexel &pt1, const Hexel &pt2) const { //2-d distance on the layer (x-y)
240  return std::sqrt(distance2(pt1, pt2));
241  }
242  double calculateLocalDensity(std::vector<KDNode> &, KDTree &, const unsigned int) const; //return max density
243  double calculateDistanceToHigher(std::vector<KDNode> &) const;
244  int findAndAssignClusters(std::vector<KDNode> &,
245  KDTree &,
246  double,
247  KDTreeBox<2> &,
248  const unsigned int,
249  std::vector<std::vector<KDNode>> &) const;
250  math::XYZPoint calculatePosition(std::vector<KDNode> &) const;
251 
252  //For keeping the density information
253  void setDensity(const std::vector<KDNode> &nd);
254 
255  // attempt to find subclusters within a given set of hexels
256  std::vector<unsigned> findLocalMaximaInCluster(const std::vector<KDNode> &);
257  math::XYZPoint calculatePositionWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
258  double calculateEnergyWithFraction(const std::vector<KDNode> &, const std::vector<double> &);
259  // outputs
260  void shareEnergy(const std::vector<KDNode> &, const std::vector<unsigned> &, std::vector<std::vector<double>> &);
261 };
262 
263 #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 >> &)
#define nullptr
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
hgcal_clustering::Density Density
Definition: HGCalCLUEAlgo.h:27
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_
#define end
Definition: vmac.h:39
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_
HGCalImagingAlgo(const edm::ParameterSet &ps)
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
double calculateDistanceToHigher(std::vector< KDNode > &) const
static int position[264][3]
Definition: ReadPGInfo.cc:289
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)
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const