CMS 3D CMS Logo

PFMultiDepthClusterizer.cc
Go to the documentation of this file.
7 #include "Math/GenVector/VectorUtil.h"
11 #include "vdt/vdtMath.h"
12 
13 #include <iterator>
14 #include <unordered_map>
15 
18 
19 public:
21 
22  ~PFMultiDepthClusterizer() override = default;
23  PFMultiDepthClusterizer(const B2DGPF&) = delete;
24  B2DGPF& operator=(const B2DGPF&) = delete;
25 
26  void update(const edm::EventSetup& es) override { _allCellsPosCalc->update(es); }
27 
29  const std::vector<bool>&,
31  const HcalPFCuts*) override;
32 
33 private:
34  std::unique_ptr<PFCPositionCalculatorBase> _allCellsPosCalc;
35  double nSigmaEta_;
36  double nSigmaPhi_;
37 
38  class ClusterLink {
39  public:
40  ClusterLink(unsigned int i, unsigned int j, double DR, int DZ, double energy) {
41  from_ = i;
42  to_ = j;
43  linkDR_ = DR;
44  linkDZ_ = DZ;
45  linkE_ = energy;
46  }
47 
48  ~ClusterLink() = default;
49 
50  unsigned int from() const { return from_; }
51  unsigned int to() const { return to_; }
52  double dR() const { return linkDR_; }
53  int dZ() const { return linkDZ_; }
54  double energy() const { return linkE_; }
55 
56  private:
57  unsigned int from_;
58  unsigned int to_;
59  double linkDR_;
60  int linkDZ_;
61  double linkE_;
62  };
63 
64  void calculateShowerShapes(const reco::PFClusterCollection&, std::vector<double>&, std::vector<double>&);
65  std::vector<ClusterLink> link(const reco::PFClusterCollection&,
66  const std::vector<double>&,
67  const std::vector<double>&);
68  std::vector<ClusterLink> prune(std::vector<ClusterLink>&, std::vector<bool>& linkedClusters);
69 
71  unsigned int point,
72  std::vector<bool>& mask,
74  const std::vector<ClusterLink>& links);
75 
77 };
78 
80 
82  : PFClusterBuilderBase(conf, cc) {
83  if (conf.exists("allCellsPositionCalc")) {
84  const edm::ParameterSet& acConf = conf.getParameterSet("allCellsPositionCalc");
85  const std::string& algoac = acConf.getParameter<std::string>("algoName");
86  _allCellsPosCalc = PFCPositionCalculatorFactory::get()->create(algoac, acConf, cc);
87  }
88 
89  nSigmaEta_ = pow(conf.getParameter<double>("nSigmaEta"), 2);
90  nSigmaPhi_ = pow(conf.getParameter<double>("nSigmaPhi"), 2);
91 }
92 
94  const std::vector<bool>& seedable,
96  const HcalPFCuts* hcalCuts) {
97  std::vector<double> etaRMS2(input.size(), 0.0);
98  std::vector<double> phiRMS2(input.size(), 0.0);
99 
100  //We need to sort the clusters for smaller to larger depth
101  // for (unsigned int i=0;i<input.size();++i)
102  // printf(" cluster%f %f \n",input[i].depth(),input[i].energy());
103 
104  //calculate cluster shapes
105  calculateShowerShapes(input, etaRMS2, phiRMS2);
106 
107  //link
108  auto&& links = link(input, etaRMS2, phiRMS2);
109  // for (const auto& link: links)
110  // printf("link %d %d %f %f\n",link.from(),link.to(),link.dR(),link.dZ());
111 
112  std::vector<bool> mask(input.size(), false);
113  std::vector<bool> linked(input.size(), false);
114 
115  //prune
116  auto&& prunedLinks = prune(links, linked);
117 
118  //printf("Pruned links\n")
119  // for (const auto& link: prunedLinks)
120  // printf("link %d %d %f %f\n",link.from(),link.to(),link.dR(),link.dZ());
121 
122  //now we need to build clusters
123  for (unsigned int i = 0; i < input.size(); ++i) {
124  //if not masked
125  if (mask[i])
126  continue;
127  //if not linked just spit it out
128  if (!linked[i]) {
129  output.push_back(input[i]);
130  // printf("Added single cluster with energy =%f \n",input[i].energy());
131  mask[i] = true;
132  continue;
133  }
134 
135  //now business: if linked and not masked gather clusters
136  reco::PFCluster cluster = input[i];
137  mask[i] = true;
138  expandCluster(cluster, i, mask, input, prunedLinks);
139  _allCellsPosCalc->calculateAndSetPosition(cluster, hcalCuts);
140  output.push_back(cluster);
141  // printf("Added linked cluster with energy =%f\n",cluster.energy());
142  }
143 }
144 
146  std::vector<double>& etaRMS2,
147  std::vector<double>& phiRMS2) {
148  //shower shapes. here do not use the fractions
149 
150  for (unsigned int i = 0; i < clusters.size(); ++i) {
151  const reco::PFCluster& cluster = clusters[i];
152  double etaSum = 0.0;
153  double phiSum = 0.0;
154  auto const& crep = cluster.positionREP();
155  for (const auto& frac : cluster.recHitFractions()) {
156  auto const& h = *frac.recHitRef();
157  auto const& rep = h.positionREP();
158  etaSum += (frac.fraction() * h.energy()) * std::abs(rep.eta() - crep.eta());
159  phiSum += (frac.fraction() * h.energy()) * std::abs(deltaPhi(rep.phi(), crep.phi()));
160  }
161  //protection for single line : assign ~ tower
162  etaRMS2[i] = std::max(etaSum / cluster.energy(), 0.1);
163  etaRMS2[i] *= etaRMS2[i];
164  phiRMS2[i] = std::max(phiSum / cluster.energy(), 0.1);
165  phiRMS2[i] *= phiRMS2[i];
166  }
167 }
168 
169 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::link(
170  const reco::PFClusterCollection& clusters, const std::vector<double>& etaRMS2, const std::vector<double>& phiRMS2) {
171  std::vector<ClusterLink> links;
172  //loop on all pairs
173  for (unsigned int i = 0; i < clusters.size(); ++i)
174  for (unsigned int j = 0; j < clusters.size(); ++j) {
175  if (i == j)
176  continue;
177 
178  const reco::PFCluster& cluster1 = clusters[i];
179  const reco::PFCluster& cluster2 = clusters[j];
180 
181  // PFCluster depth stored as double but HCAL layer clusters have integral depths only
182  auto dz = (static_cast<int>(cluster2.depth()) - static_cast<int>(cluster1.depth()));
183 
184  //Do not link at the same layer and only link inside out!
185  if (dz <= 0)
186  continue;
187 
188  auto const& crep1 = cluster1.positionREP();
189  auto const& crep2 = cluster2.positionREP();
190 
191  auto deta = crep1.eta() - crep2.eta();
192  deta = deta * deta / (etaRMS2[i] + etaRMS2[j]);
193  auto dphi = deltaPhi(crep1.phi(), crep2.phi());
194  dphi = dphi * dphi / (phiRMS2[i] + phiRMS2[j]);
195 
196  // printf("Testing Link %d -> %d (%f %f %f %f ) \n",i,j,deta,dphi,cluster1.position().Eta()-cluster2.position().Eta(),deltaPhi(cluster1.position().Phi(),cluster2.position().Phi()));
197 
198  if ((deta < nSigmaEta_) & (dphi < nSigmaPhi_))
199  links.push_back(ClusterLink(i, j, deta + dphi, std::abs(dz), cluster1.energy() + cluster2.energy()));
200  }
201 
202  return links;
203 }
204 
205 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::prune(std::vector<ClusterLink>& links,
206  std::vector<bool>& linkedClusters) {
207  std::vector<ClusterLink> goodLinks;
208  std::vector<bool> mask(links.size(), false);
209  if (links.empty())
210  return goodLinks;
211 
212  for (unsigned int i = 0; i < links.size() - 1; ++i) {
213  if (mask[i])
214  continue;
215  for (unsigned int j = i + 1; j < links.size(); ++j) {
216  if (mask[j])
217  continue;
218 
219  const ClusterLink& link1 = links[i];
220  const ClusterLink& link2 = links[j];
221 
222  if (link1.to() == link2.to()) { //found two links going to the same spot,kill one
223  //first prefer nearby layers
224  if (link1.dZ() < link2.dZ()) {
225  mask[j] = true;
226  } else if (link1.dZ() > link2.dZ()) {
227  mask[i] = true;
228  } else { //if same layer-pick based on transverse compatibility
229  if (link1.dR() < link2.dR()) {
230  mask[j] = true;
231  } else if (link1.dR() > link2.dR()) {
232  mask[i] = true;
233  } else {
234  //same distance as well -> can happen !!!!! Pick the highest SUME
235  if (link1.energy() < link2.energy())
236  mask[i] = true;
237  else
238  mask[j] = true;
239  }
240  }
241  }
242  }
243  }
244 
245  for (unsigned int i = 0; i < links.size(); ++i) {
246  if (mask[i])
247  continue;
248  goodLinks.push_back(links[i]);
249  linkedClusters[links[i].from()] = true;
250  linkedClusters[links[i].to()] = true;
251  }
252 
253  return goodLinks;
254 }
255 
257  double e1 = 0.0;
258  double e2 = 0.0;
259 
260  //find seeds
261  for (const auto& fraction : main.recHitFractions())
262  if (fraction.recHitRef()->detId() == main.seed()) {
263  e1 = fraction.recHitRef()->energy();
264  }
265 
266  for (const auto& fraction : added.recHitFractions()) {
267  main.addRecHitFraction(fraction);
268  if (fraction.recHitRef()->detId() == added.seed()) {
269  e2 = fraction.recHitRef()->energy();
270  }
271  }
272  if (e2 > e1)
273  main.setSeed(added.seed());
274 }
275 
277  unsigned int point,
278  std::vector<bool>& mask,
280  const std::vector<ClusterLink>& links) {
281  for (const auto& link : links) {
282  if (link.from() == point) {
283  //found link that starts from this guy if not masked absorb
284  if (!mask[link.from()]) {
285  absorbCluster(cluster, clusters[link.from()]);
286  mask[link.from()] = true;
287  }
288 
289  if (!mask[link.to()]) {
290  absorbCluster(cluster, clusters[link.to()]);
291  mask[link.to()] = true;
292  expandCluster(cluster, link.to(), mask, clusters, links);
293  }
294  }
295  if (link.to() == point) {
296  //found link that starts from this guy if not masked absorb
297  if (!mask[link.to()]) {
298  absorbCluster(cluster, clusters[link.to()]);
299  mask[link.to()] = true;
300  }
301 
302  if (!mask[link.from()]) {
303  absorbCluster(cluster, clusters[link.from()]);
304  mask[link.from()] = true;
305  expandCluster(cluster, link.from(), mask, clusters, links);
306  }
307  }
308  }
309 }
void calculateShowerShapes(const reco::PFClusterCollection &, std::vector< double > &, std::vector< double > &)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:218
bool exists(std::string const &parameterName) const
checks if a parameter exists
ParameterSet const & getParameterSet(std::string const &) const
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
B2DGPF & operator=(const B2DGPF &)=delete
std::vector< ClusterLink > prune(std::vector< ClusterLink > &, std::vector< bool > &linkedClusters)
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:92
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus, const HcalPFCuts *) override
static std::string const input
Definition: EdmProvDump.cc:50
double depth() const
cluster depth
Definition: PFCluster.h:82
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
void absorbCluster(reco::PFCluster &, const reco::PFCluster &)
double energy() const
cluster energy
Definition: PFCluster.h:74
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< ClusterLink > link(const reco::PFClusterCollection &, const std::vector< double > &, const std::vector< double > &)
rep
Definition: cuy.py:1189
void update(const edm::EventSetup &es) override
PFMultiDepthClusterizer(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
~PFMultiDepthClusterizer() override=default
Definition: main.py:1
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: output.py:1
#define get
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
PFMultiDepthClusterizer B2DGPF
void expandCluster(reco::PFCluster &, unsigned int point, std::vector< bool > &mask, const reco::PFClusterCollection &, const std::vector< ClusterLink > &links)