CMS 3D CMS Logo

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

#include <PFlow2DClusterizerWithTime.h>

Inheritance diagram for PFlow2DClusterizerWithTime:
PFClusterBuilderBase

Public Member Functions

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

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 11 of file PFlow2DClusterizerWithTime.h.

Member Typedef Documentation

Definition at line 12 of file PFlow2DClusterizerWithTime.h.

Constructor & Destructor Documentation

PFlow2DClusterizerWithTime::PFlow2DClusterizerWithTime ( const edm::ParameterSet conf)

Definition at line 27 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, and PFLayer::PS2.

27  :
29  _maxIterations(conf.getParameter<unsigned>("maxIterations")),
30  _stoppingTolerance(conf.getParameter<double>("stoppingTolerance")),
31  _showerSigma2(std::pow(conf.getParameter<double>("showerSigma"),2.0)),
32  _timeSigma_eb(std::pow(conf.getParameter<double>("timeSigmaEB"),2.0)),
33  _timeSigma_ee(std::pow(conf.getParameter<double>("timeSigmaEE"),2.0)),
34  _excludeOtherSeeds(conf.getParameter<bool>("excludeOtherSeeds")),
35  _minFracTot(conf.getParameter<double>("minFracTot")),
36  _maxNSigmaTime(std::pow(conf.getParameter<double>("maxNSigmaTime"),2.0)),
37  _minChi2Prob(conf.getParameter<double>("minChi2Prob")),
38  _clusterTimeResFromSeed(conf.getParameter<bool>("clusterTimeResFromSeed")),
39 
40  _layerMap({ {"PS2",(int)PFLayer::PS2},
41  {"PS1",(int)PFLayer::PS1},
42  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
43  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
44  {"NONE",(int)PFLayer::NONE},
45  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
46  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
47  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
48  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
49  {"HF_EM",(int)PFLayer::HF_EM},
50  {"HF_HAD",(int)PFLayer::HF_HAD} }) {
T getParameter(std::string const &) const
PFClusterBuilderBase(const edm::ParameterSet &conf)
const std::unordered_map< std::string, int > _layerMap
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual PFlow2DClusterizerWithTime::~PFlow2DClusterizerWithTime ( )
inlinevirtual

Definition at line 16 of file PFlow2DClusterizerWithTime.h.

16 {}
PFlow2DClusterizerWithTime::PFlow2DClusterizerWithTime ( const B2DGPF )
delete

Member Function Documentation

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

Implements PFClusterBuilderBase.

Definition at line 105 of file PFlow2DClusterizerWithTime.cc.

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

107  {
108  reco::PFClusterCollection clustersInTopo;
109  for( const auto& topocluster : input ) {
110  clustersInTopo.clear();
111  seedPFClustersFromTopo(topocluster,seedable,clustersInTopo);
112  const unsigned tolScal =
113  std::pow(std::max(1.0,clustersInTopo.size()-1.0),2.0);
114  growPFClusters(topocluster,seedable,tolScal,0,tolScal,clustersInTopo);
115  // step added by Josh Bendavid, removes low-fraction clusters
116  // did not impact position resolution with fraction cut of 1e-7
117  // decreases the size of each pf cluster considerably
118  prunePFClusters(clustersInTopo);
119  // recalculate the positions of the pruned clusters
120  if( _convergencePosCalc ) {
121  // if defined, use the special position calculation for convergence tests
122  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
123  } else {
124  if( clustersInTopo.size() == 1 && _allCellsPosCalc ) {
125  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
126  } else {
127  _positionCalc->calculateAndSetPositions(clustersInTopo);
128  }
129  }
130  for( auto& clusterout : clustersInTopo ) {
131  output.insert(output.end(),std::move(clusterout));
132  }
133  }
134 }
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:43
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
def move
Definition: eostools.py:510
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:40
void PFlow2DClusterizerWithTime::clusterTimeResolution ( reco::PFCluster cluster,
double &  res 
) const
private

Definition at line 322 of file PFlow2DClusterizerWithTime.cc.

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

Referenced by growPFClusters().

324 {
325  double sumTimeSigma2 = 0.;
326  double sumSigma2 = 0.;
327 
328  for (auto& rhf : cluster.recHitFractions())
329  {
330  const reco::PFRecHit& rh = *(rhf.recHitRef());
331  const double rhf_f = rhf.fraction();
332 
333  if (rhf_f == 0.)
334  continue;
335 
336  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 ||
337  rh.layer() == PFLayer::HCAL_BARREL2 ||
338  rh.layer() == PFLayer::ECAL_BARREL);
339  double res2 = 10000.;
340  if (isBarrel)
341  res2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
342  else
343  res2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
344 
345  sumTimeSigma2 += rhf_f * rh.time()/res2;
346  sumSigma2 += rhf_f/res2;
347  }
348  if (sumSigma2 > 0.) {
349  clusterRes2 = 1./sumSigma2;
350  cluster.setTime(sumTimeSigma2/sumSigma2);
351  } else {
352  clusterRes2 = 999999.;
353  }
354 }
bool isBarrel(GeomDetEnumerators::SubDetector m)
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:109
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
void setTime(double time)
Definition: PFCluster.h:89
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
double energy() const
rechit energy
Definition: PFRecHit.h:112
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
double time() const
timing for cleaned hits
Definition: PFRecHit.h:116
void PFlow2DClusterizerWithTime::clusterTimeResolutionFromSeed ( reco::PFCluster cluster,
double &  res 
) const
private

Definition at line 301 of file PFlow2DClusterizerWithTime.cc.

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

Referenced by growPFClusters().

303 {
304  clusterRes2 = 10000.;
305  for (auto& rhf : cluster.recHitFractions())
306  {
307  const reco::PFRecHit& rh = *(rhf.recHitRef());
308  if (rh.detId() == cluster.seed())
309  {
310  cluster.setTime(rh.time());
311  bool isBarrel = (rh.layer() == PFLayer::HCAL_BARREL1 ||
312  rh.layer() == PFLayer::HCAL_BARREL2 ||
313  rh.layer() == PFLayer::ECAL_BARREL);
314  if (isBarrel)
315  clusterRes2 = _timeResolutionCalcBarrel->timeResolution2(rh.energy());
316  else
317  clusterRes2 = _timeResolutionCalcEndcap->timeResolution2(rh.energy());
318  }
319  }
320 }
unsigned detId() const
rechit detId
Definition: PFRecHit.h:106
bool isBarrel(GeomDetEnumerators::SubDetector m)
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:109
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
void setTime(double time)
Definition: PFCluster.h:89
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
double energy() const
rechit energy
Definition: PFRecHit.h:112
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
double time() const
timing for cleaned hits
Definition: PFRecHit.h:116
double PFlow2DClusterizerWithTime::dist2Time ( const reco::PFCluster cluster,
const reco::PFRecHitRef refhit,
int  cell_layer,
double  prev_timeres2 
) const
private

Definition at line 356 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().

358 {
359  const double deltaT = cluster.time()-refhit->time();
360  const double t2 = deltaT*deltaT;
361  double res2 = 100.;
362 
363  if (cell_layer == PFLayer::HCAL_BARREL1 ||
364  cell_layer == PFLayer::HCAL_BARREL2 ||
365  cell_layer == PFLayer::ECAL_BARREL) {
367  const double resCluster2 = prev_timeres2;
368  res2 = resCluster2 + _timeResolutionCalcBarrel->timeResolution2(
369  refhit->energy());
370  }
371  else {
372  return t2/_timeSigma_eb;
373  }
374  }
375  else if (cell_layer == PFLayer::HCAL_ENDCAP ||
376  cell_layer == PFLayer::HF_EM ||
377  cell_layer == PFLayer::HF_HAD ||
378  cell_layer == PFLayer::ECAL_ENDCAP) {
380  const double resCluster2 = prev_timeres2;
381  res2 = resCluster2 + _timeResolutionCalcEndcap->timeResolution2(
382  refhit->energy());
383  }
384  else {
385  return t2/_timeSigma_ee;
386  }
387  }
388 
389  double distTime2 = t2/res2;
390  if (distTime2 > _maxNSigmaTime)
391  return 999.; // reject hit
392 
393  return distTime2;
394 }
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcBarrel
std::unique_ptr< CaloRecHitResolutionProvider > _timeResolutionCalcEndcap
double time() const
cluster time
Definition: PFCluster.h:85
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 156 of file PFlow2DClusterizerWithTime.cc.

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

Referenced by buildClusters().

161  {
162  if( iter >= _maxIterations ) {
163  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
164  <<"reached " << _maxIterations << " iterations, terminated position "
165  << "fit with diff = " << diff;
166  }
167  if( iter >= _maxIterations ||
168  diff <= _stoppingTolerance*toleranceScaling) return;
169  // reset the rechits in this cluster, keeping the previous position
170  std::vector<reco::PFCluster::REPPoint> clus_prev_pos;
171  // also calculate and keep the previous time resolution
172  std::vector<double> clus_prev_timeres2;
173 
174  for( auto& cluster : clusters) {
175  const reco::PFCluster::REPPoint& repp = cluster.positionREP();
176  clus_prev_pos.emplace_back(repp.rho(),repp.eta(),repp.phi());
177  if( _convergencePosCalc ) {
178  if( clusters.size() == 1 && _allCellsPosCalc ) {
179  _allCellsPosCalc->calculateAndSetPosition(cluster);
180  } else {
181  _positionCalc->calculateAndSetPosition(cluster);
182  }
183  }
184  double resCluster2;
186  clusterTimeResolutionFromSeed(cluster, resCluster2);
187  else
188  clusterTimeResolution(cluster, resCluster2);
189  clus_prev_timeres2.push_back(resCluster2);
190  cluster.resetHitsAndFractions();
191  }
192  // loop over topo cluster and grow current PFCluster hypothesis
193  std::vector<double> dist2, frac;
194 
195  // Store chi2 values and nhits to calculate chi2 prob. inline
196  std::vector<double> clus_chi2;
197  std::vector<size_t> clus_chi2_nhits;
198 
199  double fractot = 0, fraction = 0;
200  for( const reco::PFRecHitFraction& rhf : topo.recHitFractions() ) {
201  const reco::PFRecHitRef& refhit = rhf.recHitRef();
202  int cell_layer = (int)refhit->layer();
203  if( cell_layer == PFLayer::HCAL_BARREL2 &&
204  std::abs(refhit->positionREP().eta()) > 0.34 ) {
205  cell_layer *= 100;
206  }
207  const double recHitEnergyNorm =
208  _recHitEnergyNorms.find(cell_layer)->second;
209  const math::XYZPoint& topocellpos_xyz = refhit->position();
210  dist2.clear(); frac.clear(); fractot = 0;
211 
212  // add rechits to clusters, calculating fraction based on distance
213  // need position in vector to get cluster time resolution
214  for (size_t iCluster = 0; iCluster < clusters.size(); ++iCluster) {
215  reco::PFCluster& cluster = clusters[iCluster];
216  const math::XYZPoint& clusterpos_xyz = cluster.position();
217  fraction = 0.0;
218  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
219  double d2 = deltav.Mag2()/_showerSigma2;
220 
221  double d2time = dist2Time(cluster, refhit, cell_layer,
222  clus_prev_timeres2[iCluster]);
223  d2 += d2time;
224 
225  if (_minChi2Prob > 0. && !passChi2Prob(iCluster, d2time, clus_chi2,
226  clus_chi2_nhits))
227  d2 = 999.;
228 
229  dist2.emplace_back( d2);
230 
231  if( d2 > 100 ) {
232  LOGDRESSED("PFlow2DClusterizerWithTime:growAndStabilizePFClusters")
233  << "Warning! :: pfcluster-topocell distance is too large! d= "
234  << d2;
235  }
236  // fraction assignment logic
237  if( refhit->detId() == cluster.seed() && _excludeOtherSeeds ) {
238  fraction = 1.0;
239  } else if ( seedable[refhit.key()] && _excludeOtherSeeds ) {
240  fraction = 0.0;
241  } else {
242  fraction = cluster.energy()/recHitEnergyNorm * vdt::fast_expf( -0.5*d2 );
243  }
244  fractot += fraction;
245  frac.emplace_back(fraction);
246  }
247  for( unsigned i = 0; i < clusters.size(); ++i ) {
248  if( fractot > _minFracTot ||
249  ( refhit->detId() == clusters[i].seed() && fractot > 0.0 ) ) {
250  frac[i]/=fractot;
251  } else {
252  continue;
253  }
254  // if the fraction has been set to 0, the cell
255  // is now added to the cluster - careful ! (PJ, 19/07/08)
256  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
257  // (PJ, 15/09/08 <- similar to what existed before the
258  // previous bug fix, but keeps the close seeds inside,
259  // even if their fraction was set to zero.)
260  // Also add a protection to keep the seed in the cluster
261  // when the latter gets far from the former. These cases
262  // (about 1% of the clusters) need to be studied, as
263  // they create fake photons, in general.
264  // (PJ, 16/09/08)
265  if( dist2[i] < 100.0 || frac[i] > 0.9999 ) {
266  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit,frac[i]));
267  }
268  }
269  }
270  // recalculate positions and calculate convergence parameter
271  double diff2 = 0.0;
272  for( unsigned i = 0; i < clusters.size(); ++i ) {
273  if( _convergencePosCalc ) {
274  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
275  } else {
276  if( clusters.size() == 1 && _allCellsPosCalc ) {
277  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
278  } else {
279  _positionCalc->calculateAndSetPosition(clusters[i]);
280  }
281  }
282  const double delta2 =
283  reco::deltaR2(clusters[i].positionREP(),clus_prev_pos[i]);
284  if( delta2 > diff2 ) diff2 = delta2;
285  }
286  diff = std::sqrt(diff2);
287  dist2.clear(); frac.clear(); clus_prev_pos.clear();// avoid badness
288  clus_chi2.clear(); clus_chi2_nhits.clear(); clus_prev_timeres2.clear();
289  growPFClusters(topo,seedable,toleranceScaling,iter+1,diff,clusters);
290 }
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:126
int i
Definition: DBlmapReader.cc:9
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
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:264
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:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
double dist2Time(const reco::PFCluster &, const reco::PFRecHitRef &, int cell_layer, double prev_timeres2) const
double energy() const
cluster energy
Definition: PFCluster.h:82
DetId seed() const
return DetId of seed
Definition: CaloCluster.h:202
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
std::unordered_map< int, double > _recHitEnergyNorms
float fast_expf(float x)
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFCluster.h:54
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
void clusterTimeResolutionFromSeed(reco::PFCluster &cluster, double &res) const
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
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 396 of file PFlow2DClusterizerWithTime.cc.

References _minChi2Prob.

Referenced by growPFClusters().

398 {
399  if (iCluster >= clus_chi2.size()) { // first hit
400  clus_chi2.push_back(dist2);
401  clus_chi2_nhits.push_back(1);
402  return true;
403  }
404  else {
405  double chi2 = clus_chi2[iCluster];
406  size_t nhitsCluster = clus_chi2_nhits[iCluster];
407  chi2 += dist2;
408  if (TMath::Prob(chi2, nhitsCluster) >= _minChi2Prob) {
409  clus_chi2[iCluster] = chi2;
410  clus_chi2_nhits[iCluster] = nhitsCluster + 1;
411  return true;
412  }
413  }
414  return false;
415 }
void PFlow2DClusterizerWithTime::prunePFClusters ( reco::PFClusterCollection clusters) const
private

Definition at line 293 of file PFlow2DClusterizerWithTime.cc.

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

Referenced by buildClusters().

293  {
294  for( auto& cluster : clusters ) {
295  cluster.pruneUsing( [&](const reco::PFRecHitFraction& rhf)
296  {return rhf.fraction() > _minFractionToKeep;}
297  );
298  }
299 }
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 137 of file PFlow2DClusterizerWithTime.cc.

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

Referenced by buildClusters().

139  {
140  const auto& recHitFractions = topo.recHitFractions();
141  for( const auto& rhf : recHitFractions ) {
142  if( !seedable[rhf.recHitRef().key()] ) continue;
143  initialPFClusters.push_back(reco::PFCluster());
144  reco::PFCluster& current = initialPFClusters.back();
145  current.addRecHitFraction(rhf);
146  current.setSeed(rhf.recHitRef()->detId());
147  if( _convergencePosCalc ) {
148  _convergencePosCalc->calculateAndSetPosition(current);
149  } else {
150  _positionCalc->calculateAndSetPosition(current);
151  }
152  }
153 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
std::unique_ptr< PosCalc > _positionCalc
void setSeed(const DetId &id)
Definition: CaloCluster.h:118
void addRecHitFraction(const reco::PFRecHitFraction &frac)
add a given fraction of the rechit
Definition: PFCluster.cc:57
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
void PFlow2DClusterizerWithTime::update ( const edm::EventSetup es)
inlinevirtual

Reimplemented from PFClusterBuilderBase.

Definition at line 20 of file PFlow2DClusterizerWithTime.h.

References _allCellsPosCalc, _convergencePosCalc, and PFClusterBuilderBase::_positionCalc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

20  {
21  _positionCalc->update(es);
22  if( _allCellsPosCalc ) _allCellsPosCalc->update(es);
23  if( _convergencePosCalc ) _convergencePosCalc->update(es);
24  }
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 44 of file PFlow2DClusterizerWithTime.h.

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

const bool PFlow2DClusterizerWithTime::_clusterTimeResFromSeed
private

Definition at line 40 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters().

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

Definition at line 36 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters().

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

Definition at line 42 of file PFlow2DClusterizerWithTime.h.

Referenced by for().

const unsigned PFlow2DClusterizerWithTime::_maxIterations
private

Definition at line 31 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters().

const double PFlow2DClusterizerWithTime::_maxNSigmaTime
private

Definition at line 38 of file PFlow2DClusterizerWithTime.h.

Referenced by dist2Time().

const double PFlow2DClusterizerWithTime::_minChi2Prob
private

Definition at line 39 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters(), and passChi2Prob().

const double PFlow2DClusterizerWithTime::_minFracTot
private

Definition at line 37 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters().

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

Definition at line 43 of file PFlow2DClusterizerWithTime.h.

Referenced by for(), and growPFClusters().

const double PFlow2DClusterizerWithTime::_showerSigma2
private

Definition at line 33 of file PFlow2DClusterizerWithTime.h.

Referenced by growPFClusters().

const double PFlow2DClusterizerWithTime::_stoppingTolerance
private

Definition at line 32 of file PFlow2DClusterizerWithTime.h.

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 34 of file PFlow2DClusterizerWithTime.h.

Referenced by dist2Time().

const double PFlow2DClusterizerWithTime::_timeSigma_ee
private

Definition at line 35 of file PFlow2DClusterizerWithTime.h.

Referenced by dist2Time().