CMS 3D CMS Logo

HGCalImagingAlgo.h
Go to the documentation of this file.
1 #ifndef RecoLocalCalo_HGCalRecAlgos_HGCalImagingAlgo_h
2 #define RecoLocalCalo_HGCalRecAlgos_HGCalImagingAlgo_h
3 
13 
16 
18 
19 // C/C++ headers
20 #include <string>
21 #include <vector>
22 #include <set>
23 #include <numeric>
24 
25 #include "KDTreeLinkerAlgoT.h"
26 
27 
28 template <typename T>
29 std::vector<size_t> sorted_indices(const std::vector<T> &v) {
30 
31  // initialize original index locations
32  std::vector<size_t> idx(v.size());
33  std::iota (std::begin(idx), std::end(idx), 0);
34 
35  // sort indices based on comparing values in v
36  std::sort(idx.begin(), idx.end(),
37  [&v](size_t i1, size_t i2) {
38  return v[i1] > v[i2];
39  });
40 
41  return idx;
42 }
43 
45 {
46 
47 
48 public:
49 
50 enum VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 };
51 
53  sigma2_(1.0),
54  algoId_(reco::CaloCluster::undefined),
56 }
57 
58 HGCalImagingAlgo(const std::vector<double>& vecDeltas_in, double kappa_in, double ecut_in,
59  reco::CaloCluster::AlgoId algoId_in,
60  bool dependSensor_in,
61  const std::vector<double>& dEdXweights_in,
62  const std::vector<double>& thicknessCorrection_in,
63  const std::vector<double>& fcPerMip_in,
64  double fcPerEle_in,
65  const std::vector<double>& nonAgedNoises_in,
66  double noiseMip_in,
67  VerbosityLevel the_verbosity = pERROR) :
68  vecDeltas_(vecDeltas_in), kappa_(kappa_in),
69  ecut_(ecut_in),
70  sigma2_(1.0),
71  algoId_(algoId_in),
72  dependSensor_(dependSensor_in),
73  dEdXweights_(dEdXweights_in),
74  thicknessCorrection_(thicknessCorrection_in),
75  fcPerMip_(fcPerMip_in),
76  fcPerEle_(fcPerEle_in),
77  nonAgedNoises_(nonAgedNoises_in),
78  noiseMip_(noiseMip_in),
79  verbosity_(the_verbosity),
81  points_(2*(maxlayer+1)),
82  minpos_(2*(maxlayer+1),{
83  {0.0f,0.0f}
84  }),
85  maxpos_(2*(maxlayer+1),{ {0.0f,0.0f} })
86 {
87 }
88 
89 HGCalImagingAlgo(const std::vector<double>& vecDeltas_in, double kappa_in, double ecut_in,
90  double showerSigma,
91  reco::CaloCluster::AlgoId algoId_in,
92  bool dependSensor_in,
93  const std::vector<double>& dEdXweights_in,
94  const std::vector<double>& thicknessCorrection_in,
95  const std::vector<double>& fcPerMip_in,
96  double fcPerEle_in,
97  const std::vector<double>& nonAgedNoises_in,
98  double noiseMip_in,
99  VerbosityLevel the_verbosity = pERROR) : vecDeltas_(vecDeltas_in), kappa_(kappa_in),
100  ecut_(ecut_in),
101  sigma2_(std::pow(showerSigma,2.0)),
102  algoId_(algoId_in),
103  dependSensor_(dependSensor_in),
104  dEdXweights_(dEdXweights_in),
105  thicknessCorrection_(thicknessCorrection_in),
106  fcPerMip_(fcPerMip_in),
107  fcPerEle_(fcPerEle_in),
108  nonAgedNoises_(nonAgedNoises_in),
109  noiseMip_(noiseMip_in),
110  verbosity_(the_verbosity),
112  points_(2*(maxlayer+1)),
113  minpos_(2*(maxlayer+1),{
114  {0.0f,0.0f}
115  }),
116  maxpos_(2*(maxlayer+1),{ {0.0f,0.0f} })
117 {
118 }
119 
121 {
122 }
123 
124 void setVerbosity(VerbosityLevel the_verbosity)
125 {
126  verbosity_ = the_verbosity;
127 }
128 
129 void populate(const HGCRecHitCollection &hits);
130 // 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
131 // different hit collections (or else use reset)
132 void makeClusters();
133 // this is the method to get the cluster collection out
134 std::vector<reco::BasicCluster> getClusters(bool);
135 // needed to switch between EE and HE with the same algorithm object (to get a single cluster collection)
138 }
139 // use this if you want to reuse the same cluster object but don't want to accumulate clusters (hardly useful?)
140 void reset(){
141  clusters_v_.clear();
142  layerClustersPerLayer_.clear();
143  for( auto& it: points_)
144  {
145  it.clear();
146  std::vector<KDNode>().swap(it);
147  }
148  for(unsigned int i = 0; i < minpos_.size(); i++)
149  {
150  minpos_[i][0]=0.; minpos_[i][1]=0.;
151  maxpos_[i][0]=0.; maxpos_[i][1]=0.;
152  }
153 }
154 void computeThreshold();
155 
158 
159 //max number of layers
160 static const unsigned int maxlayer = 52;
161 
162 
163 private:
164 // last layer per subdetector
165 static const unsigned int lastLayerEE = 28;
166 static const unsigned int lastLayerFH = 40;
167 
168 // The two parameters used to identify clusters
169 std::vector<double> vecDeltas_;
170 double kappa_;
171 
172 // The hit energy cutoff
173 double ecut_;
174 
175 // for energy sharing
176 double sigma2_; // transverse shower size
177 
178 // The vector of clusters
179 std::vector<reco::BasicCluster> clusters_v_;
180 
182 
183 // The algo id
185 
186 // various parameters used for calculating the noise levels for a given sensor (and whether to use them)
188 std::vector<double> dEdXweights_;
189 std::vector<double> thicknessCorrection_;
190 std::vector<double> fcPerMip_;
191 double fcPerEle_;
192 std::vector<double> nonAgedNoises_;
193 double noiseMip_;
194 std::vector<std::vector<double> > thresholds_;
195 std::vector<std::vector<double> > v_sigmaNoise_;
196 
197 // The verbosity level
199 
200 // initialization bool
202 
203 struct Hexel {
204 
205  double x;
206  double y;
207  double z;
209  double weight;
210  double fraction;
212  double rho;
213  double delta;
215  bool isBorder;
216  bool isHalo;
218  float sigmaNoise;
219  float thickness;
221 
222  Hexel(const HGCRecHit &hit, DetId id_in, bool isHalf, float sigmaNoise_in, float thickness_in, const hgcal::RecHitTools *tools_in) :
223  isHalfCell(isHalf),
224  weight(0.), fraction(1.0), detid(id_in), rho(0.), delta(0.),
225  nearestHigher(-1), isBorder(false), isHalo(false),
226  clusterIndex(-1), sigmaNoise(sigmaNoise_in), thickness(thickness_in),
227  tools(tools_in)
228  {
229  const GlobalPoint position( tools->getPosition( detid ) );
230  weight = hit.energy();
231  x = position.x();
232  y = position.y();
233  z = position.z();
234  }
235  Hexel() :
236  x(0.),y(0.),z(0.),isHalfCell(false),
237  weight(0.), fraction(1.0), detid(), rho(0.), delta(0.),
238  nearestHigher(-1), isBorder(false), isHalo(false),
239  clusterIndex(-1),
240  sigmaNoise(0.),
241  thickness(0.),
242  tools(nullptr)
243  {
244  }
245  bool operator > (const Hexel& rhs) const {
246  return (rho > rhs.rho);
247  }
248 
249 };
250 
253 
254 
255 std::vector<std::vector<std::vector< KDNode> > > layerClustersPerLayer_;
256 
257 std::vector<size_t> sort_by_delta(const std::vector<KDNode> &v) const {
258  std::vector<size_t> idx(v.size());
259  std::iota (std::begin(idx), std::end(idx), 0);
260  sort(idx.begin(), idx.end(),
261  [&v](size_t i1, size_t i2) {
262  return v[i1].data.delta > v[i2].data.delta;
263  });
264  return idx;
265 }
266 
267 std::vector<std::vector<KDNode> > points_; //a vector of vectors of hexels, one for each layer
268 //@@EM todo: the number of layers should be obtained programmatically - the range is 1-n instead of 0-n-1...
269 
270 std::vector<std::array<float,2> > minpos_;
271 std::vector<std::array<float,2> > maxpos_;
272 
273 
274 //these functions should be in a helper class.
275 inline double distance2(const Hexel &pt1, const Hexel &pt2) const{ //distance squared
276  const double dx = pt1.x - pt2.x;
277  const double dy = pt1.y - pt2.y;
278  return (dx*dx + dy*dy);
279 } //distance squaredq
280 inline double distance(const Hexel &pt1, const Hexel &pt2) const{ //2-d distance on the layer (x-y)
281  return std::sqrt(distance2(pt1,pt2));
282 }
283 double calculateLocalDensity(std::vector<KDNode> &, KDTree &, const unsigned int) const; //return max density
284 double calculateDistanceToHigher(std::vector<KDNode> &) const;
285 int findAndAssignClusters(std::vector<KDNode> &, KDTree &, double, KDTreeBox &, const unsigned int, std::vector<std::vector<KDNode> >&) const;
286 math::XYZPoint calculatePosition(std::vector<KDNode> &) const;
287 
288 // attempt to find subclusters within a given set of hexels
289 std::vector<unsigned> findLocalMaximaInCluster(const std::vector<KDNode>&);
290 math::XYZPoint calculatePositionWithFraction(const std::vector<KDNode>&, const std::vector<double>&);
291 double calculateEnergyWithFraction(const std::vector<KDNode>&, const std::vector<double>&);
292 // outputs
293 void shareEnergy(const std::vector<KDNode>&,
294  const std::vector<unsigned>&,
295  std::vector<std::vector<double> >&);
296 };
297 
298 #endif
constexpr float energy() const
Definition: CaloRecHit.h:31
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
hgcal::RecHitTools rhtools_
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
static const unsigned int maxlayer
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int) const
void getEventSetup(const edm::EventSetup &es)
std::vector< std::vector< double > > thresholds_
std::vector< std::vector< double > > v_sigmaNoise_
math::XYZPoint Point
point in the space
#define nullptr
double distance2(const Hexel &pt1, const Hexel &pt2) const
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:73
HGCalImagingAlgo(const std::vector< double > &vecDeltas_in, double kappa_in, double ecut_in, reco::CaloCluster::AlgoId algoId_in, bool dependSensor_in, const std::vector< double > &dEdXweights_in, const std::vector< double > &thicknessCorrection_in, const std::vector< double > &fcPerMip_in, double fcPerEle_in, const std::vector< double > &nonAgedNoises_in, double noiseMip_in, VerbosityLevel the_verbosity=pERROR)
void shareEnergy(const std::vector< KDNode > &, const std::vector< unsigned > &, std::vector< std::vector< double > > &)
static const unsigned int lastLayerEE
double distance(const Hexel &pt1, const Hexel &pt2) const
T sqrt(T t)
Definition: SSEVec.h:18
static const unsigned int lastLayerFH
std::vector< reco::BasicCluster > clusters_v_
reco::CaloCluster::AlgoId algoId_
std::vector< std::vector< KDNode > > points_
std::vector< double > vecDeltas_
std::vector< double > dEdXweights_
#define end
Definition: vmac.h:39
std::vector< size_t > sorted_indices(const std::vector< T > &v)
std::vector< std::vector< std::vector< KDNode > > > layerClustersPerLayer_
std::vector< reco::BasicCluster > getClusters(bool)
std::vector< unsigned > findLocalMaximaInCluster(const std::vector< KDNode > &)
std::vector< std::array< float, 2 > > minpos_
const hgcal::RecHitTools * tools
Definition: DetId.h:18
std::vector< double > thicknessCorrection_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
HGCalImagingAlgo(const std::vector< double > &vecDeltas_in, double kappa_in, double ecut_in, double showerSigma, reco::CaloCluster::AlgoId algoId_in, bool dependSensor_in, const std::vector< double > &dEdXweights_in, const std::vector< double > &thicknessCorrection_in, const std::vector< double > &fcPerMip_in, double fcPerEle_in, const std::vector< double > &nonAgedNoises_in, double noiseMip_in, VerbosityLevel the_verbosity=pERROR)
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
std::vector< double > fcPerMip_
std::vector< double > nonAgedNoises_
void setVerbosity(VerbosityLevel the_verbosity)
KDTreeLinkerAlgo< Hexel, 2 > KDTree
fixed size matrix
#define begin
Definition: vmac.h:32
double calculateDistanceToHigher(std::vector< KDNode > &) const
static int position[264][3]
Definition: ReadPGInfo.cc:509
VerbosityLevel verbosity_
math::XYZPoint calculatePosition(std::vector< KDNode > &) const
bool operator>(const Hexel &rhs) const
virtual ~HGCalImagingAlgo()
std::vector< std::array< float, 2 > > maxpos_
Hexel(const HGCRecHit &hit, DetId id_in, bool isHalf, float sigmaNoise_in, float thickness_in, const hgcal::RecHitTools *tools_in)
KDTreeNodeInfoT< Hexel, 2 > KDNode
void populate(const HGCRecHitCollection &hits)
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox &, const unsigned int, std::vector< std::vector< KDNode > > &) const