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 Member Functions | Private Attributes
EcalCleaningAlgo Class Reference

#include <EcalCleaningAlgo.h>

Public Member Functions

EcalRecHit::Flags checkTopology (const DetId &id, const EcalRecHitCollection &rhs)
 
 EcalCleaningAlgo (const edm::ParameterSet &p)
 
void setFlags (EcalRecHitCollection &rhs)
 

Private Member Functions

float e4e1 (const DetId &id, const EcalRecHitCollection &rhs)
 yet another function to calculate swiss cross More...
 
float e6e2 (const DetId &id, const EcalRecHitCollection &rhs)
 
bool isNearCrack (const DetId &detid)
 in EB, check if we are near a crack More...
 
const std::vector< DetIdneighbours (const DetId &id)
 return the id of the 4 neighbours in the swiss cross More...
 
float recHitE (const DetId id, const EcalRecHitCollection &recHits, bool useTimingInfo)
 

Private Attributes

float cThreshold_barrel_
 
float cThreshold_double_
 
float cThreshold_endcap_
 
float e4e1_a_barrel_
 
float e4e1_a_endcap_
 
float e4e1_b_barrel_
 
float e4e1_b_endcap_
 
float e4e1Treshold_barrel_
 
float e4e1Treshold_endcap_
 
float e6e2thresh_
 
float ignoreOutOfTimeThresh_
 ignore kOutOfTime above threshold when calculating e4e1 More...
 
float tightenCrack_e1_double_
 
float tightenCrack_e1_single_
 
float tightenCrack_e4e1_single_
 
float tightenCrack_e6e2_double_
 

Detailed Description

Definition at line 21 of file EcalCleaningAlgo.h.

Constructor & Destructor Documentation

EcalCleaningAlgo::EcalCleaningAlgo ( const edm::ParameterSet p)

Definition at line 14 of file EcalCleaningAlgo.cc.

References cThreshold_barrel_, cThreshold_double_, cThreshold_endcap_, e4e1_a_barrel_, e4e1_a_endcap_, e4e1_b_barrel_, e4e1_b_endcap_, e4e1Treshold_barrel_, e4e1Treshold_endcap_, e6e2thresh_, edm::ParameterSet::getParameter(), ignoreOutOfTimeThresh_, tightenCrack_e1_double_, tightenCrack_e1_single_, tightenCrack_e4e1_single_, and tightenCrack_e6e2_double_.

14  {
15 
16  cThreshold_barrel_ = p.getParameter<double>("cThreshold_barrel");;
17  cThreshold_endcap_ = p.getParameter<double>("cThreshold_endcap");;
18  e4e1_a_barrel_ = p.getParameter<double>("e4e1_a_barrel");
19  e4e1_b_barrel_ = p.getParameter<double>("e4e1_b_barrel");
20  e4e1_a_endcap_ = p.getParameter<double>("e4e1_a_endcap");
21  e4e1_b_endcap_ = p.getParameter<double>("e4e1_b_endcap");
22  e4e1Treshold_barrel_ = p.getParameter<double>("e4e1Threshold_barrel");
23  e4e1Treshold_endcap_ = p.getParameter<double>("e4e1Threshold_endcap");
24 
25  cThreshold_double_ = p.getParameter<double>("cThreshold_double");
26 
27  ignoreOutOfTimeThresh_ =p.getParameter<double>("ignoreOutOfTimeThresh");
28  tightenCrack_e1_single_ =p.getParameter<double>("tightenCrack_e1_single");
29  tightenCrack_e4e1_single_=p.getParameter<double>("tightenCrack_e4e1_single");
30  tightenCrack_e1_double_ =p.getParameter<double>("tightenCrack_e1_double");
31  tightenCrack_e6e2_double_=p.getParameter<double>("tightenCrack_e6e2_double");
32  e6e2thresh_= p.getParameter<double>("e6e2thresh");
33 
34 
35 
36 }
T getParameter(std::string const &) const
float tightenCrack_e4e1_single_
float ignoreOutOfTimeThresh_
ignore kOutOfTime above threshold when calculating e4e1
float tightenCrack_e6e2_double_

Member Function Documentation

EcalRecHit::Flags EcalCleaningAlgo::checkTopology ( const DetId id,
const EcalRecHitCollection rhs 
)

check topology, return : kGood : not anomalous kWeird : spike kDiWeird : dispike

Flag spikey channels

Mark single spikes. Spike definition:

Barrel: e> cThreshold_barrel_ && e4e1 > e4e1_a_barrel_ * log10(e) + e4e1_b_barrel_

Near cracks: energy threshold is multiplied by tightenCrack_e1_single e4e1 threshold is divided by tightenCrack_e4e1_single

Endcap : e> cThreshold_endcap_ && e4e1> e4e1_a_endcap_ * log10(e) + e4e1_b_endcap_

Mark double spikes (barrel only) e> cThreshold_double_ && e6e2 > e6e2thresh_;

Near cracks: energy threshold multiplied by tightenCrack_e1_double e6e2 threshold divided by tightenCrack_e6e2_double

Out of time hits above e4e1_IgnoreOutOfTimeThresh_ are ignored in topological quantities

Definition at line 69 of file EcalCleaningAlgo.cc.

References a, b, cThreshold_barrel_, cThreshold_double_, cThreshold_endcap_, e4e1(), e4e1_a_barrel_, e4e1_a_endcap_, e4e1_b_barrel_, e4e1_b_endcap_, e6e2(), e6e2thresh_, EcalBarrel, EcalEndcap, relval_parameters_module::energy, isNearCrack(), EcalRecHit::kDiWeird, EcalRecHit::kGood, EcalRecHit::kWeird, recHitE(), tightenCrack_e1_double_, tightenCrack_e1_single_, tightenCrack_e4e1_single_, and tightenCrack_e6e2_double_.

Referenced by setFlags().

70  {
71 
72 
73  float a=0,b=0,e4e1thresh=0,ethresh=0;
74 
75 
76  if( id.subdetId() == EcalBarrel) {
77  a= e4e1_a_barrel_;
78  b= e4e1_b_barrel_;
79  ethresh=cThreshold_barrel_;
80 
81  }
82  else if( id.subdetId() == EcalEndcap){
83  a= e4e1_a_endcap_;
85  ethresh=cThreshold_endcap_;
86  }
87 
88 
89 
90  // for energies below threshold, we don't apply e4e1 cut
91  float energy = recHitE(id,rhs,false);
92 
93  if (energy< ethresh) return EcalRecHit::kGood;
94  if (isNearCrack(id) && energy < ethresh*tightenCrack_e1_single_)
95  return EcalRecHit::kGood;
96 
97 
98  float e4e1value = e4e1(id,rhs);
99  e4e1thresh = a* log10(energy) + b;
100 
101  // near cracks the cut is tighter by a factor
102  if (isNearCrack(id)) {
103  e4e1thresh/=tightenCrack_e4e1_single_;
104  }
105 
106  // identify spike
107  if (e4e1value < e4e1thresh) return EcalRecHit::kWeird;
108 
109 
110 
111  // now for double spikes
112 
113  // no checking for double spikes in EE
114  if( id.subdetId() == EcalEndcap) return EcalRecHit::kGood;
115 
116  float e6e2value = e6e2(id,rhs);
117  float e6e2thresh = e6e2thresh_ ;
119  return EcalRecHit::kGood;
120 
121  if (energy < cThreshold_double_) return EcalRecHit::kGood;
122 
123  // near cracks the cut is tighter by a factor
124  if (id.subdetId() == EcalBarrel && isNearCrack(id))
125  e6e2thresh/=tightenCrack_e6e2_double_;
126 
127  // identify double spike
128  if (e6e2value < e6e2thresh) return EcalRecHit::kDiWeird;
129 
130  return EcalRecHit::kGood;
131 
132 }
float tightenCrack_e4e1_single_
float recHitE(const DetId id, const EcalRecHitCollection &recHits, bool useTimingInfo)
float e4e1(const DetId &id, const EcalRecHitCollection &rhs)
yet another function to calculate swiss cross
bool isNearCrack(const DetId &detid)
in EB, check if we are near a crack
float tightenCrack_e6e2_double_
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
float e6e2(const DetId &id, const EcalRecHitCollection &rhs)
float EcalCleaningAlgo::e4e1 ( const DetId id,
const EcalRecHitCollection rhs 
)
private

yet another function to calculate swiss cross

Definition at line 139 of file EcalCleaningAlgo.cc.

References i, neighbours(), and recHitE().

Referenced by checkTopology(), and e6e2().

140  {
141 
142 
143  float s4 = 0;
144  float e1 = recHitE( id, rhs, false );
145 
146 
147  if ( e1 == 0 ) return 0;
148  const std::vector<DetId>& neighs = neighbours(id);
149  for (size_t i=0; i<neighs.size(); ++i)
150  // avoid hits out of time when making s4
151  s4+=recHitE(neighs[i],rhs, true);
152 
153  return s4 / e1;
154 
155 
156 }
int i
Definition: DBlmapReader.cc:9
const std::vector< DetId > neighbours(const DetId &id)
return the id of the 4 neighbours in the swiss cross
float recHitE(const DetId id, const EcalRecHitCollection &recHits, bool useTimingInfo)
float EcalCleaningAlgo::e6e2 ( const DetId id,
const EcalRecHitCollection rhs 
)
private

Compute e6 over e2 around xtal 1, where 2 is the most energetic in the swiss cross around 1

  | | |   
+-+-+-+-+
| |1|2| |
+-+-+-+-+
  | | |                

Definition at line 161 of file EcalCleaningAlgo.cc.

References e4e1(), i, neighbours(), and recHitE().

Referenced by checkTopology().

162  {
163 
164  float s4_1 = 0;
165  float s4_2 = 0;
166  float e1 = recHitE( id, rhs , false );
167 
168 
169  float maxene=0;
170  DetId maxid;
171 
172  if ( e1 == 0 ) return 0;
173 
174  const std::vector<DetId>& neighs = neighbours(id);
175 
176  // find the most energetic neighbour ignoring time info
177  for (size_t i=0; i<neighs.size(); ++i){
178  float ene = recHitE(neighs[i],rhs,false);
179  if (ene>maxene) {
180  maxene=ene;
181  maxid = neighs[i];
182  }
183  }
184 
185  float e2=maxene;
186 
187  s4_1 = e4e1(id,rhs)* e1;
188  s4_2 = e4e1(maxid,rhs)* e2;
189 
190  return (s4_1 + s4_2) / (e1+e2) -1. ;
191 
192 }
int i
Definition: DBlmapReader.cc:9
const std::vector< DetId > neighbours(const DetId &id)
return the id of the 4 neighbours in the swiss cross
float recHitE(const DetId id, const EcalRecHitCollection &recHits, bool useTimingInfo)
float e4e1(const DetId &id, const EcalRecHitCollection &rhs)
yet another function to calculate swiss cross
Definition: DetId.h:18
bool EcalCleaningAlgo::isNearCrack ( const DetId detid)
private

in EB, check if we are near a crack

Definition at line 255 of file EcalCleaningAlgo.cc.

References EcalEndcap, EBDetId::isNextToBoundary(), and EEDetId::isNextToRingBoundary().

Referenced by checkTopology().

255  {
256 
257  if (id.subdetId() == EcalEndcap) {
259  } else {
260  return EBDetId::isNextToBoundary(id);
261  }
262 }
static bool isNextToBoundary(EBDetId id)
Definition: EBDetId.cc:121
static bool isNextToRingBoundary(EEDetId id)
Definition: EEDetId.cc:375
const std::vector< DetId > EcalCleaningAlgo::neighbours ( const DetId id)
private

return the id of the 4 neighbours in the swiss cross

four neighbours in the swiss cross around id

Definition at line 229 of file EcalCleaningAlgo.cc.

References EcalBarrel, EcalEndcap, EBDetId::offsetBy(), EEDetId::offsetBy(), and run_regression::ret.

Referenced by e4e1(), and e6e2().

229  {
230 
231  std::vector<DetId> ret;
232 
233  if ( id.subdetId() == EcalBarrel) {
234 
235  ret.push_back( EBDetId::offsetBy( id, 1, 0 ));
236  ret.push_back( EBDetId::offsetBy( id, -1, 0 ));
237  ret.push_back( EBDetId::offsetBy( id, 0, 1 ));
238  ret.push_back( EBDetId::offsetBy( id, 0,-1 ));
239  }
240  // nobody understands what polymorphism is for, sgrunt !
241  else if (id.subdetId() == EcalEndcap) {
242  ret.push_back( EEDetId::offsetBy( id, 1, 0 ));
243  ret.push_back( EEDetId::offsetBy( id, -1, 0 ));
244  ret.push_back( EEDetId::offsetBy( id, 0, 1 ));
245  ret.push_back( EEDetId::offsetBy( id, 0,-1 ));
246 
247  }
248 
249 
250  return ret;
251 
252 }
EEDetId offsetBy(int nrStepsX, int nrStepsY) const
Definition: EEDetId.cc:474
EBDetId offsetBy(int nrStepsEta, int nrStepsPhi) const
Definition: EBDetId.cc:61
float EcalCleaningAlgo::recHitE ( const DetId  id,
const EcalRecHitCollection recHits,
bool  useTimingInfo 
)
private

Definition at line 197 of file EcalCleaningAlgo.cc.

References e4e1Treshold_barrel_, e4e1Treshold_endcap_, EcalBarrel, EcalEndcap, edm::SortedCollection< T, SORT >::end(), edm::SortedCollection< T, SORT >::find(), ignoreOutOfTimeThresh_, EcalRecHit::kOutOfTime, and dtDQMClient_cfg::threshold.

Referenced by checkTopology(), e4e1(), and e6e2().

200 {
201  if ( id.rawId() == 0 ) return 0;
202 
203 
205  if ( id.subdetId() == EcalEndcap) threshold = e4e1Treshold_endcap_;
206 
207  EcalRecHitCollection::const_iterator it = recHits.find( id );
208  if ( it != recHits.end() ){
209  float ene= (*it).energy();
210 
211  // ignore out of time in EB when making e4e1 if so configured
212  if (useTimingInfo){
213 
214  if (id.subdetId()==EcalBarrel &&
215  it->checkFlag(EcalRecHit::kOutOfTime)
216  && ene>ignoreOutOfTimeThresh_) return 0;
217  }
218 
219  // ignore hits below threshold
220  if (ene < threshold) return 0;
221 
222  // else return the energy of this hit
223  return ene;
224  }
225  return 0;
226 }
float ignoreOutOfTimeThresh_
ignore kOutOfTime above threshold when calculating e4e1
std::vector< EcalRecHit >::const_iterator const_iterator
const_iterator end() const
iterator find(key_type k)
void EcalCleaningAlgo::setFlags ( EcalRecHitCollection rhs)

Definition at line 266 of file EcalCleaningAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), checkTopology(), edm::SortedCollection< T, SORT >::end(), and EcalRecHit::kGood.

Referenced by Vispa.Plugins.EdmBrowser.EventContentView.LabelItem::__init__(), and EcalRecHitProducer::produce().

266  {
268  //changing the collection on place
269  for (rh=rhs.begin(); rh!=rhs.end(); ++rh){
270  EcalRecHit::Flags state=checkTopology(rh->id(),rhs);
271  if (state!=EcalRecHit::kGood) {
272  rh->unsetFlag(EcalRecHit::kGood);
273  rh->setFlag(state);
274  }
275  }
276 }
EcalRecHit::Flags checkTopology(const DetId &id, const EcalRecHitCollection &rhs)
std::vector< EcalRecHit >::iterator iterator
const_iterator end() const
const_iterator begin() const

Member Data Documentation

float EcalCleaningAlgo::cThreshold_barrel_
private

Definition at line 67 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::cThreshold_double_
private

Definition at line 78 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::cThreshold_endcap_
private

Definition at line 68 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::e4e1_a_barrel_
private

Definition at line 69 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::e4e1_a_endcap_
private

Definition at line 71 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::e4e1_b_barrel_
private

Definition at line 70 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::e4e1_b_endcap_
private

Definition at line 72 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::e4e1Treshold_barrel_
private

Definition at line 74 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

float EcalCleaningAlgo::e4e1Treshold_endcap_
private

Definition at line 75 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

float EcalCleaningAlgo::e6e2thresh_
private

Definition at line 81 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::ignoreOutOfTimeThresh_
private

ignore kOutOfTime above threshold when calculating e4e1

Definition at line 63 of file EcalCleaningAlgo.h.

Referenced by EcalCleaningAlgo(), and recHitE().

float EcalCleaningAlgo::tightenCrack_e1_double_
private

Definition at line 79 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::tightenCrack_e1_single_
private

Definition at line 76 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::tightenCrack_e4e1_single_
private

Definition at line 77 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().

float EcalCleaningAlgo::tightenCrack_e6e2_double_
private

Definition at line 80 of file EcalCleaningAlgo.h.

Referenced by checkTopology(), and EcalCleaningAlgo().