CMS 3D CMS Logo

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
 
virtual ~RecHitTopologicalCleanerBase ()=default
 

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 62 of file SpikeAndDoubleSpikeCleaner.cc.

References SpikeAndDoubleSpikeCleaner::spike_cleaning::_doubleSpikeS6S2, SpikeAndDoubleSpikeCleaner::spike_cleaning::_doubleSpikeThresh, SpikeAndDoubleSpikeCleaner::spike_cleaning::_eneThreshMod, SpikeAndDoubleSpikeCleaner::spike_cleaning::_fracThreshMod, _layerMap, SpikeAndDoubleSpikeCleaner::spike_cleaning::_minS4S1_a, SpikeAndDoubleSpikeCleaner::spike_cleaning::_minS4S1_b, SpikeAndDoubleSpikeCleaner::spike_cleaning::_singleSpikeThresh, _thresholds, clean(), PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, mps_splice::entry, Exception, PFLayer::HCAL_BARREL1, PFLayer::HCAL_BARREL2, PFLayer::HCAL_ENDCAP, PFLayer::HF_EM, PFLayer::HF_HAD, info(), createfilelist::int, PFLayer::NONE, PFLayer::PS1, PFLayer::PS2, muonDTDigis_cfi::pset, AlCaHLTBitMon_QueryRunRegistry::string, and particleFlowRecHitECAL_cfi::thresholds.

62  :
64  _layerMap({ {"PS2",(int)PFLayer::PS2},
65  {"PS1",(int)PFLayer::PS1},
66  {"ECAL_ENDCAP",(int)PFLayer::ECAL_ENDCAP},
67  {"ECAL_BARREL",(int)PFLayer::ECAL_BARREL},
68  {"NONE",(int)PFLayer::NONE},
69  {"HCAL_BARREL1",(int)PFLayer::HCAL_BARREL1},
70  {"HCAL_BARREL2_RING0",(int)PFLayer::HCAL_BARREL2},
71  // hack to deal with ring1 in HO
72  {"HCAL_BARREL2_RING1",100*(int)PFLayer::HCAL_BARREL2},
73  {"HCAL_ENDCAP",(int)PFLayer::HCAL_ENDCAP},
74  {"HF_EM",(int)PFLayer::HF_EM},
75  {"HF_HAD",(int)PFLayer::HF_HAD} }) {
76  const std::vector<edm::ParameterSet>& thresholds =
77  conf.getParameterSetVector("cleaningByDetector");
78  for( const auto& pset : thresholds ) {
79  spike_cleaning info;
80  const std::string& det = pset.getParameter<std::string>("detector");
81  info._minS4S1_a = pset.getParameter<double>("minS4S1_a");
82  info._minS4S1_b = pset.getParameter<double>("minS4S1_b");
83  info._doubleSpikeS6S2 = pset.getParameter<double>("doubleSpikeS6S2");
84  info._eneThreshMod = pset.getParameter<double>("energyThresholdModifier");
85  info._fracThreshMod =
86  pset.getParameter<double>("fractionThresholdModifier");
87  info._doubleSpikeThresh = pset.getParameter<double>("doubleSpikeThresh");
88  info._singleSpikeThresh = pset.getParameter<double>("singleSpikeThresh");
89  auto entry = _layerMap.find(det);
90  if( entry == _layerMap.end() ) {
91  throw cms::Exception("InvalidDetectorLayer")
92  << "Detector layer : " << det << " is not in the list of recognized"
93  << " detector layers!";
94  }
95  _thresholds.emplace(_layerMap.find(det)->second,info);
96  }
97 }
VParameterSet const & getParameterSetVector(std::string const &name) const
static const TGPicture * info(bool iBackgroundIsBlack)
RecHitTopologicalCleanerBase(const edm::ParameterSet &conf)
std::unordered_map< int, spike_cleaning > _thresholds
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 101 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(), RhoEtaPhi::eta(), PFLayer::HCAL_BARREL2, hfClusterShapes_cfi::hits, mps_fire::i, training_settings::idx, input, createfilelist::int, gen::k, reco::PFRecHit::layer(), min(), reco::PFRecHit::neighbours4(), RhoEtaPhi::phi(), and reco::PFRecHit::positionREP().

Referenced by SpikeAndDoubleSpikeCleaner().

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

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

Definition at line 30 of file SpikeAndDoubleSpikeCleaner.h.

Referenced by clean(), and SpikeAndDoubleSpikeCleaner().