CMS 3D CMS Logo

EcalRecHitWorkerRecover.cc
Go to the documentation of this file.
2 
8 
11 
13 
16 
18  : EcalRecHitWorkerBaseClass(ps, c), ecalScaleTokens_(c), tpgscaleTokens_(c) {
19  rechitMaker_ = std::make_unique<EcalRecHitSimpleAlgo>();
20  // isolated channel recovery
21  singleRecoveryMethod_ = ps.getParameter<std::string>("singleChannelRecoveryMethod");
22  singleRecoveryThreshold_ = ps.getParameter<double>("singleChannelRecoveryThreshold");
23  sum8RecoveryThreshold_ = ps.getParameter<double>("sum8ChannelRecoveryThreshold");
24  killDeadChannels_ = ps.getParameter<bool>("killDeadChannels");
25  recoverEBIsolatedChannels_ = ps.getParameter<bool>("recoverEBIsolatedChannels");
26  recoverEEIsolatedChannels_ = ps.getParameter<bool>("recoverEEIsolatedChannels");
27  recoverEBVFE_ = ps.getParameter<bool>("recoverEBVFE");
28  recoverEEVFE_ = ps.getParameter<bool>("recoverEEVFE");
29  recoverEBFE_ = ps.getParameter<bool>("recoverEBFE");
30  recoverEEFE_ = ps.getParameter<bool>("recoverEEFE");
38 
39  dbStatusToBeExcludedEE_ = ps.getParameter<std::vector<int> >("dbStatusToBeExcludedEE");
40  dbStatusToBeExcludedEB_ = ps.getParameter<std::vector<int> >("dbStatusToBeExcludedEB");
41 
42  logWarningEtThreshold_EB_FE_ = ps.getParameter<double>("logWarningEtThreshold_EB_FE");
43  logWarningEtThreshold_EE_FE_ = ps.getParameter<double>("logWarningEtThreshold_EE_FE");
44 
45  tpDigiToken_ =
46  c.consumes<EcalTrigPrimDigiCollection>(ps.getParameter<edm::InputTag>("triggerPrimitiveDigiCollection"));
47 
50 }
51 
57  // geometry...
64  recoveredDetIds_EB_.clear();
65  recoveredDetIds_EE_.clear();
66  eventSetup_ = &es;
67 }
68 
70  const EcalUncalibratedRecHit& uncalibRH,
72  DetId detId = uncalibRH.id();
73  uint32_t flags = (0xF & uncalibRH.flags());
74 
75  // get laser coefficient
76  //float lasercalib = laser->getLaserCorrection( detId, evt.time());
77 
78  // killDeadChannels_ = true, means explicitely kill dead channels even if the recovered energies are computed in the code
79  // if you don't want to store the recovered energies in the rechit you can produce LogWarnings if logWarningEtThreshold_EB(EE)_FE>0
80  // logWarningEtThreshold_EB(EE)_FE_<0 will not compute the recovered energies at all (faster)
81 
82  if (killDeadChannels_) {
87  EcalRecHit hit(detId, 0., 0., EcalRecHit::kDead);
88  hit.setFlag(EcalRecHit::kDead);
89  insertRecHit(hit, result); // insert trivial rechit with kDead flag
90  return true;
91  }
93  EcalTrigTowerDetId ttDetId(((EBDetId)detId).tower());
94  std::vector<DetId> vid = ttMap_->constituentsOf(ttDetId);
95  for (std::vector<DetId>::const_iterator dit = vid.begin(); dit != vid.end(); ++dit) {
96  EcalRecHit hit((*dit), 0., 0., EcalRecHit::kDead);
97  hit.setFlag(EcalRecHit::kDead);
98  insertRecHit(hit, result); // insert trivial rechit with kDead flag
99  }
101  return true; // if you don't want log warning just return true
102  }
104  EEDetId id(detId);
105  EcalScDetId sc(1 + (id.ix() - 1) / 5, 1 + (id.iy() - 1) / 5, id.zside());
106  std::vector<DetId> eeC;
107  for (int dx = 1; dx <= 5; ++dx) {
108  for (int dy = 1; dy <= 5; ++dy) {
109  int ix = (sc.ix() - 1) * 5 + dx;
110  int iy = (sc.iy() - 1) * 5 + dy;
111  int iz = sc.zside();
112  if (EEDetId::validDetId(ix, iy, iz)) {
113  eeC.push_back(EEDetId(ix, iy, iz));
114  }
115  }
116  }
117  for (size_t i = 0; i < eeC.size(); ++i) {
118  EcalRecHit hit(eeC[i], 0., 0., EcalRecHit::kDead);
119  hit.setFlag(EcalRecHit::kDead);
120  insertRecHit(hit, result); // insert trivial rechit with kDead flag
121  }
123  return true; // if you don't want log warning just return true
124  }
125  }
126 
128  // recover as single dead channel
130 
131  // channel recovery. Accepted new RecHit has the flag AcceptRecHit=TRUE
132  bool AcceptRecHit = true;
133  float ebEn = ebDeadChannelCorrector.correct(
135  EcalRecHit hit(detId, ebEn, 0., EcalRecHit::kDead);
136 
137  if (hit.energy() != 0 and AcceptRecHit == true) {
139  } else {
140  // recovery failed
141  hit.setFlag(EcalRecHit::kDead);
142  }
144 
146  // recover as single dead channel
148 
149  // channel recovery. Accepted new RecHit has the flag AcceptRecHit=TRUE
150  bool AcceptRecHit = true;
151  float eeEn = eeDeadChannelCorrector.correct(
153  EcalRecHit hit(detId, eeEn, 0., EcalRecHit::kDead);
154  if (hit.energy() != 0 and AcceptRecHit == true) {
156  } else {
157  // recovery failed
158  hit.setFlag(EcalRecHit::kDead);
159  }
161 
162  } else if (flags == EcalRecHitWorkerRecover::EB_VFE) {
163  // recover as dead VFE
164  EcalRecHit hit(detId, 0., 0.);
165  hit.setFlag(EcalRecHit::kDead);
166  // recovery not implemented
168  } else if (flags == EcalRecHitWorkerRecover::EB_FE) {
169  // recover as dead TT
170 
171  EcalTrigTowerDetId ttDetId(((EBDetId)detId).tower());
173  evt.getByToken(tpDigiToken_, pTPDigis);
174  const EcalTrigPrimDigiCollection* tpDigis = nullptr;
175  tpDigis = pTPDigis.product();
176 
178  // recover the whole trigger tower
179  if (tp != tpDigis->end()) {
181  //std::vector<DetId> vid = ecalMapping_->dccTowerConstituents( ecalMapping_->DCCid( ttDetId ), ecalMapping_->iTT( ttDetId ) );
182  std::vector<DetId> vid = ttMap_->constituentsOf(ttDetId);
183  float tpEt = ecalScale.getTPGInGeV(tp->compressedEt(), tp->id());
184  float tpEtThreshEB = logWarningEtThreshold_EB_FE_;
185  if (tpEt > tpEtThreshEB) {
186  edm::LogWarning("EnergyInDeadEB_FE") << "TP energy in the dead TT = " << tpEt << " at " << ttDetId;
187  }
189  // democratic energy sharing
190 
191  for (std::vector<DetId>::const_iterator dit = vid.begin(); dit != vid.end(); ++dit) {
192  if (alreadyInserted(*dit))
193  continue;
194  float theta = ebGeom_->getGeometry(*dit)->getPosition().theta();
195  float tpEt = ecalScale.getTPGInGeV(tp->compressedEt(), tp->id());
197  EcalRecHit hit(*dit, tpEt / ((float)vid.size()) / sin(theta), 0.);
199  if (tp->compressedEt() == 0xFF)
200  hit.setFlag(EcalRecHit::kTPSaturated);
201  if (tp->sFGVB())
202  hit.setFlag(EcalRecHit::kL1SpikeFlag);
204  }
205  }
206  } else {
207  // tp not found => recovery failed
208  std::vector<DetId> vid = ttMap_->constituentsOf(ttDetId);
209  for (std::vector<DetId>::const_iterator dit = vid.begin(); dit != vid.end(); ++dit) {
210  if (alreadyInserted(*dit))
211  continue;
212  EcalRecHit hit(*dit, 0., 0.);
213  hit.setFlag(EcalRecHit::kDead);
215  }
216  }
217  }
218  } else if (flags == EcalRecHitWorkerRecover::EE_FE) {
219  // Structure for recovery:
220  // ** SC --> EEDetId constituents (eeC) --> associated Trigger Towers (aTT) --> EEDetId constituents (aTTC)
221  // ** energy for a SC EEDetId = [ sum_aTT(energy) - sum_aTTC(energy) ] / N_eeC
222  // .. i.e. the total energy of the TTs covering the SC minus
223  // .. the energy of the recHits in the TTs but not in the SC
224  //std::vector<DetId> vid = ecalMapping_->dccTowerConstituents( ecalMapping_->DCCid( ttDetId ), ecalMapping_->iTT( ttDetId ) );
225  // due to lack of implementation of the EcalTrigTowerDetId ix,iy methods in EE we compute Et recovered energies (in EB we compute E)
226 
227  EEDetId eeId(detId);
228  EcalScDetId sc((eeId.ix() - 1) / 5 + 1, (eeId.iy() - 1) / 5 + 1, eeId.zside());
229  std::set<DetId> eeC;
230  for (int dx = 1; dx <= 5; ++dx) {
231  for (int dy = 1; dy <= 5; ++dy) {
232  int ix = (sc.ix() - 1) * 5 + dx;
233  int iy = (sc.iy() - 1) * 5 + dy;
234  int iz = sc.zside();
235  if (EEDetId::validDetId(ix, iy, iz)) {
236  EEDetId id(ix, iy, iz);
238  eeC.insert(id);
239  } // check status
240  }
241  }
242  }
243 
245  evt.getByToken(tpDigiToken_, pTPDigis);
246  const EcalTrigPrimDigiCollection* tpDigis = nullptr;
247  tpDigis = pTPDigis.product();
248 
249  // associated trigger towers
250  std::set<EcalTrigTowerDetId> aTT;
251  for (std::set<DetId>::const_iterator it = eeC.begin(); it != eeC.end(); ++it) {
252  aTT.insert(ttMap_->towerOf(*it));
253  }
254 
257  // associated trigger towers: total energy
258  float totE = 0;
259  // associated trigger towers: EEDetId constituents
260  std::set<DetId> aTTC;
261  bool atLeastOneTPSaturated = false;
262  for (std::set<EcalTrigTowerDetId>::const_iterator it = aTT.begin(); it != aTT.end(); ++it) {
263  // add the energy of this trigger tower
265  if (itTP != tpDigis->end()) {
266  std::vector<DetId> v = ttMap_->constituentsOf(*it);
267 
268  // from the constituents, remove dead channels
269  std::vector<DetId>::iterator ttcons = v.begin();
270  while (ttcons != v.end()) {
272  ttcons = v.erase(ttcons);
273  } else {
274  ++ttcons;
275  }
276  } // while
277 
278  if (itTP->compressedEt() == 0xFF) { // In the case of a saturated trigger tower, a fraction
279  atLeastOneTPSaturated =
280  true; //of the saturated energy is put in: number of xtals in dead region/total xtals in TT *63.75
281 
282  //Alternative recovery algorithm that I will now investigate.
283  //Estimate energy sums the energy in the working channels, then decides how much energy
284  //to put here depending on that. Duncan 20101203
285 
286  totE += estimateEnergy(itTP->id().ietaAbs(), &result, eeC, v, tpgscale);
287 
288  /*
289  These commented out lines use
290  64GeV*fraction of the TT overlapping the dead FE
291 
292  int count = 0;
293  for (std::vector<DetId>::const_iterator idsit = v.begin(); idsit != v.end(); ++ idsit){
294  std::set<DetId>::const_iterator itFind = eeC.find(*idsit);
295  if (itFind != eeC.end())
296  ++count;
297  }
298  //std::cout << count << ", " << v.size() << std::endl;
299  totE+=((float)count/(float)v.size())* ((it->ietaAbs()>26)?2*ecalScale_.getTPGInGeV( itTP->compressedEt(), itTP->id() ):ecalScale_.getTPGInGeV( itTP->compressedEt(), itTP->id() ));*/
300  } else {
301  totE += ((it->ietaAbs() > 26) ? 2 : 1) * ecalScale.getTPGInGeV(itTP->compressedEt(), itTP->id());
302  }
303 
304  // get the trigger tower constituents
305 
306  if (itTP->compressedEt() == 0) { // If there's no energy in TT, the constituents are removed from the recovery.
307  for (size_t i = 0; i < v.size(); ++i)
308  eeC.erase(v[i]);
309  } else if (itTP->compressedEt() != 0xFF) {
310  //If it's saturated the energy has already been determined, so we do not want to subtract any channels
311  for (size_t j = 0; j < v.size(); ++j) {
312  aTTC.insert(v[j]);
313  }
314  }
315  }
316  }
317  // remove crystals of dead SC
318  // (this step is not needed if sure that SC crystals are not
319  // in the recHit collection)
320 
321  for (std::set<DetId>::const_iterator it = eeC.begin(); it != eeC.end(); ++it) {
322  aTTC.erase(*it);
323  }
324  // compute the total energy for the dead SC
325  const EcalRecHitCollection* hits = &result;
326  for (std::set<DetId>::const_iterator it = aTTC.begin(); it != aTTC.end(); ++it) {
328  if (jt != hits->end()) {
329  float energy = jt->energy(); // Correct conversion to Et
330  float eta = geo_->getPosition(jt->id()).eta();
331  float pf = 1.0 / cosh(eta);
332  // use Et instead of E, consistent with the Et estimation of the associated TT
333  totE -= energy * pf;
334  }
335  }
336 
337  float scEt = totE;
338  float scEtThreshEE = logWarningEtThreshold_EE_FE_;
339  if (scEt > scEtThreshEE) {
340  edm::LogWarning("EnergyInDeadEE_FE") << "TP energy in the dead TT = " << scEt << " at " << sc;
341  }
342 
343  // assign the energy to the SC crystals
344  if (!killDeadChannels_ || recoverEEFE_) { // if eeC is empty, i.e. there are no hits
345  // in the tower, nothing is returned. No negative values from noise.
346  for (std::set<DetId>::const_iterator it = eeC.begin(); it != eeC.end(); ++it) {
347  float eta = geo_->getPosition(*it).eta(); //Convert back to E from Et for the recovered hits
348  float pf = 1.0 / cosh(eta);
349  EcalRecHit hit(*it, totE / ((float)eeC.size() * pf), 0);
350 
351  if (atLeastOneTPSaturated)
352  hit.setFlag(EcalRecHit::kTPSaturated);
355 
356  } // for
357  } // if
358  }
359  return true;
360 }
361 
364  const std::set<DetId>& sId,
365  const std::vector<DetId>& vId,
366  const EcalTPGScale& tpgscale) {
367  float xtalE = 0;
368  int count = 0;
369  for (std::vector<DetId>::const_iterator vIdit = vId.begin(); vIdit != vId.end(); ++vIdit) {
370  std::set<DetId>::const_iterator sIdit = sId.find(*vIdit);
371  if (sIdit == sId.end()) {
372  float energy = hits->find(*vIdit)->energy();
373  float eta = geo_->getPosition(*vIdit).eta();
374  float pf = 1.0 / cosh(eta);
375  xtalE += energy * pf;
376  count++;
377  }
378  }
379 
380  if (count == 0) { // If there are no overlapping crystals return saturated value.
381 
382  double etsat = tpgscale.getTPGInGeV(0xFF,
383  ttMap_->towerOf(*vId.begin())); // get saturation value for the first
384  // constituent, for the others it's the same
385 
386  return etsat / cosh(ieta) * (ieta > 26 ? 2 : 1); // account for duplicated TT in EE for ieta>26
387  } else
388  return xtalE * ((vId.size() / (float)count) - 1) * (ieta > 26 ? 2 : 1);
389 }
390 
392  // skip already inserted DetId's and raise a log warning
393  if (alreadyInserted(hit.id())) {
394  edm::LogWarning("EcalRecHitWorkerRecover") << "DetId already recovered! Skipping...";
395  return;
396  }
398  if (it == collection.end()) {
399  // insert the hit in the collection
400  collection.push_back(hit);
401  } else {
402  // overwrite existing recHit
403  *it = hit;
404  }
405  if (hit.id().subdetId() == EcalBarrel) {
406  recoveredDetIds_EB_.insert(hit.id());
407  } else if (hit.id().subdetId() == EcalEndcap) {
408  recoveredDetIds_EE_.insert(hit.id());
409  } else {
410  edm::LogError("EcalRecHitWorkerRecover::InvalidDetId") << "Invalid DetId " << hit.id().rawId();
411  }
412 }
413 
415  bool res = false;
416  if (id.subdetId() == EcalBarrel) {
417  res = (recoveredDetIds_EB_.find(id) != recoveredDetIds_EB_.end());
418  } else if (id.subdetId() == EcalEndcap) {
419  res = (recoveredDetIds_EE_.find(id) != recoveredDetIds_EE_.end());
420  } else {
421  edm::LogError("EcalRecHitWorkerRecover::InvalidDetId") << "Invalid DetId " << id.rawId();
422  }
423  return res;
424 }
425 
426 // In the future, this will be used to calibrate the TT energy. There is a dependance on
427 // eta at lower energies that can be corrected for here after more validation.
428 float EcalRecHitWorkerRecover::recCheckCalib(float eTT, int ieta) { return eTT; }
429 
430 // return false is the channel has status in the list of statusestoexclude
431 // true otherwise (channel ok)
432 // Careful: this function works on raw (encoded) channel statuses
433 bool EcalRecHitWorkerRecover::checkChannelStatus(const DetId& id, const std::vector<int>& statusestoexclude) {
434  if (!chStatus_.isValid())
435  edm::LogError("ObjectNotFound") << "Channel Status not set";
436 
438  uint16_t dbStatus = 0;
439  if (chIt != chStatus_->end()) {
440  dbStatus = chIt->getEncodedStatusCode();
441  } else {
442  edm::LogError("ObjectNotFound") << "No channel status found for xtal " << id.rawId()
443  << "! something wrong with EcalChannelStatus in your DB? ";
444  }
445 
446  for (std::vector<int>::const_iterator status = statusestoexclude.begin(); status != statusestoexclude.end();
447  ++status) {
448  if (*status == dbStatus)
449  return false;
450  }
451 
452  return true;
453 }
454 
edm::ESGetToken< CaloSubdetectorGeometry, EcalBarrelGeometryRecord > pEBGeomToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EventSetup * eventSetup_
float estimateEnergy(int ieta, EcalRecHitCollection *hits, const std::set< DetId > &sId, const std::vector< DetId > &vId, const EcalTPGScale &tpgscale)
int zside() const
Definition: EcalScDetId.h:64
bool alreadyInserted(const DetId &id)
bool run(const edm::Event &evt, const EcalUncalibratedRecHit &uncalibRH, EcalRecHitCollection &result) override
std::vector< DetId > constituentsOf(const EcalTrigTowerDetId &id) const
Get the constituent detids for this tower id.
EcalDeadChannelRecoveryAlgos< EEDetId > eeDeadChannelCorrector
edm::ESGetToken< EcalElectronicsMapping, EcalMappingRcd > pEcalMappingToken_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeometryToken_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
Definition: PV3DBase.h:73
T const * product() const
Definition: Handle.h:70
int ix() const
Definition: EEDetId.h:77
std::vector< T >::const_iterator const_iterator
std::vector< int > dbStatusToBeExcludedEE_
edm::ESGetToken< EcalLaserDbService, EcalLaserDbRecord > laserToken_
edm::ESHandle< EcalElectronicsMapping > pEcalMapping_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::ESGetToken< EcalChannelStatus, EcalChannelStatusRcd > chStatusToken_
int zside(DetId const &)
edm::ESHandle< CaloTopology > caloTopology_
EcalTPGScale::Tokens tpgscaleTokens_
Log< level::Error, false > LogError
float recCheckCalib(float energy, int ieta)
void setCaloTopology(const CaloTopology *topology)
Definition: Electron.h:6
EcalTrigTowerDetId towerOf(const DetId &id) const
Get the tower id for this det id (or null if not known)
void insertRecHit(const EcalRecHit &hit, EcalRecHitCollection &collection)
EcalTPGScale::Tokens ecalScaleTokens_
std::unique_ptr< EcalRecHitSimpleAlgo > rechitMaker_
T const * product() const
Definition: ESHandle.h:86
EcalChannelStatusMap EcalChannelStatus
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:50
edm::ESHandle< EcalTrigTowerConstituentsMap > ttMap_
void setParameters(const edm::ParameterSet &ps)
edm::ESHandle< CaloSubdetectorGeometry > pEBGeom_
int iy() const
Definition: EcalScDetId.h:76
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
const CaloSubdetectorGeometry * ebGeom_
std::vector< EcalRecHit >::iterator iterator
unsigned int id
const_iterator end() const
Definition: DetId.h:17
void set(const edm::EventSetup &es) override
edm::ESGetToken< EcalTrigTowerConstituentsMap, IdealGeometryRecord > ttMapToken_
int zside() const
Definition: EEDetId.h:71
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
edm::EDGetTokenT< EcalTrigPrimDigiCollection > tpDigiToken_
std::vector< Item >::const_iterator const_iterator
std::set< DetId > recoveredDetIds_EE_
double getTPGInGeV(const EcalTriggerPrimitiveDigi &tpDigi) const
Definition: EcalTPGScale.cc:17
iterator find(key_type k)
EcalRecHitWorkerRecover(const edm::ParameterSet &, edm::ConsumesCollector &c)
edm::ESHandle< CaloGeometry > caloGeometry_
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< int > dbStatusToBeExcludedEB_
float correct(const DetIdT id, const EcalRecHitCollection &hit_collection, std::string algo, double single8Cut, double sum8Cut, bool *accFlag)
edm::ESGetToken< CaloTopology, CaloTopologyRecord > caloTopologyToken_
Log< level::Warning, false > LogWarning
EcalDeadChannelRecoveryAlgos< EBDetId > ebDeadChannelCorrector
Geom::Theta< T > theta() const
edm::ESHandle< EcalChannelStatus > chStatus_
edm::ESHandle< EcalLaserDbService > laser
const EcalElectronicsMapping * ecalMapping_
bool checkChannelStatus(const DetId &id, const std::vector< int > &statusestoexclude)
int iy() const
Definition: EEDetId.h:83
std::set< DetId > recoveredDetIds_EB_
int ix() const
Definition: EcalScDetId.h:70