CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SpikeAndDoubleSpikeCleaner.cc
Go to the documentation of this file.
2 #include <cmath>
3 
4 namespace {
5  bool greaterByEnergy(const std::pair<unsigned,double>& a,
6  const std::pair<unsigned,double>& b) {
7  return a.second > b.second;
8  }
9  std::pair<double,double> dCrack(double phi, double eta) {
10  constexpr double oneOverCrystalSize=1.0/0.0175;
11  constexpr double pi=M_PI;
12  constexpr double twopi=2*pi;
13  // the result below is from unrolling
14  // the lazy-eval loop in PFClusterAlgo::dCrack
15  constexpr double cPhi[18] = {2.97025, 2.621184149601134, 2.272118299202268,
16  1.9230524488034024, 1.5739865984045365,
17  1.2249207480056705, 0.8758548976068048,
18  0.5267890472079388, 0.1777231968090729,
19  -0.17134265358979306, -0.520408503988659,
20  -0.8694743543875245, -1.2185402047863905,
21  -1.5676060551852569, -1.9166719055841224,
22  -2.265737755982988, -2.6148036063818543,
23  -2.9638694567807207};
24  constexpr double cEta[9] = {0.0,
25  4.44747e-01, -4.44747e-01,
26  7.92824e-01, -7.92824e-01,
27  1.14090e+00, -1.14090e+00,
28  1.47464e+00, -1.47464e+00};
29  // shift for eta < 0
30  constexpr double delta_cPhi = 0.00638;
31  // let's calculate dphi
32  double defi = 0;
33  if( eta < 0 ) phi+=delta_cPhi;
34  if( phi >= -pi && phi <= pi ) {
35  //the problem of the extrema
36  if( phi < cPhi[17] || phi >= cPhi[0] ) {
37  if( phi < 0 ) phi += 2*pi;
38  defi = std::min(std::abs(phi-cPhi[0]),std::abs(phi-cPhi[17]-twopi));
39  } else { // between these extrema
40  bool OK = false;
41  unsigned i = 16;
42  while(!OK) {
43  if( phi < cPhi[i] ) {
44  defi = std::min(std::abs(phi-cPhi[i+1]),std::abs(phi-cPhi[i]));
45  OK=true;
46  } else {
47  i -= 1;
48  }
49  }// end while
50  }
51  } else { // if there's a problem assume we're in a crack
52  defi = 0;
53  }
54  // let's calculate deta
55  double deta = 999.0;
56  for( const double etaGap : cEta ) {
57  deta = std::min(deta,std::abs(eta-etaGap));
58  }
59  defi *= oneOverCrystalSize;
60  deta *= oneOverCrystalSize;
61  return std::make_pair(defi,deta);
62  }
63 }
64 
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} }) {
80  const std::vector<edm::ParameterSet>& thresholds =
81  conf.getParameterSetVector("cleaningByDetector");
82  for( const auto& pset : thresholds ) {
83  spike_cleaning info;
84  const std::string& det = pset.getParameter<std::string>("detector");
85  info._minS4S1_a = pset.getParameter<double>("minS4S1_a");
86  info._minS4S1_b = pset.getParameter<double>("minS4S1_b");
87  info._doubleSpikeS6S2 = pset.getParameter<double>("doubleSpikeS6S2");
88  info._eneThreshMod = pset.getParameter<double>("energyThresholdModifier");
89  info._fracThreshMod =
90  pset.getParameter<double>("fractionThresholdModifier");
91  info._doubleSpikeThresh = pset.getParameter<double>("doubleSpikeThresh");
92  info._singleSpikeThresh = pset.getParameter<double>("singleSpikeThresh");
93  auto entry = _layerMap.find(det);
94  if( entry == _layerMap.end() ) {
95  throw cms::Exception("InvalidDetectorLayer")
96  << "Detector layer : " << det << " is not in the list of recognized"
97  << " detector layers!";
98  }
99  _thresholds.emplace(_layerMap.find(det)->second,info);
100  }
101 }
102 
103 
106  std::vector<bool>& mask ) {
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
VParameterSet const & getParameterSetVector(std::string const &name) const
static const TGPicture * info(bool iBackgroundIsBlack)
SpikeAndDoubleSpikeCleaner(const edm::ParameterSet &conf)
void clean(const edm::Handle< reco::PFRecHitCollection > &input, std::vector< bool > &mask)
const PFRecHitRefVector & neighbours4() const
Definition: PFRecHit.h:84
key_type key() const
Accessor for product key.
Definition: Ref.h:264
#define constexpr
static std::string const input
Definition: EdmProvDump.cc:43
const Double_t pi
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
tuple conf
Definition: dbtoconf.py:185
std::pair< int, edm::FunctionWithDict > OK
Definition: findMethod.cc:136
std::unordered_map< int, spike_cleaning > _thresholds
double b
Definition: hdecay.h:120
Geom::Phi< T > phi() const
double energy() const
rechit energy
Definition: PFRecHit.h:112
double a
Definition: hdecay.h:121
const REPPoint & positionREP() const
Definition: PFRecHit.h:133
const std::unordered_map< std::string, int > _layerMap