CMS 3D CMS Logo

EcalRecHitWorkerSimple.cc
Go to the documentation of this file.
1 
28 
30 public:
32  ~EcalRecHitWorkerSimple() override;
33 
34  void set(const edm::EventSetup& es) override;
35  bool run(const edm::Event& evt, const EcalUncalibratedRecHit& uncalibRH, EcalRecHitCollection& result) override;
36 
37 protected:
38  double EBLaserMIN_;
39  double EELaserMIN_;
40  double EBLaserMAX_;
41  double EELaserMAX_;
42 
48  std::vector<int> v_chstatus_;
56 
57  // Associate reco flagbit ( outer vector) to many db status flags (inner vector)
58  std::vector<std::vector<uint32_t> > v_DB_reco_flags_;
59 
60  uint32_t setFlagBits(const std::vector<std::vector<uint32_t> >& map, const uint32_t& status);
61 
62  uint32_t flagmask_; // do not propagate channels with these flags on
63 
67 
69 };
70 
74  v_chstatus_ = StringToEnumValue<EcalChannelStatusCode::Code>(
75  ps.getParameter<std::vector<std::string> >("ChannelStatusToBeExcluded"));
76  killDeadChannels_ = ps.getParameter<bool>("killDeadChannels");
77  laserCorrection_ = ps.getParameter<bool>("laserCorrection");
78  EBLaserMIN_ = ps.getParameter<double>("EBLaserMIN");
79  EELaserMIN_ = ps.getParameter<double>("EELaserMIN");
80  EBLaserMAX_ = ps.getParameter<double>("EBLaserMAX");
81  EELaserMAX_ = ps.getParameter<double>("EELaserMAX");
82 
83  skipTimeCalib_ = ps.getParameter<bool>("skipTimeCalib");
84 
86  if (!skipTimeCalib_) {
89  }
92  if (laserCorrection_)
94 
95  // Traslate string representation of flagsMapDBReco into enum values
96  const edm::ParameterSet& p = ps.getParameter<edm::ParameterSet>("flagsMapDBReco");
97  std::vector<std::string> recoflagbitsStrings = p.getParameterNames();
98  v_DB_reco_flags_.resize(32);
99 
100  for (unsigned int i = 0; i != recoflagbitsStrings.size(); ++i) {
101  EcalRecHit::Flags recoflagbit = (EcalRecHit::Flags)StringToEnumValue<EcalRecHit::Flags>(recoflagbitsStrings[i]);
102  std::vector<std::string> dbstatus_s = p.getParameter<std::vector<std::string> >(recoflagbitsStrings[i]);
103  std::vector<uint32_t> dbstatuses;
104  for (unsigned int j = 0; j != dbstatus_s.size(); ++j) {
105  EcalChannelStatusCode::Code dbstatus =
106  (EcalChannelStatusCode::Code)StringToEnumValue<EcalChannelStatusCode::Code>(dbstatus_s[j]);
107  dbstatuses.push_back(dbstatus);
108  }
109 
110  v_DB_reco_flags_[recoflagbit] = dbstatuses;
111  }
112 
113  flagmask_ = 0;
116  flagmask_ |= 0x1 << EcalRecHit::kDead;
117  flagmask_ |= 0x1 << EcalRecHit::kKilled;
120 }
121 
123  ical = es.getHandle(icalToken_);
124 
125  if (!skipTimeCalib_) {
128  }
129 
130  agc = es.getHandle(agcToken_);
132  if (laserCorrection_)
134 }
135 
137  const EcalUncalibratedRecHit& uncalibRH,
139  DetId detid = uncalibRH.id();
140 
142  EcalChannelStatusCode::Code dbstatus = chit->getStatusCode();
143 
144  // check for channels to be excluded from reconstruction
145  if (!v_chstatus_.empty()) {
146  std::vector<int>::const_iterator res = std::find(v_chstatus_.begin(), v_chstatus_.end(), dbstatus);
147  if (res != v_chstatus_.end())
148  return false;
149  }
150 
151  uint32_t flagBits = setFlagBits(v_DB_reco_flags_, dbstatus);
152 
153  float offsetTime = 0; // the global time phase
154  const EcalIntercalibConstantMap& icalMap = ical->getMap();
155  if (detid.subdetId() == EcalEndcap) {
157  if (!skipTimeCalib_)
158  offsetTime = offtime->getEEValue();
159  } else {
161  if (!skipTimeCalib_)
162  offsetTime = offtime->getEBValue();
163  }
164 
165  // first intercalibration constants
166  EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(detid);
167  EcalIntercalibConstant icalconst = 1;
168  if (icalit != icalMap.end()) {
169  icalconst = (*icalit);
170  } else {
171  edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << detid.rawId()
172  << "! something wrong with EcalIntercalibConstants in your DB? ";
173  }
174 
175  // get laser coefficient
176  float lasercalib = 1.;
177  if (laserCorrection_)
178  lasercalib = laser->getLaserCorrection(detid, evt.time());
179 
180  // get time calibration coefficient
181  EcalTimeCalibConstant itimeconst = 0;
182 
183  if (!skipTimeCalib_) {
184  const EcalTimeCalibConstantMap& itimeMap = itime->getMap();
186 
187  if (itime != itimeMap.end()) {
188  itimeconst = (*itime);
189  } else {
190  edm::LogError("EcalRecHitError") << "No time calib const found for xtal " << detid.rawId()
191  << "! something wrong with EcalTimeCalibConstants in your DB? ";
192  }
193  }
194 
195  // make the rechit and put in the output collection, unless recovery has to take care of it
196  if (!(flagmask_ & flagBits) || !killDeadChannels_) {
197  EcalRecHit myrechit(rechitMaker_->makeRecHit(uncalibRH,
198  icalconst * lasercalib,
199  (itimeconst + offsetTime),
200  /*recoflags_ 0*/
201  flagBits));
202 
203  if (detid.subdetId() == EcalBarrel && (lasercalib < EBLaserMIN_ || lasercalib > EBLaserMAX_))
205  if (detid.subdetId() == EcalEndcap && (lasercalib < EELaserMIN_ || lasercalib > EELaserMAX_))
206  myrechit.setFlag(EcalRecHit::kPoorCalib);
207  result.push_back(myrechit);
208  }
209 
210  return true;
211 }
212 
213 // Take our association map of dbstatuses-> recHit flagbits and return the apporpriate flagbit word
214 uint32_t EcalRecHitWorkerSimple::setFlagBits(const std::vector<std::vector<uint32_t> >& map, const uint32_t& status) {
215  for (unsigned int i = 0; i != map.size(); ++i) {
216  if (std::find(map[i].begin(), map[i].end(), status) != map[i].end())
217  return 0x1 << i;
218  }
219 
220  return 0;
221 }
222 
224 
edm::ESGetToken< EcalTimeCalibConstants, EcalTimeCalibConstantsRcd > itimeToken_
edm::ESGetToken< EcalLaserDbService, EcalLaserDbRecord > laserToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
EcalRecHitWorkerSimple(const edm::ParameterSet &, edm::ConsumesCollector &c)
edm::ESHandle< EcalADCToGeVConstant > agc
std::vector< std::vector< uint32_t > > v_DB_reco_flags_
edm::ESGetToken< EcalTimeOffsetConstant, EcalTimeOffsetConstantRcd > offtimeToken_
void setFlag(int flag)
set the flags (from Flags or ESFlags)
Definition: EcalRecHit.h:184
edm::ESHandle< EcalChannelStatus > chStatus
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
Definition: Electron.h:6
edm::ESGetToken< EcalIntercalibConstants, EcalIntercalibConstantsRcd > icalToken_
edm::ESHandle< EcalTimeCalibConstants > itime
EcalChannelStatusMap EcalChannelStatus
edm::Timestamp time() const
Definition: EventBase.h:64
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
EcalRecHitSimpleAlgo * rechitMaker_
float getLaserCorrection(DetId const &xid, edm::Timestamp const &iTime) const
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::vector< int > v_chstatus_
const_iterator find(uint32_t rawId) const
void setADCToGeVConstant(const float &value) override
make rechits from dataframes
Definition: DetId.h:17
EcalTimeCalibConstantMap EcalTimeCalibConstants
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< Item >::const_iterator const_iterator
uint32_t setFlagBits(const std::vector< std::vector< uint32_t > > &map, const uint32_t &status)
EcalRecHit makeRecHit(const EcalUncalibratedRecHit &uncalibRH, const float &intercalibConstant, const float &timeIntercalib=0, const uint32_t &flags=0) const override
Compute parameters.
edm::ESHandle< EcalTimeOffsetConstant > offtime
edm::ESGetToken< EcalADCToGeVConstant, EcalADCToGeVConstantRcd > agcToken_
edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > chStatusToken_
EcalIntercalibConstantMap EcalIntercalibConstants
float EcalTimeCalibConstant
#define DEFINE_EDM_PLUGIN(factory, type, name)
const_iterator end() const
bool run(const edm::Event &evt, const EcalUncalibratedRecHit &uncalibRH, EcalRecHitCollection &result) override
void set(const edm::EventSetup &es) override
edm::ESHandle< EcalIntercalibConstants > ical
float EcalIntercalibConstant
edm::ESHandle< EcalLaserDbService > laser