CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
PFlow2DClusterizerWithTime Class Reference
Inheritance diagram for PFlow2DClusterizerWithTime:
PFClusterBuilderBase

Public Member Functions

void buildClusters (const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus) override
 
B2DGPFoperator= (const B2DGPF &)=delete
 
 PFlow2DClusterizerWithTime (const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
 
 PFlow2DClusterizerWithTime (const B2DGPF &)=delete
 
void update (const edm::EventSetup &es) override
 
 ~PFlow2DClusterizerWithTime () override=default
 
- Public Member Functions inherited from PFClusterBuilderBase
std::ostream & operator<< (std::ostream &o) const
 
PFCBBoperator= (const PFCBB &)=delete
 
 PFClusterBuilderBase (const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
 
 PFClusterBuilderBase (const PFCBB &)=delete
 
void reset ()
 
virtual ~PFClusterBuilderBase ()=default
 

Private Types

typedef PFlow2DClusterizerWithTime B2DGPF
 

Private Member Functions

void clusterTimeResolution (reco::PFCluster &cluster, double &res) const
 
void clusterTimeResolutionFromSeed (reco::PFCluster &cluster, double &res) const
 
double dist2Time (const reco::PFCluster &, const reco::PFRecHitRef &, int cell_layer, double prev_timeres2) const
 
void growPFClusters (const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
 
bool passChi2Prob (size_t iCluster, double dist2, std::vector< double > &clus_chi2, std::vector< size_t > &clus_chi2_nhits) const
 
void prunePFClusters (reco::PFClusterCollection &) const
 
void seedPFClustersFromTopo (const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
 

Private Attributes

std::unique_ptr
< PFCPositionCalculatorBase
_allCellsPosCalc
 
const bool _clusterTimeResFromSeed
 
std::unique_ptr
< PFCPositionCalculatorBase
_convergencePosCalc
 
const bool _excludeOtherSeeds
 
const std::unordered_map
< std::string, int > 
_layerMap
 
const unsigned _maxIterations
 
const double _maxNSigmaTime
 
const double _minChi2Prob
 
const double _minFracTot
 
std::unordered_map< int, double > _recHitEnergyNorms
 
const double _showerSigma2
 
const double _stoppingTolerance
 
std::unique_ptr
< CaloRecHitResolutionProvider
_timeResolutionCalcBarrel
 
std::unique_ptr
< CaloRecHitResolutionProvider
_timeResolutionCalcEndcap
 
const double _timeSigma_eb
 
const double _timeSigma_ee
 

Additional Inherited Members

- Public Types inherited from PFClusterBuilderBase
typedef PFCPositionCalculatorBase PosCalc
 
- Protected Attributes inherited from PFClusterBuilderBase
const float _minFractionToKeep
 
unsigned _nClustersFound
 
unsigned _nSeeds
 
std::unique_ptr< PosCalc_positionCalc
 

Detailed Description

Definition at line 16 of file PFlow2DClusterizerWithTime.cc.

Member Typedef Documentation

Definition at line 17 of file PFlow2DClusterizerWithTime.cc.

Constructor & Destructor Documentation

PFlow2DClusterizerWithTime::PFlow2DClusterizerWithTime ( const edm::ParameterSet conf,
edm::ConsumesCollector cc 
)

Definition at line 90 of file PFlow2DClusterizerWithTime.cc.

References PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, edm::ParameterSet::getParameterSetVector(), PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PFLayer::HCAL_ENDCAP, PFLayer::HF_EM, PFLayer::HF_HAD, PFLayer::NONE, PFLayer::PS1, PFLayer::PS2, and lowPtGsfElectronSeeds_cfi::thresholds().

91  : PFClusterBuilderBase(conf, cc),
92  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
93  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
94  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"), 2.0)),
95  _timeSigma_eb(std::pow(conf.getParameter<double>("timeSigmaEB"), 2.0)),
96  _timeSigma_ee(std::pow(conf.getParameter<double>("timeSigmaEE"), 2.0)),
97  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
98  _minFracTot(conf.getParameter<double>("minFracTot")),
99  _maxNSigmaTime(std::pow(conf.getParameter<double>("maxNSigmaTime"), 2.0)),
100  _minChi2Prob(conf.getParameter<double>("minChi2Prob")),
101  _clusterTimeResFromSeed(conf.getParameter<bool>("clusterTimeResFromSeed")),
102 
103  _layerMap({{"PS2", (int)PFLayer::PS2},
104  {"PS1", (int)PFLayer::PS1},
105  {"ECAL_ENDCAP", (int)PFLayer::ECAL_ENDCAP},
106  {"ECAL_BARREL", (int)PFLayer::ECAL_BARREL},
107  {"NONE", (int)PFLayer::NONE},
108  {"HCAL_BARREL1", (int)PFLayer::HCAL_BARREL1},
109  {"HCAL_BARREL2_RING0", (int)PFLayer::HCAL_BARREL2},
110  {"HCAL_BARREL2_RING1", 100 * (int)PFLayer::HCAL_BARREL2},
111  {"HCAL_ENDCAP", (int)PFLayer::HCAL_ENDCAP},
112  {"HF_EM", (int)PFLayer::HF_EM},
113  {"HF_HAD", (int)PFLayer::HF_HAD}}) {
const std::unordered_map< std::string, int > _layerMap
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PFClusterBuilderBase(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
PFlow2DClusterizerWithTime::~PFlow2DClusterizerWithTime ( )
overridedefault
PFlow2DClusterizerWithTime::PFlow2DClusterizerWithTime ( const B2DGPF )
delete

Member Function Documentation

void PFlow2DClusterizerWithTime::buildClusters ( const reco::PFClusterCollection input,
const std::vector< bool > &  seedable,
reco::PFClusterCollection outclus 
)
overridevirtual

Implements PFClusterBuilderBase.

Definition at line 147 of file PFlow2DClusterizerWithTime.cc.

References _allCellsPosCalc, _convergencePosCalc, PFClusterBuilderBase::_positionCalc, growPFClusters(), SiStripPI::max, eostools::move(), funct::pow(), prunePFClusters(), and seedPFClustersFromTopo().

149  {
150  reco::PFClusterCollection clustersInTopo;
151  for (const auto& topocluster : input) {
152  clustersInTopo.clear();
153  seedPFClustersFromTopo(topocluster, seedable, clustersInTopo);
154  const unsigned tolScal = std::pow(std::max(1.0, clustersInTopo.size() - 1.0), 2.0);
155  growPFClusters(topocluster, seedable, tolScal, 0, tolScal, clustersInTopo);
156  // step added by Josh Bendavid, removes low-fraction clusters
157  // did not impact position resolution with fraction cut of 1e-7
158  // decreases the size of each pf cluster considerably
159  prunePFClusters(clustersInTopo);
160  // recalculate the positions of the pruned clusters
161  if (_convergencePosCalc) {
162  // if defined, use the special position calculation for convergence tests
163  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
164  } else {
165  if (clustersInTopo.size() == 1 && _allCellsPosCalc) {
166  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
167  } else {
168  _positionCalc->calculateAndSetPositions(clustersInTopo);
169  }
170  }
171  for (auto& clusterout : clustersInTopo) {
172  output.insert(output.end(), std::move(clusterout));
173  }
174  }
175 }
void prunePFClusters(reco::PFClusterCollection &) const
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
std::unique_ptr< PosCalc > _positionCalc
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
static std::string const input
Definition: EdmProvDump.cc:47
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
def move
Definition: eostools.py:511
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
void PFlow2DClusterizerWithTime::clusterTimeResolution ( reco::PFCluster cluster,
double &  res 
) const
private

Definition at line 355 of file PFlow2DClusterizerWithTime.cc.

References _timeResolutionCalcBarrel, _timeResolutionCalcEndcap, PFLayer::ECAL_BARREL, reco::PFRecHit::energy(), PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PixelPluginsPhase0_cfi::isBarrel, reco::PFRecHit::layer(), reco::PFCluster::recHitFractions(), reco::PFCluster::setTime(), reco::PFCluster::setTimeError(), mathSSE::sqrt(), and reco::PFRecHit::time().

Referenced by growPFClusters().

355  {
356  double sumTimeSigma2 = 0.;
357  double sumSigma2 = 0.;
358 
359  for (auto& rhf : cluster.recHitFractions()) {
360  const reco::PFRecHit& rh = *(rhf.recHitRef());
361  const double rhf_f = rhf.fraction();
362 
363  if (rhf_f == 0.)
364  continue;
365 
366  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 || rh.layer() == PFLayer::HCAL_BARREL2 ||
367  rh.layer() == PFLayer::ECAL_BARREL);
368  double res2 = 10000.;
369  if (isBarrel)
370  res2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
371  else
372  res2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
373 
374  sumTimeSigma2 += rhf_f * rh.time() / res2;
375  sumSigma2 += rhf_f / res2;
376  }
377  if (sumSigma2 > 0.) {
378  clusterRes2 = 1. / sumSigma2;
379  cluster.setTime(sumTimeSigma2 / sumSigma2);
380  cluster.setTimeError(std::sqrt(float(clusterRes2)));
381  } else {
382  clusterRes2 = 999999.;
383  }
384 }
float time() const
timing for cleaned hits
Definition: PFRecHit.h:102
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
T sqrt(T t)
Definition: SSEVec.h:19
float energy() const
rechit energy
Definition: PFRecHit.h:99
void setTimeError(float timeError)
Definition: PFCluster.h:88
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
void PFlow2DClusterizerWithTime::clusterTimeResolutionFromSeed ( reco::PFCluster cluster,
double &  res 
) const
private

Definition at line 338 of file PFlow2DClusterizerWithTime.cc.

References _timeResolutionCalcBarrel, _timeResolutionCalcEndcap, reco::PFRecHit::detId(), PFLayer::ECAL_BARREL, reco::PFRecHit::energy(), PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PixelPluginsPhase0_cfi::isBarrel, reco::PFRecHit::layer(), reco::PFCluster::recHitFractions(), reco::CaloCluster::seed(), reco::PFCluster::setTime(), reco::PFCluster::setTimeError(), mathSSE::sqrt(), and reco::PFRecHit::time().

Referenced by growPFClusters().

338  {
339  clusterRes2 = 10000.;
340  for (auto& rhf : cluster.recHitFractions()) {
341  const reco::PFRecHit& rh = *(rhf.recHitRef());
342  if (rh.detId() == cluster.seed()) {
343  cluster.setTime(rh.time());
344  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 || rh.layer() == PFLayer::HCAL_BARREL2 ||
345  rh.layer() == PFLayer::ECAL_BARREL);
346  if (isBarrel)
347  clusterRes2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
348  else
349  clusterRes2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
350  cluster.setTimeError(std::sqrt(float(clusterRes2)));
351  }
352  }
353 }
float time() const
timing for cleaned hits
Definition: PFRecHit.h:102
unsigned detId() const
rechit detId
Definition: PFRecHit.h:93
void setTime(float time, float timeError=0)
Definition: PFCluster.h:84
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:96
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
T sqrt(T t)
Definition: SSEVec.h:19
float energy() const
rechit energy
Definition: PFRecHit.h:99
void setTimeError(float timeError)
Definition: PFCluster.h:88
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
double PFlow2DClusterizerWithTime::dist2Time ( const reco::PFCluster cluster,
const reco::PFRecHitRef refhit,
int  cell_layer,
double  prev_timeres2 
) const
private

Definition at line 386 of file PFlow2DClusterizerWithTime.cc.

References _maxNSigmaTime, _timeResolutionCalcBarrel, _timeResolutionCalcEndcap, _timeSigma_eb, _timeSigma_ee, PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PFLayer::HCAL_ENDCAP, PFLayer::HF_EM, PFLayer::HF_HAD, and reco::PFCluster::time().

Referenced by growPFClusters().

389  {
390  const double deltaT = cluster.time() - refhit->time();
391  const double t2 = deltaT * deltaT;
392  double res2 = 100.;
393 
394  if (cell_layer == PFLayer::HCAL_BARREL1 || cell_layer == PFLayer::HCAL_BARREL2 ||
395  cell_layer == PFLayer::ECAL_BARREL) {
397  const double resCluster2 = prev_timeres2;
398  res2 = resCluster2 + _timeResolutionCalcBarrel->timeResolution2(refhit->energy());
399  } else {
400  return t2 / _timeSigma_eb;
401  }
402  } else if (cell_layer == PFLayer::HCAL_ENDCAP || cell_layer == PFLayer::HF_EM || cell_layer == PFLayer::HF_HAD ||
403  cell_layer == PFLayer::ECAL_ENDCAP) {
405  const double resCluster2 = prev_timeres2;
406  res2 = resCluster2 + _timeResolutionCalcEndcap->timeResolution2(refhit->energy());
407  } else {
408  return t2 / _timeSigma_ee;
409  }
410  }
411 
412  double distTime2 = t2 / res2;
413  if (distTime2 > _maxNSigmaTime)
414  return 999.; // reject hit
415 
416  return distTime2;
417 }
float time() const
Definition: PFCluster.h:77
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
void PFlow2DClusterizerWithTime::growPFClusters ( const reco::PFCluster topo,
const std::vector< bool > &  seedable,
const unsigned  toleranceScaling,
const unsigned  iter,
double  dist,
reco::PFClusterCollection clusters 
) const
private

Definition at line 196 of file PFlow2DClusterizerWithTime.cc.

References _allCellsPosCalc, _clusterTimeResFromSeed, _convergencePosCalc, _excludeOtherSeeds, _maxIterations, _minChi2Prob, _minFracTot, PFClusterBuilderBase::_positionCalc, _recHitEnergyNorms, _showerSigma2, _stoppingTolerance, funct::abs(), clusterTimeResolution(), clusterTimeResolutionFromSeed(), reco::deltaR2(), change_name::diff, dist2Time(), reco::PFCluster::energy(), myMath::fast_expf(), DivergingColor::frac, HLT_FULL_cff::fraction, PFLayer::HCAL_BARREL2, mps_fire::i, edm::Ref< C, T, F >::key(), LOGDRESSED, passChi2Prob(), reco::CaloCluster::position(), reco::PFCluster::recHitFractions(), reco::CaloCluster::seed(), and mathSSE::sqrt().

Referenced by buildClusters().

201  {
202  if (iter >= _maxIterations) {
203  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
204  << "reached " << _maxIterations << " iterations, terminated position "
205  << "fit with diff = " << diff;
206  }
207  if (iter >= _maxIterations || diff <= _stoppingTolerance * toleranceScaling)
208  return;
209  // reset the rechits in this cluster, keeping the previous position
210  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
211  // also calculate and keep the previous time resolution
212  std::vector<double> clus_prev_timeres2;
213 
214  for (auto& cluster : clusters) {
215  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
216  clus_prev_pos.emplace_back(repp.rho(), repp.eta(), repp.phi());
217  if (_convergencePosCalc) {
218  if (clusters.size() == 1 && _allCellsPosCalc) {
219  _allCellsPosCalc->calculateAndSetPosition(cluster);
220  } else {
221  _positionCalc->calculateAndSetPosition(cluster);
222  }
223  }
224  double resCluster2;
226  clusterTimeResolutionFromSeed(cluster, resCluster2);
227  else
228  clusterTimeResolution(cluster, resCluster2);
229  clus_prev_timeres2.push_back(resCluster2);
230  cluster.resetHitsAndFractions();
231  }
232  // loop over topo cluster and grow current PFCluster hypothesis
233  std::vector<double> dist2, frac;
234 
235  // Store chi2 values and nhits to calculate chi2 prob. inline
236  std::vector<double> clus_chi2;
237  std::vector<size_t> clus_chi2_nhits;
238 
239  double fractot = 0;
240  for (const reco::PFRecHitFraction& rhf : topo.recHitFractions()) {
241  const reco::PFRecHitRef& refhit = rhf.recHitRef();
242  int cell_layer = (int)refhit->layer();
243  if (cell_layer == PFLayer::HCAL_BARREL2 && std::abs(refhit->positionREP().eta()) > 0.34) {
244  cell_layer *= 100;
245  }
246  const double recHitEnergyNorm = _recHitEnergyNorms.find(cell_layer)->second;
247  math::XYZPoint topocellpos_xyz(refhit->position());
248  dist2.clear();
249  frac.clear();
250  fractot = 0;
251 
252  // add rechits to clusters, calculating fraction based on distance
253  // need position in vector to get cluster time resolution
254  for (size_t iCluster = 0; iCluster < clusters.size(); ++iCluster) {
255  reco::PFCluster& cluster = clusters[iCluster];
256  const math::XYZPoint& clusterpos_xyz = cluster.position();
257  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
258  double d2 = deltav.Mag2() / _showerSigma2;
259 
260  double d2time = dist2Time(cluster, refhit, cell_layer, clus_prev_timeres2[iCluster]);
261  d2 += d2time;
262 
263  if (_minChi2Prob > 0. && !passChi2Prob(iCluster, d2time, clus_chi2, clus_chi2_nhits))
264  d2 = 999.;
265 
266  dist2.emplace_back(d2);
267 
268  if (d2 > 100) {
269  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
270  << "Warning! :: pfcluster-topocell distance is too large! d= " << d2;
271  }
272  // fraction assignment logic
273  double fraction;
274  if (refhit->detId() == cluster.seed() && _excludeOtherSeeds) {
275  fraction = 1.0;
276  } else if (seedable[refhit.key()] && _excludeOtherSeeds) {
277  fraction = 0.0;
278  } else {
279  fraction = cluster.energy() / recHitEnergyNorm * vdt::fast_expf(-0.5 * d2);
280  }
281  fractot += fraction;
282  frac.emplace_back(fraction);
283  }
284  for (unsigned i = 0; i < clusters.size(); ++i) {
285  if (fractot > _minFracTot || (refhit->detId() == clusters[i].seed() && fractot > 0.0)) {
286  frac[i] /= fractot;
287  } else {
288  continue;
289  }
290  // if the fraction has been set to 0, the cell
291  // is now added to the cluster - careful ! (PJ, 19/07/08)
292  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
293  // (PJ, 15/09/08 <- similar to what existed before the
294  // previous bug fix, but keeps the close seeds inside,
295  // even if their fraction was set to zero.)
296  // Also add a protection to keep the seed in the cluster
297  // when the latter gets far from the former. These cases
298  // (about 1% of the clusters) need to be studied, as
299  // they create fake photons, in general.
300  // (PJ, 16/09/08)
301  if (dist2[i] < 100.0 || frac[i] > 0.9999) {
302  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit, frac[i]));
303  }
304  }
305  }
306  // recalculate positions and calculate convergence parameter
307  double diff2 = 0.0;
308  for (unsigned i = 0; i < clusters.size(); ++i) {
309  if (_convergencePosCalc) {
310  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
311  } else {
312  if (clusters.size() == 1 && _allCellsPosCalc) {
313  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
314  } else {
315  _positionCalc->calculateAndSetPosition(clusters[i]);
316  }
317  }
318  const double delta2 = reco::deltaR2(clusters[i].positionREP(), clus_prev_pos[i]);
319  if (delta2 > diff2)
320  diff2 = delta2;
321  }
322  diff = std::sqrt(diff2);
323  dist2.clear();
324  frac.clear();
325  clus_prev_pos.clear(); // avoid badness
326  clus_chi2.clear();
327  clus_chi2_nhits.clear();
328  clus_prev_timeres2.clear();
329  growPFClusters(topo, seedable, toleranceScaling, iter + 1, diff, clusters);
330 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
bool passChi2Prob(size_t iCluster, double dist2, std::vector< double > &clus_chi2, std::vector< size_t > &clus_chi2_nhits) const
std::unique_ptr< PosCalc > _positionCalc
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
#define LOGDRESSED(x)
key_type key() const
Accessor for product key.
Definition: Ref.h:250
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
void clusterTimeResolution(reco::PFCluster &cluster, double &res) const
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double dist2Time(const reco::PFCluster &, const reco::PFRecHitRef &, int cell_layer, double prev_timeres2) const
double energy() const
cluster energy
Definition: PFCluster.h:74
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:219
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
float fast_expf(float x)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:48
void clusterTimeResolutionFromSeed(reco::PFCluster &cluster, double &res) const
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
std::unordered_map< int, double > _recHitEnergyNorms
B2DGPF& PFlow2DClusterizerWithTime::operator= ( const B2DGPF )
delete
bool PFlow2DClusterizerWithTime::passChi2Prob ( size_t  iCluster,
double  dist2,
std::vector< double > &  clus_chi2,
std::vector< size_t > &  clus_chi2_nhits 
) const
private

Definition at line 419 of file PFlow2DClusterizerWithTime.cc.

References _minChi2Prob, and HLT_FULL_cff::chi2.

Referenced by growPFClusters().

422  {
423  if (iCluster >= clus_chi2.size()) { // first hit
424  clus_chi2.push_back(dist2);
425  clus_chi2_nhits.push_back(1);
426  return true;
427  } else {
428  double chi2 = clus_chi2[iCluster];
429  size_t nhitsCluster = clus_chi2_nhits[iCluster];
430  chi2 += dist2;
431  if (TMath::Prob(chi2, nhitsCluster) >= _minChi2Prob) {
432  clus_chi2[iCluster] = chi2;
433  clus_chi2_nhits[iCluster] = nhitsCluster + 1;
434  return true;
435  }
436  }
437  return false;
438 }
void PFlow2DClusterizerWithTime::prunePFClusters ( reco::PFClusterCollection clusters) const
private

Definition at line 332 of file PFlow2DClusterizerWithTime.cc.

References PFClusterBuilderBase::_minFractionToKeep, and reco::PFRecHitFraction::fraction().

Referenced by buildClusters().

332  {
333  for (auto& cluster : clusters) {
334  cluster.pruneUsing([&](const reco::PFRecHitFraction& rhf) { return rhf.fraction() > _minFractionToKeep; });
335  }
336 }
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
double fraction() const
void PFlow2DClusterizerWithTime::seedPFClustersFromTopo ( const reco::PFCluster topo,
const std::vector< bool > &  seedable,
reco::PFClusterCollection initialPFClusters 
) const
private

Definition at line 177 of file PFlow2DClusterizerWithTime.cc.

References _convergencePosCalc, PFClusterBuilderBase::_positionCalc, reco::PFCluster::addRecHitFraction(), reco::PFCluster::recHitFractions(), and reco::CaloCluster::setSeed().

Referenced by buildClusters().

179  {
180  const auto& recHitFractions = topo.recHitFractions();
181  for (const auto& rhf : recHitFractions) {
182  if (!seedable[rhf.recHitRef().key()])
183  continue;
184  initialPFClusters.push_back(reco::PFCluster());
185  reco::PFCluster& current = initialPFClusters.back();
186  current.addRecHitFraction(rhf);
187  current.setSeed(rhf.recHitRef()->detId());
188  if (_convergencePosCalc) {
189  _convergencePosCalc->calculateAndSetPosition(current);
190  } else {
191  _positionCalc->calculateAndSetPosition(current);
192  }
193  }
194 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
std::unique_ptr< PosCalc > _positionCalc
void setSeed(const DetId &id)
Definition: CaloCluster.h:146
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:65
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:33
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
void PFlow2DClusterizerWithTime::update ( const edm::EventSetup es)
inlineoverridevirtual

Reimplemented from PFClusterBuilderBase.

Definition at line 26 of file PFlow2DClusterizerWithTime.cc.

References _allCellsPosCalc, _convergencePosCalc, and PFClusterBuilderBase::_positionCalc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

26  {
27  _positionCalc->update(es);
28  if (_allCellsPosCalc)
29  _allCellsPosCalc->update(es);
31  _convergencePosCalc->update(es);
32  }
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
std::unique_ptr< PosCalc > _positionCalc
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc

Member Data Documentation

std::unique_ptr<PFCPositionCalculatorBase> PFlow2DClusterizerWithTime::_allCellsPosCalc
private

Definition at line 52 of file PFlow2DClusterizerWithTime.cc.

Referenced by buildClusters(), growPFClusters(), if(), and update().

const bool PFlow2DClusterizerWithTime::_clusterTimeResFromSeed
private

Definition at line 48 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

std::unique_ptr<PFCPositionCalculatorBase> PFlow2DClusterizerWithTime::_convergencePosCalc
private
const bool PFlow2DClusterizerWithTime::_excludeOtherSeeds
private

Definition at line 44 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

const std::unordered_map<std::string, int> PFlow2DClusterizerWithTime::_layerMap
private

Definition at line 50 of file PFlow2DClusterizerWithTime.cc.

Referenced by for().

const unsigned PFlow2DClusterizerWithTime::_maxIterations
private

Definition at line 39 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

const double PFlow2DClusterizerWithTime::_maxNSigmaTime
private

Definition at line 46 of file PFlow2DClusterizerWithTime.cc.

Referenced by dist2Time().

const double PFlow2DClusterizerWithTime::_minChi2Prob
private

Definition at line 47 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters(), and passChi2Prob().

const double PFlow2DClusterizerWithTime::_minFracTot
private

Definition at line 45 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

std::unordered_map<int, double> PFlow2DClusterizerWithTime::_recHitEnergyNorms
private

Definition at line 51 of file PFlow2DClusterizerWithTime.cc.

Referenced by for(), and growPFClusters().

const double PFlow2DClusterizerWithTime::_showerSigma2
private

Definition at line 41 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

const double PFlow2DClusterizerWithTime::_stoppingTolerance
private

Definition at line 40 of file PFlow2DClusterizerWithTime.cc.

Referenced by growPFClusters().

std::unique_ptr<CaloRecHitResolutionProvider> PFlow2DClusterizerWithTime::_timeResolutionCalcBarrel
private
std::unique_ptr<CaloRecHitResolutionProvider> PFlow2DClusterizerWithTime::_timeResolutionCalcEndcap
private
const double PFlow2DClusterizerWithTime::_timeSigma_eb
private

Definition at line 42 of file PFlow2DClusterizerWithTime.cc.

Referenced by dist2Time().

const double PFlow2DClusterizerWithTime::_timeSigma_ee
private

Definition at line 43 of file PFlow2DClusterizerWithTime.cc.

Referenced by dist2Time().