CMS 3D CMS Logo

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 
48  std::vector<int> depths;
49  std::vector<double> rhE_norm;
50 
51  if (det==std::string("HCAL_BARREL1") || det==std::string("HCAL_ENDCAP")) {
52  depths= pset.getParameter<std::vector<int> >("depths");
53  rhE_norm = pset.getParameter<std::vector<double> >("recHitEnergyNorm");
54  } else {
55  depths.push_back(0);
56  rhE_norm.push_back(pset.getParameter<double>("recHitEnergyNorm"));
57  }
58 
59  if( rhE_norm.size()!=depths.size() ) {
60  throw cms::Exception("InvalidPFRecHitThreshold")
61  << "PFlowClusterizerThreshold mismatch with the numbers of depths";
62  }
63 
64 
65  auto entry = _layerMap.find(det);
66  if( entry == _layerMap.end() ) {
67  throw cms::Exception("InvalidDetectorLayer")
68  << "Detector layer : " << det << " is not in the list of recognized"
69  << " detector layers!";
70  }
71  _recHitEnergyNorms.emplace(_layerMap.find(det)->second,std::make_pair(depths,rhE_norm));
72  }
73 
74  _allCellsPosCalc.reset(nullptr);
75  if( conf.exists("allCellsPositionCalc") ) {
76  const edm::ParameterSet& acConf =
77  conf.getParameterSet("allCellsPositionCalc");
78  const std::string& algoac =
79  acConf.getParameter<std::string>("algoName");
80  PosCalc* accalc =
81  PFCPositionCalculatorFactory::get()->create(algoac, acConf);
82  _allCellsPosCalc.reset(accalc);
83  }
84  // if necessary a third pos calc for convergence testing
85  _convergencePosCalc.reset(nullptr);
86  if( conf.exists("positionCalcForConvergence") ) {
87  const edm::ParameterSet& convConf =
88  conf.getParameterSet("positionCalcForConvergence");
89  const std::string& algoconv =
90  convConf.getParameter<std::string>("algoName");
91  PosCalc* convcalc =
92  PFCPositionCalculatorFactory::get()->create(algoconv, convConf);
93  _convergencePosCalc.reset(convcalc);
94  }
95 }
96 
99  const std::vector<bool>& seedable,
101  reco::PFClusterCollection clustersInTopo;
102  for( const auto& topocluster : input ) {
103  clustersInTopo.clear();
104  seedPFClustersFromTopo(topocluster,seedable,clustersInTopo);
105  const unsigned tolScal =
106  std::pow(std::max(1.0,clustersInTopo.size()-1.0),2.0);
107  growPFClusters(topocluster,seedable,tolScal,0,tolScal,clustersInTopo);
108  // step added by Josh Bendavid, removes low-fraction clusters
109  // did not impact position resolution with fraction cut of 1e-7
110  // decreases the size of each pf cluster considerably
111  prunePFClusters(clustersInTopo);
112  // recalculate the positions of the pruned clusters
113  if( _convergencePosCalc ) {
114  // if defined, use the special position calculation for convergence tests
115  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
116  } else {
117  if( clustersInTopo.size() == 1 && _allCellsPosCalc ) {
118  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
119  } else {
120  _positionCalc->calculateAndSetPositions(clustersInTopo);
121  }
122  }
123  for( auto& clusterout : clustersInTopo ) {
124  output.insert(output.end(),std::move(clusterout));
125  }
126  }
127 }
128 
131  const std::vector<bool>& seedable,
132  reco::PFClusterCollection& initialPFClusters) const {
133  const auto& recHitFractions = topo.recHitFractions();
134  for( const auto& rhf : recHitFractions ) {
135  if( !seedable[rhf.recHitRef().key()] ) continue;
136  initialPFClusters.push_back(reco::PFCluster());
137  reco::PFCluster& current = initialPFClusters.back();
138  current.addRecHitFraction(rhf);
139  current.setSeed(rhf.recHitRef()->detId());
140  if( _convergencePosCalc ) {
141  _convergencePosCalc->calculateAndSetPosition(current);
142  } else {
143  _positionCalc->calculateAndSetPosition(current);
144  }
145  }
146 }
147 
150  const std::vector<bool>& seedable,
151  const unsigned toleranceScaling,
152  const unsigned iter,
153  double diff,
155  if( iter >= _maxIterations ) {
156  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
157  <<"reached " << _maxIterations << " iterations, terminated position "
158  << "fit with diff = " << diff;
159  }
160  if( iter >= _maxIterations ||
161  diff <= _stoppingTolerance*toleranceScaling) return;
162  // reset the rechits in this cluster, keeping the previous position
163  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
164  for( auto& cluster : clusters) {
165  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
166  clus_prev_pos.emplace_back(repp.rho(),repp.eta(),repp.phi());
167  if( _convergencePosCalc ) {
168  if( clusters.size() == 1 && _allCellsPosCalc ) {
169  _allCellsPosCalc->calculateAndSetPosition(cluster);
170  } else {
171  _positionCalc->calculateAndSetPosition(cluster);
172  }
173  }
174  cluster.resetHitsAndFractions();
175  }
176  // loop over topo cluster and grow current PFCluster hypothesis
177  std::vector<double> dist2, frac;
178  double fractot = 0;
179  for( const reco::PFRecHitFraction& rhf : topo.recHitFractions() ) {
180  const reco::PFRecHitRef& refhit = rhf.recHitRef();
181  int cell_layer = (int)refhit->layer();
182  if( cell_layer == PFLayer::HCAL_BARREL2 &&
183  std::abs(refhit->positionREP().eta()) > 0.34 ) {
184  cell_layer *= 100;
185  }
186 
187  math::XYZPoint topocellpos_xyz(refhit->position());
188  dist2.clear(); frac.clear(); fractot = 0;
189 
190  double recHitEnergyNorm=0.;
191  auto const& recHitEnergyNormDepthPair = _recHitEnergyNorms.find(cell_layer)->second;
192 
193  for (unsigned int j=0; j<recHitEnergyNormDepthPair.second.size(); ++j) {
194  int depth=recHitEnergyNormDepthPair.first[j];
195 
196  if( ( cell_layer == PFLayer::HCAL_BARREL1 && refhit->depth()== depth)
197  || ( cell_layer == PFLayer::HCAL_ENDCAP && refhit->depth()== depth)
198  || ( cell_layer != PFLayer::HCAL_ENDCAP && cell_layer != PFLayer::HCAL_BARREL1)
199  ) recHitEnergyNorm = recHitEnergyNormDepthPair.second[j];
200  }
201 
202  // add rechits to clusters, calculating fraction based on distance
203  for( auto& cluster : clusters ) {
204  const math::XYZPoint& clusterpos_xyz = cluster.position();
205  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
206  const double d2 = deltav.Mag2()/_showerSigma2;
207  dist2.emplace_back( d2 );
208  if( d2 > 100 ) {
209  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
210  << "Warning! :: pfcluster-topocell distance is too large! d= "
211  << d2;
212  }
213 
214  // fraction assignment logic
215  double fraction;
216  if( refhit->detId() == cluster.seed() && _excludeOtherSeeds ) {
217  fraction = 1.0;
218  } else if ( seedable[refhit.key()] && _excludeOtherSeeds ) {
219  fraction = 0.0;
220  } else {
221  fraction = cluster.energy()/recHitEnergyNorm * vdt::fast_expf( -0.5*d2 );
222  }
223  fractot += fraction;
224  frac.emplace_back(fraction);
225  }
226  for( unsigned i = 0; i < clusters.size(); ++i ) {
227  if( fractot > _minFracTot ||
228  ( refhit->detId() == clusters[i].seed() && fractot > 0.0 ) ) {
229  frac[i]/=fractot;
230  } else {
231  continue;
232  }
233  // if the fraction has been set to 0, the cell
234  // is now added to the cluster - careful ! (PJ, 19/07/08)
235  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
236  // (PJ, 15/09/08 <- similar to what existed before the
237  // previous bug fix, but keeps the close seeds inside,
238  // even if their fraction was set to zero.)
239  // Also add a protection to keep the seed in the cluster
240  // when the latter gets far from the former. These cases
241  // (about 1% of the clusters) need to be studied, as
242  // they create fake photons, in general.
243  // (PJ, 16/09/08)
244  if( dist2[i] < 100.0 || frac[i] > 0.9999 ) {
245  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit,frac[i]));
246  }
247  }
248  }
249  // recalculate positions and calculate convergence parameter
250  double diff2 = 0.0;
251  for( unsigned i = 0; i < clusters.size(); ++i ) {
252  if( _convergencePosCalc ) {
253  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
254  } else {
255  if( clusters.size() == 1 && _allCellsPosCalc ) {
256  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
257  } else {
258  _positionCalc->calculateAndSetPosition(clusters[i]);
259  }
260  }
261  const double delta2 =
262  reco::deltaR2(clusters[i].positionREP(),clus_prev_pos[i]);
263  if( delta2 > diff2 ) diff2 = delta2;
264  }
265  diff = std::sqrt(diff2);
266  dist2.clear(); frac.clear(); clus_prev_pos.clear();// avoid badness
267  growPFClusters(topo,seedable,toleranceScaling,iter+1,diff,clusters);
268 }
269 
272  for( auto& cluster : clusters ) {
273  cluster.pruneUsing( [&](const reco::PFRecHitFraction& rhf)
274  {return rhf.fraction() > _minFractionToKeep;}
275  );
276  }
277 }
278 
279 
T getParameter(std::string const &) const
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus) override
std::unordered_map< int, std::pair< std::vector< int >, std::vector< double > > > _recHitEnergyNorms
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
std::unique_ptr< PosCalc > _positionCalc
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
key_type key() const
Accessor for product key.
Definition: Ref.h:265
#define LOGDRESSED(x)
double fraction() const
static std::string const input
Definition: EdmProvDump.cc:44
void setSeed(const DetId &id)
Definition: CaloCluster.h:123
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
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:99
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
def move(src, dest)
Definition: eostools.py:510
T get(const Candidate &c)
Definition: component.h:55
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
Basic2DGenericPFlowClusterizer(const edm::ParameterSet &conf)