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
Basic2DGenericPFlowClusterizer Class Reference

#include <Basic2DGenericPFlowClusterizer.h>

Inheritance diagram for Basic2DGenericPFlowClusterizer:
PFClusterBuilderBase

Public Member Functions

 Basic2DGenericPFlowClusterizer (const edm::ParameterSet &conf)
 
 Basic2DGenericPFlowClusterizer (const B2DGPF &)=delete
 
void buildClusters (const reco::PFClusterCollection &, const std::vector< bool > &, reco::PFClusterCollection &outclus)
 
B2DGPFoperator= (const B2DGPF &)=delete
 
void update (const edm::EventSetup &es)
 
virtual ~Basic2DGenericPFlowClusterizer ()
 
- 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
Basic2DGenericPFlowClusterizer 
B2DGPF
 

Private Member Functions

void growPFClusters (const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) 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
 
std::unique_ptr
< PFCPositionCalculatorBase
_convergencePosCalc
 
const bool _excludeOtherSeeds
 
const std::unordered_map
< std::string, int > 
_layerMap
 
const unsigned _maxIterations
 
const double _minFracTot
 
std::unordered_map< int, double > _recHitEnergyNorms
 
const double _showerSigma2
 
const double _stoppingTolerance
 

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 9 of file Basic2DGenericPFlowClusterizer.h.

Member Typedef Documentation

Definition at line 10 of file Basic2DGenericPFlowClusterizer.h.

Constructor & Destructor Documentation

Basic2DGenericPFlowClusterizer::Basic2DGenericPFlowClusterizer ( const edm::ParameterSet conf)

Definition at line 25 of file Basic2DGenericPFlowClusterizer.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.

25  :
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} }) {
T getParameter(std::string const &) const
const std::unordered_map< std::string, int > _layerMap
PFClusterBuilderBase(const edm::ParameterSet &conf)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
virtual Basic2DGenericPFlowClusterizer::~Basic2DGenericPFlowClusterizer ( )
inlinevirtual

Definition at line 14 of file Basic2DGenericPFlowClusterizer.h.

14 {}
Basic2DGenericPFlowClusterizer::Basic2DGenericPFlowClusterizer ( const B2DGPF )
delete

Member Function Documentation

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

Implements PFClusterBuilderBase.

Definition at line 81 of file Basic2DGenericPFlowClusterizer.cc.

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

83  {
84  reco::PFClusterCollection clustersInTopo;
85  for( const auto& topocluster : input ) {
86  clustersInTopo.clear();
87  seedPFClustersFromTopo(topocluster,seedable,clustersInTopo);
88  const unsigned tolScal =
89  std::pow(std::max(1.0,clustersInTopo.size()-1.0),2.0);
90  growPFClusters(topocluster,seedable,tolScal,0,tolScal,clustersInTopo);
91  // step added by Josh Bendavid, removes low-fraction clusters
92  // did not impact position resolution with fraction cut of 1e-7
93  // decreases the size of each pf cluster considerably
94  prunePFClusters(clustersInTopo);
95  // recalculate the positions of the pruned clusters
96  if( _convergencePosCalc ) {
97  // if defined, use the special position calculation for convergence tests
98  _convergencePosCalc->calculateAndSetPositions(clustersInTopo);
99  } else {
100  if( clustersInTopo.size() == 1 && _allCellsPosCalc ) {
101  _allCellsPosCalc->calculateAndSetPosition(clustersInTopo.back());
102  } else {
103  _positionCalc->calculateAndSetPositions(clustersInTopo);
104  }
105  }
106  for( auto& clusterout : clustersInTopo ) {
107  output.insert(output.end(),std::move(clusterout));
108  }
109  }
110 }
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
std::unique_ptr< PosCalc > _positionCalc
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
static std::string const input
Definition: EdmProvDump.cc:43
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
def move
Definition: eostools.py:510
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
void prunePFClusters(reco::PFClusterCollection &) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void seedPFClustersFromTopo(const reco::PFCluster &, const std::vector< bool > &, reco::PFClusterCollection &) const
void Basic2DGenericPFlowClusterizer::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 132 of file Basic2DGenericPFlowClusterizer.cc.

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

Referenced by buildClusters().

137  {
138  if( iter >= _maxIterations ) {
139  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
140  <<"reached " << _maxIterations << " iterations, terminated position "
141  << "fit with diff = " << diff;
142  }
143  if( iter >= _maxIterations ||
144  diff <= _stoppingTolerance*toleranceScaling) 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, fraction = 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 &&
166  std::abs(refhit->positionREP().eta()) > 0.34 ) {
167  cell_layer *= 100;
168  }
169  const double recHitEnergyNorm =
170  _recHitEnergyNorms.find(cell_layer)->second;
171  const math::XYZPoint& topocellpos_xyz = refhit->position();
172  dist2.clear(); frac.clear(); fractot = 0;
173  // add rechits to clusters, calculating fraction based on distance
174  for( auto& cluster : clusters ) {
175  const math::XYZPoint& clusterpos_xyz = cluster.position();
176  fraction = 0.0;
177  const math::XYZVector deltav = clusterpos_xyz - topocellpos_xyz;
178  const double d2 = deltav.Mag2()/_showerSigma2;
179  dist2.emplace_back( d2 );
180  if( d2 > 100 ) {
181  LOGDRESSED("Basic2DGenericPFlowClusterizer:growAndStabilizePFClusters")
182  << "Warning! :: pfcluster-topocell distance is too large! d= "
183  << d2;
184  }
185  // fraction assignment logic
186  if( refhit->detId() == cluster.seed() && _excludeOtherSeeds ) {
187  fraction = 1.0;
188  } else if ( seedable[refhit.key()] && _excludeOtherSeeds ) {
189  fraction = 0.0;
190  } else {
191  fraction = cluster.energy()/recHitEnergyNorm * vdt::fast_expf( -0.5*d2 );
192  }
193  fractot += fraction;
194  frac.emplace_back(fraction);
195  }
196  for( unsigned i = 0; i < clusters.size(); ++i ) {
197  if( fractot > _minFracTot ||
198  ( refhit->detId() == clusters[i].seed() && fractot > 0.0 ) ) {
199  frac[i]/=fractot;
200  } else {
201  continue;
202  }
203  // if the fraction has been set to 0, the cell
204  // is now added to the cluster - careful ! (PJ, 19/07/08)
205  // BUT KEEP ONLY CLOSE CELLS OTHERWISE MEMORY JUST EXPLOSES
206  // (PJ, 15/09/08 <- similar to what existed before the
207  // previous bug fix, but keeps the close seeds inside,
208  // even if their fraction was set to zero.)
209  // Also add a protection to keep the seed in the cluster
210  // when the latter gets far from the former. These cases
211  // (about 1% of the clusters) need to be studied, as
212  // they create fake photons, in general.
213  // (PJ, 16/09/08)
214  if( dist2[i] < 100.0 || frac[i] > 0.9999 ) {
215  clusters[i].addRecHitFraction(reco::PFRecHitFraction(refhit,frac[i]));
216  }
217  }
218  }
219  // recalculate positions and calculate convergence parameter
220  double diff2 = 0.0;
221  for( unsigned i = 0; i < clusters.size(); ++i ) {
222  if( _convergencePosCalc ) {
223  _convergencePosCalc->calculateAndSetPosition(clusters[i]);
224  } else {
225  if( clusters.size() == 1 && _allCellsPosCalc ) {
226  _allCellsPosCalc->calculateAndSetPosition(clusters[i]);
227  } else {
228  _positionCalc->calculateAndSetPosition(clusters[i]);
229  }
230  }
231  const double delta2 =
232  reco::deltaR2(clusters[i].positionREP(),clus_prev_pos[i]);
233  if( delta2 > diff2 ) diff2 = delta2;
234  }
235  diff = std::sqrt(diff2);
236  dist2.clear(); frac.clear(); clus_prev_pos.clear();// avoid badness
237  growPFClusters(topo,seedable,toleranceScaling,iter+1,diff,clusters);
238 }
int i
Definition: DBlmapReader.cc:9
void growPFClusters(const reco::PFCluster &, const std::vector< bool > &, const unsigned toleranceScaling, const unsigned iter, double dist, reco::PFClusterCollection &) const
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:264
#define LOGDRESSED(x)
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc
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
std::unordered_map< int, double > _recHitEnergyNorms
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
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
B2DGPF& Basic2DGenericPFlowClusterizer::operator= ( const B2DGPF )
delete
void Basic2DGenericPFlowClusterizer::prunePFClusters ( reco::PFClusterCollection clusters) const
private

Definition at line 241 of file Basic2DGenericPFlowClusterizer.cc.

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

Referenced by buildClusters().

241  {
242  for( auto& cluster : clusters ) {
243  cluster.pruneUsing( [&](const reco::PFRecHitFraction& rhf)
244  {return rhf.fraction() > _minFractionToKeep;}
245  );
246  }
247 }
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
double fraction() const
void Basic2DGenericPFlowClusterizer::seedPFClustersFromTopo ( const reco::PFCluster topo,
const std::vector< bool > &  seedable,
reco::PFClusterCollection initialPFClusters 
) const
private

Definition at line 113 of file Basic2DGenericPFlowClusterizer.cc.

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

Referenced by buildClusters().

115  {
116  const auto& recHitFractions = topo.recHitFractions();
117  for( const auto& rhf : recHitFractions ) {
118  if( !seedable[rhf.recHitRef().key()] ) continue;
119  initialPFClusters.push_back(reco::PFCluster());
120  reco::PFCluster& current = initialPFClusters.back();
121  current.addRecHitFraction(rhf);
122  current.setSeed(rhf.recHitRef()->detId());
123  if( _convergencePosCalc ) {
124  _convergencePosCalc->calculateAndSetPosition(current);
125  } else {
126  _positionCalc->calculateAndSetPosition(current);
127  }
128  }
129 }
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
std::unique_ptr< PosCalc > _positionCalc
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
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
void Basic2DGenericPFlowClusterizer::update ( const edm::EventSetup es)
inlinevirtual

Reimplemented from PFClusterBuilderBase.

Definition at line 18 of file Basic2DGenericPFlowClusterizer.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().

18  {
19  _positionCalc->update(es);
20  if( _allCellsPosCalc ) _allCellsPosCalc->update(es);
21  if( _convergencePosCalc ) _convergencePosCalc->update(es);
22  }
std::unique_ptr< PosCalc > _positionCalc
std::unique_ptr< PFCPositionCalculatorBase > _convergencePosCalc
std::unique_ptr< PFCPositionCalculatorBase > _allCellsPosCalc

Member Data Documentation

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

Definition at line 36 of file Basic2DGenericPFlowClusterizer.h.

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

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

Definition at line 32 of file Basic2DGenericPFlowClusterizer.h.

Referenced by growPFClusters().

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

Definition at line 34 of file Basic2DGenericPFlowClusterizer.h.

Referenced by for().

const unsigned Basic2DGenericPFlowClusterizer::_maxIterations
private

Definition at line 29 of file Basic2DGenericPFlowClusterizer.h.

Referenced by growPFClusters().

const double Basic2DGenericPFlowClusterizer::_minFracTot
private

Definition at line 33 of file Basic2DGenericPFlowClusterizer.h.

Referenced by growPFClusters().

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

Definition at line 35 of file Basic2DGenericPFlowClusterizer.h.

Referenced by for(), and growPFClusters().

const double Basic2DGenericPFlowClusterizer::_showerSigma2
private

Definition at line 31 of file Basic2DGenericPFlowClusterizer.h.

Referenced by growPFClusters().

const double Basic2DGenericPFlowClusterizer::_stoppingTolerance
private

Definition at line 30 of file Basic2DGenericPFlowClusterizer.h.

Referenced by growPFClusters().