ATLAS Offline Software
Loading...
Searching...
No Matches
TileRawChannelBuilder.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// Tile includes
12
13// Atlas includes
17
18// Gaudi includes
19
20
21static const InterfaceID IID_ITileRawChannelBuilder("TileRawChannelBuilder", 1, 0);
22
26
27
32
36
40
44
49 , const std::string& name, const IInterface* parent)
50 : AthAlgTool(type, name, parent)
52 , m_rChUnit(TileRawChannelUnit::ADCcounts)
53 , m_bsflags(0)
54 , m_tileID(nullptr)
55 , m_tileHWID(nullptr)
56 , m_trigType(0)
57 , m_idophys(false)
58 , m_idolas(false)
59 , m_idoped(false)
60 , m_idocis(false)
62 , m_capdaq(100)
63 , m_evtCounter(0)
64 , m_chCounter(0)
65 , m_nChL(0)
66 , m_nChH(0)
67 , m_RChSumL(0.0)
68 , m_RChSumH(0.0)
70 , m_tileInfo(nullptr)
71{
73 memset(m_error, 0, sizeof(m_error));
74
75 declareProperty("calibrateEnergy", m_calibrateEnergy = false);
76 declareProperty("correctTime", m_correctTime = false);
77 declareProperty("AmpMinForAmpCorrection", m_ampMinThresh = 15.0);
78 declareProperty("TimeMinForAmpCorrection", m_timeMinThresh = -12.5);
79 declareProperty("TimeMaxForAmpCorrection", m_timeMaxThresh = 12.5);
80 declareProperty("RunType", m_runType = 0);
81 declareProperty("DataPoolSize", m_dataPoollSize = -1);
82 declareProperty("UseDSPCorrection", m_useDSP = true);
83 declareProperty("TileInfoName", m_infoName = "TileInfo");
84 declareProperty("FirstSample",m_firstSample = 0);
85
86}
87
93
98
99 ATH_MSG_INFO( "TileRawChannelBuilder::initialize()" );
100
102 m_idophys = ((m_trigType == 0) || (m_trigType == 1));
103 m_idolas = ((m_trigType == 2) || (m_trigType == 3));
104 m_idoped = ((m_trigType == 4) || (m_trigType == 5));
105 m_idocis = ((m_trigType == 8) || (m_trigType == 9));
106 m_cischan = -1;
107 m_capdaq = 0;
108 m_chCounter = 0;
109 m_evtCounter = 0;
110 m_rawChannelCnt = nullptr;
111 m_nChL = m_nChH = 0;
112 m_RChSumL = m_RChSumH = 0.0;
113 m_evtCounter = -1;
114 // retrieve TileID helpers and TileIfno from det store
115 ATH_CHECK( detStore()->retrieve(m_tileID, "TileID") );
116 ATH_CHECK( detStore()->retrieve(m_tileHWID, "TileHWID") );
117
118 ATH_CHECK( detStore()->retrieve(m_tileInfo, m_infoName) );
119 m_i_ADCmax = m_tileInfo->ADCmax();
124 m_ADCmaskValueMinusEps = m_tileInfo->ADCmaskValue() - 0.01; // indicates channels which were masked in background dataset
125
126 // access tools and store them
127 ATH_CHECK( m_noiseFilterTools.retrieve() );
128 ATH_MSG_DEBUG( "Successfully retrieve NoiseFilterTools: " << m_noiseFilterTools );
129
131 ATH_MSG_DEBUG( "Obsolete calibrateEnergy flag is set to True in jobOptions - disabling it" );
132 m_calibrateEnergy = false;
133 }
134
135 // check if we want to keep ADC counts or convert them to pCb
138
139 // if unit is not pCb, but already MeV one can use method TileRawChannelContainer::set_unit()
140 // later to declare that
141
142 // 8 upper bits of bsflags:
143 // UUPPSTTT
144 // 31,30 - units
145 // 29,28 - pulse type = 3 for offline reco
146 // 27 - 7(=0) or 9(=1) samples
147 // 24,25,26 - TileFragHash::TYPE - OF algorithm type
148 int nsamp = (m_tileInfo->NdigitSamples() > 7) ? 1 : 0;
149 m_bsflags = (m_rChUnit << 30) | (3 << 28) | (nsamp << 27) | (m_rChType << 24);
150
151 // bits 12-15 - various options
152 if (m_correctTime) m_bsflags |= 0x1000;
153
154 if (msgLvl(MSG::DEBUG)) {
155 msg(MSG::DEBUG) << "TileRawChannelBuilder created, storing rc in '"
156 << m_rawChannelContainerKey.key() << "'" << endmsg;
157 msg(MSG::DEBUG) << " calibrate energy = " << m_calibrateEnergy << endmsg;
158 msg(MSG::DEBUG) << " correct time = " << m_correctTime << endmsg;
159 msg(MSG::DEBUG) << " run type = " << m_runType << endmsg;
160 }
161
162 if (m_dataPoollSize < 0) m_dataPoollSize = m_tileHWID->channel_hash_max();
163
164 ATH_CHECK( m_cablingSvc.retrieve());
165
166 m_cabling = m_cablingSvc->cablingService();
167 if (!m_cabling) {
168 ATH_MSG_ERROR( "Unable to retrieve TileCablingService" );
169 return StatusCode::FAILURE;
170 }
171
173
174 int runPeriod = m_cabling->runPeriod();
175 std::ostringstream os;
176 if (runPeriod==3) {
177 if ( m_demoFragIDs.empty() ) {
178 std::vector<int> v = { 0x10d }; // LBA14 is demonstrator in RUN3
179 m_demoFragIDs = v;
180 }
181 os << " in RUN3";
182 }
183
184 if ( !m_demoFragIDs.empty() ) {
186 os << " (frag IDs):";
187 for (int fragID : m_demoFragIDs) {
188 if (fragID>0)
189 os << " 0x" << std::hex << fragID << std::dec;
190 else
191 os << " " << fragID;
192 }
193 ATH_MSG_INFO("Enable special treatment for demonstrator modules" << os.str());
194 }
195
196 if (m_calibrateEnergy) {
197 ATH_CHECK( m_tileToolEmscale.retrieve() );
198 } else {
199 m_tileToolEmscale.disable();
200 }
201
202 if (m_correctTime) {
203 ATH_CHECK( m_tileToolTiming.retrieve() );
204 } else {
205 m_tileToolTiming.disable();
206 }
207
209 ATH_CHECK( m_tileIdTransforms.retrieve() );
210 } else {
211 m_tileIdTransforms.disable();
212 }
213
216
217 if (m_useDSP && !m_DSPContainerKey.key().empty()) {
218 ATH_CHECK( m_DSPContainerKey.initialize() );
219 }
220 else {
222 }
223
224 return StatusCode::SUCCESS;
225}
226
228 ATH_MSG_INFO( "Finalizing" );
229 return StatusCode::SUCCESS;
230}
231
232StatusCode TileRawChannelBuilder::createContainer(const EventContext& ctx) {
233 initLog(ctx);
234
235 // create TRC container
236 m_rawChannelCnt = std::make_unique<TileMutableRawChannelContainer>(true, m_rChType, m_rChUnit, SG::VIEW_ELEMENTS);
237 ATH_CHECK( m_rawChannelCnt->status() );
238 m_rawChannelCnt->set_bsflags(m_bsflags);
239
240 ATH_MSG_DEBUG( "Created TileRawChannelContainer '" << m_rawChannelContainerKey.key() << "'" );
241
242 return StatusCode::SUCCESS;
243}
244
245void TileRawChannelBuilder::initLog(const EventContext& ctx) {
246
247 const TileDQstatus* DQstatus = SG::makeHandle (m_DQstatusKey, ctx).get();
248
249 // update only if there is new event
250 if (m_evtCounter != ctx.evt()) {
251
252 m_evtCounter = ctx.evt();
253 if (m_runType != 0) m_trigType = m_runType;
254 else m_trigType = DQstatus->trigType();
255
256 if (0 == m_trigType) {
257 m_idophys = (DQstatus->calibMode() == 0);
258 m_idolas = false;
259 m_idoped = false;
260 m_idocis = (DQstatus->calibMode() == 1);
261 } else {
262 m_idophys = (m_trigType <= 1);
263 m_idolas = ((m_trigType == 2) || (m_trigType == 3));
264 m_idoped = ((m_trigType == 4) || (m_trigType == 5));
265 m_idocis = ((m_trigType == 8) || (m_trigType == 9));
266 }
267
268 const unsigned int *cispar = DQstatus->cispar();
269 if (0 == cispar[7]) { // if capdaq not set, it can't be CIS event
270 if (m_idocis) { // cis flag was set incorrectly, change to ped
271 m_idoped = true;
272 m_idocis = false;
273 }
274 m_capdaq = 0.0;
275 } else {
276 m_capdaq = (cispar[7] < 10) ? 5.2 : 100.0;
277 }
278 m_cischan = cispar[8] - 1; // channel where CIS is fired (-1 = all channels)
279
280 ATH_MSG_DEBUG( "Trig type is " << m_trigType
281 << "; dophys is " << ((m_idophys) ? "true" : "false")
282 << "; dolas is " << ((m_idolas) ? "true" : "false")
283 << "; doped is " << ((m_idoped) ? "true" : "false")
284 << "; docis is " << ((m_idocis) ? "true" : "false")
285 << "; capacitor is " << m_capdaq
286 << "; cis chan is " << m_cischan );
287 }
288}
289
290TileRawChannel* TileRawChannelBuilder::rawChannel(const TileDigits* digits, const EventContext& /*ctx*/) {
291 ++m_chCounter;
292 ATH_MSG_WARNING( "Default constructor for rawChannel!" );
293 DataPool<TileRawChannel> tileRchPool (100);
294 TileRawChannel *rawCh = tileRchPool.nextElementPtr();
295 rawCh->assign (digits->adc_HWID(), 0.0, 0.0, 0.0);
296 return rawCh;
297}
298
299void TileRawChannelBuilder::fill_drawer_errors(const EventContext& ctx,
300 const TileDigitsCollection* coll)
301{
302 const TileDQstatus* DQstatus = SG::makeHandle (m_DQstatusKey, ctx).get();
303
304 int frag = coll->identify();
305 int ros = (frag >> 8);
306 int drawer = (frag & 0xff);
307
308 m_lastDrawer = frag;
309
310 memset(m_error, 0, sizeof(m_error));
311 int dmuerr[MAX_DMUS] = {0};
312 int nch = 0;
313 bool bigain = DQstatus->isBiGain();
314 if (!bigain) { // in bigain runs we don't have DQ status fragment
315 for (int ch = 0; ch < MAX_CHANNELS; ch += 3) {
316 if (!DQstatus->isAdcDQgood(ros, drawer, ch, 0)) {
317 m_error[ch + 2] = m_error[ch + 1] = m_error[ch] = -3;
318 dmuerr[ch / 3] = 3;
319 nch += 3;
320 }
321 }
322 }
323 if (nch == MAX_CHANNELS) { // all bad - nothing to do
324 m_badDrawer = true;
325 ATH_MSG_VERBOSE( "Drawer 0x" << MSG::hex << frag << MSG::dec
326 << " is bad - skipping bad patterns check " );
327 return;
328 } else {
329 m_badDrawer = false;
330 ATH_MSG_VERBOSE( "Drawer 0x" << MSG::hex << frag << MSG::dec
331 << " looking for bad patterns in digits" );
332 }
333
334 float mindig, maxdig;
335 int nchbad[2] = { 0, 0 };
336
337 // Iterate over all digits in this collection
339 TileDigitsCollection::const_iterator lastDigit = coll->end();
340
341 for (; digitItr != lastDigit; ++digitItr) {
342 const TileDigits * pDigits = (*digitItr);
343 HWIdentifier adcId = pDigits->adc_HWID();
344 int channel = m_tileHWID->channel(adcId);
345 int gain = m_tileHWID->adc(adcId);
346
347 if (m_error[channel]) {
348 ATH_MSG_VERBOSE( "BadCh " << ros
349 << "/" << drawer
350 << "/" << channel
351 << "/" << gain << " BAD DQ STATUS ");
352
353 } else {
354
355 int err = CorruptedData(ros, drawer, channel, gain, pDigits->samples(), mindig, maxdig, m_ADCmaxMinusEps, m_ADCmaskValueMinusEps);
356
357 if (err) {
358
359 m_error[channel] = err;
360 if (err > -5) {
361 ++dmuerr[channel / 3];
362 ++nchbad[channel / 24];
363 }
364
365 if (msgLvl(MSG::VERBOSE)) {
366
367 msg(MSG::VERBOSE) << "BadCh " << ros
368 << "/" << drawer
369 << "/" << channel
370 << "/" << gain;
371 if (err < -5) msg(MSG::VERBOSE) << " Warning " << err;
372 else msg(MSG::VERBOSE) << " Error " << err;
373 if (mindig > m_ADCmaskValueMinusEps) msg(MSG::VERBOSE) << " BADDQ";
374 if (maxdig > m_ADCmaxMinusEps) msg(MSG::VERBOSE) << " Overflow";
375 if (mindig < 0.1) msg(MSG::VERBOSE) << " Underflow";
376 if (err < 0) msg(MSG::VERBOSE) << " Const";
377
378 msg(MSG::VERBOSE) << " samp=";
379 std::vector<float> digits = pDigits->samples();
380 for (unsigned int i = 0; i < digits.size(); ++i) {
381 msg(MSG::VERBOSE) << " " << digits[i];
382 }
383 msg(MSG::VERBOSE) << endmsg;
384 }
385
386 } else {
387 if (mindig < 0.01) err += 1;
388 if (maxdig > m_ADCmaxMinusEps) err += 2;
389 if (err) m_error[channel] = err - 10;
390 }
391 }
392 }
393
394 // check if we want to mask half a drawer
395 // in this case set error = -4 for channels which were good before
396
397 int ndmubad[2] = { 0, 0 };
398 int dmu = 0;
399 for (; dmu < MAX_DMUS / 2; ++dmu) { // first half
400 if (dmuerr[dmu] > 1)
401 ++ndmubad[0]; // count DMUs with at least two bad channels
402 }
403 for (; dmu < MAX_DMUS; ++dmu) { // second half
404 if (dmuerr[dmu] > 1)
405 ++ndmubad[1]; // count DMUs with at least two bad channels
406 }
407
408 int ndmulimit[2] = { 3, 3 }; // max number of bad DMUs when half-drawer is not yet masked
409 // if 4 DMUs will be bad - mask whole half-drawer
410 if (frag > 0x2ff) { // if extended barrel
411 if (frag == 0x30e || frag == 0x411)
412 ndmulimit[0] = 4; // in EB special one DMU is always bad (missing)
413 ndmulimit[1] = 5; // in second half of EB 4 DMUs ara always bad (missing)
414 // only if 7 DMUs are bad, mask whole half-drawer
415 }
416
417 bool printall = true;
418 for (int p = 0; p < 2; ++p) {
419 if (ndmubad[p] > ndmulimit[p] && nchbad[p] > 0) {
420 if (msgLvl(MSG::VERBOSE)) {
421 msg(MSG::VERBOSE) << "Drawer 0x" << MSG::hex << frag << MSG::dec
422 << " masking whole " << ((p) ? "second" : "first")
423 << " half" << endmsg;
424 if (printall) {
425 msg(MSG::VERBOSE) << "nDMuErr ";
426 for (int d = 0; d < MAX_DMUS; ++d) {
427 msg(MSG::VERBOSE) << " " << dmuerr[d];
428 }
429 msg(MSG::VERBOSE) << " total " << ndmubad[p] << " errors" << endmsg;
430
431 msg(MSG::VERBOSE) << "ChErr ";
432 int ch = 0;
433 while (ch < MAX_CHANNELS-2) {
434 msg(MSG::VERBOSE) << " " << m_error[ch++];
435 msg(MSG::VERBOSE) << " " << m_error[ch++];
436 msg(MSG::VERBOSE) << " " << m_error[ch++];
437 msg(MSG::VERBOSE) << " ";
438 }
439
440 msg(MSG::VERBOSE) << " total " << nchbad[p]
441 << " bad patterns" << endmsg;
442 printall = false;
443 }
444 }
445 int ch = (p) ? MAX_CHANNELS / 2 : 0;
446 int chmax = (p) ? MAX_CHANNELS : MAX_CHANNELS / 2;
447 for (; ch < chmax; ++ch) {
448 if (m_error[ch] == 0 || m_error[ch] < -5) { // channel was good before
449 m_error[ch] = -4;
450 }
451 }
452 }
453 }
454
455}
456
458 static const char * const errname[26] = {
459 "-10 - good signal",
460 "-9 - underflow",
461 "-8 - overflow",
462 "-7 - underflow and overflow",
463 "-6 - constant signal",
464 "-5 - disconnected channel",
465 "-4 - half a drawer masked",
466 "-3 - bad DQ status",
467 "-2 - underflow in all samples",
468 "-1 - overflow in all samples",
469 "0 - unknown error",
470 "1 - jump from zero to saturation",
471 "2 - samples with zeros",
472 "3 - at least two saturated. others - close to pedestal",
473 "4 - two distinct levels with at least 2 samples each",
474 "5 - pedestal with jump up in one sample",
475 "6 - pedestal with jump down in one sample",
476 "7 - signal with jump up in one sample",
477 "8 - signal with jump down in one sample",
478 "9 - base line above threshold in low gain",
479 "10 - jump down in first sample in low gain",
480 "11 - jump down in last sample in low gain",
481 "12 - jump up in one sample above const",
482 "13 - jump down in one sample below const",
483 "14 - unrecoverable timing jump",
484 "15 - unknown error"
485 };
486
487 return errname[std::min(25, std::max(0, int((ped + 500) * 1e-4)))];
488}
489
490
491StatusCode TileRawChannelBuilder::build(const TileDigitsCollection* coll, const EventContext& ctx)
492{
493
494 int frag = coll->identify();
495
496 // make sure that error array is up-to-date
497 if (frag != m_lastDrawer && m_notUpgradeCabling) {
498 fill_drawer_errors(ctx, coll);
499 }
500
501 // Iterate over all digits in this collection
503 TileDigitsCollection::const_iterator lastDigit = coll->end();
504
505 for (; digitItr != lastDigit; ++digitItr) {
506
507 TileRawChannel* rch = rawChannel((*digitItr), ctx);
508
510
511 int err = m_error[m_tileHWID->channel(rch->adc_HWID())];
512
513 if (err) {
514 if (err == -8 || err == -7) m_overflows.push_back(std::make_pair(rch, (*digitItr)));
515 float ped = rch->pedestal() + 100000 + 10000 * err;
516 rch->setPedestal(ped);
517 if (msgLvl(MSG::VERBOSE) && !m_badDrawer) {
518 if (err < -5) {
519 msg(MSG::VERBOSE) << "BadCh " << m_tileHWID->to_string(rch->adc_HWID())
520 << " warning = " << BadPatternName(ped) << endmsg;
521 } else {
522 msg(MSG::VERBOSE) << "BadCh " << m_tileHWID->to_string(rch->adc_HWID())
523 << " error = " << BadPatternName(ped) << endmsg;
524 }
525 }
526 }
527
528 }
529
530 ATH_CHECK( m_rawChannelCnt->push_back (rch) );
531 }
532
533 IdentifierHash hash = m_rawChannelCnt->hashFunc().hash(coll->identify());
534 TileRawChannelCollection* rawChannelCollection = m_rawChannelCnt->indexFindPtr(hash);
535 rawChannelCollection->setLvl1Id(coll->getLvl1Id());
536 rawChannelCollection->setLvl1Type(coll->getLvl1Type());
537 rawChannelCollection->setDetEvType(coll->getDetEvType());
538 rawChannelCollection->setRODBCID(coll->getRODBCID());
539
540 return StatusCode::SUCCESS;
541}
542
543StatusCode TileRawChannelBuilder::commitContainer(const EventContext& ctx)
544{
545
546 const TileDQstatus* DQstatus = SG::makeHandle (m_DQstatusKey, ctx).get();
547
548 ToolHandleArray<ITileRawChannelTool>::iterator itrTool = m_noiseFilterTools.begin();
549 ToolHandleArray<ITileRawChannelTool>::iterator endTool = m_noiseFilterTools.end();
550
551 if ( m_useDSP && !m_DSPContainerKey.key().empty() &&
552 (DQstatus->incompleteDigits() || m_chCounter<12288) && itrTool!=endTool )
553 {
554 const TileRawChannelContainer * dspCnt = SG::makeHandle (m_DSPContainerKey, ctx).get();
555 ATH_MSG_DEBUG( "Incomplete container - use noise filter corrections from DSP container" );
556
557 uint32_t bsFlags = dspCnt->get_bsflags();
558 std::vector<IdentifierHash> hashes = m_rawChannelCnt->GetAllCurrentHashes();
559 std::vector<IdentifierHash> dspHashes = dspCnt->GetAllCurrentHashes();
560 if (bsFlags == 0) {
561 ATH_MSG_WARNING("Problem in applying noise corrections: DSP container ("
562 << m_DSPContainerKey.key() << ") seems to be emtpy!");
563 } else if (hashes != dspHashes) {
564 ATH_MSG_ERROR( " Error in applying noise corrections; "
565 "hash vectors do not match.");
566 } else {
567 // Go through all TileRawChannelCollections
568 for (IdentifierHash hash : hashes) {
569 TileRawChannelCollection* coll = m_rawChannelCnt->indexFindPtr (hash);
570 const TileRawChannelCollection* dcoll = dspCnt->indexFindPtr (hash);
571
572 if (coll->identify() != dcoll->identify()) {
573
574 ATH_MSG_ERROR( " Error in applying noise corrections " << MSG::hex
575 << " collection IDs 0x" << coll->identify() << " and 0x" << dcoll->identify()
576 << " do not match " << MSG::dec );
577 break;
578 }
579
580 // iterate over all channels in a collection
583
584 for (TileRawChannel* rch : *coll) {
585 HWIdentifier adc_id = rch->adc_HWID();
586 while (dspItr != dspLast && adc_id != (*dspItr)->adc_HWID()) {
587 ++dspItr;
588 }
589 if (dspItr != dspLast) {
590 float corr = (*dspItr)->pedestal();
591 ATH_MSG_VERBOSE( "Ch "<<m_tileHWID->to_string(adc_id)
592 <<" amp " << rch->amplitude() << " ped " << rch->pedestal()
593 << " corr " << corr );
594 if (corr<10000.) {
595 rch->setAmplitude (rch->amplitude() - corr); // just baseline shift
596 rch->setPedestal (rch->pedestal() + corr); // just baseline shift
597 } else {
598 float ped = rch->pedestal();
599 if (corr > ped) {
600 rch->setPedestal (fmod(ped,10000.) + int(corr)/10000 * 10000); // changing error status
601 ATH_MSG_VERBOSE( "New error status in ped "<<rch->pedestal());
602 }
603 }
604 } else {
605 ATH_MSG_WARNING(" Problem in applying noise corrections "
606 << " can not find channel in DSP container with HWID "
607 << m_tileHWID->to_string(adc_id) );
608 dspItr = dcoll->begin();
609 }
610 }
611 }
612 }
613
614 } else {
615
616 for (ToolHandle<ITileRawChannelTool>& noiseFilterTool : m_noiseFilterTools) {
617 if (noiseFilterTool->process(*m_rawChannelCnt.get(), ctx).isFailure()) {
618 ATH_MSG_ERROR( " Error status returned from noise filter " );
619 } else {
620 ATH_MSG_DEBUG( "Noise filter applied to the container" );
621 }
622 }
623
624 }
625
626 ATH_MSG_DEBUG( " nCh=" << m_chCounter
627 << " nChH/L=" << m_nChH << "/" << m_nChL
628 << " RChSumH/L=" << m_RChSumH << "/" << m_RChSumL );
629
631 ATH_CHECK( rawChannelsContainer.record(std::move(m_rawChannelCnt)) );
632
633 endLog();
634
635 return StatusCode::SUCCESS;
636}
637
639 m_chCounter = 0;
640 m_nChL = m_nChH = 0;
641 m_RChSumL = m_RChSumH = 0.0;
642
643}
644
645double TileRawChannelBuilder::correctAmp(double phase, bool of2) {
646
647 double corr = 1.0;
648 if (of2) {
649 // estimation from Belen for rel 14.0.0
650 /*double a,b,c;
651 if(fabs(phase)<5.){
652 a=0.137; b=0.0877; c=0.0865;
653 }else{
654 a=0.565; b=0.116; c=0.0751;
655 }
656 corr=(1+(a+b*phase+c*phase*phase)/100.);
657 */
658
659 // estimation from Vakhtang for rel 14.4.0
660 /*double k = (phase < 0.0 ? 0.0009400 : 0.0010160);
661 corr = (1.0 + k * phase * phase);
662 */
663
664 // Parabolic correction from Tigran
665 double a1,a2,b,c;
666 a1 = phase < 0.0 ? 0.000940774 : 0.00102111;
667 a2 = phase < 0.0 ? 0.000759051 : 0.000689625;
668 b = phase < 0.0 ? -2.0 * 7.0 * (a1 - a2) : 2.0 * 12.5 * (a1 - a2);
669 c = phase < 0.0 ? 1.0 - 7.0 * 7.0 * (a1-a2) : 1.0 - 12.5 * 12.5 * (a1-a2);
670 if (phase < 12.5 && phase > -7.0) corr = a1 * phase * phase + 1.0;
671 else corr = phase * ( a2 * phase + b) + c;
672
673
674 } else {
675 /*double a,b,c;
676 if(phase<0){
677 a=1.0002942; b=0.0003528; c=0.0005241;
678 }else{
679 a=1.0001841; b=-0.0004182; c=0.0006167;
680 }
681 corr = a + phase * ( b + c * phase);
682 */
683
684 /*double k = (phase < 0.0 ? 0.0005241 : 0.0006167);
685 corr = (1.0 + k * phase * phase);
686 */
687
688 // 4th degree polynomial correction from Tigran
689 double k1 = (phase < 0.0 ? -0.0000326707:0.000380336);
690 double k2 = (phase < 0.0 ? -0.000560962:-0.000670487);
691 double k3 = (phase < 0.0 ? -0.00000807869:0.00000501773);
692 double k4 = (phase < 0.0 ? -0.000000145008:0.0000000584647);
693
694 corr = 1.0 / (1.0 + (k1 + (k2 + (k3 + k4 *phase)*phase)*phase)*phase);
695
696
697 }
698
699 return corr;
700}
701
702
703// Time correction for shifted pulses by Tigran
704double TileRawChannelBuilder::correctTime(double phase, bool of2) {
705
706 double correction = 0.0;
707
708 if (of2) {
709 if(phase < 0) {
710 correction = (-0.00695743 + (0.0020673 - (0.0002976 + 0.00000361305 * phase) * phase) * phase) * phase;
711 } else {
712 correction = (0.0130013 + (0.00128769 + (-0.000550218 + 0.00000755344 * phase) * phase) * phase) * phase;
713 }
714 }
715 // OF1 does not need correction
716
717 return correction;
718}
719
720
721
722int TileRawChannelBuilder::CorruptedData(int ros, int drawer, int channel, int gain,
723 const std::vector<float> & digits, float &dmin, float &dmax, float ADCmaxMinusEps, float ADCmaskValueMinusEps) {
724 bool eb = (ros > 2);
725 bool ebsp = ((ros == 3 && drawer == 14) || (ros == 4 && drawer == 17));
726 bool empty = ((eb && ((channel > 23 && channel < 30) || channel > 41)) || (ebsp && channel < 3));
727 bool not_gap = !(empty || (eb && (channel == 0 || channel == 1 || channel == 12 || channel == 13))
728 || (ebsp && (channel == 18 || channel == 19)));
729
730 const float epsilon = 4.1; // allow +/- 2 counts fluctuations around const value
731 const float delta[4] = { 29.9, 29.9, 49.9, 99.9 }; // jump levels between constLG, constHG, non-constLG, non-constHG
732 const float level1 = 99.9; // jump from this level to m_i_ADCmax is bad
733 const float level2 = 149.9; // base line at this level in low gain is bad
734 const float narrowLevel[2] = { 29.9, 49.9 }; // minimal amplitude for narrow pulses
735 const float delt = std::min(std::min(std::min(delta[0], delta[1]), std::min(delta[2], delta[3])),
736 std::min(narrowLevel[0], narrowLevel[1]));
737 const float secondMaxLevel = 0.3;
738
739 int error = 0;
740
741 unsigned int nSamp = digits.size();
742 if (nSamp) {
743 dmin = dmax = digits[0];
744 unsigned int pmin = 0;
745 unsigned int pmax = 0;
746 unsigned int nzero = (dmin < 0.01) ? 1 : 0;
747 unsigned int nover = (dmax > ADCmaxMinusEps) ? 1 : 0;
748
749 for (unsigned int i = 1; i < nSamp; ++i) {
750 float dig = digits[i];
751 if (dig > dmax) {
752 dmax = dig;
753 pmax = i;
754 } else if (dig < dmin) {
755 dmin = dig;
756 pmin = i;
757 }
758 if (dig < 0.01) ++nzero;
759 else if (dig > ADCmaxMinusEps) ++nover;
760 }
761
762 float dmaxmin = dmax - dmin;
763 //std::cout << " ros " << ros << " drawer " << drawer << " channel " << channel << " not_gap " << not_gap << " nzero " << nzero << " nover " << nover << std::endl;
764
765 if (dmin > ADCmaxMinusEps) { // overflow in all samples
766 error = (dmin > ADCmaskValueMinusEps) ? -3 : -1; // dmin=m_tileInfo->ADCmaskValue() - masking in overlay job (set in TileDigitsMaker)
767
768 } else if (dmax < 0.01) { // underflow in all samples
769 error = (empty) ? -5 : -2; // set different type of errors for exsiting and non-existing channels
770
771 } else if (dmaxmin < 0.01) { // constant value in all samples
772 error = -6;
773
774 } else if (nzero && nover) { // jump from zero to saturation
775 error = 1;
776
777 } else if ((nzero && (not_gap || empty)) || nzero > 1) { // one sample at zero in normal channel
778 error = 2; // or 2 samples at zero in gap/crack/MBTS
779
780 } else if (gain == 0 && dmin > level2) { // baseline above threshold in low gain is bad
781 error = 9;
782
783 } else if (dmaxmin > delt) { // check that max-min is above minimal allowed jump
784
785 float abovemin = dmax;
786 float belowmax = dmin;
787 unsigned int nmin = 0;
788 unsigned int nmax = 0;
789 for (unsigned int i = 0; i < nSamp; ++i) {
790 float smp = digits[i];
791 if (smp - dmin < epsilon) {
792 ++nmin;
793 }
794 if (dmax - smp < epsilon) {
795 ++nmax;
796 }
797 if (smp < abovemin && smp > dmin) {
798 abovemin = smp;
799 }
800 if (smp > belowmax && smp < dmax) {
801 belowmax = smp;
802 }
803 }
804 // more than two different values - shift index by 2, i.e. use thresholds for non-const levels
805 int gainInd = (abovemin != dmax || belowmax != dmin) ? gain + 2 : gain;
806 bool big_jump = (dmaxmin > delta[gainInd]);
807 bool max_in_middle = (pmax > 0 && pmax < nSamp - 1);
808 bool min_in_middle = (pmin > 0 && pmin < nSamp - 1);
809
810 if (nover > 1 && belowmax < level1) { // at least two saturated. others - close to pedestal
811 error = 3;
812 } else if (nmax + nmin == nSamp && big_jump) {
813 if (nmax > 1 && nmin > 1) { // at least 2 samples at two distinct levels
814 error = 4;
815 } else if (nmax == 1) {
816 if (max_in_middle) { // jump up in one sample, but not at the edge
817 error = 5;
818 }
819 } else if (nmin == 1) { // jump down in one sample
820 error = 6;
821 }
822 }
823 if (error == 0 && dmaxmin > narrowLevel[gain]) {
824 float secondMax = dmaxmin * secondMaxLevel;
825 float dminPlus = dmin + secondMax;
826 float dmaxMinus = dmax - secondMax;
827 if (not_gap) { // jumps above two (or one) neighbour samples
828 if (max_in_middle && std::max(digits[pmax - 1], digits[pmax + 1]) < dminPlus) {
829 error = 7; // jump up in one sample in the middle, which is much higher than two neighbours
830 } else if (min_in_middle && std::min(digits[pmin - 1], digits[pmin + 1]) > dmaxMinus) {
831 error = 8; // jump down in one sample, which is much lower than two neighbours
832 } else if (big_jump && gain == 0) { // check first and last sample only in low gain
833 if (pmin == 0 && digits[1] > dmax - secondMax) {
834 error = 10; // jump down in first sample. which is much lower than next one
835 } else if (pmin == nSamp - 1 && digits[pmin - 1] > dmax - secondMax) {
836 error = 11; // jump down in last sample. which is much lower than previous one
837 }
838 }
839 }
840 if (!error && big_jump) { // jumps above all samples
841 if ((max_in_middle || gain == 0) && nmax == 1 && belowmax < dminPlus) {
842 error = 12; // jump up in one sample in the middle, which is much higher than all others
843 } else if ((min_in_middle || gain == 0) && nmin == 1 && abovemin > dmaxMinus) {
844 error = 13; // jump down in one sample, which is much lower than all others (
845 }
846 }
847 }
848 }
849
850 } else {
851 dmin = dmax = 0.0;
852 }
853
854 return error;
855}
#define endmsg
#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)
Helpers for checking error return status codes and reporting errors.
const int nmax(200)
Handle class for recording to StoreGate.
static const InterfaceID IID_ITileRawChannelBuilder("TileRawChannelBuilder", 1, 0)
std::vector< std::pair< TileRawChannel *, const TileDigits * > > Overflows_t
static const Attributes_t empty
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
a typed memory pool that saves time spent allocation small object.
Definition DataPool.h:63
pointer nextElementPtr()
obtain the next available element in pool by pointer pool is resized if its limit has been reached On...
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
virtual std::vector< IdentifierHash > GetAllCurrentHashes() const override final
Returns a collection of all hashes availiable in this IDC.
virtual const T * indexFindPtr(IdentifierHash hashId) const override final
return pointer on the found entry or null if out of range using hashed index - fast version,...
This is a "hash" representation of an Identifier.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
Class that holds Data Quality fragment information and provides functions to extract the data quality...
bool isAdcDQgood(int partition, int drawer, int ch, int gain) const
returns status of single ADC returns False if there are any errors
int trigType() const
Trigger type.
uint32_t calibMode() const
Calibration mode.
bool incompleteDigits() const
A few extra items (from TileBeamInfoProvider).
const uint32_t * cispar() const
CIS parameters.
bool isBiGain() const
returns gain mode of run
const std::vector< float > & samples() const
Definition TileDigits.h:58
Hash table for Tile fragments (==drawers ==collections in StoreGate).
virtual StatusCode finalize()
ToolHandle< TileCondToolTiming > m_tileToolTiming
Overflows_t & getOverflowedChannels(void)
static const char * BadPatternName(float ped)
TileRawChannelUnit::UNIT m_rChUnit
StatusCode build(const TileDigitsCollection *collection, const EventContext &ctx)
static const InterfaceID & interfaceID()
AlgTool InterfaceID.
float m_timeMinThresh
correct amplitude is time is above time min threshold
std::unique_ptr< TileMutableRawChannelContainer > m_rawChannelCnt
SG::ReadHandleKey< TileDQstatus > m_DQstatusKey
static double correctTime(double phase, bool of2=true)
Time correction factor.
virtual ~TileRawChannelBuilder()
Destructor.
float m_ADCmaskValueMinusEps
indicates channels which were masked in background dataset
SG::ReadHandleKey< TileRawChannelContainer > m_DSPContainerKey
float m_ampMinThresh
correct amplitude if it's above amplitude threshold (in ADC counts)
virtual StatusCode createContainer(const EventContext &ctx)
Create container in SG with name given by parameter (m_rawChannelContainerKey).
ServiceHandle< TileCablingSvc > m_cablingSvc
Name of Tile cabling service.
float m_timeMaxThresh
correct amplitude is time is below time max threshold
static int CorruptedData(int ros, int drawer, int channel, int gain, const std::vector< float > &digits, float &dmin, float &dmax, float ADCmaxMinusEps, float ADCmaskValueMinusEps)
ToolHandle< TileCondIdTransforms > m_tileIdTransforms
void fill_drawer_errors(const EventContext &ctx, const TileDigitsCollection *collection)
std::string getTileRawChannelContainerID(void)
virtual StatusCode commitContainer(const EventContext &ctx)
Commit RawChannelContiner in SG and make const.
static double correctAmp(double phase, bool of2=true)
Amplitude correction factor according to the time when using weights for tau=0 without iterations.
const TileCablingService * m_cabling
TileCabling instance.
virtual TileRawChannel * rawChannel(const TileDigits *digits, const EventContext &ctx)
Builder virtual method to be implemented by subclasses.
void initLog(const EventContext &ctx)
ToolHandleArray< ITileRawChannelTool > m_noiseFilterTools
virtual StatusCode initialize()
Initializer.
ToolHandle< TileCondToolEmscale > m_tileToolEmscale
SG::WriteHandleKey< TileRawChannelContainer > m_rawChannelContainerKey
Gaudi::Property< std::vector< int > > m_demoFragIDs
TileRawChannelBuilder(const std::string &type, const std::string &name, const IInterface *parent)
Constructor.
void assign(const HWIdentifier &id, float amplitude, float time, float quality, float ped=0.0)
float pedestal(void) const
void setPedestal(float ped)
uint32_t getLvl1Type() const
Getter for level1 type.
uint32_t getDetEvType() const
Getter for detector event type.
void setLvl1Type(uint32_t lvl1Type)
Setter for level1 type from ROD header.
void setLvl1Id(uint32_t lvl1Id)
Setter for level1 id from ROD header.
void setRODBCID(uint32_t rodBCID)
Setter for BCID from ROD header.
void setDetEvType(uint32_t detEvType)
Setter for detector event type from ROD header.
uint32_t getLvl1Id() const
Getter for level1 id.
uint32_t getRODBCID() const
Getter for BCID from.
uint32_t get_bsflags() const
HWIdentifier adc_HWID(void) const
Definition TileRawData.h:53
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.