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 
25  : PFClusterBuilderBase(conf),
26  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
27  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
28  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"), 2.0)),
29  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
30  _minFracTot(conf.getParameter<double>("minFracTot")),
31  _layerMap({{"PS2", (int)PFLayer::PS2},
32  {"PS1", (int)PFLayer::PS1},
33  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
34  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
35  {"NONE", (int)PFLayer::NONE},
36  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
37  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
38  {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
39  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
40  {"HF_EM", (int)PFLayer::HF_EM},
41  {"HF_HAD", (int)PFLayer::HF_HAD}}) {
42  const std::vector<edm::ParameterSet>& thresholds = conf.getParameterSetVector("recHitEnergyNorms");
43  for (const auto& pset : thresholds) {
44  const std::string& det = pset.getParameter<std::string>("detector");
45 
46  std::vector<int> depths;
47  std::vector<double> rhE_norm;
48 
49  if (det == std::string("HCAL_BARREL1") || det == std::string("HCAL_ENDCAP")) {
50  depths = pset.getParameter<std::vector<int> >("depths");
51  rhE_norm = pset.getParameter<std::vector<double> >("recHitEnergyNorm");
52  } else {
53  depths.push_back(0);
54  rhE_norm.push_back(pset.getParameter<double>("recHitEnergyNorm"));
55  }
56 
57  if (rhE_norm.size() != depths.size()) {
58  throw cms::Exception("InvalidPFRecHitThreshold")
59  << "PFlowClusterizerThreshold mismatch with the numbers of depths";
60  }
61 
62  auto entry = _layerMap.find(det);
63  if (entry == _layerMap.end()) {
64  throw cms::Exception("InvalidDetectorLayer") << "Detector layer : " << det << " is not in the list of recognized"
65  << " detector layers!";
66  }
67  _recHitEnergyNorms.emplace(_layerMap.find(det)->second, std::make_pair(depths, rhE_norm));
68  }
69 
70  if (conf.exists("allCellsPositionCalc")) {
71  const edm::ParameterSet& acConf = conf.getParameterSet("allCellsPositionCalc");
72  const std::string& algoac = acConf.getParameter<std::string>("algoName");
73  _allCellsPosCalc = PFCPositionCalculatorFactory::get()->create(algoac, acConf);
74  }
75  // if necessary a third pos calc for convergence testing
76  if (conf.exists("positionCalcForConvergence")) {
77  const edm::ParameterSet& convConf = conf.getParameterSet("positionCalcForConvergence");
78  const std::string& algoconv = convConf.getParameter<std::string>("algoName");
79  _convergencePosCalc = PFCPositionCalculatorFactory::get()->create(algoconv, convConf);
80  }
81 }
82 
84  const std::vector<bool>& seedable,
86  reco::PFClusterCollection clustersInTopo;
87  for (const auto& topocluster : input) {
88  clustersInTopo.clear();
89  seedPFClustersFromTopo(topocluster, seedable, clustersInTopo);
90  const unsigned tolScal = std::pow(std::max(1.0, clustersInTopo.size() - 1.0), 2.0);
91  growPFClusters(topocluster, seedable, tolScal, 0, tolScal, clustersInTopo);
92  // step added by Josh Bendavid, removes low-fraction clusters
93  // did not impact position resolution with fraction cut of 1e-7
94  // decreases the size of each pf cluster considerably
95  prunePFClusters(clustersInTopo);
96  // recalculate the positions of the pruned clusters
97  if (_convergencePosCalc) {
98  // if defined, use the special position calculation for convergence tests
99  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
100  } else {
101  if (clustersInTopo.size() == 1 && _allCellsPosCalc) {
102  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
103  } else {
104  _positionCalc->calculateAndSetPositions(clustersInTopo);
105  }
106  }
107  for (auto& clusterout : clustersInTopo) {
108  output.insert(output.end(), std::move(clusterout));
109  }
110  }
111 }
112 
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()])
119  continue;
120  initialPFClusters.push_back(reco::PFCluster());
121  reco::PFCluster& current = initialPFClusters.back();
122  current.addRecHitFraction(rhf);
123  current.setSeed(rhf.recHitRef()->detId());
124  if (_convergencePosCalc) {
125  _convergencePosCalc->calculateAndSetPosition(current);
126  } else {
127  _positionCalc->calculateAndSetPosition(current);
128  }
129  }
130 }
131 
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 || diff <= _stoppingTolerance * toleranceScaling)
144  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;
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 && std::abs(refhit->positionREP().eta()) > 0.34) {
166  cell_layer *= 100;
167  }
168 
169  math::XYZPoint topocellpos_xyz(refhit->position());
170  dist2.clear();
171  frac.clear();
172  fractot = 0;
173 
174  double recHitEnergyNorm = 0.;
175  auto const& recHitEnergyNormDepthPair = _recHitEnergyNorms.find(cell_layer)->second;
176 
177  for (unsigned int j = 0; j < recHitEnergyNormDepthPair.second.size(); ++j) {
178  int depth = recHitEnergyNormDepthPair.first[j];
179 
180  if ((cell_layer == PFLayer::HCAL_BARREL1 && refhit->depth() == depth) ||
181  (cell_layer == PFLayer::HCAL_ENDCAP && refhit->depth() == depth) ||
182  (cell_layer != PFLayer::HCAL_ENDCAP && cell_layer != PFLayer::HCAL_BARREL1))
183  recHitEnergyNorm = recHitEnergyNormDepthPair.second[j];
184  }
185 
186  // add rechits to clusters, calculating fraction based on distance
187  for (auto& cluster : clusters) {
188  const math::XYZPoint& clusterpos_xyz = cluster.position();
189  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
190  const double d2 = deltav.Mag2() / _showerSigma2;
191  dist2.emplace_back(d2);
192  if (d2 > 100) {
193  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
194  << "Warning! :: pfcluster-topocell distance is too large! d= " << d2;
195  }
196 
197  // fraction assignment logic
198  double fraction;
199  if (refhit->detId() == cluster.seed() && _excludeOtherSeeds) {
200  fraction = 1.0;
201  } else if (seedable[refhit.key()] && _excludeOtherSeeds) {
202  fraction = 0.0;
203  } else {
204  fraction = cluster.energy() / recHitEnergyNorm * vdt::fast_expf(-0.5 * d2);
205  }
206  fractot += fraction;
207  frac.emplace_back(fraction);
208  }
209  for (unsigned i = 0; i < clusters.size(); ++i) {
210  if (fractot > _minFracTot || (refhit->detId() == clusters[i].seed() && fractot > 0.0)) {
211  frac[i] /= fractot;
212  } else {
213  continue;
214  }
215  // if the fraction has been set to 0, the cell
216  // is now added to the cluster - careful ! (PJ, 19/07/08)
217  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
218  // (PJ, 15/09/08 <- similar to what existed before the
219  // previous bug fix, but keeps the close seeds inside,
220  // even if their fraction was set to zero.)
221  // Also add a protection to keep the seed in the cluster
222  // when the latter gets far from the former. These cases
223  // (about 1% of the clusters) need to be studied, as
224  // they create fake photons, in general.
225  // (PJ, 16/09/08)
226  if (dist2[i] < 100.0 || frac[i] > 0.9999) {
227  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit, frac[i]));
228  }
229  }
230  }
231  // recalculate positions and calculate convergence parameter
232  double diff2 = 0.0;
233  for (unsigned i = 0; i < clusters.size(); ++i) {
234  if (_convergencePosCalc) {
235  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
236  } else {
237  if (clusters.size() == 1 && _allCellsPosCalc) {
238  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
239  } else {
240  _positionCalc->calculateAndSetPosition(clusters[i]);
241  }
242  }
243  const double delta2 = reco::deltaR2(clusters[i].positionREP(), clus_prev_pos[i]);
244  if (delta2 > diff2)
245  diff2 = delta2;
246  }
247  diff = std::sqrt(diff2);
248  dist2.clear();
249  frac.clear();
250  clus_prev_pos.clear(); // avoid badness
251  growPFClusters(topo, seedable, toleranceScaling, iter + 1, diff, clusters);
252 }
253 
255  for (auto& cluster : clusters) {
256  cluster.pruneUsing([&](const reco::PFRecHitFraction& rhf) { return rhf.fraction() > _minFractionToKeep; });
257  }
258 }
reco::PFClusterCollection
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
change_name.diff
diff
Definition: change_name.py:13
reco::CaloCluster::setSeed
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
reco::PFRecHitFraction::fraction
double fraction() const
Definition: PFRecHitFraction.h:30
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
Basic2DGenericPFlowClusterizer::_minFracTot
const double _minFracTot
Definition: Basic2DGenericPFlowClusterizer.h:36
Basic2DGenericPFlowClusterizer::_convergencePosCalc
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
Definition: Basic2DGenericPFlowClusterizer.h:41
input
static const std::string input
Definition: EdmProvDump.cc:48
MessageLogger.h
particleFlowZeroSuppressionECAL_cff.thresholds
thresholds
Definition: particleFlowZeroSuppressionECAL_cff.py:31
reco::PFRecHitFraction
Fraction of a PFRecHit (rechits can be shared between several PFCluster's)
Definition: PFRecHitFraction.h:18
Basic2DGenericPFlowClusterizer::seedPFClustersFromTopo
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
Definition: Basic2DGenericPFlowClusterizer.cc:113
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
LOGDRESSED
#define LOGDRESSED(x)
Definition: Basic2DGenericPFlowClusterizer.cc:21
mps_splice.entry
entry
Definition: mps_splice.py:68
Basic2DGenericPFlowClusterizer::growPFClusters
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
Definition: Basic2DGenericPFlowClusterizer.cc:132
PFLayer::HCAL_ENDCAP
Definition: PFLayer.h:37
reco::PFCluster::recHitFractions
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
Basic2DGenericPFlowClusterizer::Basic2DGenericPFlowClusterizer
Basic2DGenericPFlowClusterizer(const edm::ParameterSet &conf)
Definition: Basic2DGenericPFlowClusterizer.cc:24
Basic2DGenericPFlowClusterizer::_recHitEnergyNorms
std::unordered_map< int, std::pair< std::vector< int >, std::vector< double > > > _recHitEnergyNorms
Definition: Basic2DGenericPFlowClusterizer.h:39
Basic2DGenericPFlowClusterizer::_excludeOtherSeeds
const bool _excludeOtherSeeds
Definition: Basic2DGenericPFlowClusterizer.h:35
Basic2DGenericPFlowClusterizer::_maxIterations
const unsigned _maxIterations
Definition: Basic2DGenericPFlowClusterizer.h:32
PFRecHit.h
edm::Ref< PFRecHitCollection >
PFLayer::ECAL_BARREL
Definition: PFLayer.h:33
myMath::fast_expf
float fast_expf(float x)
Definition: EcalUncalibRecHitRatioMethodAlgo.h:27
DivergingColor.frac
float frac
Definition: DivergingColor.py:175
PFLayer::PS1
Definition: PFLayer.h:31
PFLayer::HCAL_BARREL2
Definition: PFLayer.h:36
PFLayer::HF_EM
Definition: PFLayer.h:38
Basic2DGenericPFlowClusterizer.h
HLT_FULL_cff.fraction
fraction
Definition: HLT_FULL_cff.py:52795
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
PFLayer::HCAL_BARREL1
Definition: PFLayer.h:35
HLT_FULL_cff.recHitEnergyNorm
recHitEnergyNorm
Definition: HLT_FULL_cff.py:13279
PFLayer::NONE
Definition: PFLayer.h:34
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
PFLayer::HF_HAD
Definition: PFLayer.h:39
edm::ParameterSet
Definition: ParameterSet.h:47
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
deltaR.h
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
Basic2DGenericPFlowClusterizer::buildClusters
void buildClusters(const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus) override
Definition: Basic2DGenericPFlowClusterizer.cc:83
math::XYZVector
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
createfilelist.int
int
Definition: createfilelist.py:10
get
#define get
HLT_FULL_cff.depths
depths
Definition: HLT_FULL_cff.py:13412
reco::PFCluster::addRecHitFraction
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
Basic2DGenericPFlowClusterizer::_stoppingTolerance
const double _stoppingTolerance
Definition: Basic2DGenericPFlowClusterizer.h:33
Exception
Definition: hltDiff.cc:246
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PFClusterBuilderBase
Definition: PFClusterBuilderBase.h:18
edm::Ref::key
key_type key() const
Accessor for product key.
Definition: Ref.h:250
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::PFCluster::REPPoint
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
PFClusterBuilderBase::_positionCalc
std::unique_ptr< PosCalc > _positionCalc
Definition: PFClusterBuilderBase.h:57
PFClusterBuilderBase::_minFractionToKeep
const float _minFractionToKeep
Definition: PFClusterBuilderBase.h:56
Basic2DGenericPFlowClusterizer::_allCellsPosCalc
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
Definition: Basic2DGenericPFlowClusterizer.h:40
PFLayer::ECAL_ENDCAP
Definition: PFLayer.h:32
PFLayer::PS2
Definition: PFLayer.h:30
Basic2DGenericPFlowClusterizer::_showerSigma2
const double _showerSigma2
Definition: Basic2DGenericPFlowClusterizer.h:34
edm::ParameterSet::getParameterSet
ParameterSet const & getParameterSet(std::string const &) const
Definition: ParameterSet.cc:2128
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
Basic2DGenericPFlowClusterizer::prunePFClusters
void prunePFClusters(reco::PFClusterCollection &) const
Definition: Basic2DGenericPFlowClusterizer.cc:254