CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFRecHitQTests.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFClusterProducer_PFEcalRecHitQTests_h
2 #define RecoParticleFlow_PFClusterProducer_PFEcalRecHitQTests_h
3 
4 #include <memory>
7 
8 
9 
10 //
11 // Quality test that checks threshold
12 //
14  public:
16 
17  }
18 
20  PFRecHitQTestBase(iConfig)
21  {
22  threshold_ = iConfig.getParameter<double>("threshold");
23 
24  }
25 
26  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
27  }
28 
29  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
30  return pass(hit);
31  }
32  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
33  return pass(hit);
34  }
35 
36  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
37  return pass(hit);
38  }
39  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
40  return pass(hit);
41  }
42 
43  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
44  return pass(hit);
45  }
46 
47  protected:
48  double threshold_;
49 
50  bool pass(const reco::PFRecHit& hit){
51  if (hit.energy()>threshold_) return true;
52 
53  return false;
54  }
55 };
56 
57 
58 
59 
60 //
61 // Quality test that checks kHCAL Severity
62 //
64 
65  public:
66 
67 
69 
70  }
71 
73  PFRecHitQTestBase(iConfig)
74  {
75  thresholds_ = iConfig.getParameter<std::vector<int> >("maxSeverities");
76  cleanThresholds_ = iConfig.getParameter<std::vector<double> >("cleaningThresholds");
77  std::vector<std::string> flags = iConfig.getParameter<std::vector<std::string> >("flags");
78  for (unsigned int i=0;i<flags.size();++i) {
79  if (flags[i] =="Standard") {
80  flags_.push_back(-1);
81  depths_.push_back(-1);
82 
83  }
84  else if (flags[i] =="HFInTime") {
86  depths_.push_back(-1);
87  }
88  else if (flags[i] =="HFDigi") {
90  depths_.push_back(-1);
91 
92  }
93  else if (flags[i] =="HFLong") {
95  depths_.push_back(1);
96 
97  }
98  else if (flags[i] =="HFShort") {
100  depths_.push_back(2);
101 
102  }
103  else {
104  flags_.push_back(-1);
105  depths_.push_back(-1);
106 
107  }
108  }
109 
110  }
111 
112  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
114  iSetup.get<HcalRecNumberingRecord>().get(topo);
115  edm::ESHandle<HcalChannelQuality> hcalChStatus;
116  iSetup.get<HcalChannelQualityRcd>().get( "withTopo", hcalChStatus );
117  theHcalChStatus_ = hcalChStatus.product();
118  edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
119  iSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
120  hcalSevLvlComputer_ = hcalSevLvlComputerHndl.product();
121  }
122 
123  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
124  return true;
125  }
126  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
127  return test(rh.detid(),rh.energy(),rh.flags(),clean);
128  }
129 
130  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
131  return test(rh.detid(),rh.energy(),rh.flags(),clean);
132  }
133  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
134  return test(rh.detid(),rh.energy(),rh.flags(),clean);
135  }
136 
137  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
138  return true;
139 
140  }
141 
142  protected:
143  std::vector<int> thresholds_;
144  std::vector<double> cleanThresholds_;
145  std::vector<int> flags_;
146  std::vector<int> depths_;
149 
150  bool test(unsigned aDETID,double energy,int flags,bool& clean){
151  const HcalDetId& detid = (HcalDetId)aDETID;
152  const HcalChannelStatus* theStatus = theHcalChStatus_->getValues(detid);
153  unsigned theStatusValue = theStatus->getValue();
154  // Now get severity of problems for the given detID, based on the rechit flag word and the channel quality status value
155  for (unsigned int i=0;i<thresholds_.size();++i) {
156  int hitSeverity =0;
157  if (energy < cleanThresholds_[i])
158  continue;
159 
160  if(flags_[i]<0) {
161  hitSeverity=hcalSevLvlComputer_->getSeverityLevel(detid, flags,theStatusValue);
162  }
163  else {
164  hitSeverity=hcalSevLvlComputer_->getSeverityLevel(detid, flags & flags_[i],theStatusValue);
165  }
166 
167  if (hitSeverity>thresholds_[i] && ((depths_[i]<0 || (depths_[i]==detid.depth())))) {
168  clean=true;
169  return false;
170  }
171  }
172  return true;
173  }
174 
175 };
176 
177 //
178 // Quality test that applies threshold and timing as a function of depth
179 //
181  public:
183 
184  }
185 
187  PFRecHitQTestBase(iConfig)
188  {
189  std::vector<edm::ParameterSet> psets = iConfig.getParameter<std::vector<edm::ParameterSet> >("cuts");
190  for (unsigned int i=0;i<psets.size();++i) {
191  depths_.push_back(psets[i].getParameter<int>("depth"));
192  minTimes_.push_back(psets[i].getParameter<double>("minTime"));
193  maxTimes_.push_back(psets[i].getParameter<double>("maxTime"));
194  thresholds_.push_back(psets[i].getParameter<double>("threshold"));
195  }
196  }
197 
198  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
199  }
200 
201  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
202  return true;
203  }
204  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
205  return test(rh.detid(),rh.energy(),rh.time(),clean);
206  }
207 
208  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
209  return test(rh.detid(),rh.energy(),rh.time(),clean);
210  }
211  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean) {
212  return test(rh.detid(),rh.energy(),rh.time(),clean);
213  }
214 
215  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
216  return true;
217  }
218 
219  protected:
220  std::vector<int> depths_;
221  std::vector<double> minTimes_;
222  std::vector<double> maxTimes_;
223  std::vector<double> thresholds_;
224 
225  bool test(unsigned aDETID,double energy,double time,bool& clean){
226  HcalDetId detid(aDETID);
227  for (unsigned int i=0;i<depths_.size();++i) {
228  if (detid.depth() == depths_[i]) {
229  if ((time <minTimes_[i] || time >maxTimes_[i] ) && energy>thresholds_[i])
230  {
231  clean=true;
232  return false;
233  }
234  break;
235  }
236  }
237  return true;
238  }
239 };
240 
241 
242 
243 
244 //
245 // Quality test that applies threshold as a function of depth
246 //
248  public:
250 
251  }
252 
254  PFRecHitQTestBase(iConfig)
255  {
256  std::vector<edm::ParameterSet> psets = iConfig.getParameter<std::vector<edm::ParameterSet> >("cuts");
257  for (unsigned int i=0;i<psets.size();++i) {
258  depths_.push_back(psets[i].getParameter<int>("depth"));
259  thresholds_.push_back(psets[i].getParameter<double>("threshold"));
260  }
261  }
262 
263  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
264  }
265 
266  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
267  return true;
268  }
269  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
270  return test(rh.detid(),rh.energy(),rh.time(),clean);
271  }
272 
273  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
274  return test(rh.detid(),rh.energy(),rh.time(),clean);
275  }
276  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
277  return test(rh.detid(),rh.energy(),rh.time(),clean);
278  }
279 
280  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
281  return true;
282  }
283 
284  protected:
285  std::vector<int> depths_;
286  std::vector<double> thresholds_;
287 
288  bool test(unsigned aDETID,double energy,double time,bool& clean){
289  HcalDetId detid(aDETID);
290  for (unsigned int i=0;i<depths_.size();++i) {
291  if (detid.depth() == depths_[i]) {
292  if ( energy<thresholds_[i])
293  {
294  clean=false;
295  return false;
296  }
297  break;
298  }
299  }
300  return true;
301  }
302 };
303 
304 
305 
306 
307 
308 //
309 // Quality test that checks HO threshold applying different threshold in rings
310 //
312  public:
314 
315  }
316 
318  PFRecHitQTestBase(iConfig)
319  {
320  threshold0_ = iConfig.getParameter<double>("threshold_ring0");
321  threshold12_ = iConfig.getParameter<double>("threshold_ring12");
322  }
323 
324  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
325  }
326 
327  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
328  return true;
329  }
330  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
331  return true;
332  }
333 
334  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
335  return true;
336  }
337  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
338  HcalDetId detid(rh.detid());
339  if (abs(detid.ieta())<=4 && hit.energy()>threshold0_)
340  return true;
341  if (abs(detid.ieta())>4 && hit.energy()>threshold12_)
342  return true;
343 
344  return false;
345  }
346 
347  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
348  return true;
349  }
350 
351  protected:
352  double threshold0_;
353  double threshold12_;
354 
355 };
356 
357 //
358 // Quality test that checks ecal quality cuts
359 //
361  public:
363 
364  }
365 
367  PFRecHitQTestBase(iConfig)
368  {
369  thresholdCleaning_ = iConfig.getParameter<double>("cleaningThreshold");
370  timingCleaning_ = iConfig.getParameter<bool>("timingCleaning");
371  topologicalCleaning_ = iConfig.getParameter<bool>("topologicalCleaning");
372  skipTTRecoveredHits_ = iConfig.getParameter<bool>("skipTTRecoveredHits");
373 
374  }
375 
376  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
377  }
378 
379  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
381  {
382  clean=true;
383  return false;
384  }
385  if ( timingCleaning_ && rh.energy() > thresholdCleaning_ &&
387  clean=true;
388  return false;
389  }
390 
391  if ( topologicalCleaning_ &&
392  ( rh.checkFlag(EcalRecHit::kWeird) ||
394  clean=true;
395  return false;
396  }
397 
398  return true;
399  }
400 
401  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
402  return true;
403  }
404 
405  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
406  return true;
407 
408  }
409 
410  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
411  return true;
412  }
413 
414  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
415  return true;
416 
417  }
418 
419 
420  protected:
425 
426 };
427 
428 
429 
430 
431 
432 //
433 // Quality test that calibrates tower 29 of HCAL
434 //
436  public:
438 
439  }
440 
442  PFRecHitQTestBase(iConfig)
443  {
444  calibFactor_ =iConfig.getParameter<double>("calibFactor");
445  }
446 
447  void beginEvent(const edm::Event& event,const edm::EventSetup& iSetup) {
448  }
449 
450  bool test(reco::PFRecHit& hit,const EcalRecHit& rh,bool& clean){
451  return true;
452  }
453  bool test(reco::PFRecHit& hit,const HBHERecHit& rh,bool& clean){
454  HcalDetId detId(hit.detId());
455  if (abs(detId.ieta())==29)
456  hit.setEnergy(hit.energy()*calibFactor_);
457  return true;
458 
459  }
460 
461  bool test(reco::PFRecHit& hit,const HFRecHit& rh,bool& clean){
462  return true;
463 
464  }
465  bool test(reco::PFRecHit& hit,const HORecHit& rh,bool& clean){
466  return true;
467  }
468 
469  bool test(reco::PFRecHit& hit,const CaloTower& rh,bool& clean){
470  CaloTowerDetId detId(hit.detId());
471  if (detId.ietaAbs()==29)
472  hit.setEnergy(hit.energy()*calibFactor_);
473  return true;
474 
475  }
476 
477  protected:
479 };
480 
481 
482 
483 #endif
T getParameter(std::string const &) const
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
int i
Definition: DBlmapReader.cc:9
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
const HcalChannelQuality * theHcalChStatus_
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
std::vector< double > minTimes_
bool pass(const reco::PFRecHit &hit)
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
const DetId & detid() const
Definition: CaloRecHit.h:20
unsigned detId() const
rechit detId
Definition: PFRecHit.h:106
PFRecHitQTestThreshold(const edm::ParameterSet &iConfig)
std::vector< int > depths_
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
const HcalSeverityLevelComputer * hcalSevLvlComputer_
std::vector< double > maxTimes_
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
PFRecHitQTestECAL(const edm::ParameterSet &iConfig)
const Item * getValues(DetId fId, bool throwOnFail=true) const
bool test(unsigned aDETID, double energy, double time, bool &clean)
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
float time() const
Definition: CaloRecHit.h:19
PFRecHitQTestHCALChannel(const edm::ParameterSet &iConfig)
std::vector< double > cleanThresholds_
std::vector< double > thresholds_
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
PFRecHitQTestHCALThresholdVsDepth(const edm::ParameterSet &iConfig)
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
std::vector< double > thresholds_
bool test(unsigned aDETID, double energy, double time, bool &clean)
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
void setEnergy(double energy)
Definition: PFRecHit.h:72
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
int depth() const
get the tower depth
Definition: HcalDetId.h:55
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
bool checkFlag(int flag) const
check if the flag is true
Definition: EcalRecHit.h:172
float energy() const
Definition: CaloRecHit.h:17
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
PFRecHitQTestHCALCalib29(const edm::ParameterSet &iConfig)
uint32_t flags() const
Definition: CaloRecHit.h:21
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
float energy() const
Definition: EcalRecHit.h:68
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
PFRecHitQTestHCALTimeVsDepth(const edm::ParameterSet &iConfig)
std::vector< T * > clean
Definition: MVATrainer.cc:156
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
PFRecHitQTestHOThreshold(const edm::ParameterSet &iConfig)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const HFRecHit &rh, bool &clean)
const T & get() const
Definition: EventSetup.h:56
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
T const * product() const
Definition: ESHandle.h:86
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
std::vector< int > depths_
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
double energy() const
rechit energy
Definition: PFRecHit.h:112
std::vector< int > thresholds_
bool test(reco::PFRecHit &hit, const CaloTower &rh, bool &clean)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const EcalRecHit &rh, bool &clean)
bool test(unsigned aDETID, double energy, int flags, bool &clean)
uint32_t getValue() const
bool test(reco::PFRecHit &hit, const HORecHit &rh, bool &clean)
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
void beginEvent(const edm::Event &event, const edm::EventSetup &iSetup)
std::vector< int > flags_
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)
bool test(reco::PFRecHit &hit, const HBHERecHit &rh, bool &clean)