CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Basic2DGenericPFlowClusterizer.cc
Go to the documentation of this file.
5 
6 #include "Math/GenVector/VectorUtil.h"
7 
8 #include "vdt/vdtMath.h"
9 
10 #include <iterator>
11 
12 #ifdef PFLOW_DEBUG
13 #define LOGVERB(x) edm::LogVerbatim(x)
14 #define LOGWARN(x) edm::LogWarning(x)
15 #define LOGERR(x) edm::LogError(x)
16 #define LOGDRESSED(x) edm::LogInfo(x)
17 #else
18 #define LOGVERB(x) LogTrace(x)
19 #define LOGWARN(x) edm::LogWarning(x)
20 #define LOGERR(x) edm::LogError(x)
21 #define LOGDRESSED(x) LogDebug(x)
22 #endif
23 
27  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
28  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
29  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"),2.0)),
30  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
31  _minFracTot(conf.getParameter<double>("minFracTot")),
32  _layerMap({ {"PS2",(int)PFLayer::PS2},
33  {"PS1",(int)PFLayer::PS1},
34  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
35  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
36  {"NONE",(int)PFLayer::NONE},
37  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
38  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
39  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
40  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
41  {"HF_EM",(int)PFLayer::HF_EM},
42  {"HF_HAD",(int)PFLayer::HF_HAD} }) {
43  const std::vector<edm::ParameterSet>& thresholds =
44  conf.getParameterSetVector("recHitEnergyNorms");
45  for( const auto& pset : thresholds ) {
46  const std::string& det = pset.getParameter<std::string>("detector");
47  const double& rhE_norm = pset.getParameter<double>("recHitEnergyNorm");
48  auto entry = _layerMap.find(det);
49  if( entry == _layerMap.end() ) {
50  throw cms::Exception("InvalidDetectorLayer")
51  << "Detector layer : " << det << " is not in the list of recognized"
52  << " detector layers!";
53  }
54  _recHitEnergyNorms.emplace(_layerMap.find(det)->second,rhE_norm);
55  }
56 
57  _allCellsPosCalc.reset(NULL);
58  if( conf.exists("allCellsPositionCalc") ) {
59  const edm::ParameterSet& acConf =
60  conf.getParameterSet("allCellsPositionCalc");
61  const std::string& algoac =
62  acConf.getParameter<std::string>("algoName");
63  PosCalc* accalc =
64  PFCPositionCalculatorFactory::get()->create(algoac, acConf);
65  _allCellsPosCalc.reset(accalc);
66  }
67  // if necessary a third pos calc for convergence testing
69  if( conf.exists("positionCalcForConvergence") ) {
70  const edm::ParameterSet& convConf =
71  conf.getParameterSet("positionCalcForConvergence");
72  const std::string& algoconv =
73  convConf.getParameter<std::string>("algoName");
74  PosCalc* convcalc =
75  PFCPositionCalculatorFactory::get()->create(algoconv, convConf);
76  _convergencePosCalc.reset(convcalc);
77  }
78 }
79 
82  const std::vector<bool>& seedable,
84  reco::PFClusterCollection clustersInTopo;
85  for( const auto& topocluster : input ) {
86  clustersInTopo.clear();
87  seedPFClustersFromTopo(topocluster,seedable,clustersInTopo);
88  const unsigned tolScal =
89  std::pow(std::max(1.0,clustersInTopo.size()-1.0),2.0);
90  growPFClusters(topocluster,seedable,tolScal,0,tolScal,clustersInTopo);
91  // step added by Josh Bendavid, removes low-fraction clusters
92  // did not impact position resolution with fraction cut of 1e-7
93  // decreases the size of each pf cluster considerably
94  prunePFClusters(clustersInTopo);
95  // recalculate the positions of the pruned clusters
96  if( _convergencePosCalc ) {
97  // if defined, use the special position calculation for convergence tests
98  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
99  } else {
100  if( clustersInTopo.size() == 1 && _allCellsPosCalc ) {
101  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
102  } else {
103  _positionCalc->calculateAndSetPositions(clustersInTopo);
104  }
105  }
106  for( auto& clusterout : clustersInTopo ) {
107  output.insert(output.end(),std::move(clusterout));
108  }
109  }
110 }
111 
114  const std::vector<bool>& seedable,
115  reco::PFClusterCollection& initialPFClusters) const {
116  const auto& recHitFractions = topo.recHitFractions();
117  for( const auto& rhf : recHitFractions ) {
118  if( !seedable[rhf.recHitRef().key()] ) continue;
119  initialPFClusters.push_back(reco::PFCluster());
120  reco::PFCluster& current = initialPFClusters.back();
121  current.addRecHitFraction(rhf);
122  current.setSeed(rhf.recHitRef()->detId());
123  if( _convergencePosCalc ) {
124  _convergencePosCalc->calculateAndSetPosition(current);
125  } else {
126  _positionCalc->calculateAndSetPosition(current);
127  }
128  }
129 }
130 
133  const std::vector<bool>& seedable,
134  const unsigned toleranceScaling,
135  const unsigned iter,
136  double diff,
138  if( iter >= _maxIterations ) {
139  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
140  <<"reached " << _maxIterations << " iterations, terminated position "
141  << "fit with diff = " << diff;
142  }
143  if( iter >= _maxIterations ||
144  diff <= _stoppingTolerance*toleranceScaling) return;
145  // reset the rechits in this cluster, keeping the previous position
146  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
147  for( auto& cluster : clusters) {
148  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
149  clus_prev_pos.emplace_back(repp.rho(),repp.eta(),repp.phi());
150  if( _convergencePosCalc ) {
151  if( clusters.size() == 1 && _allCellsPosCalc ) {
152  _allCellsPosCalc->calculateAndSetPosition(cluster);
153  } else {
154  _positionCalc->calculateAndSetPosition(cluster);
155  }
156  }
157  cluster.resetHitsAndFractions();
158  }
159  // loop over topo cluster and grow current PFCluster hypothesis
160  std::vector<double> dist2, frac;
161  double fractot = 0, fraction = 0;
162  for( const reco::PFRecHitFraction& rhf : topo.recHitFractions() ) {
163  const reco::PFRecHitRef& refhit = rhf.recHitRef();
164  int cell_layer = (int)refhit->layer();
165  if( cell_layer == PFLayer::HCAL_BARREL2 &&
166  std::abs(refhit->positionREP().eta()) > 0.34 ) {
167  cell_layer *= 100;
168  }
169  const double recHitEnergyNorm =
170  _recHitEnergyNorms.find(cell_layer)->second;
171  const math::XYZPoint& topocellpos_xyz = refhit->position();
172  dist2.clear(); frac.clear(); fractot = 0;
173  // add rechits to clusters, calculating fraction based on distance
174  for( auto& cluster : clusters ) {
175  const math::XYZPoint& clusterpos_xyz = cluster.position();
176  fraction = 0.0;
177  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
178  const double d2 = deltav.Mag2()/_showerSigma2;
179  dist2.emplace_back( d2 );
180  if( d2 > 100 ) {
181  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
182  << "Warning! :: pfcluster-topocell distance is too large! d= "
183  << d2;
184  }
185  // fraction assignment logic
186  if( refhit->detId() == cluster.seed() && _excludeOtherSeeds ) {
187  fraction = 1.0;
188  } else if ( seedable[refhit.key()] && _excludeOtherSeeds ) {
189  fraction = 0.0;
190  } else {
191  fraction = cluster.energy()/recHitEnergyNorm * vdt::fast_expf( -0.5*d2 );
192  }
193  fractot += fraction;
194  frac.emplace_back(fraction);
195  }
196  for( unsigned i = 0; i < clusters.size(); ++i ) {
197  if( fractot > _minFracTot ||
198  ( refhit->detId() == clusters[i].seed() && fractot > 0.0 ) ) {
199  frac[i]/=fractot;
200  } else {
201  continue;
202  }
203  // if the fraction has been set to 0, the cell
204  // is now added to the cluster - careful ! (PJ, 19/07/08)
205  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
206  // (PJ, 15/09/08 <- similar to what existed before the
207  // previous bug fix, but keeps the close seeds inside,
208  // even if their fraction was set to zero.)
209  // Also add a protection to keep the seed in the cluster
210  // when the latter gets far from the former. These cases
211  // (about 1% of the clusters) need to be studied, as
212  // they create fake photons, in general.
213  // (PJ, 16/09/08)
214  if( dist2[i] < 100.0 || frac[i] > 0.9999 ) {
215  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit,frac[i]));
216  }
217  }
218  }
219  // recalculate positions and calculate convergence parameter
220  double diff2 = 0.0;
221  for( unsigned i = 0; i < clusters.size(); ++i ) {
222  if( _convergencePosCalc ) {
223  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
224  } else {
225  if( clusters.size() == 1 && _allCellsPosCalc ) {
226  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
227  } else {
228  _positionCalc->calculateAndSetPosition(clusters[i]);
229  }
230  }
231  const double delta2 =
232  reco::deltaR2(clusters[i].positionREP(),clus_prev_pos[i]);
233  if( delta2 > diff2 ) diff2 = delta2;
234  }
235  diff = std::sqrt(diff2);
236  dist2.clear(); frac.clear(); clus_prev_pos.clear();// avoid badness
237  growPFClusters(topo,seedable,toleranceScaling,iter+1,diff,clusters);
238 }
239 
242  for( auto& cluster : clusters ) {
243  cluster.pruneUsing( [&](const reco::PFRecHitFraction& rhf)
244  {return rhf.fraction() > _minFractionToKeep;}
245  );
246  }
247 }
248 
249 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
VParameterSet const & getParameterSetVector(std::string const &name) const
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
std::unique_ptr< PosCalc > _positionCalc
bool exists(std::string const &parameterName) const
checks if a parameter exists
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
#define NULL
Definition: scimark2.h:8
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
key_type key() const
Accessor for product key.
Definition: Ref.h:264
#define LOGDRESSED(x)
double fraction() const
static std::string const input
Definition: EdmProvDump.cc:43
void setSeed(const DetId &id)
Definition: CaloCluster.h:118
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
std::unordered_map< int, double > _recHitEnergyNorms
tuple conf
Definition: dbtoconf.py:185
const std::unordered_map< std::string, int > _layerMap
ParameterSet const & getParameterSet(std::string const &) const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:57
float fast_expf(float x)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:54
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 prunePFClusters(reco::PFClusterCollection &) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
T get(const Candidate &c)
Definition: component.h:55
PFCPositionCalculatorBase PosCalc
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus)
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
Basic2DGenericPFlowClusterizer(const edm::ParameterSet &conf)