CMS 3D CMS Logo

Basic2DGenericPFlowClusterizer.cc
Go to the documentation of this file.
9 
10 #include "Math/GenVector/VectorUtil.h"
11 #include "vdt/vdtMath.h"
12 
13 #include <iterator>
14 #include <unordered_map>
15 
18 
19 public:
21 
22  ~Basic2DGenericPFlowClusterizer() override = default;
23  Basic2DGenericPFlowClusterizer(const B2DGPF&) = delete;
24  B2DGPF& operator=(const B2DGPF&) = delete;
25 
26  void update(const edm::EventSetup& es) override {
27  _positionCalc->update(es);
28  if (_allCellsPosCalc)
29  _allCellsPosCalc->update(es);
31  _convergencePosCalc->update(es);
32  }
33 
35  const std::vector<bool>&,
37  const HcalPFCuts*) override;
38 
39 private:
40  const unsigned _maxIterations;
41  const double _stoppingTolerance;
42  const double _showerSigma2;
43  const bool _excludeOtherSeeds;
44  const double _minFracTot;
45  const std::unordered_map<std::string, int> _layerMap;
46 
47  std::unordered_map<int, std::pair<std::vector<int>, std::vector<double> > > _recHitEnergyNorms;
48  std::unique_ptr<PFCPositionCalculatorBase> _allCellsPosCalc;
49  std::unique_ptr<PFCPositionCalculatorBase> _convergencePosCalc;
50 
52  const std::vector<bool>&,
54  const HcalPFCuts*) const;
55 
56  void growPFClusters(const reco::PFCluster&,
57  const std::vector<bool>&,
58  const unsigned toleranceScaling,
59  const unsigned iter,
60  double dist,
62  const HcalPFCuts*) const;
63 
65 };
66 
68 
69 #ifdef PFLOW_DEBUG
70 #define LOGVERB(x) edm::LogVerbatim(x)
71 #define LOGWARN(x) edm::LogWarning(x)
72 #define LOGERR(x) edm::LogError(x)
73 #define LOGDRESSED(x) edm::LogInfo(x)
74 #else
75 #define LOGVERB(x) LogTrace(x)
76 #define LOGWARN(x) edm::LogWarning(x)
77 #define LOGERR(x) edm::LogError(x)
78 #define LOGDRESSED(x) LogDebug(x)
79 #endif
80 
83  : PFClusterBuilderBase(conf, cc),
84  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
85  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
86  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"), 2.0)),
87  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
88  _minFracTot(conf.getParameter<double>("minFracTot")),
89  _layerMap({{"PS2", (int)PFLayer::PS2},
90  {"PS1", (int)PFLayer::PS1},
91  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
92  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
93  {"NONE", (int)PFLayer::NONE},
94  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
95  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
96  {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
97  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
98  {"HF_EM", (int)PFLayer::HF_EM},
99  {"HF_HAD", (int)PFLayer::HF_HAD}}) {
100  const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("recHitEnergyNorms");
101  for (const auto& pset : thresholds) {
102  const std::string& det = pset.getParameter<std::string>("detector");
103 
104  std::vector<int> depths;
105  std::vector<double> rhE_norm;
106 
107  if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
108  depths = pset.getParameter<std::vector<int> >("depths");
109  rhE_norm = pset.getParameter<std::vector<double> >("recHitEnergyNorm");
110  } else {
111  depths.push_back(0);
112  rhE_norm.push_back(pset.getParameter<double>("recHitEnergyNorm"));
113  }
114 
115  if (rhE_norm.size() != depths.size()) {
116  throw cms::Exception("InvalidPFRecHitThreshold")
117  << "PFlowClusterizerThreshold mismatch with the numbers of depths";
118  }
119 
120  auto entry = _layerMap.find(det);
121  if (entry == _layerMap.end()) {
122  throw cms::Exception("InvalidDetectorLayer") << "Detector layer : " << det << " is not in the list of recognized"
123  << " detector layers!";
124  }
125  _recHitEnergyNorms.emplace(_layerMap.find(det)->second, std::make_pair(depths, rhE_norm));
126  }
127 
128  const auto& acConf = conf.getParameterSet("allCellsPositionCalc");
129  if (!acConf.empty()) {
130  const std::string& algoac = acConf.getParameter<std::string>("algoName");
131  if (!algoac.empty())
132  _allCellsPosCalc = PFCPositionCalculatorFactory::get()->create(algoac, acConf, cc);
133  }
134  // if necessary a third pos calc for convergence testing
135  const auto& convConf = conf.getParameterSet("positionCalcForConvergence");
136  if (!convConf.empty()) {
137  const std::string& algoconv = convConf.getParameter<std::string>("algoName");
138  if (!algoconv.empty())
139  _convergencePosCalc = PFCPositionCalculatorFactory::get()->create(algoconv, convConf, cc);
140  }
141 }
142 
144  const std::vector<bool>& seedable,
146  const HcalPFCuts* hcalCuts) {
147  reco::PFClusterCollection clustersInTopo;
148  for (const auto& topocluster : input) {
149  clustersInTopo.clear();
150  seedPFClustersFromTopo(topocluster, seedable, clustersInTopo, hcalCuts);
151  const unsigned tolScal = std::pow(std::max(1.0, clustersInTopo.size() - 1.0), 2.0);
152  growPFClusters(topocluster, seedable, tolScal, 0, tolScal, clustersInTopo, hcalCuts);
153  // step added by Josh Bendavid, removes low-fraction clusters
154  // did not impact position resolution with fraction cut of 1e-7
155  // decreases the size of each pf cluster considerably
156  prunePFClusters(clustersInTopo);
157  // recalculate the positions of the pruned clusters
158  if (_convergencePosCalc) {
159  // if defined, use the special position calculation for convergence tests
160  _convergencePosCalc->calculateAndSetPositions(clustersInTopo, hcalCuts);
161  } else {
162  if (clustersInTopo.size() == 1 && _allCellsPosCalc) {
163  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back(), hcalCuts);
164  } else {
165  _positionCalc->calculateAndSetPositions(clustersInTopo, hcalCuts);
166  }
167  }
168  for (auto& clusterout : clustersInTopo) {
169  output.insert(output.end(), std::move(clusterout));
170  }
171  }
172 }
173 
175  const std::vector<bool>& seedable,
176  reco::PFClusterCollection& initialPFClusters,
177  const HcalPFCuts* hcalCuts) const {
178  const auto& recHitFractions = topo.recHitFractions();
179  for (const auto& rhf : recHitFractions) {
180  if (!seedable[rhf.recHitRef().key()])
181  continue;
182  initialPFClusters.push_back(reco::PFCluster());
183  reco::PFCluster& current = initialPFClusters.back();
184  current.addRecHitFraction(rhf);
185  current.setSeed(rhf.recHitRef()->detId());
186  if (_convergencePosCalc) {
187  _convergencePosCalc->calculateAndSetPosition(current, hcalCuts);
188  } else {
189  _positionCalc->calculateAndSetPosition(current, hcalCuts);
190  }
191  }
192 }
193 
195  const std::vector<bool>& seedable,
196  const unsigned toleranceScaling,
197  const unsigned iter,
198  double diff,
200  const HcalPFCuts* hcalCuts) const {
201  if (iter >= _maxIterations) {
202  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
203  << "reached " << _maxIterations << " iterations, terminated position "
204  << "fit with diff = " << diff;
205  }
206  if (iter >= _maxIterations || diff <= _stoppingTolerance * toleranceScaling)
207  return;
208  // reset the rechits in this cluster, keeping the previous position
209  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
210  for (auto& cluster : clusters) {
211  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
212  clus_prev_pos.emplace_back(repp.rho(), repp.eta(), repp.phi());
213  if (_convergencePosCalc) {
214  if (clusters.size() == 1 && _allCellsPosCalc) {
215  _allCellsPosCalc->calculateAndSetPosition(cluster, hcalCuts);
216  } else {
217  _positionCalc->calculateAndSetPosition(cluster, hcalCuts);
218  }
219  }
220  cluster.resetHitsAndFractions();
221  }
222  // loop over topo cluster and grow current PFCluster hypothesis
223  std::vector<double> dist2, frac;
224  double fractot = 0;
225  for (const reco::PFRecHitFraction& rhf : topo.recHitFractions()) {
226  const reco::PFRecHitRef& refhit = rhf.recHitRef();
227  int cell_layer = (int)refhit->layer();
228  if (cell_layer == PFLayer::HCAL_BARREL2 && std::abs(refhit->positionREP().eta()) > 0.34) {
229  cell_layer *= 100;
230  }
231 
232  math::XYZPoint topocellpos_xyz(refhit->position());
233  dist2.clear();
234  frac.clear();
235  fractot = 0;
236 
237  double recHitEnergyNorm = 0.;
238  auto const& recHitEnergyNormDepthPair = _recHitEnergyNorms.find(cell_layer)->second;
239 
240  if (hcalCuts != nullptr && // this means, cutsFromDB is set to True in PFClusterProducer.cc
241  (cell_layer == PFLayer::HCAL_BARREL1 || cell_layer == PFLayer::HCAL_ENDCAP)) {
242  HcalDetId thisId = refhit->detId();
243  const HcalPFCut* item = hcalCuts->getValues(thisId.rawId());
244  recHitEnergyNorm = item->noiseThreshold();
245  } else {
246  for (unsigned int j = 0; j < recHitEnergyNormDepthPair.second.size(); ++j) {
247  int depth = recHitEnergyNormDepthPair.first[j];
248  if ((cell_layer == PFLayer::HCAL_BARREL1 && refhit->depth() == depth) ||
249  (cell_layer == PFLayer::HCAL_ENDCAP && refhit->depth() == depth) ||
250  (cell_layer != PFLayer::HCAL_ENDCAP && cell_layer != PFLayer::HCAL_BARREL1))
251  recHitEnergyNorm = recHitEnergyNormDepthPair.second[j];
252  }
253  }
254 
255  // add rechits to clusters, calculating fraction based on distance
256  for (auto& cluster : clusters) {
257  const math::XYZPoint& clusterpos_xyz = cluster.position();
258  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
259  const double d2 = deltav.Mag2() / _showerSigma2;
260  dist2.emplace_back(d2);
261  if (d2 > 100) {
262  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
263  << "Warning! :: pfcluster-topocell distance is too large! d= " << d2;
264  }
265 
266  // fraction assignment logic
267  double fraction;
268  if (refhit->detId() == cluster.seed() && _excludeOtherSeeds) {
269  fraction = 1.0;
270  } else if (seedable[refhit.key()] && _excludeOtherSeeds) {
271  fraction = 0.0;
272  } else {
273  fraction = cluster.energy() / recHitEnergyNorm * vdt::fast_expf(-0.5 * d2);
274  }
275  fractot += fraction;
276  frac.emplace_back(fraction);
277  }
278  for (unsigned i = 0; i < clusters.size(); ++i) {
279  if (fractot > _minFracTot || (refhit->detId() == clusters[i].seed() && fractot > 0.0)) {
280  frac[i] /= fractot;
281  } else {
282  continue;
283  }
284  // if the fraction has been set to 0, the cell
285  // is now added to the cluster - careful ! (PJ, 19/07/08)
286  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
287  // (PJ, 15/09/08 <- similar to what existed before the
288  // previous bug fix, but keeps the close seeds inside,
289  // even if their fraction was set to zero.)
290  // Also add a protection to keep the seed in the cluster
291  // when the latter gets far from the former. These cases
292  // (about 1% of the clusters) need to be studied, as
293  // they create fake photons, in general.
294  // (PJ, 16/09/08)
295  if (dist2[i] < 100.0 || frac[i] > 0.9999) {
296  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit, frac[i]));
297  }
298  }
299  }
300  // recalculate positions and calculate convergence parameter
301  double diff2 = 0.0;
302  for (unsigned i = 0; i < clusters.size(); ++i) {
303  if (_convergencePosCalc) {
304  _convergencePosCalc->calculateAndSetPosition(clusters[i], hcalCuts);
305  } else {
306  if (clusters.size() == 1 && _allCellsPosCalc) {
307  _allCellsPosCalc->calculateAndSetPosition(clusters[i], hcalCuts);
308  } else {
309  _positionCalc->calculateAndSetPosition(clusters[i], hcalCuts);
310  }
311  }
312  const double delta2 = reco::deltaR2(clusters[i].positionREP(), clus_prev_pos[i]);
313  if (delta2 > diff2)
314  diff2 = delta2;
315  }
316  diff = std::sqrt(diff2);
317  dist2.clear();
318  frac.clear();
319  clus_prev_pos.clear(); // avoid badness
320  growPFClusters(topo, seedable, toleranceScaling, iter + 1, diff, clusters, hcalCuts);
321 }
322 
324  for (auto& cluster : clusters) {
325  cluster.pruneUsing([&](const reco::PFRecHitFraction& rhf) { return rhf.fraction() > _minFractionToKeep; });
326  }
327 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::unique_ptr< PosCalc > _positionCalc
const std::unordered_map< std::string, int > _layerMap
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
#define LOGDRESSED(x)
key_type key() const
Accessor for product key.
Definition: Ref.h:244
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &, const HcalPFCuts *) const
const Item * getValues(DetId fId, bool throwOnFail=true) const
static std::string const input
Definition: EdmProvDump.cc:50
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
Basic2DGenericPFlowClusterizer(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
void prunePFClusters(reco::PFClusterCollection &) const
T sqrt(T t)
Definition: SSEVec.h:23
~Basic2DGenericPFlowClusterizer() override=default
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
B2DGPF & operator=(const B2DGPF &)=delete
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::unordered_map< int, std::pair< std::vector< int >, std::vector< double > > > _recHitEnergyNorms
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
float fast_expf(float x)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: output.py:1
#define get
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &, const HcalPFCuts *) const
void update(const edm::EventSetup &es) override
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus, const HcalPFCuts *) override
def move(src, dest)
Definition: eostools.py:511