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