CMS 3D CMS Logo

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

#include <SpikeAndDoubleSpikeCleaner.h>

Inheritance diagram for SpikeAndDoubleSpikeCleaner:
RecHitTopologicalCleanerBase

Classes

struct  spike_cleaning
 

Public Member Functions

void clean (const edm::Handle< reco::PFRecHitCollection > &input, std::vector< bool > &mask)
 
SpikeAndDoubleSpikeCleaneroperator= (const SpikeAndDoubleSpikeCleaner &)=delete
 
 SpikeAndDoubleSpikeCleaner (const edm::ParameterSet &conf)
 
 SpikeAndDoubleSpikeCleaner (const SpikeAndDoubleSpikeCleaner &)=delete
 
- Public Member Functions inherited from RecHitTopologicalCleanerBase
const std::string & name () const
 
RecHitTopologicalCleanerBaseoperator= (const RecHitTopologicalCleanerBase &)=delete
 
 RecHitTopologicalCleanerBase (const edm::ParameterSet &conf)
 
 RecHitTopologicalCleanerBase (const RecHitTopologicalCleanerBase &)=delete
 

Private Attributes

const std::unordered_map
< std::string, int > 
_layerMap
 
std::unordered_map< int,
spike_cleaning
_thresholds
 

Detailed Description

Definition at line 8 of file SpikeAndDoubleSpikeCleaner.h.

Constructor & Destructor Documentation

SpikeAndDoubleSpikeCleaner::SpikeAndDoubleSpikeCleaner ( const edm::ParameterSet conf)

Definition at line 66 of file SpikeAndDoubleSpikeCleaner.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 fff_deletion::thresholds.

66  :
68  _layerMap({ {"PS2",(int)PFLayer::PS2},
69  {"PS1",(int)PFLayer::PS1},
70  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
71  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
72  {"NONE",(int)PFLayer::NONE},
73  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
74  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
75  // hack to deal with ring1 in HO
76  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
77  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
78  {"HF_EM",(int)PFLayer::HF_EM},
79  {"HF_HAD",(int)PFLayer::HF_HAD} }) {
RecHitTopologicalCleanerBase(const edm::ParameterSet &conf)
const std::unordered_map< std::string, int > _layerMap
SpikeAndDoubleSpikeCleaner::SpikeAndDoubleSpikeCleaner ( const SpikeAndDoubleSpikeCleaner )
delete

Member Function Documentation

void SpikeAndDoubleSpikeCleaner::clean ( const edm::Handle< reco::PFRecHitCollection > &  input,
std::vector< bool > &  mask 
)
virtual

Implements RecHitTopologicalCleanerBase.

Definition at line 105 of file SpikeAndDoubleSpikeCleaner.cc.

References SpikeAndDoubleSpikeCleaner::spike_cleaning::_doubleSpikeS6S2, SpikeAndDoubleSpikeCleaner::spike_cleaning::_doubleSpikeThresh, SpikeAndDoubleSpikeCleaner::spike_cleaning::_eneThreshMod, SpikeAndDoubleSpikeCleaner::spike_cleaning::_fracThreshMod, SpikeAndDoubleSpikeCleaner::spike_cleaning::_minS4S1_a, SpikeAndDoubleSpikeCleaner::spike_cleaning::_minS4S1_b, SpikeAndDoubleSpikeCleaner::spike_cleaning::_singleSpikeThresh, _thresholds, funct::abs(), PFLayer::ECAL_BARREL, reco::PFRecHit::energy(), PFLayer::HCAL_BARREL2, i, edm::Ref< C, T, F >::key(), reco::PFRecHit::layer(), bookConverter::min, reco::PFRecHit::neighbours4(), and reco::PFRecHit::positionREP().

106  {
107  //need to run over energy sorted rechits
108  std::vector<std::pair<unsigned,double> > ordered_hits;
109  for( unsigned i = 0; i < input->size(); ++i ) {
110  std::pair<unsigned,double> val = std::make_pair(i,input->at(i).energy());
111  auto pos = std::upper_bound(ordered_hits.begin(),ordered_hits.end(),
112  val, greaterByEnergy);
113  ordered_hits.insert(pos,val);
114  }
115 
116  for( const auto& idx_e : ordered_hits ) {
117  const unsigned i = idx_e.first;
118  if( !mask[i] ) continue; // don't need to re-mask things :-)
119  const reco::PFRecHit& rechit = input->at(i);
120  int hitlayer = (int)rechit.layer();
121  if( hitlayer == PFLayer::HCAL_BARREL2 &&
122  std::abs(rechit.positionREP().eta()) > 0.34 ) {
123  hitlayer *= 100;
124  }
125  const spike_cleaning& clean = _thresholds.find(hitlayer)->second;
126  if( rechit.energy() < clean._singleSpikeThresh ) continue;
127  const double rhenergy = rechit.energy();
128  // single spike cleaning
129  const reco::PFRecHitRefVector& neighbours4 = rechit.neighbours4();
130  double surroundingEnergy = rechit.energy();
131  double neighbourEnergy = 0.0;
132  double layerEnergy = 0.0;
133  for( const reco::PFRecHitRef& neighbour : neighbours4 ) {
134  if( !mask[neighbour.key()] ) continue;
135  const double sum = neighbour->energy(); //energyUp is just rechit energy?
136  surroundingEnergy += sum;
137  neighbourEnergy += sum;
138  layerEnergy += neighbour->energy();
139  }
140  // wannaBeSeed.energyUp()/wannaBeSeed.energy() : 1.;
141  // Fraction 1 is the balance between the hit and its neighbours
142  // from both layers
143  const double fraction1 = surroundingEnergy/rhenergy;
144  // removed spurious comments from old pfcluster algo...
145  // look there if you want more history
146  const double f1Cut = ( clean._minS4S1_a*std::log10(rechit.energy()) +
147  clean._minS4S1_b );
148  if( fraction1 < f1Cut ) {
149  const double eta = rechit.positionREP().eta();
150  const double aeta = std::abs(eta);
151  const double phi = rechit.positionREP().phi();
152  std::pair<double,double> dcr = dCrack(phi,eta);
153  const double dcrmin = ( rechit.layer() == PFLayer::ECAL_BARREL ?
154  std::min(dcr.first,dcr.second):
155  dcr.second );
156  if( aeta < 5.0 &&
157  ( (aeta < 2.85 && dcrmin > 1.0) ||
158  (rhenergy > clean._eneThreshMod*clean._singleSpikeThresh &&
159  fraction1 < f1Cut/clean._fracThreshMod ) ) ) {
160  mask[i] = false;
161  }
162  }//if initial fraction cut (single spike)
163  // double spike removal
164  if( mask[i] && rhenergy > clean._doubleSpikeThresh ) {
165  //Determine energy surrounding the seed and the most energetic neighbour
166  double surroundingEnergyi = 0.0;
167  double enmax = -999.0;
168  reco::PFRecHitRef mostEnergeticNeighbour;
169  const reco::PFRecHitRefVector& neighbours4i = rechit.neighbours4();
170  for( reco::PFRecHitRef neighbour : neighbours4i ) {
171  if( !mask[neighbour.key()] ) continue;
172  const double nenergy = neighbour->energy();
173  surroundingEnergyi += nenergy;
174  if( nenergy > enmax ) {
175  enmax = nenergy;
176  mostEnergeticNeighbour = neighbour;
177  }
178  }
179  // is there an energetic neighbour
180  if( enmax > 0.0 ) {
181  double surroundingEnergyj = 0.0;
182  const reco::PFRecHitRefVector& neighbours4j =
183  mostEnergeticNeighbour->neighbours4();
184  for( const reco::PFRecHitRef& neighbour : neighbours4j ) {
185  //if( !mask[neighbour] && neighbour != i) continue; // leave out?
186  surroundingEnergyj += neighbour->energy();
187  }
188  // The energy surrounding the double spike candidate
189  const double surroundingEnergyFraction =
190  (surroundingEnergyi+surroundingEnergyj) /
191  (rechit.energy()+mostEnergeticNeighbour->energy()) - 1.;
192  if ( surroundingEnergyFraction < clean._doubleSpikeS6S2 ) {
193  const double eta = rechit.positionREP().eta();
194  const double aeta = std::abs(eta);
195  const double phi = rechit.positionREP().phi();
196  std::pair<double,double> dcr = dCrack(phi,eta);
197  const double dcrmin = ( rechit.layer() == PFLayer::ECAL_BARREL ?
198  std::min(dcr.first,dcr.second):
199  dcr.second );
200  if( aeta < 5.0 &&
201  ( (aeta < 2.85 && dcrmin > 1.0) ||
202  (rhenergy > clean._eneThreshMod*clean._doubleSpikeThresh &&
203  surroundingEnergyFraction < clean._doubleSpikeS6S2/clean._fracThreshMod )
204  ) ) {
205  mask[i] = false;
206  mask[mostEnergeticNeighbour.key()] = false;
207  }
208  }
209  } // was there an energetic neighbour ?
210  }// if double spike thresh
211  } // rechit loop
212 }
int i
Definition: DBlmapReader.cc:9
void clean(const edm::Handle< reco::PFRecHitCollection > &input, std::vector< bool > &mask)
const PFRecHitRefVector & neighbours4() const
Definition: PFRecHit.h:79
T eta() const
PFLayer::Layer layer() const
rechit layer
Definition: PFRecHit.h:104
U second(std::pair< T, U > const &p)
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unordered_map< int, spike_cleaning > _thresholds
key_type key() const
Accessor for product key.
Definition: Ref.h:266
double energy() const
rechit energy
Definition: PFRecHit.h:107
const REPPoint & positionREP() const
Definition: PFRecHit.h:125
Definition: DDAxes.h:10
SpikeAndDoubleSpikeCleaner& SpikeAndDoubleSpikeCleaner::operator= ( const SpikeAndDoubleSpikeCleaner )
delete

Member Data Documentation

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

Definition at line 29 of file SpikeAndDoubleSpikeCleaner.h.

Referenced by for().

std::unordered_map<int,spike_cleaning> SpikeAndDoubleSpikeCleaner::_thresholds
private

Definition at line 30 of file SpikeAndDoubleSpikeCleaner.h.

Referenced by clean(), and for().