ATLAS Offline Software
Loading...
Searching...
No Matches
LArRODMonAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5// ********************************************************************
6//
7// NAME: LArRODMonAlg.cxx
8// PACKAGE: LArMonitoring
9//
10// ********************************************************************
11
12#include "LArRODMonAlg.h"
13
19
22
24
26#include "GaudiKernel/ConcurrencyFlags.h"
27#include <cmath>
28
29
30/*---------------------------------------------------------*/
32
33/*---------------------------------------------------------*/
34StatusCode
36 ATH_MSG_VERBOSE( "In LArRODMonAlg::initialize() ");
37
38 ATH_CHECK(detStore()->retrieve(m_LArOnlineIDHelper, "LArOnlineID"));
39
42 ATH_CHECK(m_digitContainerKey.initialize());
43 ATH_CHECK(m_headerContainerKey.initialize());
44
45 ATH_CHECK(m_keyOFC.initialize());
46 ATH_CHECK(m_keyShape.initialize());
47 ATH_CHECK(m_keyHVScaleCorr.initialize());
48 ATH_CHECK(m_keyPedestal.initialize());
49
50 ATH_CHECK(m_adc2mevKey.initialize());
51
52 ATH_CHECK(m_noiseCDOKey.initialize());
53 ATH_CHECK(m_cablingKey.initialize());
54
55 ATH_CHECK(m_bcContKey.initialize());
56 ATH_CHECK(m_bcMask.buildBitMask(m_problemsToMask,msg()));
57
58 ATH_CHECK( m_eventInfoKey.initialize() );
59
60
61 //Check properties ...
62 if (Gaudi::Concurrency::ConcurrencyFlags::numThreads() > 1) { //MT environment
63 if (m_doDspTestDump) {
64 ATH_MSG_ERROR("Property 'DoDspTestDump' must not be true if nThreads>1");
65 return StatusCode::FAILURE;
66 }
67 if (m_doCellsDump) {
68 ATH_MSG_ERROR("Property 'DoCellsDump' must not be true if nThreads>1");
69 return StatusCode::FAILURE;
70 }
71 }
72
73
74 auto pairCmp = [](const std::pair<int,int>& p1, const std::pair<int,int>& p2) {return (p1.first<p2.first);};
75 if (!std::is_sorted(m_E_precision.begin(),m_E_precision.end(),pairCmp)) {
76 ATH_MSG_ERROR("Configuration problem: Energy ranges not in ascending order!");
77 return StatusCode::FAILURE;
78 }
79
80 if (!std::is_sorted(m_T_precision.begin(),m_T_precision.end(),pairCmp)) {
81 ATH_MSG_ERROR("Configuration problem: Time ranges not in ascending order!");
82 return StatusCode::FAILURE;
83 }
84
85 if (!std::is_sorted(m_Q_precision.begin(),m_Q_precision.end(),pairCmp)) {
86 ATH_MSG_ERROR("Configuration problem: Quality ranges not in ascending order!");
87 return StatusCode::FAILURE;
88 }
89
90
91 // Open output files for DspTest
92 if (m_doDspTestDump) {
93 m_fai.open(m_AiFileName,std::ios::trunc);
94 m_fdig.open(m_DigitsFileName,std::ios::trunc);
95 m_fen.open(m_EnergyFileName,std::ios::trunc);
96 }
97
98 // Output file
99 if (m_doCellsDump) {
100 m_fdump.open(m_DumpCellsFileName);
101 m_fdump<<"febid channel CellID slot FT barrel_ec posneg partition E_off E_on T_off T_on Q_off Q_on event "<<std::endl;
102 }
103
104
105 ATH_MSG_DEBUG("Setting an offset time of " << m_timeOffset << " BC, i.e. " << m_timeOffset*m_BC << " ps");
106
108
110}
111
112StatusCode LArRODMonAlg::fillHistograms(const EventContext& ctx) const {
113 ATH_MSG_VERBOSE( "In LArRODMonAlg::fillHistograms()");
114
115 //monitored variables
116 auto weight_e = Monitored::Scalar<float>("weight_e",1.);
117 auto weight_q = Monitored::Scalar<float>("weight_q",1.);
118 auto weight_t = Monitored::Scalar<float>("weight_t",1.);
119 auto numE = Monitored::Scalar<int>("numE",1.);
120 auto numQ = Monitored::Scalar<int>("numQ",1.);
121 auto numT = Monitored::Scalar<int>("numT",1.);
122 auto gain = Monitored::Scalar<int>("gain",-1);
123 auto partition = Monitored::Scalar<int>("partition",-1);
124 auto partitionI = Monitored::Scalar<int>("partitionI",-1);
125 auto lb = Monitored::Scalar<int>("LBN",0);
126 auto sweetc = Monitored::Scalar<float>("Sweetc",1.);
127
129
131 const ILArPedestal* pedestals=*pedestalHdl;
132
133 //retrieve BadChannel info:
135 const LArBadChannelCont* bcCont{*bcContHdl};
136
137 bool isEventFlaggedByLArNoisyROAlg = false; // keep default as false
138 bool isEventFlaggedByLArNoisyROAlgInTimeW = false; // keep deault as false
139
140 if ( thisEventInfo->isEventFlagBitSet(xAOD::EventInfo::LAr,0) ) {
141 isEventFlaggedByLArNoisyROAlg = true;
142 ATH_MSG_DEBUG( "Event flagged as Noisy" );
143 }
144
145 if ( thisEventInfo->isEventFlagBitSet(xAOD::EventInfo::LAr,3) ) {
146 isEventFlaggedByLArNoisyROAlgInTimeW = true;
147 ATH_MSG_DEBUG( " !!! Noisy event found by LArNoisyROAlg in Time window of 500ms!!!" );
148 }
149
150 // Noise bursts cleaning (LArNoisyRO_Std or TimeWindowVeto) added by B.Trocme - 19/7/12
151 if (m_removeNoiseBursts && (isEventFlaggedByLArNoisyROAlg || isEventFlaggedByLArNoisyROAlgInTimeW)) return StatusCode::SUCCESS;
152
153 // Retrieve stream info
154 int nStreams = m_streams.size();
155 // if ((nStreams == 1 && m_streams[0] == "all") || nStreams <= 0) selectstreams = false;
156
157 lb = thisEventInfo->lumiBlock();
158
159 //Fixme: Use LArTrigStreamMatching also here.
160 const int streamsize = nStreams + 1;
161 std::vector<int> hasStream(streamsize,0);
162 // for (int str = 0; str < nStreams + 1; str++) hasStream[str] = 0;
163
164 bool hasstrlist = false;
165 const std::vector< xAOD::EventInfo::StreamTag >& evtStreamTags=thisEventInfo->streamTags();
166 for (const auto& evtStreamTag : evtStreamTags) {
167 std::vector<std::string>::const_iterator evtStreamTagIt=std::find(m_streams.begin(),m_streams.end(),evtStreamTag.name());
168 if (evtStreamTagIt!=m_streams.end()) {
169 const unsigned str=evtStreamTagIt-m_streams.begin();
170 ATH_MSG_VERBOSE( "Keeping Stream Tag: " << evtStreamTag.type() << "_" << evtStreamTag.name());
171 hasStream[str] = 1;
172 hasstrlist = true;
173 }
174 }
175 if (! hasstrlist) hasStream[nStreams] = 1;
176
177
179
181
183
184 std::set<HWIdentifier> ignoreFEBs;
185
188 if (!febCont) {
189 ATH_MSG_WARNING( "No LArFEB container found in TDS" );
190 } else {
191 for (const auto* febH : *febCont) {
192 if (((m_doCheckSum && febH->ChecksumVerification()==false)) ||
193 (m_doRodStatus && febH->RodStatus()!=0))
194 ignoreFEBs.insert(febH->FEBId());
195 }
196 }
197 ATH_MSG_DEBUG("Found " << ignoreFEBs.size() << " FEBs with checksum errors or statatus errors. Will ignore these FEBs.");
198 }
199
200 std::vector<ERRCOUNTER> errcounters(N_PARTITIONS);
201 ERRCOUNTER allEC;
202
203 std::vector<unsigned> errsPerFEB;
204 errsPerFEB.resize(m_LArOnlineIDHelper->febHashMax(),0);
205
206 const bool ignoreFebs=(ignoreFEBs.size()>0);
207 std::set<HWIdentifier>::const_iterator ignoreFebsEnd=ignoreFEBs.end();
208
209 unsigned count_gain[(int)CaloGain::LARNGAIN] = {0};
210
211 //Build an association of channels in the two LArRawChannelContainers.
212 //The LArRawChannelContainers are unordered collections of LArRawChannels.
213 //But we know that they have the same order because they were built from the same source (namely the LArDigits and RawChannels in the Bytestream)
214 //and we know that the LArRawChannels built offline are a subset of the LArRawChannelContainers read from Bytestream.
215 //Therfore we can search much more efficiently
216 LArRawChannelContainer::const_iterator rcDigIt=rawColl_fromDigits->begin();
217 LArRawChannelContainer::const_iterator rcDigIt_e=rawColl_fromDigits->end();
218 LArRawChannelContainer::const_iterator rcBSIt=rawColl_fromBytestream->begin();
219 LArRawChannelContainer::const_iterator rcBSIt_e=rawColl_fromBytestream->end();
220
221 LArDigitContainer::const_iterator digIt=pLArDigitContainer->begin();
222 LArDigitContainer::const_iterator digIt_e=pLArDigitContainer->end();
223
224
225 //Loop over indices in LArRawChannelContainer built offline (the small one)
226 ATH_MSG_DEBUG( "Entering the LArRawChannel loop." );
227
228 for (;rcDigIt!=rcDigIt_e;++rcDigIt) {
229 const HWIdentifier idDig=rcDigIt->hardwareID();
230 const HWIdentifier febId=m_LArOnlineIDHelper->feb_Id(idDig);
231 // Check if this FEB should be ignored
232 if (ignoreFebs) {
233 if (ignoreFEBs.find(febId)!=ignoreFebsEnd) continue;
234 }
235 //Check if this is a bad channel
236 if (m_skipKnownProblematicChannels && m_bcMask.cellShouldBeMasked(bcCont,idDig)) continue;
237
238 const CaloGain::CaloGain gain = rcDigIt->gain();
239 //Check pedestal if needed
240 if (m_skipNullPed) {
241 const float ped = pedestals->pedestal(idDig,gain);
242 if(ped <= (1.0+LArElecCalib::ERRORCODE)) continue;
243 }
244
245 //Now look for corresponding channel in the LArRawChannelContainer read from Bytestream (the big one)
246 LArRawChannelContainer::const_iterator currIt=rcBSIt; //Remember current position in container
247 for (;rcBSIt!=rcBSIt_e && rcBSIt->hardwareID() != idDig; ++rcBSIt);
248 if (rcBSIt==rcBSIt_e) {
249 ATH_MSG_WARNING( "LArDigitContainer not in the expected order. Change of LArByteStream format?" );
250 //Wrap-around
251 for (rcBSIt=rawColl_fromBytestream->begin();rcBSIt!=currIt && rcBSIt->hardwareID() != idDig; ++rcBSIt);
252 if (rcBSIt==currIt) {
253 ATH_MSG_ERROR( "Channel " << m_LArOnlineIDHelper->channel_name(idDig) << " not found." );
254 return StatusCode::FAILURE;
255 }
256 }
257
258 //Now look for corresponding channel in the LArDigitContainer read from Bytestream
259 //Should be in almost in sync with the RawChannelContainer we are iterating over,
260 //but contains disconnected channels that are not part of the LArRawChannelContainer
261 LArDigitContainer::const_iterator currDigIt=digIt; //Remember current position in digit-container
262 for (;digIt!=digIt_e && (*digIt)->hardwareID() != idDig; ++digIt);
263 if (digIt==digIt_e) {
264 ATH_MSG_WARNING( "LArRawChannelContainer not in the expected order. Change of LArRawChannelBuilder behavior?" );
265 //Wrap-around
266 for (digIt=pLArDigitContainer->begin();digIt!=currDigIt && (*digIt)->hardwareID() != idDig; ++digIt);
267 if (digIt==currDigIt) {
268 ATH_MSG_ERROR( "Channel " << m_LArOnlineIDHelper->channel_name(idDig) << " not found in LArDigitContainer." );
269 return StatusCode::FAILURE;
270 }
271 }
272 const LArDigit* dig=*digIt;
273 const std::vector<short>& samples=dig->samples();
274 const auto [minSamplesIt, maxSamplesIt] = std::minmax_element(samples.begin(),samples.end());
275 if (m_adc_th<=0 || (*maxSamplesIt-*minSamplesIt)>m_adc_th) {
276 const diff_t compRes=compareChannel(*rcDigIt,*rcBSIt);
277 if (compRes.e_on!=compRes.e_off || compRes.t_on!=compRes.t_off || compRes.q_on!=compRes.q_off) {
279 detailedOutput(compRes,*dig,ctx);
280 ++m_ndump;
281 }
282 if (m_doCellsDump) {
283 dumpCellInfo(idDig,gain,ctx,compRes);
284 }
285 }
286 //Count errors:
287 const PARTITION p=getPartition(idDig);
288 if (compRes.e_on!=compRes.e_off) {
289 ++(errcounters[p].errors_E[gain]);
290 ++(allEC.errors_E[gain]);
291
292 IdentifierHash febHash=m_LArOnlineIDHelper->feb_Hash(febId);
293 ++(errsPerFEB[febHash]);
294 }
295 if (compRes.t_on!=compRes.t_off) {//Time-error
296 ++(errcounters[p].errors_T[gain]);
297 ++(allEC.errors_T[gain]);
298 }
299
300 if (compRes.q_on!=compRes.q_off) {//Quality-error
301 ++(errcounters[p].errors_Q[gain]);
302 ++(allEC.errors_Q[gain]);
303 }
304
305 }
306 else {
307 ATH_MSG_DEBUG( "Samples : "<< *maxSamplesIt << " " << *minSamplesIt );
308 }
309 }//end loop over rawColl_fromDigits
310
311 ATH_MSG_DEBUG( "End of rawChannels loop" );
312
313 for (unsigned i=0;i<m_LArOnlineIDHelper->febHashMax();++i) {
314 const HWIdentifier febid=m_LArOnlineIDHelper->feb_Id(i);
315 const std::string pn=getPartitionName(febid);
316 sweetc = errsPerFEB[i];
317 fill(m_tools[m_histoGroups.at(pn)],sweetc);
318 }
319
320
321 ATH_MSG_VERBOSE( "*Number of errors in Energy Computation : " );
322 ATH_MSG_VERBOSE( "* Low Gain : " << allEC.errors_E[2] << " / " << count_gain[CaloGain::LARLOWGAIN] );
323 ATH_MSG_VERBOSE( "* Medium Gain : " << allEC.errors_E[1] << " / " << count_gain[CaloGain::LARMEDIUMGAIN] );
324 ATH_MSG_VERBOSE( "* High Gain : " << allEC.errors_E[0] << " / " << count_gain[CaloGain::LARHIGHGAIN] );
325 ATH_MSG_VERBOSE( "*Number of errors in Time Computation : " );
326 ATH_MSG_VERBOSE( "* Low Gain : " << allEC.errors_T[2] << " / " << count_gain[CaloGain::LARLOWGAIN] );
327 ATH_MSG_VERBOSE( "* Medium Gain : " << allEC.errors_T[1] << " / " << count_gain[CaloGain::LARMEDIUMGAIN] );
328 ATH_MSG_VERBOSE( "* High Gain : " << allEC.errors_T[0] << " / " << count_gain[CaloGain::LARHIGHGAIN] );
329 ATH_MSG_VERBOSE( "*Number of errors in Quality Computation : " );
330 ATH_MSG_VERBOSE( "* Low Gain : " << allEC.errors_Q[2] << " / " << count_gain[CaloGain::LARLOWGAIN] );
331 ATH_MSG_VERBOSE( "* Medium Gain : " << allEC.errors_Q[1] << " / " << count_gain[CaloGain::LARMEDIUMGAIN] );
332 ATH_MSG_VERBOSE( "* High Gain : " << allEC.errors_Q[0] << " / " << count_gain[CaloGain::LARHIGHGAIN] );
333
334 for (unsigned p=0;p<N_PARTITIONS;++p) {
335 unsigned allErrsPartE=0;
336 unsigned allErrsPartT=0;
337 unsigned allErrsPartQ=0;
338 partition = p;
339 for (unsigned g=0;g<3;++g) {
340 gain = g;
341 weight_e = (float)errcounters[p].errors_E[g];
342 weight_q = (float)errcounters[p].errors_Q[g];
343 weight_t = (float)errcounters[p].errors_T[g];
344 fill(m_MonGroupName, partition, gain, weight_e, weight_q, weight_t);
345
346 allErrsPartE+=errcounters[p].errors_E[g];
347 allErrsPartT+=errcounters[p].errors_T[g];
348 allErrsPartQ+=errcounters[p].errors_Q[g];
349 }
350 partitionI = p;
351 numE = (float)allErrsPartE;
352 numT = (float)allErrsPartT;
353 numQ = (float)allErrsPartQ;
354 fill(m_MonGroupName, lb, partitionI, numE, numT, numQ);
355 }//end loop over partitions
356
357 /*
358 for(int str = 0; str < nStreams + 1; str++) {
359 if (hasStream[str] == 1) {
360
361 m_hEErrors_LB_stream->Fill((float)m_curr_lb,(float)str,(float)allErr_E);
362 m_hTErrors_LB_stream->Fill((float)m_curr_lb,(float)str,(float)allErr_T);
363 m_hQErrors_LB_stream->Fill((float)m_curr_lb,(float)str,(float)allErr_Q);
364 FIXME
365 }
366 }
367 */
368
369 return StatusCode::SUCCESS;
370}
371
372
373
375 const LArRawChannel& rcBS) const {
376 diff_t result;
377 const HWIdentifier chid=rcDig.channelID();
378 const int slot_fD = m_LArOnlineIDHelper->slot(chid);
379 const int feedthrough_fD = m_LArOnlineIDHelper->feedthrough(chid);
380 const float timeOffline = rcDig.time() - m_timeOffset*m_BC;
381
382 const float& en_fB=rcBS.energy();
383 const float& en_fD=rcDig.energy();
384
385 ATH_MSG_VERBOSE(chid.getString()<<" | "<<timeOffline<<" "<<rcBS.time()<<" | "<<en_fB<<" "<<en_fD);
386
387 if (fabs(timeOffline) > m_peakTime_cut*1000.){
388 ATH_MSG_DEBUG( " timeOffline too large " << timeOffline );
389 return result;
390 }
391
392 const std::string & hg=getPartitionName(chid);
393 auto slot = Monitored::Scalar<int>("slot",slot_fD);
394 auto ft = Monitored::Scalar<int>("FT",feedthrough_fD);
395
396 auto pairValueCmp = [](const int& a, const std::pair<int,int>& b){return a<b.first;};
397 //Energy check:
398 //Find expected precision given the cell-energy:
399 auto e_Precision=std::upper_bound(m_E_precision.begin(),m_E_precision.end(),std::abs(en_fB),pairValueCmp);
400 const size_t energyrange=e_Precision-m_E_precision.begin();
401 auto erange = Monitored::Scalar<int>("Erange",energyrange);
402 auto DiffE = Monitored::Scalar<float>("Ediff",en_fD - en_fB);
403 fill(m_MonGroupName,DiffE,erange);
404
405 auto Eon = Monitored::Scalar<float>("Eon",en_fB);
406 auto Eoff = Monitored::Scalar<float>("Eoff",en_fD);
407 fill(m_tools[m_histoGroups.at(hg)],DiffE,erange,Eon,Eoff);
408
409 if (std::abs(en_fD-en_fB) > e_Precision->second) {
410 //Fill results object for error counting, and dumping (if needed)compRes.e_on!=compRes.e_off)
411 result.e_off=en_fD;
412 result.e_on=en_fB;
413 auto weight_e = Monitored::Scalar<float>("weight_e",1.);
414 fill(m_tools[m_histoGroups.at(hg)],slot,ft,weight_e);
415 }
416
417 const float q_fB=rcBS.quality();
418 const float q_fD=rcDig.quality();
419 const float t_fB=rcBS.time();
420
421 if ((rcDig.provenance() & 0x2000) == 0 || q_fD==0 || t_fB==0 || q_fB==0 || timeOffline==0) {
422 // Skip time/Quality comparison if
423 // * provenance bits indicate the LArRawChannel builder didn't calculate these quantities
424 // * the offline time is zero (may happen if the OFC amplitude < 0.1 )
425 // * online value are zero
426 // * offline quality is 0 (avoid div-by-zero later on)
427 ATH_MSG_VERBOSE("Skip time/Quality comparison, not computed either online or offline");
428 return result;
429 }
430
431 auto DiffT = Monitored::Scalar<float>("Tdiff",timeOffline - t_fB);
432 //Find expected precision given the cell-time:
433 auto t_Precision=std::upper_bound(m_T_precision.begin(),m_T_precision.end(),std::abs(t_fB),pairValueCmp);
434
435 auto Ton = Monitored::Scalar<float>("Ton",t_fB);
436 auto Toff = Monitored::Scalar<float>("Toff",timeOffline);
437 fill(m_tools[m_histoGroups.at(hg)],DiffT,Ton,Toff);
438 fill(m_MonGroupName,DiffT);
439 if (fabs(DiffT) > t_Precision->second) {
440 auto weight_t = Monitored::Scalar<float>("weight_t",1.);
441 fill(m_tools[m_histoGroups.at(hg)],slot,ft,weight_t);
442 result.t_on=t_fB;
443 result.t_off=timeOffline;
444 }
445
446
447 //Quality check:
448 float qdiff = 65535.0; // max possible quality
449 if (q_fD > 0.) qdiff = (q_fD - q_fB)/std::sqrt(q_fD);
450
451 //Find expected precision given the cell-quality:
452 auto q_Precision=std::upper_bound(m_Q_precision.begin(),m_Q_precision.end(),std::abs(q_fB),pairValueCmp);
453
454 auto DiffQ = Monitored::Scalar<float>("Qdiff", qdiff);
455 auto Qon = Monitored::Scalar<float>("Qon",q_fB);
456 auto Qoff = Monitored::Scalar<float>("Qoff",q_fD);
457 fill(m_tools[m_histoGroups.at(hg)],DiffQ,Qon,Qoff);
458 fill(m_MonGroupName,DiffQ);
459
460 if (fabs(DiffQ) > q_Precision->second) {
461 auto weight_q = Monitored::Scalar<float>("weight_q",1.);
462 fill(m_tools[m_histoGroups.at(hg)],slot,ft,weight_q);
463 result.q_on=q_fB;
464 result.q_off=q_fD;
465 }
466 return result;
467}
468
470 const LArDigit& dig,
471 const EventContext& ctx) const{
472
473
474 const HWIdentifier chid=dig.channelID();
475 const auto gain=dig.gain();
476
478 const CaloNoise *noisep = *noiseHdl;
479
481 const ILArPedestal* pedestals=*pedestalHdl;
482
484 const ILArOFC* ofcs=*ofcHdl;
485 const ILArOFC::OFCRef_t& OFCa = ofcs->OFC_a(chid,gain);
486 const ILArOFC::OFCRef_t& OFCb = ofcs->OFC_b(chid,gain);
488
489 const ILArShape* shapes=*shapeHdl;
490 const ILArShape::ShapeRef_t& shape = shapes->Shape(chid,gain);
491 const ILArShape::ShapeRef_t& shapeDer = shapes->ShapeDer(chid,gain);
492
494 const ILArHVScaleCorr* hvScaleCorrs=*hvScaleCorrHdl;
495
497 const LArADC2MeV* adc2mev=*adc2MeVHdl;
498 const auto& polynom_adc2mev=adc2mev->ADC2MEV(chid,gain);
499 const float escale = (polynom_adc2mev)[1];
500 float ramp0 = (polynom_adc2mev)[0];
502 const LArOnOffIdMapping* cabling=*cablingHdl;
503
504 const float hvscale = hvScaleCorrs->HVScaleCorr(chid);
505 const int channel=m_LArOnlineIDHelper->channel(chid);
506 const HWIdentifier febid=m_LArOnlineIDHelper->feb_Id(chid);
507 const std::vector<short>& samples=dig.samples();
508 if (gain == 0) ramp0 = 0.; // no ramp intercepts in HG
509 const float ped = pedestals->pedestal(chid,dig.gain());
510
511
512 //Log-file output if requested ...
514 ATH_MSG_WARNING("Channel " << m_LArOnlineIDHelper->channel_name(chid) << ", gain " << gain
515 <<", run " << ctx.eventID().run_number() << ", evt " << ctx.eventID().event_number());
516 if (cmp.e_on!=cmp.e_off) {
517 ATH_MSG_WARNING("DSP Energy Error:");
518 ATH_MSG_WARNING( " Eonl = " << cmp.e_on << " , Eoff = " << cmp.e_off
519 << " , Eoff - Eonl = " << cmp.e_off-cmp.e_on);
520 }
521 else {
522 ATH_MSG_WARNING("Eonl=Eofl="<< cmp.e_on);
523 }
524 if(cmp.t_off!=cmp.t_on ) {
525 ATH_MSG_WARNING( "DSP Time Error:");
526 ATH_MSG_WARNING( " Tonl = " << cmp.t_on << " , Toff = " << cmp.t_off
527 << " , Toff - Tonl = " << cmp.t_off - cmp.t_on);
528 }
529 if (cmp.q_off!=cmp.q_on) {
530 ATH_MSG_WARNING( "DSP Quality Error");
531 ATH_MSG_WARNING( " Qonl = " << cmp.q_on << " , Qoff = " << cmp.q_off
532 << " (Qoff - Qnl)/sqrt(Qoff) = " << (cmp.q_off - cmp.q_on)/std::sqrt(cmp.q_off));
533 }
534 }
535 if (m_ndump<m_max_dump) {
536 std::string output;
537 output="Digits : ";
538 for (const short s : samples) {output+=std::to_string(s)+ " ";}
539 ATH_MSG_INFO(output);
540
541 output="OFCa : ";
542 for (const auto o : OFCa) {output+=std::to_string(o)+" ";}
543 ATH_MSG_INFO(output);
544
545 output="OFCb : ";
546 for (const auto o : OFCb) {output+=std::to_string(o)+" ";}
547 ATH_MSG_INFO(output);
548
549 output="Shape : ";
550 for (const auto s : shape) {output+=std::to_string(s)+" ";}
551 ATH_MSG_INFO(output);
552
553 output="ShapeDer : ";
554 for (const auto s : shapeDer) {output+=std::to_string(s)+" ";}
555 ATH_MSG_INFO( output );
556
557 ATH_MSG_INFO( "Escale: "<<escale<<" intercept: "<<ramp0<<" pedestal: "<<ped<<" gain: "<<dig.gain() );
558 const Identifier cellid=cabling->cnvToIdentifier(chid);
559 const float noise=noisep->getNoise(cellid,gain);
560 ATH_MSG_INFO( "Noise: "<<noise);
561 ATH_MSG_INFO( "HVScaleCorr: "<<hvscale);
562 double emon=0.;
563 const unsigned nOFCSamp=std::min(samples.size(),OFCa.size());
564 for (unsigned k=0; k<nOFCSamp; ++k) emon += (samples.at(k)-ped)*OFCa.at(k);
565 emon *= escale;
566 emon += ramp0;
567 ATH_MSG_INFO( "intercept + Escale*Sum[(sample-ped)*OFCa] "<<emon);
568 }//end log-file dump
569
570 if(m_doCellsDump) {
571 float sumai = 0.;
572 for (const float a : OFCa) {sumai += a;}
573 float pedplusoffset=0;
574 if (escale*sumai != 0) pedplusoffset = ped - ramp0/(escale*sumai);
575 else pedplusoffset = 0;
576 const float inv_Escale = 1. / escale;
577
578 m_fai << channel<<"\t"<< ped<<"\t"<< pedplusoffset<<"\t"
579 << OFCa[0]*escale<<"\t"<< OFCa[1]*escale<<"\t"<< OFCa[2]*escale<<"\t"<< OFCa[3]*escale<<"\t"<< OFCa[4]*escale<<"\t"
580 << OFCb[0]*escale<<"\t"<< OFCb[1]*escale<<"\t"<< OFCb[2]*escale<<"\t"<< OFCb[3]*escale<<"\t"<< OFCb[4]*escale<<"\t"
581 << shape[0]*inv_Escale<<"\t"<< shape[1]*inv_Escale<<"\t"<< shape[2]*inv_Escale<<"\t"<< shape[3]*inv_Escale<<"\t"<< shape[4]*inv_Escale<<"\t"
582 << shapeDer[0]*inv_Escale<<"\t"<< shapeDer[1]*inv_Escale<<"\t"<< shapeDer[2]*inv_Escale<<"\t"<< shape[3]*inv_Escale<<"\t"<< shapeDer[4]*inv_Escale << std::endl;
583
584 // Dumping file with offline energy and online energ
585 m_fen << m_ndump << " " << cmp.e_off << " " << cmp.e_on ;
586 m_fen << " // FEB " << febid.get_identifier32().get_compact() << " ( channel " << channel << " ), event " << ctx.eventID().event_number() << std::endl;
587
588 // Dumping file with digits
589 m_fdig << channel << " ";
590 for (const short d : samples) {
591 m_fdig << d << " ";
592 }
593 m_fdig << " // FEB " << febid.get_identifier32().get_compact() << " ( channel " << channel << " ), event " << ctx.eventID().event_number() << std::endl;
594 } //end if m_doCellsDump
595 return;
596}
597
599 const int gain,
600 const EventContext& ctx,
601 const diff_t & cmp) const {
602
603 // Dumping cell information in a text file so that it can be compared to unhappy channels lists for instance ...
604
605 const int barrel_ec=m_LArOnlineIDHelper->barrel_ec(chid);
606 const int posneg=m_LArOnlineIDHelper->pos_neg(chid);
607 const int slot = m_LArOnlineIDHelper->slot(chid);
608 const int FT = m_LArOnlineIDHelper->feedthrough(chid);
609 const int channel = m_LArOnlineIDHelper->channel(chid);
610 const HWIdentifier febid= m_LArOnlineIDHelper->feb_Id(chid);
611 m_fdump << "0x" << std::hex << std::setw(8)<<febid.get_identifier32().get_compact() << " " << std::dec << std::setw(3) << std::right << channel << " 0x" << std::hex << std::setw(8)<<chid.get_identifier32().get_compact()
612 <<std::dec << std::setw(3) << std::right << slot << std::setw(3) << std::right << FT << std::setw(3) << std::right << barrel_ec << std::setw(3) << std::right<< posneg << std::setw(6) << std::right << getPartitionName(chid)
613 << " " << gain << " " << " " << cmp.e_off << " "<< cmp.e_on << " "<<cmp.t_off << " "<<cmp.t_on <<" "<<cmp.q_off << " "<<cmp.q_on <<ctx.eventID().event_number()<<std::endl;
614 return;
615}
616
618 errors_E.fill(0);
619 errors_T.fill(0);
620 errors_Q.fill(0);
621 return;
622}
623
625 m_fai.close();
626 m_fdig.close();
627 m_fen.close();
628 m_fdump.close();
629 return StatusCode::SUCCESS;
630}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
LArBadXCont< LArBadChannel > LArBadChannelCont
static Double_t a
const ServiceHandle< StoreGateSvc > & detStore() const
virtual StatusCode initialize() override
initialize
SG::ReadHandle< xAOD::EventInfo > GetEventInfo(const EventContext &) const
Return a ReadHandle for an EventInfo object (get run/event numbers, etc.)
ToolHandleArray< GenericMonitoringTool > m_tools
Array of Generic Monitoring Tools.
float getNoise(const IdentifierHash h, const int gain) const
Accessor by IdentifierHash and gain.
Definition CaloNoise.h:34
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
virtual const float & HVScaleCorr(const HWIdentifier &id) const =0
virtual OFCRef_t OFC_b(const HWIdentifier &id, int gain, int tbin=0) const =0
virtual OFCRef_t OFC_a(const HWIdentifier &id, int gain, int tbin=0) const =0
access to OFCs by online ID, gain, and tbin (!=0 for testbeam)
LArVectorProxy OFCRef_t
This class defines the interface for accessing Optimal Filtering coefficients for each channel provid...
Definition ILArOFC.h:26
virtual float pedestal(const HWIdentifier &id, int gain) const =0
LArVectorProxy ShapeRef_t
This class defines the interface for accessing Shape (Nsample variable, Dt = 25 ns fixed) @stereotype...
Definition ILArShape.h:26
virtual ShapeRef_t Shape(const HWIdentifier &id, int gain, int tbin=0, int mode=0) const =0
virtual ShapeRef_t ShapeDer(const HWIdentifier &id, int gain, int tbin=0, int mode=0) const =0
value_type get_compact() const
Get the compact id.
This is a "hash" representation of an Identifier.
std::string getString() const
Provide a string form of the identifier - hexadecimal.
Identifier32 get_identifier32() const
Get the 32-bit version Identifier, will be invalid if >32 bits needed.
const LArVectorProxy ADC2MEV(const HWIdentifier &id, int gain) const
Definition LArADC2MeV.h:32
Liquid Argon digit base class.
Definition LArDigit.h:25
CaloGain::CaloGain gain() const
Definition LArDigit.h:72
const std::vector< short > & samples() const
Definition LArDigit.h:78
const HWIdentifier & channelID() const
Definition LArDigit.h:69
Container class for LArFebHeader.
std::array< unsigned, 3 > errors_T
std::array< unsigned, 3 > errors_Q
std::array< unsigned, 3 > errors_E
Gaudi::Property< std::vector< std::pair< int, int > > > m_T_precision
Gaudi::Property< std::vector< std::pair< int, int > > > m_E_precision
const char * getPartitionName(const HWIdentifier chid) const
virtual StatusCode initialize() override final
initialize
Gaudi::Property< float > m_timeOffset
Gaudi::Property< bool > m_printEnergyErrors
Gaudi::Property< bool > m_doRodStatus
Gaudi::Property< std::string > m_DigitsFileName
Gaudi::Property< std::string > m_AiFileName
SG::ReadCondHandleKey< CaloNoise > m_noiseCDOKey
SG::ReadCondHandleKey< ILArPedestal > m_keyPedestal
Gaudi::Property< std::vector< std::pair< int, int > > > m_Q_precision
void dumpCellInfo(const HWIdentifier chid, const int gain, const EventContext &ctx, const diff_t &comp) const
Dump a cell's information and calculated energies into a txt file.
Gaudi::Property< bool > m_skipNullPed
const LArOnlineID * m_LArOnlineIDHelper
SG::ReadHandleKey< LArFebHeaderContainer > m_headerContainerKey
Gaudi::Property< std::vector< std::string > > m_streams
SG::ReadCondHandleKey< LArADC2MeV > m_adc2mevKey
virtual StatusCode fillHistograms(const EventContext &ctx) const override final
adds event to the monitoring histograms
Gaudi::Property< bool > m_removeNoiseBursts
virtual ~LArRODMonAlg()
Default destructor.
Gaudi::Property< float > m_peakTime_cut
PARTITION getPartition(const HWIdentifier chid) const
Gaudi::Property< bool > m_skipKnownProblematicChannels
SG::ReadCondHandleKey< ILArOFC > m_keyOFC
Gaudi::Property< std::string > m_DumpCellsFileName
SG::ReadDecorHandleKey< xAOD::EventInfo > m_eventInfoKey
Gaudi::Property< std::vector< std::string > > m_problemsToMask
Gaudi::Property< bool > m_doCellsDump
virtual StatusCode finalize() override final
diff_t compareChannel(const LArRawChannel &rcDig, const LArRawChannel &rcBS) const
SG::ReadHandleKey< LArRawChannelContainer > m_channelKey_fromBytestream
SG::ReadCondHandleKey< LArBadChannelCont > m_bcContKey
LArBadChannelMask m_bcMask
Gaudi::Property< std::string > m_EnergyFileName
SG::ReadCondHandleKey< ILArShape > m_keyShape
Gaudi::Property< unsigned > m_max_dump
SG::ReadCondHandleKey< LArOnOffIdMapping > m_cablingKey
const float m_BC
std::atomic< unsigned > m_ndump
Gaudi::Property< bool > m_doCheckSum
Gaudi::Property< bool > m_doDspTestDump
Gaudi::Property< std::string > m_MonGroupName
Gaudi::Property< std::vector< std::string > > m_partitions
Gaudi::Property< short > m_adc_th
std::map< std::string, int > m_histoGroups
SG::ReadHandleKey< LArDigitContainer > m_digitContainerKey
void detailedOutput(const LArRODMonAlg::diff_t &, const LArDigit &dig, const EventContext &ctx) const
SG::ReadHandleKey< LArRawChannelContainer > m_channelKey_fromDigits
SG::ReadCondHandleKey< ILArHVScaleCorr > m_keyHVScaleCorr
Liquid Argon ROD output object base class.
uint16_t provenance() const
int time() const
uint16_t quality() const
HWIdentifier channelID() const
int energy() const
Declare a monitored scalar variable.
@ LAr
The LAr calorimeter.
int lb
Definition globals.cxx:23
@ LARMEDIUMGAIN
Definition CaloGain.h:18
@ LARLOWGAIN
Definition CaloGain.h:18
@ LARNGAIN
Definition CaloGain.h:19
@ LARHIGHGAIN
Definition CaloGain.h:18
std::vector< V > buildToolMap(const ToolHandleArray< GenericMonitoringTool > &tools, const std::string &baseName, int nHist)
Builds an array of indices (base case)
const T * get(const ReadCondHandleKey< T > &key, const EventContext &ctx)
Convenience function to retrieve an object given a ReadCondHandleKey.
void fill(H5::Group &out_file, size_t iterations)