CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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(NULL);
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> etaRMS(input.size(),0.0);
41  std::vector<double> phiRMS(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,etaRMS,phiRMS);
51 
52  //link
53  std::vector<ClusterLink> links = link(input,etaRMS,phiRMS);
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  std::vector<ClusterLink> prunedLinks;
64  if (links.size())
65  prunedLinks = prune(links,linked);
66 
67  //printf("Pruned links\n")
68  // for (const auto& link: prunedLinks)
69  // printf("link %d %d %f %f\n",link.from(),link.to(),link.dR(),link.dZ());
70 
71 
72 
73  //now we need to build clusters
74  for (unsigned int i=0;i<input.size();++i) {
75  //if not masked
76  if( mask[i])
77  continue;
78  //if not linked just spit it out
79  if (!linked[i]) {
80  output.push_back(input[i]);
81  // printf("Added single cluster with energy =%f \n",input[i].energy());
82  mask[i] = true;
83  continue;
84  }
85 
86 
87  //now business: if linked and not masked gather clusters
88  reco::PFCluster cluster = input[i];
89  mask[i] = true;
90  expandCluster(cluster,i, mask,input,prunedLinks);
91  _allCellsPosCalc->calculateAndSetPosition(cluster);
92  output.push_back(cluster);
93  // printf("Added linked cluster with energy =%f\n",cluster.energy());
94 
95  }
96 }
97 
98 
100 calculateShowerShapes(const reco::PFClusterCollection& clusters,std::vector<double>& etaRMS,std::vector<double>& phiRMS) {
101  float etaSum;
102  float phiSum;
103 
104  //shower shapes. here do not use the fractions
105 
106  for( unsigned int i=0;i<clusters.size();++i ) {
107  const reco::PFCluster& cluster = clusters[i];
108  etaSum=0.0;phiSum=0.0;
109  auto const & crep = cluster.positionREP();
110  for (const auto& frac : cluster.recHitFractions()) {
111  auto const & h = *frac.recHitRef();
112  auto const & rep = h.positionREP();
113  etaSum += (frac.fraction()*h.energy())*std::abs(rep.eta()-crep.eta());
114  phiSum += (frac.fraction()*h.energy())*std::abs(deltaPhi(rep.phi(),crep.phi()));
115  }
116  //protection for single line : assign ~ tower
117  etaRMS[i] = std::max(etaSum/cluster.energy(),0.1);
118  phiRMS[i] = std::max(phiSum/cluster.energy(),0.1);
119 
120  }
121 
122 }
123 
124 
125 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::
126 link(const reco::PFClusterCollection& clusters ,const std::vector<double>& etaRMS,const std::vector<double>& phiRMS) {
127 
128  std::vector<ClusterLink> links;
129  //loop on all pairs
130  for (unsigned int i=0;i<clusters.size();++i)
131  for (unsigned int j=0;j<clusters.size();++j) {
132  if (i==j )
133  continue;
134 
135  const reco::PFCluster& cluster1 = clusters[i];
136  const reco::PFCluster& cluster2 = clusters[j];
137 
138  float dz = (cluster2.depth() - cluster1.depth());
139 
140  //Do not link at the same layer and only link inside out!
141  if (dz<0.0 || fabs(dz)<0.2)
142  continue;
143  auto const & crep1 = cluster1.positionREP();
144  auto const & crep2 = cluster2.positionREP();
145 
146  float deta =(crep1.eta()-crep2.eta())*(crep1.eta()-crep2.eta())/(etaRMS[i]*etaRMS[i]+etaRMS[j]*etaRMS[j]);
147  float dphi = deltaPhi(crep1.phi(),crep2.phi())*deltaPhi(crep1.phi(),crep2.phi())/(phiRMS[i]*phiRMS[i]+phiRMS[j]*phiRMS[j]);
148 
149  // 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()));
150 
151  if (deta<nSigmaEta_ && dphi<nSigmaPhi_ )
152  links.push_back(ClusterLink(i,j,deta+dphi,fabs(dz),cluster1.energy()+cluster2.energy()));
153  }
154 
155  return links;
156 }
157 
158 std::vector<PFMultiDepthClusterizer::ClusterLink> PFMultiDepthClusterizer::
159 prune(std::vector<ClusterLink>& links,std::vector<bool>& linkedClusters) {
160  std::vector<ClusterLink> goodLinks ;
161  std::vector<bool> mask(links.size(),false);
162 
163  for (unsigned int i=0;i<links.size()-1;++i) {
164  if (mask[i])
165  continue;
166  for (unsigned int j=i+1;j<links.size();++j) {
167 
168  if (mask[j])
169  continue;
170 
171  const ClusterLink& link1 = links[i];
172  const ClusterLink& link2 = links[j];
173 
174  if (link1.to() == link2.to()) {//found two links going to the same spot,kill one
175  //first prefer nearby layers
176  if (link1.dZ() < link2.dZ()) {
177  mask[j]=true;
178  }
179  else if (link1.dZ() > link2.dZ()) {
180  mask[i] = true;
181  }
182  else if (fabs(link1.dZ()-link2.dZ())<0.2) { //if same layer-pick based on transverse compatibility
183  if (link1.dR()<link2.dR()) {
184  mask[j]=true;
185  }
186  else if (link1.dR()>link2.dR()) {
187  mask[i] = true;
188  }
189  else {
190  //same distance as well -> can happen !!!!! Pick the highest SUME
191  if (link1.energy()<link2.energy())
192  mask[i] = true;
193  else
194  mask[j] = true;
195 
196  }
197  }
198  }
199  }
200  }
201 
202  for (unsigned int i=0;i<links.size();++i) {
203  if (mask[i])
204  continue;
205  goodLinks.push_back(links[i]);
206  linkedClusters[links[i].from()]=true;
207  linkedClusters[links[i].to()]=true;
208  }
209 
210 
211  return goodLinks;
212 }
213 
214 
215 void
218 
219  double e1=0.0;
220  double e2=0.0;
221 
222  //find seeds
223  for ( const auto& fraction : main.recHitFractions())
224  if (fraction.recHitRef()->detId() == main.seed()) {
225  e1=fraction.recHitRef()->energy();
226  }
227 
228  for ( const auto& fraction : added.recHitFractions()) {
230  if (fraction.recHitRef()->detId() == added.seed()) {
231  e2=fraction.recHitRef()->energy();
232  }
233  }
234  if (e2>e1)
235  main.setSeed(added.seed());
236 
237 }
238 
239 
240 void
242 expandCluster(reco::PFCluster& cluster,unsigned int point,std::vector<bool>& mask,const reco::PFClusterCollection& clusters, const std::vector<ClusterLink>& links) {
243  for (const auto& link : links) {
244  if (link.from() == point) {
245  //found link that starts from this guy if not masked absorb
246  if (!mask[link.from()]) {
247  absorbCluster(cluster,clusters[link.from()]);
248  mask[link.from()] = true;
249  }
250 
251  if (!mask[link.to()]) {
252  absorbCluster(cluster,clusters[link.to()]);
253  mask[link.to()]=true;
254  expandCluster(cluster,link.to(),mask,clusters,links);
255  }
256  }
257  if (link.to() == point) {
258  //found link that starts from this guy if not masked absorb
259  if (!mask[link.to()]) {
260  absorbCluster(cluster,clusters[link.to()]);
261  mask[link.to()] = true;
262  }
263 
264  if (!mask[link.from()]) {
265  absorbCluster(cluster,clusters[link.from()]);
266  mask[link.from()]=true;
267  expandCluster(cluster,link.from(),mask,clusters,links);
268  }
269  }
270 
271 
272 
273  }
274 
275 }
T getParameter(std::string const &) const
void calculateShowerShapes(const reco::PFClusterCollection &, std::vector< double > &, std::vector< double > &)
int i
Definition: DBlmapReader.cc:9
string rep
Definition: cuy.py:1188
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
#define NULL
Definition: scimark2.h:8
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus)
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:118
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:93
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
std::vector< ClusterLink > link(const reco::PFClusterCollection &, const std::vector< double > &, const std::vector< double > &)
double energy() const
cluster energy
Definition: PFCluster.h:82
Float e1
Definition: deltaR.h:20
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
ParameterSet const & getParameterSet(std::string const &) const
Float e2
Definition: deltaR.h:21
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:57
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
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:87