CMS 3D CMS Logo

PFMultiDepthClusterizer.cc
Go to the documentation of this file.
6 #include "Math/GenVector/VectorUtil.h"
7 
8 #include "vdt/vdtMath.h"
9 
10 #include <iterator>
11 
13  if (conf.exists("allCellsPositionCalc")) {
14  const edm::ParameterSet& acConf = conf.getParameterSet("allCellsPositionCalc");
15  const std::string& algoac = acConf.getParameter<std::string>("algoName");
16  _allCellsPosCalc = PFCPositionCalculatorFactory::get()->create(algoac, acConf);
17  }
18 
19  nSigmaEta_ = pow(conf.getParameter<double>("nSigmaEta"), 2);
20  nSigmaPhi_ = pow(conf.getParameter<double>("nSigmaPhi"), 2);
21 }
22 
24  const std::vector<bool>& seedable,
26  std::vector<double> etaRMS2(input.size(), 0.0);
27  std::vector<double> phiRMS2(input.size(), 0.0);
28 
29  //We need to sort the clusters for smaller to larger depth
30  // for (unsigned int i=0;i<input.size();++i)
31  // printf(" cluster%f %f \n",input[i].depth(),input[i].energy());
32 
33  //calculate cluster shapes
34  calculateShowerShapes(input, etaRMS2, phiRMS2);
35 
36  //link
37  auto&& links = link(input, etaRMS2, phiRMS2);
38  // for (const auto& link: links)
39  // printf("link %d %d %f %f\n",link.from(),link.to(),link.dR(),link.dZ());
40 
41  std::vector<bool> mask(input.size(), false);
42  std::vector<bool> linked(input.size(), false);
43 
44  //prune
45  auto&& prunedLinks = prune(links, linked);
46 
47  //printf("Pruned links\n")
48  // for (const auto& link: prunedLinks)
49  // printf("link %d %d %f %f\n",link.from(),link.to(),link.dR(),link.dZ());
50 
51  //now we need to build clusters
52  for (unsigned int i = 0; i < input.size(); ++i) {
53  //if not masked
54  if (mask[i])
55  continue;
56  //if not linked just spit it out
57  if (!linked[i]) {
58  output.push_back(input[i]);
59  // printf("Added single cluster with energy =%f \n",input[i].energy());
60  mask[i] = true;
61  continue;
62  }
63 
64  //now business: if linked and not masked gather clusters
65  reco::PFCluster cluster = input[i];
66  mask[i] = true;
67  expandCluster(cluster, i, mask, input, prunedLinks);
68  _allCellsPosCalc->calculateAndSetPosition(cluster);
69  output.push_back(cluster);
70  // printf("Added linked cluster with energy =%f\n",cluster.energy());
71  }
72 }
73 
75  std::vector<double>& etaRMS2,
76  std::vector<double>& phiRMS2) {
77  //shower shapes. here do not use the fractions
78 
79  for (unsigned int i = 0; i < clusters.size(); ++i) {
80  const reco::PFCluster& cluster = clusters[i];
81  double etaSum = 0.0;
82  double phiSum = 0.0;
83  auto const& crep = cluster.positionREP();
84  for (const auto& frac : cluster.recHitFractions()) {
85  auto const& h = *frac.recHitRef();
86  auto const& rep = h.positionREP();
87  etaSum += (frac.fraction() * h.energy()) * std::abs(rep.eta() - crep.eta());
88  phiSum += (frac.fraction() * h.energy()) * std::abs(deltaPhi(rep.phi(), crep.phi()));
89  }
90  //protection for single line : assign ~ tower
91  etaRMS2[i] = std::max(etaSum / cluster.energy(), 0.1);
92  etaRMS2[i] *= etaRMS2[i];
93  phiRMS2[i] = std::max(phiSum / cluster.energy(), 0.1);
94  phiRMS2[i] *= phiRMS2[i];
95  }
96 }
97 
98 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::link(
99  const reco::PFClusterCollection& clusters, const std::vector<double>& etaRMS2, const std::vector<double>& phiRMS2) {
100  std::vector<ClusterLink> links;
101  //loop on all pairs
102  for (unsigned int i = 0; i < clusters.size(); ++i)
103  for (unsigned int j = 0; j < clusters.size(); ++j) {
104  if (i == j)
105  continue;
106 
107  const reco::PFCluster& cluster1 = clusters[i];
108  const reco::PFCluster& cluster2 = clusters[j];
109 
110  auto dz = (cluster2.depth() - cluster1.depth());
111 
112  //Do not link at the same layer and only link inside out!
113  if (dz < 0.0f || std::abs(dz) < 0.2f)
114  continue;
115 
116  auto const& crep1 = cluster1.positionREP();
117  auto const& crep2 = cluster2.positionREP();
118 
119  auto deta = crep1.eta() - crep2.eta();
120  deta = deta * deta / (etaRMS2[i] + etaRMS2[j]);
121  auto dphi = deltaPhi(crep1.phi(), crep2.phi());
122  dphi = dphi * dphi / (phiRMS2[i] + phiRMS2[j]);
123 
124  // 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()));
125 
126  if ((deta < nSigmaEta_) & (dphi < nSigmaPhi_))
127  links.push_back(ClusterLink(i, j, deta + dphi, std::abs(dz), cluster1.energy() + cluster2.energy()));
128  }
129 
130  return links;
131 }
132 
133 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::prune(std::vector<ClusterLink>& links,
134  std::vector<bool>& linkedClusters) {
135  std::vector<ClusterLink> goodLinks;
136  std::vector<bool> mask(links.size(), false);
137  if (links.empty())
138  return goodLinks;
139 
140  for (unsigned int i = 0; i < links.size() - 1; ++i) {
141  if (mask[i])
142  continue;
143  for (unsigned int j = i + 1; j < links.size(); ++j) {
144  if (mask[j])
145  continue;
146 
147  const ClusterLink& link1 = links[i];
148  const ClusterLink& link2 = links[j];
149 
150  if (link1.to() == link2.to()) { //found two links going to the same spot,kill one
151  //first prefer nearby layers
152  if (link1.dZ() < link2.dZ()) {
153  mask[j] = true;
154  } else if (link1.dZ() > link2.dZ()) {
155  mask[i] = true;
156  } else if (fabs(link1.dZ() - link2.dZ()) < 0.2) { //if same layer-pick based on transverse compatibility
157  if (link1.dR() < link2.dR()) {
158  mask[j] = true;
159  } else if (link1.dR() > link2.dR()) {
160  mask[i] = true;
161  } else {
162  //same distance as well -> can happen !!!!! Pick the highest SUME
163  if (link1.energy() < link2.energy())
164  mask[i] = true;
165  else
166  mask[j] = true;
167  }
168  }
169  }
170  }
171  }
172 
173  for (unsigned int i = 0; i < links.size(); ++i) {
174  if (mask[i])
175  continue;
176  goodLinks.push_back(links[i]);
177  linkedClusters[links[i].from()] = true;
178  linkedClusters[links[i].to()] = true;
179  }
180 
181  return goodLinks;
182 }
183 
185  double e1 = 0.0;
186  double e2 = 0.0;
187 
188  //find seeds
189  for (const auto& fraction : main.recHitFractions())
190  if (fraction.recHitRef()->detId() == main.seed()) {
191  e1 = fraction.recHitRef()->energy();
192  }
193 
194  for (const auto& fraction : added.recHitFractions()) {
195  main.addRecHitFraction(fraction);
196  if (fraction.recHitRef()->detId() == added.seed()) {
197  e2 = fraction.recHitRef()->energy();
198  }
199  }
200  if (e2 > e1)
201  main.setSeed(added.seed());
202 }
203 
205  unsigned int point,
206  std::vector<bool>& mask,
208  const std::vector<ClusterLink>& links) {
209  for (const auto& link : links) {
210  if (link.from() == point) {
211  //found link that starts from this guy if not masked absorb
212  if (!mask[link.from()]) {
213  absorbCluster(cluster, clusters[link.from()]);
214  mask[link.from()] = true;
215  }
216 
217  if (!mask[link.to()]) {
218  absorbCluster(cluster, clusters[link.to()]);
219  mask[link.to()] = true;
220  expandCluster(cluster, link.to(), mask, clusters, links);
221  }
222  }
223  if (link.to() == point) {
224  //found link that starts from this guy if not masked absorb
225  if (!mask[link.to()]) {
226  absorbCluster(cluster, clusters[link.to()]);
227  mask[link.to()] = true;
228  }
229 
230  if (!mask[link.from()]) {
231  absorbCluster(cluster, clusters[link.from()]);
232  mask[link.from()] = true;
233  expandCluster(cluster, link.from(), mask, clusters, links);
234  }
235  }
236  }
237 }
reco::PFClusterCollection
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
mps_fire.i
i
Definition: mps_fire.py:428
input
static const std::string input
Definition: EdmProvDump.cc:48
MessageLogger.h
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
deltaPhi.h
PFMultiDepthClusterizer::calculateShowerShapes
void calculateShowerShapes(const reco::PFClusterCollection &, std::vector< double > &, std::vector< double > &)
Definition: PFMultiDepthClusterizer.cc:74
PFMultiDepthClusterizer::expandCluster
void expandCluster(reco::PFCluster &, unsigned int point, std::vector< bool > &mask, const reco::PFClusterCollection &, const std::vector< ClusterLink > &links)
Definition: PFMultiDepthClusterizer.cc:204
PFMultiDepthClusterizer::nSigmaEta_
double nSigmaEta_
Definition: PFMultiDepthClusterizer.h:27
reco::PFCluster::recHitFractions
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
PFRecHit.h
PFMultiDepthClusterizer::PFMultiDepthClusterizer
PFMultiDepthClusterizer(const edm::ParameterSet &conf)
Definition: PFMultiDepthClusterizer.cc:12
PFMultiDepthClusterizer::_allCellsPosCalc
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
Definition: PFMultiDepthClusterizer.h:26
reco::PFCluster::energy
double energy() const
cluster energy
Definition: PFCluster.h:74
PFMultiDepthClusterizer::prune
std::vector< ClusterLink > prune(std::vector< ClusterLink > &, std::vector< bool > &linkedClusters)
Definition: PFMultiDepthClusterizer.cc:133
DivergingColor.frac
float frac
Definition: DivergingColor.py:175
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
HLT_FULL_cff.fraction
fraction
Definition: HLT_FULL_cff.py:52795
reco::PFCluster::depth
double depth() const
cluster depth
Definition: PFCluster.h:82
h
reco::PFCluster::positionREP
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:92
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
PFMultiDepthClusterizer::absorbCluster
void absorbCluster(reco::PFCluster &, const reco::PFCluster &)
Definition: PFMultiDepthClusterizer.cc:184
edm::ParameterSet
Definition: ParameterSet.h:47
deltaR.h
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
PFMultiDepthClusterizer::nSigmaPhi_
double nSigmaPhi_
Definition: PFMultiDepthClusterizer.h:28
PFMultiDepthClusterizer::link
std::vector< ClusterLink > link(const reco::PFClusterCollection &, const std::vector< double > &, const std::vector< double > &)
Definition: PFMultiDepthClusterizer.cc:98
StorageManager_cfg.e1
e1
Definition: StorageManager_cfg.py:16
PFMultiDepthClusterizer.h
cuy.rep
rep
Definition: cuy.py:1190
get
#define get
reco::CaloCluster::seed
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
main
Definition: main.py:1
PVValHelper::dz
Definition: PVValidationHelpers.h:50
electronStore.links
links
Definition: electronStore.py:149
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PFMultiDepthClusterizer::buildClusters
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus) override
Definition: PFMultiDepthClusterizer.cc:23
PFClusterBuilderBase
Definition: PFClusterBuilderBase.h:18
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
point
*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
edm::ParameterSet::getParameterSet
ParameterSet const & getParameterSet(std::string const &) const
Definition: ParameterSet.cc:2128