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  cluster_offset(0),
71  sigma2(1.0),
72  algoId(algoId_in),
73  dependSensor(dependSensor_in),
74  dEdXweights(dEdXweights_in),
75  thicknessCorrection(thicknessCorrection_in),
76  fcPerMip(fcPerMip_in),
77  fcPerEle(fcPerEle_in),
78  nonAgedNoises(nonAgedNoises_in),
79  noiseMip(noiseMip_in),
80  verbosity(the_verbosity),
82  points(2*(maxlayer+1)),
83  minpos(2*(maxlayer+1),{
84  {0.0f,0.0f}
85  }),
86  maxpos(2*(maxlayer+1),{ {0.0f,0.0f} })
87 {
88 }
89 
90 HGCalImagingAlgo(const std::vector<double>& vecDeltas_in, double kappa_in, double ecut_in,
91  double showerSigma,
92  reco::CaloCluster::AlgoId algoId_in,
93  bool dependSensor_in,
94  const std::vector<double>& dEdXweights_in,
95  const std::vector<double>& thicknessCorrection_in,
96  const std::vector<double>& fcPerMip_in,
97  double fcPerEle_in,
98  const std::vector<double>& nonAgedNoises_in,
99  double noiseMip_in,
100  VerbosityLevel the_verbosity = pERROR) : vecDeltas(vecDeltas_in), kappa(kappa_in),
101  ecut(ecut_in),
102  cluster_offset(0),
103  sigma2(std::pow(showerSigma,2.0)),
104  algoId(algoId_in),
105  dependSensor(dependSensor_in),
106  dEdXweights(dEdXweights_in),
107  thicknessCorrection(thicknessCorrection_in),
108  fcPerMip(fcPerMip_in),
109  fcPerEle(fcPerEle_in),
110  nonAgedNoises(nonAgedNoises_in),
111  noiseMip(noiseMip_in),
112  verbosity(the_verbosity),
114  points(2*(maxlayer+1)),
115  minpos(2*(maxlayer+1),{
116  {0.0f,0.0f}
117  }),
118  maxpos(2*(maxlayer+1),{ {0.0f,0.0f} })
119 {
120 }
121 
123 {
124 }
125 
126 void setVerbosity(VerbosityLevel the_verbosity)
127 {
128  verbosity = the_verbosity;
129 }
130 
131 void populate(const HGCRecHitCollection &hits);
132 // 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
133 // different hit collections (or else use reset)
134 void makeClusters();
135 // this is the method to get the cluster collection out
136 std::vector<reco::BasicCluster> getClusters(bool);
137 // needed to switch between EE and HE with the same algorithm object (to get a single cluster collection)
140 }
141 // use this if you want to reuse the same cluster object but don't want to accumulate clusters (hardly useful?)
142 void reset(){
143  current_v.clear();
144  clusters_v.clear();
145  cluster_offset = 0;
146  for( auto& it: points)
147  {
148  it.clear();
149  std::vector<KDNode>().swap(it);
150  }
151  for(unsigned int i = 0; i < minpos.size(); i++)
152  {
153  minpos[i][0]=0.; minpos[i][1]=0.;
154  maxpos[i][0]=0.; maxpos[i][1]=0.;
155  }
156 }
157 void computeThreshold();
158 
161 
162 //max number of layers
163 static const unsigned int maxlayer = 52;
164 
165 
166 private:
167 // last layer per subdetector
168 static const unsigned int lastLayerEE = 28;
169 static const unsigned int lastLayerFH = 40;
170 // maximum number of wafers per Layer: 666 (V7), 794 (V8)
171 static const unsigned int maxNumberOfWafersPerLayer = 796;
172 
173 // The two parameters used to identify clusters
174 std::vector<double> vecDeltas;
175 double kappa;
176 
177 // The hit energy cutoff
178 double ecut;
179 
180 // The current offset into the temporary cluster structure
181 unsigned int cluster_offset;
182 
183 // for energy sharing
184 double sigma2; // transverse shower size
185 
186 // The vector of clusters
187 std::vector<reco::BasicCluster> clusters_v;
188 
190 
191 // The algo id
193 
194 // various parameters used for calculating the noise levels for a given sensor (and whether to use them)
196 std::vector<double> dEdXweights;
197 std::vector<double> thicknessCorrection;
198 std::vector<double> fcPerMip;
199 double fcPerEle;
200 std::vector<double> nonAgedNoises;
201 double noiseMip;
202 std::vector<std::vector<double> > thresholds;
203 std::vector<std::vector<double> > v_sigmaNoise;
204 
205 // The verbosity level
207 
208 // initialization bool
210 
211 struct Hexel {
212 
213  double x;
214  double y;
215  double z;
217  double weight;
218  double fraction;
220  double rho;
221  double delta;
223  bool isBorder;
224  bool isHalo;
226  float sigmaNoise;
227  float thickness;
229 
230  Hexel(const HGCRecHit &hit, DetId id_in, bool isHalf, float sigmaNoise_in, float thickness_in, const hgcal::RecHitTools *tools_in) :
231  isHalfCell(isHalf),
232  weight(0.), fraction(1.0), detid(id_in), rho(0.), delta(0.),
233  nearestHigher(-1), isBorder(false), isHalo(false),
234  clusterIndex(-1), sigmaNoise(sigmaNoise_in), thickness(thickness_in),
235  tools(tools_in)
236  {
237  const GlobalPoint position( std::move( tools->getPosition( detid ) ) );
238  weight = hit.energy();
239  x = position.x();
240  y = position.y();
241  z = position.z();
242  }
243  Hexel() :
244  x(0.),y(0.),z(0.),isHalfCell(false),
245  weight(0.), fraction(1.0), detid(), rho(0.), delta(0.),
246  nearestHigher(-1), isBorder(false), isHalo(false),
247  clusterIndex(-1),
248  sigmaNoise(0.),
249  thickness(0.),
250  tools(0)
251  {
252  }
253  bool operator > (const Hexel& rhs) const {
254  return (rho > rhs.rho);
255  }
256 
257 };
258 
261 
262 
263 // A vector of vectors of KDNodes holding an Hexel in the clusters - to be used to build CaloClusters of DetIds
264 std::vector< std::vector<KDNode> > current_v;
265 
266 std::vector<size_t> sort_by_delta(const std::vector<KDNode> &v){
267  std::vector<size_t> idx(v.size());
268  std::iota (std::begin(idx), std::end(idx), 0);
269  sort(idx.begin(), idx.end(),
270  [&v](size_t i1, size_t i2) {
271  return v[i1].data.delta > v[i2].data.delta;
272  });
273  return idx;
274 }
275 
276 std::vector<std::vector<KDNode> > points; //a vector of vectors of hexels, one for each layer
277 //@@EM todo: the number of layers should be obtained programmatically - the range is 1-n instead of 0-n-1...
278 
279 std::vector<std::array<float,2> > minpos;
280 std::vector<std::array<float,2> > maxpos;
281 
282 
283 //these functions should be in a helper class.
284 inline double distance2(const Hexel &pt1, const Hexel &pt2) { //distance squared
285  const double dx = pt1.x - pt2.x;
286  const double dy = pt1.y - pt2.y;
287  return (dx*dx + dy*dy);
288 } //distance squaredq
289 inline double distance(const Hexel &pt1, const Hexel &pt2) { //2-d distance on the layer (x-y)
290  return std::sqrt(distance2(pt1,pt2));
291 }
292 double calculateLocalDensity(std::vector<KDNode> &, KDTree &, const unsigned int); //return max density
293 double calculateDistanceToHigher(std::vector<KDNode> &);
294 int findAndAssignClusters(std::vector<KDNode> &, KDTree &, double, KDTreeBox &, const unsigned int);
295 math::XYZPoint calculatePosition(std::vector<KDNode> &);
296 
297 // attempt to find subclusters within a given set of hexels
298 std::vector<unsigned> findLocalMaximaInCluster(const std::vector<KDNode>&);
299 math::XYZPoint calculatePositionWithFraction(const std::vector<KDNode>&, const std::vector<double>&);
300 double calculateEnergyWithFraction(const std::vector<KDNode>&, const std::vector<double>&);
301 // outputs
302 void shareEnergy(const std::vector<KDNode>&,
303  const std::vector<unsigned>&,
304  std::vector<std::vector<double> >&);
305 };
306 
307 #endif
double calculateEnergyWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
hgcal::RecHitTools rhtools_
std::vector< std::vector< KDNode > > current_v
math::XYZPoint calculatePositionWithFraction(const std::vector< KDNode > &, const std::vector< double > &)
double distance(const Hexel &pt1, const Hexel &pt2)
std::vector< double > dEdXweights
reco::CaloCluster::AlgoId algoId
double calculateDistanceToHigher(std::vector< KDNode > &)
static const unsigned int maxlayer
std::vector< std::array< float, 2 > > minpos
void getEventSetup(const edm::EventSetup &es)
static const unsigned int maxNumberOfWafersPerLayer
std::vector< double > vecDeltas
std::vector< std::vector< double > > thresholds
std::vector< std::array< float, 2 > > maxpos
T y() const
Definition: PV3DBase.h:63
math::XYZPoint Point
point in the space
std::vector< double > thicknessCorrection
std::vector< double > nonAgedNoises
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:66
std::vector< size_t > sort_by_delta(const std::vector< KDNode > &v)
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
float energy() const
Definition: CaloRecHit.h:17
T sqrt(T t)
Definition: SSEVec.h:18
static const unsigned int lastLayerFH
unsigned int cluster_offset
T z() const
Definition: PV3DBase.h:64
VerbosityLevel verbosity
#define end
Definition: vmac.h:37
std::vector< size_t > sorted_indices(const std::vector< T > &v)
std::vector< reco::BasicCluster > getClusters(bool)
std::vector< unsigned > findLocalMaximaInCluster(const std::vector< KDNode > &)
double calculateLocalDensity(std::vector< KDNode > &, KDTree &, const unsigned int)
const hgcal::RecHitTools * tools
Definition: DetId.h:18
int findAndAssignClusters(std::vector< KDNode > &, KDTree &, double, KDTreeBox &, const unsigned int)
std::vector< std::vector< KDNode > > points
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:73
math::XYZPoint calculatePosition(std::vector< KDNode > &)
std::vector< std::vector< double > > v_sigmaNoise
void setVerbosity(VerbosityLevel the_verbosity)
KDTreeLinkerAlgo< Hexel, 2 > KDTree
fixed size matrix
#define begin
Definition: vmac.h:30
double distance2(const Hexel &pt1, const Hexel &pt2)
static int position[264][3]
Definition: ReadPGInfo.cc:509
bool operator>(const Hexel &rhs) const
virtual ~HGCalImagingAlgo()
T x() const
Definition: PV3DBase.h:62
Hexel(const HGCRecHit &hit, DetId id_in, bool isHalf, float sigmaNoise_in, float thickness_in, const hgcal::RecHitTools *tools_in)
KDTreeNodeInfoT< Hexel, 2 > KDNode
std::vector< double > fcPerMip
void populate(const HGCRecHitCollection &hits)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
def move(src, dest)
Definition: eostools.py:510
std::vector< reco::BasicCluster > clusters_v