ATLAS Offline Software
LArCalibPatchingAlg.icc
Go to the documentation of this file.
1 //Dear emacs, this is -*-c++-*-
2 /*
3  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
4 */
5 
6 #include <tuple>
7 
8 template<class CONDITIONSCONTAINER>
9 LArCalibPatchingAlg<CONDITIONSCONTAINER>::LArCalibPatchingAlg (const std::string& name, ISvcLocator* pSvcLocator) :
10  AthAlgorithm(name,pSvcLocator),
11  m_onlineHelper(0),
12  m_caloId(0),
13  m_contIn(0),
14  m_contOut(0),
15  m_patchMethod(PhiAverage) { }
16 
17 
18 template<class CONDITIONSCONTAINER>
19 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::initialize() {
20 
21  ATH_CHECK( m_BCKey.initialize() );
22  ATH_CHECK( m_cablingKey.initialize() );
23  ATH_CHECK( m_CLKey.initialize() );
24 
25 
26  if(m_patchMethodProp=="FEBNeighbor") {
27  m_patchMethod=FEBNeighbor;
28  return StatusCode::SUCCESS;
29  }
30  else if (m_patchMethodProp=="PhiNeighbor") {
31  m_patchMethod=PhiNeighbor;
32  return StatusCode::SUCCESS;
33  }
34  else if (m_patchMethodProp=="PhiAverage") {
35  if (typeid(CONDITIONSCONTAINER)==typeid(LArAutoCorrComplete)) {
36  ATH_MSG_ERROR ( "PhiAverage not implemented for LArAutoCorrComlete."
37  << "Please choose other patching strategy" );
38  return StatusCode::FAILURE;
39  }
40  m_patchMethod=PhiAverage;
41  return StatusCode::SUCCESS;
42  }
43  else if (m_patchMethodProp=="FEBAverage") {
44  if (typeid(CONDITIONSCONTAINER)==typeid(LArCaliWaveContainer)) {
45  ATH_MSG_ERROR ( "FEBAverage not implemented for CaliWaveContainer."
46  << "Please choose other patching strategy" );
47  return StatusCode::FAILURE;
48  }
49  m_patchMethod=FEBAverage;
50  return StatusCode::SUCCESS;
51  } else if (m_patchMethodProp=="SetZero") {
52  if (typeid(CONDITIONSCONTAINER)==typeid(LArCaliWaveContainer)) {
53  ATH_MSG_ERROR ( "SetZero not implemented for CaliWaveContainer."
54  << "Please choose other patching strategy" );
55  return StatusCode::FAILURE;
56  }
57  m_patchMethod=SetZero;
58  return StatusCode::SUCCESS;
59  }
60 
61  ATH_MSG_ERROR ( "Unknown patching method: " << m_patchMethodProp );
62  ATH_MSG_ERROR ( "Allowed values: [Empty, FEBNeighbor, PhiNeighbor, PhiAverage, SetZero]" );
63 
64  return StatusCode::FAILURE;
65 }
66 
67 template<class CONDITIONSCONTAINER>
68 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::stop() {
69  ATH_MSG_INFO ( "Entering LArCalibPatchingAlg" );
70  if(m_isSC) {
71  const LArOnline_SuperCellID *onlHelper=nullptr;
72  ATH_CHECK( detStore()->retrieve(onlHelper, "LArOnline_SuperCellID") );
73  m_onlineHelper=dynamic_cast<const LArOnlineID_Base*>(onlHelper);
74  const CaloCell_SuperCell_ID *caloid=nullptr;
75  ATH_CHECK( detStore()->retrieve(caloid, "CaloCell_SuperCell_ID") );
76  m_caloId = dynamic_cast<const CaloCell_Base_ID*>(caloid);
77  } else {
78  const LArOnlineID *onlHelper=nullptr;
79  ATH_CHECK( detStore()->retrieve(onlHelper, "LArOnlineID") );
80  m_onlineHelper=dynamic_cast<const LArOnlineID_Base*>(onlHelper);
81  const CaloCell_ID *caloid=nullptr;
82  ATH_CHECK( detStore()->retrieve(caloid, "CaloCell_ID") );
83  m_caloId = dynamic_cast<const CaloCell_Base_ID*>(caloid);
84  }
85 
86  if(m_isSC) m_bcMask.setSC();
87  ATH_CHECK(m_bcMask.buildBitMask(m_problemsToPatch,msg()));
88 
89  const EventContext& ctx = Gaudi::Hive::currentContext();
90 
91  SG::ReadCondHandle<LArBadChannelCont> readHandle{m_BCKey, ctx};
92  const LArBadChannelCont *bcCont {*readHandle};
93  if(!bcCont) {
94  ATH_MSG_ERROR( "Do not have Bad chan container !!!" );
95  return StatusCode::FAILURE;
96  }
97 
98  SG::ReadCondHandle<LArOnOffIdMapping> cablingHdl{m_cablingKey, ctx};
99  const LArOnOffIdMapping* cabling = *cablingHdl;
100  if(!cabling) {
101  ATH_MSG_ERROR( "Do not have OnOff Id mapping !!!" );
102  return StatusCode::FAILURE;
103  }
104 
105  SG::ReadCondHandle<LArCalibLineMapping> clHdl{m_CLKey, ctx};
106  const LArCalibLineMapping *clCont {*clHdl};
107  if(!clCont) {
108  ATH_MSG_ERROR( "Do not have calib line mapping !!!" );
109  return StatusCode::FAILURE;
110  }
111 
112  if (m_newContainerKey.size()) {
113  //New container key give -> different containers for reading and writing
114  ATH_CHECK( detStore()->retrieve(m_contIn,m_containerKey) ); //const-retrieve
115 
116  m_contOut=new CONDITIONSCONTAINER();
117  m_contOut->setGroupingType((LArConditionsContainerBase::GroupingType)m_contIn->groupingType());
118  ATH_CHECK( m_contOut->initialize() );
119  ATH_CHECK( detStore()->record(m_contOut,m_newContainerKey) );
120  ATH_CHECK( setSymlink(m_contOut) );
121  ATH_MSG_INFO ( "Loaded input container " << m_containerKey
122  << ", write to new container " << m_newContainerKey );
123  }
124  else { //Same container for reading and writing
125  if (m_unlock) {
126  ATH_CHECK( detStore()->retrieve(m_contIn,m_containerKey) ); //const-retrieve
127  m_contOut=const_cast<CONDITIONSCONTAINER*>(m_contIn);
128  }
129  else{
130  ATH_CHECK( detStore()->retrieve(m_contOut,m_containerKey) ); //non-const retrieve
131  m_contIn=const_cast<const CONDITIONSCONTAINER*>(m_contOut);
132  ATH_MSG_INFO ( "Work on container '" << m_containerKey << "'" );
133  }
134  }
135 
136  LArBadChanBitPacking packing;
137  LArBadChanSCBitPacking scpacking;
138 
139  unsigned maxgain;
140  if(m_isSC) maxgain=CaloGain::LARMEDIUMGAIN; else maxgain=CaloGain::LARNGAIN;
141  for (unsigned igain=CaloGain::LARHIGHGAIN;
142  igain<maxgain ; ++igain ) {
143  CONTIT it=m_contIn->begin(igain);
144  CONTIT it_e=m_contIn->end(igain);
145  for (;it!=it_e;it++) {
146  const HWIdentifier chid = it.channelId();
147  if (!cabling->isOnlineConnected(chid)) continue; //Don't care about disconnected channels
148  if (m_bcMask.cellShouldBeMasked(bcCont,chid)) {
149  const std::string bcType = m_isSC ? scpacking.stringStatus(bcCont->status(chid)) : packing.stringStatus(bcCont->status(chid));
150  ATH_MSG_INFO ( "Found problematic channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec << " [" << bcType << "]"
151  <<" Gain:" << igain << " " << m_onlineHelper->channel_name(chid) << ". Trying to patch." );
152  if (patch(chid,igain, bcCont, cabling, clCont)) {
153  ATH_MSG_INFO ( "Sucessfully patched channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec <<" Gain:" << igain );
154  } else {
155  ATH_MSG_WARNING ( "Failed to patch channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec <<" Gain:" << igain );
156  }
157  }//end if channel is in bad-channel database
158  else
159  if (it->isEmpty()) { //check if data-object is empty (eg the default instance
160  ATH_MSG_ERROR ( "The channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec
161  <<" Gain:" << igain << " " << m_onlineHelper->channel_name(chid)
162  << " has no calibration but is not (yet?) flagged in the bad-channel database" );
163  if (m_patchAllMissing) {
164  ATH_MSG_INFO ( "Will try to patch anyway." );
165  if (patch(chid,igain, bcCont, cabling, clCont)) {
166  ATH_MSG_INFO ( "Sucessfully patched channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec <<" Gain:" << igain );
167  } else {
168  ATH_MSG_WARNING ( "Failed to patch channel 0x" << MSG::hex << chid.get_identifier32().get_compact() << MSG::dec <<" Gain:" << igain );
169  }
170  }//end if m_patchAllMissing
171  else
172  ATH_MSG_ERROR ( "Channel remains un-patched!" );
173  }//end if isEmpty
174  }//end loop over all channels
175  }//end loop over all gains
176 
177  std::vector<unsigned> completedChans = m_contOut->completeCorrectionChannels();
178  if (completedChans.size()>0 && m_useCorrChannel) {
179  msg() << MSG::INFO << "Artificially inserted correction subsets in COOL channels";
180  for(size_t j=0;j<completedChans.size();++j)
181  msg() << MSG::INFO << " " << completedChans[j];
182  msg() << MSG::INFO << endmsg;
183  }
184 
185 
186  ATH_MSG_INFO ( "Done with LArCalibPatchingAlg" );
187  ATH_MSG_DEBUG ( detStore()->dump() );
188  return StatusCode::SUCCESS;
189 }
190 
191 template<class CONDITIONSCONTAINER>
192 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::patch(const HWIdentifier chid, const int gain, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, const LArCalibLineMapping *clCont) {
193 
194  if(m_doNotPatchCBs.size()>0) { // check if this channel is not on excluded CalibBoard
195  const std::vector<HWIdentifier>& cLids=clCont->calibSlotLine(chid);
196  for(unsigned cl=0; cl<cLids.size(); ++cl) {
197  const HWIdentifier calibModuleID = m_onlineHelper->calib_module_Id(cLids[cl]);
198  if (std::find(m_doNotPatchCBs.begin(),m_doNotPatchCBs.end(),calibModuleID.get_identifier32().get_compact()) != m_doNotPatchCBs.end()) { // we should not patch this channel
199  return false;
200  }
201  } // over CLs
202  }
203 
204  if (m_patchMethod==FEBNeighbor){
205  const HWIdentifier febId=m_onlineHelper->feb_Id(chid);
206  const int febChan=m_onlineHelper->channel(chid);
207  HWIdentifier chid_patch(0);
208  if (febChan>0) {
209  //try lower channel
210  chid_patch=m_onlineHelper->channel_Id(febId, febChan-1);
211  if (cabling->isOnlineConnected(chid_patch) && bcCont->status(chid_patch).good()) {
212  const LArCondObj patch=m_contIn->get(chid_patch,gain); //Should be the const-get method
213  if (!patch.isEmpty()) {
214  StatusCode sc=m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel);
215  if (sc.isFailure())
216  ATH_MSG_ERROR ( "Failed to insert correction for channel 0x" << MSG::hex << chid.get_compact()
217  << MSG::dec <<", gain " << gain << "." );
218  else
219  ATH_MSG_INFO ( "Replaced channel 0x" << MSG::hex << chid.get_compact()
220  << " by it's left FEB neighbor 0x" << chid_patch.get_compact() << MSG::dec );
221  return true;
222  }
223  }//end if patch connected and good
224  }//end if febChan
225  if (febChan<(m_onlineHelper->channelInSlotMax(febId)-1)) {
226  chid_patch=m_onlineHelper->channel_Id(febId, febChan+1);
227  if (cabling->isOnlineConnected(chid_patch) && bcCont->status(chid_patch).good()) {
228  const LArCondObj patch=m_contIn->get(chid_patch,gain); //Should be the const-get method
229  if (!patch.isEmpty()) {
230  StatusCode sc=m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel);
231  if (sc.isFailure()) {
232  ATH_MSG_ERROR ( "Failed to insert correction for channel 0x" << MSG::hex << chid.get_compact()
233  << MSG::dec << ", gain " << gain << "." );
234  }
235  else {
236  ATH_MSG_INFO ( "Replaced channel 0x" << MSG::hex << chid.get_compact()
237  << " by it's right FEB neighbor 0x" << chid_patch.get_compact() << MSG::dec );
238  return true;
239  }
240  }//end if patch connected and good
241  }//end if chan<max
242  }
243  ATH_MSG_ERROR ( "None of the FEB neighbors is good!" );
244  return false;
245  } else if (m_patchMethod==PhiNeighbor) {
246  // (*m_log) << MSG::ERROR << "Patching Method 'Phi-neighbor' not yet implemented." << endmsg;
247  try {
248  const Identifier id=cabling->cnvToIdentifier(chid);
249  int eta, phi, phi_min, phi_max, phi_range;
250  Identifier regionID;
251  regionID=m_caloId->region_id(id);
252  eta = m_caloId->eta(id);
253  phi = m_caloId->phi(id);
254  phi_min =m_caloId->phi_min(regionID);
255  phi_max =m_caloId->phi_max(regionID);
256  phi_range=phi_max-phi_min+1;
257  if (eta==CaloCell_ID::NOT_VALID || phi==CaloCell_ID::NOT_VALID || phi_min==CaloCell_ID::NOT_VALID || phi_max==CaloCell_ID::NOT_VALID) {
258  ATH_MSG_ERROR ( "CaloCell_ID returned NOT_VALID for offline id 0x"<< std::hex << id.get_compact()
259  <<", online id 0x" << chid.get_compact() << std::dec );
260  return false;
261  }
262  ATH_MSG_VERBOSE ( "Problem channel has phi="<< phi << " eta=" << eta );
263  //Try both phi-neighbors
264  int phi_list[4]={phi-1,phi+1,phi-2,phi+2};
265  for (unsigned i=0;i<4;i++) {
266  int phi_patch=phi_list[i];
267  if (phi_patch<m_caloId->phi_min(regionID)) phi_patch=phi_patch+phi_range;
268  if (phi_patch>m_caloId->phi_max(regionID)) phi_patch=phi_patch-phi_range+phi_min;
269  ATH_MSG_VERBOSE ( "Iteration " << i << " Using cell with phi="
270  << phi_patch << " eta=" << eta );
271  //std::cout << "i=" << i << " Using cell with phi="
272  // << phi_patch << " eta=" << eta << std::endl;
273  Identifier patch_id=m_caloId->cell_id(regionID,eta,phi_patch);
274  HWIdentifier chid_patch=cabling->createSignalChannelID(patch_id);
275  if (bcCont->status(chid_patch).good()) {
276  const LArCondObj patch=m_contIn->get(chid_patch,gain);
277  if (!patch.isEmpty()) {
278  StatusCode sc=m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel);
279  if (sc.isFailure())
280  ATH_MSG_ERROR ( "Failed to insert correction for channel 0x" << MSG::hex << chid.get_compact()
281  << MSG::dec << ", gain " << gain << "." );
282  else
283  ATH_MSG_INFO ( "Replaced channel 0x" << MSG::hex << chid.get_compact()
284  << " by neighbor 0x" << chid_patch.get_compact() << MSG::dec << " (phi " << phi << " to " << phi_patch <<")" );
285  return true;
286  }//end if isEmtpy()
287  }//end if neighbor is good
288  }// end for loop
289  ATH_MSG_ERROR ( "All phi-neighbors of channel 0x" << MSG::hex << chid.get_compact() << MSG::dec
290  << " are either absent or bad." );
291  return false;
292  }catch(LArID_Exception& except) {
293  ATH_MSG_ERROR ( "LArID_Exception caught!" );
294  return false;
295  }
296  }
297  else if (m_patchMethod==PhiAverage) {
298  LArCondObj patch;
299  if (!getAverage(chid,gain,patch,bcCont,cabling)){
300  ATH_MSG_ERROR ( "Failed get phi-average!" );
301  return false;
302  } else {
303  ATH_MSG_DEBUG ( "Got a phi-average..." );
304  }
305  ATH_CHECK( m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel), false );
306  return true;
307  }
308  else if (m_patchMethod==FEBAverage) {
309  LArCondObj patch;
310  if (!getAverage(chid,gain,patch,bcCont, cabling, false)){
311  ATH_MSG_ERROR ( "Failed get FEB-average!" );
312  return false;
313  } else {
314  ATH_MSG_DEBUG ( "Got a FEB-average..." );
315  }
316  ATH_CHECK( m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel), false );
317  return true;
318  }
319  else if (m_patchMethod==SetZero) {
320  LArCondObj patch;
321  if (!setZero(chid, gain, patch)){
322  ATH_MSG_ERROR ( "Failed set Zero!" );
323  return false;
324  } else {
325  ATH_MSG_DEBUG ( "Set zero ..." );
326  }
327  ATH_CHECK( m_contOut->insertCorrection(chid,patch,gain,m_useCorrChannel), false );
328  return true;
329  }
330  else //failed...
331  ATH_MSG_ERROR ( "Unknown correction method." );
332  return false;
333 }
334 
335 
336 template<class CONDITIONSCONTAINER>
337 std::vector<HWIdentifier>& LArCalibPatchingAlg<CONDITIONSCONTAINER>::getPhiRing(const HWIdentifier chid, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, unsigned distance) {
338  if (distance==0) distance=1;
339  m_idList.clear();
340  try {
341  const Identifier id=cabling->cnvToIdentifier(chid);
342  int eta,phi, phi_min, phi_max;
343  Identifier regionID;
344  regionID=m_caloId->region_id(id);
345  eta = m_caloId->eta(id);
346  phi = m_caloId->phi(id);
347  phi_min =m_caloId->phi_min(regionID);
348  phi_max =m_caloId->phi_max(regionID);
349  ATH_MSG_VERBOSE ( "Assembling phi-ring for eta=" << eta << " phi=" << phi );
350  //std::cout << "Assembling phi-ring for eta=" << eta << " phi=" << phi << std::endl;
351  if (eta==CaloCell_ID::NOT_VALID || phi==CaloCell_ID::NOT_VALID || phi_min==CaloCell_ID::NOT_VALID || phi_max==CaloCell_ID::NOT_VALID) {
352  ATH_MSG_ERROR ( "CaloCell_ID returned NOT_VALID for offline id 0x"<< std::hex << id.get_compact()
353  <<", online id 0x" << chid.get_compact() << std::dec );
354  return m_idList;
355  }
356 
357  if ((phi_max-phi_min)%distance) {
358  ATH_MSG_ERROR ( "Can't divide " << (phi_min-phi_max) << " by " << distance );
359  return m_idList;
360  }
361 
362  int nSteps=(phi_max-phi_min)/distance;
363 
364  for (int i=1;i<=nSteps;i++) {
365  int phi_patch=phi+i*distance;
366  if (phi_patch>phi_max) phi_patch=phi_patch-phi_max+phi_min-1;
367  ATH_MSG_VERBOSE ( "i=" << i << " Adding cell with phi="
368  << phi_patch << " eta=" << eta );
369  //std::cout << "i=" << i << " Adding cell with phi="
370  // << phi_patch << " eta=" << eta << std::endl;
371  Identifier patch_id=m_caloId->cell_id(regionID,eta,phi_patch);
372  HWIdentifier chid_patch=cabling->createSignalChannelID(patch_id);
373  if (bcCont->status(chid_patch).good()) {
374  m_idList.push_back(chid_patch);
375  }
376  else
377  ATH_MSG_VERBOSE ( "This cell is bad as well. Ignored." );
378  }//end loop over phi-steps
379 
380  }catch(LArID_Exception& except) {
381  ATH_MSG_ERROR ( "LArID_Exception caught!" );
382  }
383  return m_idList;
384 }
385 
386 template<class CONDITIONSCONTAINER>
387 std::vector<HWIdentifier>& LArCalibPatchingAlg<CONDITIONSCONTAINER>::getFEBChans(const HWIdentifier chid, const LArBadChannelCont* bcCont) {
388  m_idList.clear();
389  HWIdentifier febid = m_onlineHelper->feb_Id(chid);
390  ATH_MSG_VERBOSE ( "Assembling list of channels for FEB=" << febid );
391  for (int i=0;i<128;++i) {
392  HWIdentifier fchan = m_onlineHelper->channel_Id(febid,i);
393  if(fchan == chid) continue;
394  ATH_MSG_VERBOSE ( " Adding channel =" << i );
395  if (bcCont->status(fchan).good()) {
396  m_idList.push_back(fchan);
397  } else {
398  ATH_MSG_VERBOSE ( "This channel is bad as well. Ignored." );
399  }
400  }//end loop over chans
401 
402  return m_idList;
403 }
404 
405 
406 template<class CONDITIONSCONTAINER>
407 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArRampP1& patch, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, bool isphi) {
408 
409  std::vector<HWIdentifier>& symCells =
410  isphi ? getPhiRing(chid, bcCont, cabling) : getFEBChans(chid, bcCont);
411  if (symCells.empty()) {
412  ATH_MSG_ERROR ( "No symmetry cells found!" );
413  return false;
414  }
415  size_t s=m_contIn->get(symCells[0],gain).m_vRamp.size();
416  patch.m_vRamp.clear();
417  patch.m_vRamp.resize(s);
418  unsigned nCells=0;
419  for (HWIdentifier hwid : symCells) {
420  const LArRampP1& ramp=m_contIn->get(hwid,gain);
421  if (ramp.m_vRamp.size()==0) continue; //This one is empty...
422  if (ramp.m_vRamp.size()!=s) {
423  if(isphi) {
424  ATH_MSG_WARNING ("Cell with same phi but different size of ramp polynom found!" );
425  } else {
426  ATH_MSG_WARNING ("Cell with same FEB but different size of ramp polynom found!" );
427  }
428  continue;
429  }
430  msg() << MSG::DEBUG << "Adding cell 0x"<< std::hex << hwid.get_compact() << std::dec << " Ramp:";
431  for (size_t i=0;i<s;i++) {
432  patch.m_vRamp[i]+=ramp.m_vRamp[i];
433  msg() << MSG::DEBUG << ramp.m_vRamp[i] << " ";
434  }
435  msg() << MSG::DEBUG << endmsg;
436  nCells++;
437  }
438  if (nCells==0) {
439  if(isphi) {
440  ATH_MSG_ERROR ( "No good ramp with same phi found!" );
441  } else {
442  ATH_MSG_ERROR ( "No good ramp with same FEB found!" );
443  }
444  return false;
445  }
446  for (size_t i=0;i<s;i++)
447  patch.m_vRamp[i]=patch.m_vRamp[i]/nCells;
448 
449  //FIXME: We should somehow watch the rms....
450  msg() << MSG::INFO << "Patched Ramp (based on " << nCells << " channels):" ;
451  for (size_t i=0;i<s;i++)
452  msg() << MSG::INFO << " " << patch.m_vRamp[i];
453  msg() << MSG::INFO << endmsg;
454  return true;
455 }
456 
457 
458 template<class CONDITIONSCONTAINER>
459 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArOFCP1& patch, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, bool /*isphi*/) {
460 
461  std::vector<HWIdentifier>& symCells=getPhiRing(chid, bcCont, cabling);
462 
463  if (symCells.empty()) {
464  ATH_MSG_ERROR ( "No symmetry cells found!" );
465  return false;
466  }
467 
468  const size_t nPhases=m_contIn->get(symCells[0],gain).m_vOFC_a.size();
469  if (!nPhases) {
470  ATH_MSG_ERROR ( "OFC of neighbor nPhase=0!" );
471  return false;
472  }
473  const size_t nSamples=m_contIn->get(symCells[0],gain).mvOFC_a[0].size();
474  if (!nSamples) {
475  ATH_MSG_ERROR ( "OFC of neighbor nSamples=0!" );
476  return false;
477  }
478 
479  float timeOffset = m_contIn->get(symCells[0],gain).m_timeOffset;
480  float timeBinWidth = m_contIn->get(symCells[0],gain).m_timeBinWidth;
481  std::vector<std::vector<float> > ofc_a;
482  std::vector<std::vector<float> > ofc_b;
483  ofc_a.resize(nPhases,std::vector<float>(nSamples));
484  ofc_b.resize(nPhases,std::vector<float>(nSamples));
485 
486  unsigned nCells=0;
487  for (HWIdentifier hwid : symCells) {
488  const LArOFCP1& ofc=m_contIn->get(hwid,gain);
489  if (ofc.OFC_aSize()==0 || ofc.OFC_bSize()==0) continue; //This one is empty...
490  if (ofc.OFC_aSize()!=nPhases) {
491  ATH_MSG_WARNING ("Cell with same phi but different nPhases found! Ignored" );
492  continue;
493  }
494  if (ofc.timeOffset()!=timeOffset) {
495  ATH_MSG_WARNING ("Cell with same phi but different time-offset found! Ignored" );
496  continue;
497  }
498  if (ofc.timeBinWidth()!=timeBinWidth) {
499  ATH_MSG_WARNING ("Cell with same phi but different time-offset found! Ignored" );
500  continue;
501  }
502 
503  ATH_MSG_DEBUG ( "Adding cell 0x"<< std::hex << hwid.get_compact() << std::dec );
504  for (size_t iPhase=0;iPhase<nPhases;++iPhase) {
505  //Check size of vector?
506  for (size_t iSample=0;iSample<nSamples;++iSample) {
507  ofc_a[iPhase][iSample]+=ofc.OFC_a(iPhase)[iSample];
508  ofc_b[iPhase][iSample]+=ofc.OFC_b(iPhase)[iSample];
509  }
510  }
511  nCells++;
512  }
513  if (nCells==0) {
514  ATH_MSG_ERROR ( "No good OFC set with same phi found!" );
515  return false;
516  }
517 
518  for (size_t iPhase=0;iPhase<nPhases;++iPhase) {
519  //Check size of vector?
520  for (size_t iSample=0;iSample<nSamples;++iSample) {
521  ofc_a[iPhase][iSample]/=nCells;
522  ofc_b[iPhase][iSample]/=nCells;
523  }
524  }
525  //FIXME: We should somehow watch the rms....
526 
527  LArOFCP1 tmp (timeOffset,
528  timeBinWidth,
529  ofc_a, ofc_b);
530  patch.setFrom (tmp);
531 
532  return true;
533 }
534 
535 template<class CONDITIONSCONTAINER>
536 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArCaliWaveVec& patch,const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, bool /*isphi*/) {
537 
538  const std::vector<HWIdentifier>& symCells=getPhiRing(chid, bcCont, cabling);
539 
540 
541  if (symCells.empty()) {
542  ATH_MSG_ERROR ( "No symmetry cells found!" );
543  return false;
544  }
545  patch.clear();
546 
547  struct perDAC_t {
548  std::vector<std::tuple<HWIdentifier,double,const LArCaliWave*> > tmax_wave; //Identifier, tmax and ptr to wave, all neighbors
549  double tmaxAvg=0; //average peak time of all neighbors
550  };
551 
552  LArWaveHelper wHelper;
553 
554  // Organize waves per DAC and calculate tmax (peaking time)
555  std::map<int,perDAC_t> neighbors_per_dac; //In most cases, there is only one DAC, eg one entry in this map
556 
557  ATH_MSG_DEBUG(" symCells size "<<symCells.size());
558  for (HWIdentifier hwid : symCells) {
559  for (const LArCaliWave& cwave : m_contIn->get(hwid,gain)) {
560  ATH_MSG_VERBOSE("found wave in channel " << m_onlineHelper->channel_name(hwid) << ", gain " << gain);
561  const int ourDAC = cwave.getDAC();
562  double tmax=-1;
563  wHelper.getDMax(cwave, tmax);
564  if (tmax > 0 ) {
565  auto& dacWave=neighbors_per_dac[ourDAC];
566  dacWave.tmax_wave.emplace_back(hwid,tmax,&cwave);
567  dacWave.tmaxAvg+=tmax;
568  }
569  else {
570  ATH_MSG_WARNING("Ignoring wave with peak-time=" << tmax << ", DAC=" << ourDAC << " found in channel " << m_onlineHelper->channel_name(hwid) << ", gain " << gain);
571  }
572  }//end loop over DACs
573  }// end loop over sym-cells
574 
575  ATH_MSG_DEBUG(" neighbors_per_dac size "<<neighbors_per_dac.size());
576 
577  for (auto& dacWave : neighbors_per_dac) {
578  dacWave.second.tmaxAvg/=dacWave.second.tmax_wave.size();
579  ATH_MSG_DEBUG("Average tmax computed " << dacWave.second.tmaxAvg << " for DAC=" << dacWave.first << ", to patch channel " << m_onlineHelper->channel_name(chid) << ", gain " << gain);
580  }
581 
582  //Align waves accoring to their peak-time and sum them
583  std::vector<std::pair<LArCaliWave,unsigned> > alignedWaveSums;
584 
585  for (const auto& dacWave : neighbors_per_dac) {//loop over DAC values
586  const int& dac=dacWave.first;
587  const perDAC_t& neighbors=dacWave.second;
588  const double& tmaxAvg=neighbors.tmaxAvg;
589  bool first=true;
590 
591  for (const auto& [hwid, tmax, wave] : neighbors.tmax_wave) {
592  if (std::fabs(tmaxAvg-tmax) > wave->getWave().size()*wave->getDt()/2.0) {//Don't allow very long time shifts
593  ATH_MSG_WARNING("Peaking time of symmetric channel " << m_onlineHelper->channel_name(hwid)
594  << ", DAC=" << dac << ", gain=" << gain << " is too far off the average (" << std::fabs(tmaxAvg-tmax) << "). Ignoring.");
595  continue;
596  }
597  if (first) {
598  //First (useable) wave for this DAC value. Initialize a LArCaliWave object for it
599  alignedWaveSums.emplace_back(LArCaliWave(wHelper.Dtranslate(*wave, tmaxAvg - tmax).getWave(),wave->getDt(),wave->getDAC(),wave->getIsPulsedInt(),wave->getFlag()),1);
600  ATH_MSG_DEBUG("Adding wave of symmetric channel " << m_onlineHelper->channel_name(hwid)
601  << ", DAC=" << dac << ", gain=" << gain << " shifted by " << tmaxAvg-tmax);
602  first=false;
603  }
604  else {
605  //Subsequent waves fro this DAC value: Add them.
606  alignedWaveSums.back().first+=wHelper.Dtranslate(*wave, tmaxAvg - tmax);
607  alignedWaveSums.back().second++;
608  ATH_MSG_DEBUG("Summing wave of symmetric channel " << m_onlineHelper->channel_name(hwid)
609  << ", DAC=" << dac << ", gain=" << gain << " shifted by " << tmaxAvg-tmax);
610  }
611  }//end loop over symmetic waves for one dac value
612  if (first) {
613  ATH_MSG_ERROR("Cannot calculate patch wave for channel " << m_onlineHelper->channel_name(chid) << ", gain="<< gain << ", DAC=" << dac << ", no usable neighbors found");
614  return false;
615  }
616 
617  }//end loop over DAC values
618 
619  //Final loop to divide wave by number-of-waves to get the average, then push-back into return container
620  for (auto& waveSum : alignedWaveSums) { //Loop over DAC values
621  LArCaliWave& wave=waveSum.first;
622  wave*=(1./waveSum.second);
623  ATH_MSG_INFO("Calculated average wave for channel " << m_onlineHelper->channel_name(chid) << ", gain="<< gain << ", DAC=" << waveSum.first.getDAC()
624  << ", based on " << waveSum.second << " neighbors");
625  patch.push_back(wave);
626  }
627  return true;
628 }
629 
630 template<class CONDITIONSCONTAINER>
631 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArAutoCorrP1& patch, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* /*cabling*/, bool isphi) {
632 
633  ATH_MSG_DEBUG ( "Starting getAverage for LArAutoCorr" );
634  if(isphi) {
635  ATH_MSG_ERROR ( "No phi-average for AutoCorr!" );
636  return false;
637  }
638  std::vector<HWIdentifier>& symCells = getFEBChans(chid, bcCont);
639  if (symCells.empty()) {
640  ATH_MSG_ERROR ( "No symmetry cells found!" );
641  return false;
642  }
643  size_t s=m_contIn->get(symCells[0],gain).m_vAutoCorr.size();
644  patch.m_vAutoCorr.clear();
645  patch.m_vAutoCorr.resize(s);
646  unsigned nCells=0;
647  for (HWIdentifier hwid : symCells) {
648  const LArAutoCorrP1& ac=m_contIn->get(hwid,gain);
649  if (ac.m_vAutoCorr.size()==0) continue; //This one is empty...
650  if (ac.m_vAutoCorr.size()!=s) {
651  ATH_MSG_WARNING ("Cell with same FEB but different size of autocorr found!" );
652  continue;
653  }
654  msg() << MSG::DEBUG << "Adding cell 0x"<< std::hex << hwid.get_compact() << std::dec << " AC:";
655  for (size_t i=0;i<s;i++) {
656  patch.m_vAutoCorr[i]+=ac.m_vAutoCorr[i];
657  msg() << MSG::DEBUG << ac.m_vAutoCorr[i] << " ";
658  }
659  msg() << MSG::DEBUG << endmsg;
660  nCells++;
661  }
662  if (nCells==0) {
663  ATH_MSG_ERROR ( "No good autocorr with same FEB found!" );
664  return false;
665  }
666  for (size_t i=0;i<s;i++)
667  patch.m_vAutoCorr[i]=patch.m_vAutoCorr[i]/nCells;
668 
669  //FIXME: We should somehow watch the rms....
670  msg() << MSG::INFO << "Patched autocorr (based on " << nCells << " channels):" ;
671  for (size_t i=0;i<s;i++)
672  msg() << MSG::INFO << " " << patch.m_vAutoCorr[i];
673  msg() << MSG::INFO << endmsg;
674  return true;
675 }
676 
677 #ifdef LARRAWCONDITIONS_LARMPHYSOVERMCALP
678 template<class CONDITIONSCONTAINER>
679 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArMphysOverMcalP1& patch, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, bool isphi) {
680 
681  std::vector<HWIdentifier>& symCells=getPhiRing(chid, bcCont, cabling);
682 
683  patch.m_MphysOverMcal=0;
684 
685  unsigned nCells=0;
686  for (HWIdentifier hwid : symCells) {
687  const float mPmC=m_contIn->get(hwid,gain);
688  if (mPmC>0) {
689  patch.m_MphysOverMcal+=mPmC;
690  nCells++;
691  }
692  }
693 
694  if (nCells==0) {
695  ATH_MSG_ERROR ( "No good symmetry cells found!" );
696  return false;
697  }
698  patch.m_MphysOverMcal/=nCells;
699  return true;
700 }
701 #endif
702 
703 
704 
705 #ifdef LARRAWCONDITIONS_LARSINGLEFLOATP
706 template<class CONDITIONSCONTAINER>
707 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::getAverage(const HWIdentifier chid, const int gain, LArSingleFloatP& patch, const LArBadChannelCont* bcCont, const LArOnOffIdMapping* cabling, bool /*isphi*/) {
708 
709  std::vector<HWIdentifier>& symCells=getPhiRing(chid, bcCont, cabling);
710 
711  patch.m_data=0;
712 
713  unsigned nCells=0;
714  for (HWIdentifier hwid : symCells) {
715  const LArSingleFloatP& sf=m_contIn->get(hwid,gain);
716  if (!sf.isEmpty()) {
717  patch.m_data+=sf.m_data;
718  nCells++;
719  }
720  }
721 
722  if (nCells==0) {
723  ATH_MSG_ERROR ( "No good symmetry cells found!" );
724  return false;
725  }
726  patch.m_data/=nCells;
727  return true;
728 }
729 #endif
730 
731 template<class CONDITIONSCONTAINER>
732 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::setZero(const HWIdentifier chid, const int gain, LArRampP1& patch) {
733 
734  size_t s=m_contIn->get(chid,gain).m_vRamp.size();
735  patch.m_vRamp.clear();
736  patch.m_vRamp.resize(s);
737  for (size_t i=0; i<s; ++i) patch.m_vRamp[i] = 0.;
738  return true;
739 }
740 
741 template<class CONDITIONSCONTAINER>
742 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::setZero(const HWIdentifier chid, const int gain, LArOFCP1& patch) {
743 
744  const size_t nPhases=m_contIn->get(chid,gain).m_vOFC_a.size();
745  if (!nPhases) {
746  ATH_MSG_ERROR ( "OFC of nPhase=0 !" );
747  return false;
748  }
749  const size_t nSamples=m_contIn->get(chid,gain).mvOFC_a[0].size();
750  if (!nSamples) {
751  ATH_MSG_ERROR ( "OFC of nSamples=0 !" );
752  return false;
753  }
754 
755  LArOFCP1 newvals (patch.timeOffset(), patch.timeBinWidth(),
756  std::vector<std::vector<float> > (nPhases, std::vector<float> (nSamples, 0)),
757  std::vector<std::vector<float> > (nPhases, std::vector<float> (nSamples, 0)));
758  patch.setFrom (newvals);
759  return true;
760 }
761 
762 template<class CONDITIONSCONTAINER>
763 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::setZero(const HWIdentifier /*chid*/, const int /*gain*/, LArCaliWaveVec& /*patch*/) {
764  ATH_MSG_ERROR ( "Not implemented, should not come here !!!");
765  return false;
766 }
767 
768 template<class CONDITIONSCONTAINER>
769 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::setZero(const HWIdentifier chid, const int gain, LArAutoCorrP1& patch) {
770 
771  size_t s=m_contIn->get(chid,gain).m_vAutoCorr.size();
772  patch.m_vAutoCorr.clear();
773  patch.m_vAutoCorr.resize(s);
774  for (size_t i=0; i<s; ++i) patch.m_vAutoCorr[i] = 0.;
775  return true;
776 }
777 
778 #ifdef LARRAWCONDITIONS_LARSINGLEFLOATP
779 template<class CONDITIONSCONTAINER>
780 bool LArCalibPatchingAlg<CONDITIONSCONTAINER>::setZero(LArSingleFloatP& patch) {
781  patch.m_data=0.;
782  return true;
783 }
784 #endif
785 
786 
787 template<class CONDITIONSCONTAINER>
788 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::setSymlink(const LArRampComplete* ramp) const {
789  ATH_CHECK( detStore()->symLink(ramp, (ILArRamp*)ramp) );
790  return StatusCode::SUCCESS;
791 }
792 
793 
794 template<class CONDITIONSCONTAINER>
795 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::setSymlink(const LArOFCComplete* ofc) const {
796  ATH_CHECK( detStore()->symLink(ofc, (ILArOFC*)ofc) );
797  return StatusCode::SUCCESS;
798 }
799 
800 
801 template<class CONDITIONSCONTAINER>
802 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::setSymlink(const LArMphysOverMcalComplete* ramp) const {
803  ATH_CHECK( detStore()->symLink(ramp, (ILArMphysOverMcal*)ramp) );
804  return StatusCode::SUCCESS;
805 }
806 
807 template<class CONDITIONSCONTAINER>
808 StatusCode LArCalibPatchingAlg<CONDITIONSCONTAINER>::setSymlink(const LArAutoCorrComplete* ac) const {
809  ATH_CHECK( detStore()->symLink(ac, (ILArAutoCorr*)ac) );
810  return StatusCode::SUCCESS;
811 }
812 
813