CMS 3D CMS Logo

EcalSeverityLevelAlgo.cc
Go to the documentation of this file.
1 
14 
18 
20 
22 
24 
25 
26  timeThresh_ = p.getParameter< double> ("timeThresh");
27  chStatus_ =nullptr;
28 
29  const edm::ParameterSet & ps=p.getParameter< edm::ParameterSet >("flagMask");
30  std::vector<std::string> severities = ps.getParameterNames();
31  std::vector<std::string> flags;
32 
33  flagMask_.resize(severities.size());
34 
35  // read configuration of severities
36 
37  for (unsigned int is=0;is!=severities.size();++is){
38 
40  (EcalSeverityLevel::SeverityLevel) StringToEnumValue<EcalSeverityLevel::SeverityLevel>(severities[is]);
41  flags=ps.getParameter<std::vector<std::string> >(severities[is]);
42  uint32_t mask=0;
43  for (unsigned int ifi=0;ifi!=flags.size();++ifi){
45  (EcalRecHit::Flags)StringToEnumValue<EcalRecHit::Flags>(flags[ifi]);
46  //manipulate the mask
47  mask|=(0x1<<f);
48  }
49 
50  flagMask_[snum]=mask;
51  }
52  // read configuration of dbstatus
53 
54  const edm::ParameterSet & dbps=
55  p.getParameter< edm::ParameterSet >("dbstatusMask");
56  std::vector<std::string> dbseverities = dbps.getParameterNames();
57  std::vector<std::string> dbflags;
58 
59  dbstatusMask_.resize(dbseverities.size());
60 
61  for (unsigned int is=0;is!=dbseverities.size();++is){
62 
64  (EcalSeverityLevel::SeverityLevel) StringToEnumValue<EcalSeverityLevel::SeverityLevel>(severities[is]);
65 
66  dbflags=dbps.getParameter<std::vector<std::string> >(severities[is]);
67  uint32_t mask=0;
68  for (unsigned int ifi=0;ifi!=dbflags.size();++ifi){
70  (EcalChannelStatusCode::Code)StringToEnumValue<EcalChannelStatusCode::Code>(dbflags[ifi]);
71 
72  //manipulate the mask
73  mask|=(0x1<<f);
74  }
75 
76  dbstatusMask_[snum]=mask;
77  }
78 
79 
80 }
81 
82 
85  const EcalRecHitCollection& rhs) const{
86 
87  using namespace EcalSeverityLevel;
88 
89  // if the detid is within our rechits, evaluate from flag
91  if ( rh != rhs.end() )
92  return severityLevel(*rh);
93 
94 
95  // else evaluate from dbstatus
96  return severityLevel(id);
97 
98 }
99 
102 
103  using namespace EcalSeverityLevel;
104 
105 
107 
108  uint16_t dbStatus = chIt->getStatusCode();
109 
110  // kGood==0 we know!
111  if (0==dbStatus) return kGood;
112 
113  // check if the bit corresponding to that dbStatus is set in the mask
114  // This implementation implies that the statuses have a priority
115  for (size_t i=0; i< dbstatusMask_.size();++i){
116  uint32_t tmp = 0x1<<dbStatus;
117  if (dbstatusMask_[i] & tmp) return SeverityLevel(i);
118  }
119 
120  // no matching
121  LogDebug("EcalSeverityLevelAlgo")<<
122  "Unmatched DB status, returning kGood";
123  return kGood;
124 }
125 
126 
129 
130  using namespace EcalSeverityLevel;
131 
132  //if marked good, do not do any further test
133  if (rh.checkFlag(kGood)) return kGood;
134 
135  // check if the bit corresponding to that flag is set in the mask
136  // This implementation implies that severities have a priority...
137  for (int sev=kBad;sev>=0;--sev){
138  if(sev==kTime && rh.energy() < timeThresh_ ) continue;
139  if (rh.checkFlagMask(flagMask_[sev])) return SeverityLevel(sev);
140  }
141 
142  // no matching
143  LogDebug("EcalSeverityLevelAlgo")<< "Unmatched Flag , returning kGood";
144  return kGood;
145 
146 }
147 
148 
149 
150 
151 // Configure (x)emacs for this file ...
152 // Local Variables:
153 // mode:c++
154 // compile-command: "scram b -k"
155 // End:
#define LogDebug(id)
T getParameter(std::string const &) const
EcalSeverityLevelAlgo(const edm::ParameterSet &p)
EcalSeverityLevel::SeverityLevel severityLevel(const DetId &id) const
Evaluate status from id use channelStatus from DB.
float timeThresh_
Return kTime only if the rechit is flagged kOutOfTime and E>timeThresh_.
std::vector< EcalRecHit >::const_iterator const_iterator
std::vector< Variable::Flags > flags
Definition: MVATrainer.cc:135
bool checkFlag(int flag) const
check if the flag is true
Definition: EcalRecHit.h:189
float energy() const
Definition: EcalRecHit.h:68
double f[11][100]
std::vector< std::string > getParameterNames() const
std::vector< uint32_t > flagMask_
Configure which EcalRecHit::Flag is mapped into which EcalSeverityLevel.
const_iterator end() const
Definition: DetId.h:18
std::vector< Item >::const_iterator const_iterator
const EcalChannelStatus * chStatus_
bool checkFlagMask(uint32_t mask) const
apply a bitmask to our flags. Experts only
Definition: EcalRecHit.h:203
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
iterator find(key_type k)
SeverityLevel
const_iterator find(uint32_t rawId) const
std::vector< uint32_t > dbstatusMask_
Configure which DBStatus::Flag is mapped into which EcalSeverityLevel.