ATLAS Offline Software
Loading...
Searching...
No Matches
TileCisDefaultCalibTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6
7// Tile includes
12
13// Gaudi includes
14#include "GaudiKernel/IToolSvc.h"
15#include "GaudiKernel/MsgStream.h"
16#include "GaudiKernel/Service.h"
17#include "GaudiKernel/ThreadLocalContext.h"
18
19
20// Athena includes
23
24#include "TFile.h"
25#include "TTree.h"
26#include "TObjString.h"
27#include "TF1.h"
28#include "TGraphErrors.h"
29#include "TMap.h"
30#include "TMath.h"
31#include <cmath>
32
33TileCisDefaultCalibTool::TileCisDefaultCalibTool(const std::string& type, const std::string& name,
34 const IInterface* pParent)
35 : AthAlgTool(type, name, pParent)
36 , m_tileHWID(nullptr)
37 , m_cabling(nullptr)
38 , m_cablingSvc("TileCablingSvc", name)
39 , m_scanMap(nullptr)
40 , m_scanMapRMS(nullptr)
41 , m_tileInfo(nullptr)
42{
43 declareInterface<ITileCalibTool>(this);
44
45 declareProperty("NtupleID", m_ntupleID = "h3000");
46
47 declareProperty("removePed", m_removePed = true);
48 declareProperty("useSmallCap", m_useSmallCap = false);
49 declareProperty("phaseMin", m_phaseMin = -10);
50 declareProperty("phaseMax", m_phaseMax = 300);
51
52 declareProperty("maxPed", m_maxPed = 45.);
53
54 declareProperty("chargeMinHi", m_chargeMinHi = 0.5);
55 declareProperty("chargeMaxHi", m_chargeMaxHi = 12.5);
56 declareProperty("chargeMinLo", m_chargeMinLo = 31.5);
57 declareProperty("chargeMaxLo", m_chargeMaxLo = 800.0);
58
59 declareProperty("linfitMinHi", m_linfitMinHi = 3.0);
60 declareProperty("linfitMaxHi", m_linfitMaxHi = 10.0);
61 declareProperty("linfitMinLo", m_linfitMinLo = 300.0);
62 declareProperty("linfitMaxLo", m_linfitMaxLo = 700.0);
63
64 declareProperty("linfitMinHiDemo", m_linfitMinHiDemo = 6.0);
65 declareProperty("linfitMaxHiDemo", m_linfitMaxHiDemo = 20.0);
66 declareProperty("linfitMinLoDemo", m_linfitMinLoDemo = 300.0);
67 declareProperty("linfitMaxLoDemo", m_linfitMaxLoDemo = 700.0);
68
69 declareProperty("doSampleChecking", m_doSampleChecking = true); // do sample checking by default
70 declareProperty("TileDQstatus", m_dqStatusKey = "TileDQstatus");
71 declareProperty("TileInfoName", m_infoName = "TileInfo");
72
73 declareProperty("FragIDsDemonstrators", m_fragIDsDemonstrators, "List of Tile frag IDs of demonstrators, which have different CIS circuits than the legacy ones");
74
75 // Initialize arrays for results
81
82 // Initialize sample check arrays
85
89}
90
92
93 delete[] m_calib;
94 delete[] m_qflag;
95 delete[] m_nDAC;
96 delete[] m_nDigitalErrors;
97 delete[] m_chi2;
98 delete[] m_edgeSample;
99 delete[] m_nextToEdgeSample;
100 delete[] m_sampleBit;
101 delete[] m_bitStatus;
102 delete[] m_numSamp;
103
104}
105
107 ATH_MSG_INFO( "initialize()" );
108
109 // get TileHWID helper
110 CHECK( detStore()->retrieve(m_tileHWID) );
111
112 // get TileCabling Service
113 CHECK( m_cablingSvc.retrieve() );
114 m_cabling = m_cablingSvc->cablingService();
115 int runPeriod = m_cabling->runPeriod();
116
117 if (runPeriod==3) {
118 if ( m_fragIDsDemonstrators.empty()) {
119 m_fragIDsDemonstrators.push_back (0x10d); // LBA14 is demonstrator in RUN3
120 }
121 }
122
123 if ( not m_fragIDsDemonstrators.empty() ) {
124
126
127 std::ostringstream os;
128 for (int fragID : m_fragIDsDemonstrators) {
129 os << " 0x" << std::hex << fragID << std::dec;
130 }
131
132 ATH_MSG_INFO("Special settings in histograms for demonstrator modules (frag IDs):" << os.str());
133 }
134
135 CHECK( m_dqStatusKey.initialize() );
136
137 // get TileInfo
138 CHECK( detStore()->retrieve(m_tileInfo, m_infoName) );
139
140 // set important constants
141 m_dac2Charge[0] = 100.* 2.0 * 4.096 / double(m_tileInfo->ADCmax()); // 100 pF * 2 for legacy
142 m_dac2Charge[1] = 5.2 * 2.0 * 4.096 / double(m_tileInfo->ADCmax()); // effective 5.2 pF * 2 for 5 pF capacitor
143 m_dac2Charge[2] = 200.* 4.096 / double(m_tileInfo->ADCmax()); // 200 pF for demonstrator - similar to legacy
144 m_dac2Charge[3] = 5.2 * 4.096 / double(m_tileInfo->ADCmax()); // effective value of small capacitor is twice smaller for demonstrator
145
146 // the same overflow limit for the moment, can be changed later
147 m_maxAmp[0] = m_tileInfo->ADCmax() - m_maxPed;
148 m_maxAmp[1] = m_tileInfo->ADCmax() - m_maxPed;
149 m_maxAmp[2] = m_tileInfo->ADCmax() - m_maxPed;
150 m_maxAmp[3] = m_tileInfo->ADCmax() - m_maxPed;
151
152 m_defaultCalib[0] = 1.29;
153 m_defaultCalib[1] = 81.8;
155 m_defaultCalib[3] = m_defaultCalib[1] / 2;
156
157 if (!m_useSmallCap) {
161 m_chargeMin[3] = m_chargeMinHi * 2;
162 }
163
167 m_chargeMax[3] = m_chargeMaxHi * 2; // high gain in demo is up to 25 pC istead of 12.5 pC
168
173
178
179 ATH_CHECK( m_rawChannelContainerKey.initialize() );
181
182 return StatusCode::SUCCESS;
183}
184
185StatusCode TileCisDefaultCalibTool::initNtuple(int runNumber, int runType, TFile * rootFile) {
186 ATH_MSG_INFO( "initialize(" << runNumber << "," << runType << "," << rootFile << ")" );
187
188 return StatusCode::SUCCESS;
189}
190
192
193 ATH_MSG_DEBUG( "execute()" );
194
195 // Get the DQ digital check information
196 const EventContext& ctx = Gaudi::Hive::currentContext();
197 const TileDQstatus* theDQstatus = SG::makeHandle (m_dqStatusKey, ctx).get();
198
199 // Get event's CIS parameters
200 const uint32_t *cispar = theDQstatus->cispar();
201 uint32_t dac = cispar[6];
202 uint32_t phase = cispar[5];
203 uint32_t cap = cispar[7];
204 int cap_ind = (cap > 10) ? 0 : 1; // 100 pF or 5 pF
205 double chargeAll[4];
206 for (int i=0; i<4; ++i)
207 chargeAll[i] = dac * m_dac2Charge[i];
208
209 // Check if event should be used in calibration
210 bool pass = true;
211 if (cap == 100 && m_useSmallCap)
212 pass = false;
213 else if (cap == 5 && !m_useSmallCap) pass = false;
214 if (phase > m_phaseMax) pass = false;
215 if (phase < m_phaseMin) pass = false;
216 if (cispar[6] == 120) { // Reject garbage events at the beginning of files. This DAQ
217 pass = false; // setting isn't used during a normal CIS scan FYI.
218 }
219
220 // Get TileRawChannelContainer
222 ATH_CHECK( container.isValid() );
223
224 // Create iterator over RawChannelContainer
225 TileRawChannelContainer::const_iterator itColl = (*container).begin();
226 TileRawChannelContainer::const_iterator itCollEnd = (*container).end();
228
229 if (pass) {
230
231 // Go through all TileRawChannelCollections
232 for (; itColl != itCollEnd; ++itColl) {
233
234 int fragId = (*itColl)->identify();
235 bool demonstrator = (std::binary_search(m_fragIDsDemonstrators.begin(), m_fragIDsDemonstrators.end(), fragId));
236 int gain_offset = (demonstrator) ? 2 : 0;
237 double charge = chargeAll[cap_ind+gain_offset];
238
239 // go through all TileRawChannels in collection
240 it = (*itColl)->begin();
241 itEnd = (*itColl)->end();
242
243 for (; it != itEnd; ++it) {
244
245 // get hardware id to identify adc channel
246 HWIdentifier hwid = (*it)->adc_HWID();
247 int ros = m_tileHWID->ros(hwid); // LBA=1 LBC=2 EBA=3 EBC=4
248 int drawer = m_tileHWID->drawer(hwid); // 0 to 63
249 int chan = m_tileHWID->channel(hwid); // 0 to 47 channel not PMT
250 int gain = m_tileHWID->adc(hwid); // low=0 high=1
251
252 // check if channel is connected
253 // if( !chanIsConnected(ros,chan) ) continue;
254
255 // Is channel empty? DQ version
256 if (theDQstatus->isChEmpty(ros, drawer, chan)) continue;
257
258 // find dac maps for adc channel
259 TDACIntMap *NEvtDacMap = (m_NEvtMap)[hwid];
260 TDACIntMap *NDigitalErrorsDacMap = (m_NDigitalErrorsMap)[hwid];
261 TDACDoubleMap *MeanDacMap = (m_MeanMap)[hwid];
262 TDACDoubleMap *MeanSqDacMap = (m_MeanSqMap)[hwid];
263
264 // create new dac maps if they don't exist
265 if (NEvtDacMap == nullptr) {
266 NEvtDacMap = (m_NEvtMap)[hwid] = new TDACIntMap;
267 NDigitalErrorsDacMap = (m_NDigitalErrorsMap)[hwid] = new TDACIntMap;
268 MeanDacMap = (m_MeanMap)[hwid] = new TDACDoubleMap;
269 MeanSqDacMap = (m_MeanSqMap)[hwid] = new TDACDoubleMap;
270 }
271
272 // check that charge is less than chargeMax (depends on gain)
273 if (charge > m_chargeMin[gain+gain_offset] && charge < m_chargeMax[gain+gain_offset]) {
274
275 double amp = (*it)->amplitude();
276
277 // Hack to get rid of pedestal events: need more sophisticated method!
278 if ( (amp > 40.) || (gain == 0 && amp > 6.) || (!m_removePed) ) {
279
280 // Digital error check
281 if (!(theDQstatus->isAdcDQgood(ros, drawer, chan, gain))) {
282 ATH_MSG_DEBUG( "Skipping Module: " << ros << drawer + 1
283 << " channel: " << chan
284 << " ADC: " << gain
285 << " due to DQ error found." );
286
287 (*NDigitalErrorsDacMap)[dac] += 1;
288 continue;
289
290 } else {
291
292 // increment entries for current dac value
293 (*NEvtDacMap)[dac] += 1;
294 (*MeanDacMap)[dac] += amp;
295 (*MeanSqDacMap)[dac] += amp*amp;
296 }
297
298 }
299
300 } // end if (min < charge < max)
301 }
302 }
303 } // end if pass
304
305 if (m_doSampleChecking) {
306 // Get TileDigitsContainer
308 ATH_CHECK( digContainer.isValid() );
309
310 // Create iterator over RawDigitsContainer
311 TileDigitsContainer::const_iterator digItColl = digContainer->begin();
312 TileDigitsContainer::const_iterator digItCollEnd = digContainer->end();
313
314 for (; digItColl != digItCollEnd; ++digItColl) {
315
316 TileDigitsCollection::const_iterator digIt = (*digItColl)->begin();
317 TileDigitsCollection::const_iterator digItEnd = (*digItColl)->end();
318
319 if (digIt != digItEnd) {
320
321 int fragId = (*digItColl)->identify();
322 bool demonstrator = (std::binary_search(m_fragIDsDemonstrators.begin(), m_fragIDsDemonstrators.end(), fragId));
323 int gain_offset = (demonstrator) ? 2 : 0;
324 double charge = chargeAll[cap_ind+gain_offset];
325
326 HWIdentifier adc_id = (*digIt)->adc_HWID();
327 int ros = m_tileHWID->ros(adc_id); // LBA=1 LBC=2 EBA=3 EBC=4
328 int drawer = m_tileHWID->drawer(adc_id); // 0 to 63
329
330 // not clear how to handle this. if not 7 get off? MM - 4 June 2009
331 int numSamples = (*digIt)->NtimeSamples();
332 if (numSamples != 7) {
333 m_doSampleChecking = false;
334 break;
335 }
336
337 for (; digIt != digItEnd; ++digIt) {
338
339 adc_id = (*digIt)->adc_HWID();
340 int chan = m_tileHWID->channel(adc_id); // 0 to 47 channel not PMT
341 int gain = m_tileHWID->adc(adc_id); // low=0 high=1
342
343 std::vector<float> theDigits = (*digIt)->samples();
344
345 //MM - skip channels with digital errors
346 if (!(theDQstatus->isAdcDQgood(ros, drawer, chan, gain))) {
347 continue;
348 }
349
350 // Loop over samples for bit analysis
351 // We don't need to use the same cuts as the "Edge Sample" analysis
352 for(unsigned int sampNum = 0; sampNum < theDigits.size(); sampNum++) {
353
354 // Count the total number of samples taken by an ADC
355 m_numSamp[ros][drawer][chan][gain] += 1;
356 int k = 0;
357 int quotient = theDigits[sampNum];
358
359 // convert sample to binary number
360 while(quotient!=0) {
361 if((quotient % 2) == 1) {
362 // If the bit is one, store info in the array
363 m_sampleBit[ros][drawer][chan][gain][k] += 1;
364 }
365
366 quotient = quotient / 2;
367 k += 1;
368 } // end binary conversion
369 } //end sample loop
370
371 if (pass && charge > m_linfitMin[gain+gain_offset] && charge < m_linfitMax[gain+gain_offset]) {
372
373 std::vector<float> theDigits = (*digIt)->samples();
374 float maxSampVal = -1.;
375 int maxSampNum = -1;
376
377 for (unsigned int sampNum = 0; sampNum < theDigits.size(); sampNum++) {
378 if (theDigits[sampNum] > maxSampVal) {
379 maxSampVal = theDigits[sampNum];
380 maxSampNum = sampNum + 1;
381 }
382 }
383
384 if (maxSampNum == 1 || maxSampNum == 7) {
385 m_edgeSample[ros][drawer][chan][gain] = 1;
386 } else if (maxSampNum == 2 || maxSampNum == 6) {
387 m_nextToEdgeSample[ros][drawer][chan][gain] = 1;
388 }
389
390 } // end digits iterator
391 }
392 } // end digits collections
393
394 } // end if pass
395
396 } // end m_doSampleChecking
397
398 return StatusCode::SUCCESS;
399}
400
402
403 ATH_MSG_INFO( "finalizeCalculations()" );
404
405 // The values of 81.454 and 1.295 were derived from runs:
406 // 72652 73305 72653 72661 79259 78023 79781 78026
407 // to calibrate the detector, and looking at the mean good calibration value
408 double meanCalib[4];
409 meanCalib[0] = 1.295;
410 meanCalib[1] = 81.454;
411 // values for demonstrator were taken from run 416970
412 meanCalib[2] = 1.26282;
413 meanCalib[3] = 40.9303;
414
415 // hardware id (key to maps)
416 HWIdentifier hwid;
417
418 // dac maps
419 TDACDoubleMap* MeanDacMap;
420 TDACDoubleMap* MeanSqDacMap;
421 TDACIntMap* NEvtDacMap;
422 TDACIntMap* NDigitalErrorsDacMap;
423
424 // count number of points in dac map (for TGraph)
425 int npt, pt;
426
427 // temporary objects for loop
428 TGraphErrors* gr;
429 TGraphErrors* grrms;
430 uint32_t dac;
431 double charge, mean, mean2, meansq, rms, err, ratio;
432 int nevt, ndigerr = 0;
433 int badPts;
434 double maxRMS;
435 double maxPointInFitRange;
436
437 // linear fit for the calibration factor
438 TF1 *fslope = new TF1("fslope", "[0]*x", 0, 1000);
439
440 m_scanMap = new TMap(20000, 1);
441 m_scanMapRMS = new TMap(20000, 1);
442
443 // iterators over adc maps
444 TAdcDoubleMapIter adcIter(m_MeanMap.begin());
445 TAdcDoubleMapIter adcIterE(m_MeanMap.end());
446
447 // loop over all adcs
448 for (; adcIter != adcIterE; ++adcIter) {
449 hwid = (adcIter)->first;
450 MeanDacMap = (adcIter)->second;
451 MeanSqDacMap = m_MeanSqMap[hwid];
452 NEvtDacMap = m_NEvtMap[hwid];
453 NDigitalErrorsDacMap = m_NDigitalErrorsMap[hwid];
454
455 int ros = m_tileHWID->ros(hwid); // LBA=1 LBC=2 EBA=3 EBC=4
456 int drawer = m_tileHWID->drawer(hwid); // 0 to 63
457 int chan = m_tileHWID->channel(hwid); // 0 to 47 channel not PMT
458 int gain = m_tileHWID->adc(hwid); // low=0 high=1
459
460 int fragId = (ros << 8) | drawer;
461 bool demonstrator = (std::binary_search(m_fragIDsDemonstrators.begin(), m_fragIDsDemonstrators.end(), fragId));
462 int gain_offset = (demonstrator) ? 2 : 0;
463 int gain_ind = gain + gain_offset;
464 int cap_ind = (m_useSmallCap) ? 1 : 0;
465 double dac2Charge = m_dac2Charge[cap_ind+gain_offset];
466 double prevMean = 0.0, deltaMean = 0.0, averSlope = 0.0;
467 double prevCharge = 0.0, deltaCharge = 0.0, deltaChargeMin = 999999.;
468 double prevSlope = 0.0, slope = 0.0, maxGoodAmp=0.0, minGoodAmp = 999999.;
469 double maxGoodCharge = m_linfitMax[gain_ind], minGoodCharge = 999999.;
470 double maxCharge = 0.0, minCharge = 999999., maxAmp = 0.0, minAmp = 999999.;
471 bool signalInRange = true;
472 int nptGood = 0;
473
474 // find number of points in graph for this adc
475 npt = MeanDacMap->size();
476 m_nDAC[ros][drawer][chan][gain] = npt;
477 gr = new TGraphErrors(npt);
478 grrms = new TGraphErrors(npt);
479
480 if (npt == 0) {
481 m_calib[ros][drawer][chan][gain] = 0;
482 m_chi2[ros][drawer][chan][gain] = 0.0;
483 ATH_MSG_DEBUG( "npt==0 for adc channel "
484 << ros << "/" << drawer << "/" << chan << "/" << gain );
485 } else {
486
487 // update quality flag: adc channel is included in run
488 setBit(includedBit, m_qflag[ros][drawer][chan][gain]);
489
490 // iterator over dacs
491 TDACDoubleMapIter dacIter((*MeanDacMap).begin());
492 TDACDoubleMapIter dacIterE((*MeanDacMap).end());
493
494 // initialize current point
495 pt = 0;
496 badPts = 0;
497 maxPointInFitRange = 0.0;
498 maxRMS = 0.0;
499 for (; dacIter != dacIterE; ++dacIter) {
500 dac = (dacIter)->first;
501 mean = (dacIter)->second;
502 meansq = (*MeanSqDacMap)[dac];
503 nevt = (*NEvtDacMap)[dac];
504 ndigerr = (*NDigitalErrorsDacMap)[dac];
505
506 mean = mean / nevt;
507 mean2 = mean * mean;
508 meansq = meansq / nevt;
509
510 rms = (meansq <= mean2) ? 0. : sqrt(meansq - mean2);
511 err = sqrt(rms * rms / nevt + 0.5 * 0.5); // 0.5 is the absolute systematic uncertainty on the measurement
512
513 // find charge for this dac
514 charge = dac * dac2Charge;
515
516 if (mean < m_maxAmp[gain_ind] ) {
517 if (charge>maxCharge) maxCharge = charge;
518 if (charge<minCharge) minCharge = charge;
519 if (mean>maxAmp) maxAmp = mean;
520 if (mean<minAmp) minAmp = mean;
521 }
522
523 deltaMean = mean - prevMean;
524 deltaCharge = charge - prevCharge;
525 if (deltaCharge != 0) {
526 slope = deltaMean / deltaCharge;
527 if (deltaCharge<deltaChargeMin) deltaChargeMin = deltaCharge;
528 }
529
530 // check for problems in calibration range
531 if (charge > m_linfitMin[gain_ind] && charge < m_linfitMax[gain_ind]) {
532 if (mean > m_maxAmp[gain_ind] ) {
533 signalInRange = false;
534 ATH_MSG_DEBUG( "Too high amp in "
535 << ros << "/" << drawer
536 << "/" << chan << "/" << gain
537 << " charge " << charge
538 << " amp " << mean
539 << " => Ignoring this point ");
540 } else {
541 if (prevSlope != 0) {
542 double R = std::abs(slope/prevSlope);
543 double R1 = (nptGood>1) ? (slope / averSlope) : 1.;
544 if (R<0.025 || (nptGood>2 && R1>2.5 && R>2.5) || (nptGood==2 && R1>4.1) ) {
545 // amplitude not changed (within +/- 2.5%) or changed too much, ignore this point
546 ATH_MSG_WARNING( "Wrong amp in "
547 << ros << "/" << drawer
548 << "/" << chan << "/" << gain
549 << " charge " << charge
550 << " amp " << mean
551 << " dC " << deltaCharge
552 << " dA " << deltaMean
553 << " expected dA " << prevSlope*deltaCharge
554 << " aver dA " << averSlope*deltaCharge
555 << " rms " << rms
556 << " => Removing this point ");
557 continue;
558 }
559 }
560 if (signalInRange) {
561 maxGoodCharge = charge;
562 if (charge<minGoodCharge) minGoodCharge = charge;
563 if (mean>maxGoodAmp) maxGoodAmp = mean;
564 if (mean<minGoodAmp) minGoodAmp = mean;
565 }
566 prevMean = mean;
567 prevCharge = charge;
568 if (slope > 0) {
569 prevSlope = slope;
570 }
571 ++nptGood;
572 averSlope += (slope - averSlope) / nptGood;
573 if (slope < 0) {
574 --nptGood;
575 }
576 if (averSlope < 0) {
577 nptGood = 0;
578 averSlope = 0.;
579 }
580 }
581
582 if (rms < 0.01) badPts++;
583 if (mean > maxPointInFitRange) maxPointInFitRange = mean;
584 if (rms > maxRMS) maxRMS = rms;
585
586 } else if (charge <= m_linfitMin[gain_ind] && charge > 0) {
587
588 slope = mean/charge;
589 if (prevSlope !=0 ) {
590 double R = slope / prevSlope;
591 if (R>0.4 && R<2.5) {
592 prevMean = mean;
593 prevCharge = charge;
594 prevSlope = (prevSlope + slope)/2.;
595 }
596 } else {
597 double R = slope / m_defaultCalib[gain_ind];
598 if (R>0.4 && R<2.5) {
599 prevSlope = slope;
600 }
601 }
602 } else if (nptGood < 3 && signalInRange && charge >= m_linfitMax[gain_ind] && mean < m_maxAmp[gain_ind] && prevSlope != 0) {
603 double R = slope/prevSlope;
604 if (R>0.9 && R<1.1) {
605 maxGoodCharge = charge;
606 if (charge<minGoodCharge) minGoodCharge = charge;
607 ++nptGood;
608 }
609 }
610
611 // set point and errors in tgraph
612 gr->SetPoint(pt, charge, mean);
613 gr->SetPointError(pt, 0.0, err);
614
615 grrms->SetPoint(pt, charge, mean);
616 grrms->SetPointError(pt, 0.0, rms);
617
618 pt++;
619 } // end of for all DAC values
620
621 // remove empty points at the end
622 for (int i=npt-1; i>=pt; --i)
623 gr->RemovePoint(i);
624
625 slope = (prevMean > 0 && prevCharge > 0) ? prevMean/prevCharge : m_defaultCalib[gain_ind];
626 fslope->SetParameter(0, slope);
627 if (maxGoodCharge > m_linfitMax[gain_ind]) {
628 ATH_MSG_WARNING( "Extending fit range for "
629 << ros << "/" << drawer
630 << "/" << chan << "/" << gain
631 << " up to " << maxGoodCharge
632 << " pC; N good points " << nptGood);
633 }
634
635 if (deltaChargeMin>999.) deltaChargeMin = 0.1;
636 else deltaChargeMin *= 0.5;
637 gr->Fit("fslope", "q", "", minGoodCharge-deltaChargeMin, maxGoodCharge+deltaChargeMin);
638 ATH_MSG_VERBOSE( "Fit for "
639 << ros << "/" << drawer
640 << "/" << chan << "/" << gain
641 << " minC " << minGoodCharge
642 << " maxC " << maxGoodCharge
643 << " deltaC " << deltaChargeMin*2
644 << " minA " << minGoodAmp
645 << " maxA " << maxGoodAmp
646 << " nptG " << nptGood
647 << " slope " << slope
648 << " fit " << fslope->GetParameter(0));
649 if (nptGood < 3) {
650 const char * bms[3] = {"No good points to fit ",
651 "Only one point to fit ",
652 "Only two points to fit "};
653 const char * ems[3] = {" => Put zero slope",
654 "",
655 ""};
656 ATH_MSG_WARNING( bms[nptGood]
657 << ros << "/" << drawer
658 << "/" << chan << "/" << gain
659 << " charge " << charge
660 << " amp " << mean
661 << " slope " << slope
662 << " fit " << fslope->GetParameter(0)
663 << ems[nptGood]);
664 }
665 slope = (nptGood > 0) ? fslope->GetParameter(0) : 0.;
666
667 if (slope > 0.) {
668 if (maxGoodCharge > minGoodCharge) {
669 averSlope = (maxGoodAmp-minGoodAmp)/(maxGoodCharge-minGoodCharge);
670 if (slope > 2.5 * averSlope) {
671 ATH_MSG_WARNING( "Average slope and fit slope do not match "
672 << ros << "/" << drawer
673 << "/" << chan << "/" << gain
674 << " charge " << charge
675 << " amp " << mean
676 << " avslope " << averSlope
677 << " fit " << slope
678 << " => Put zero slope");
679 slope = 0.;
680 }
681 } else if (maxCharge > minCharge) {
682 averSlope = (maxAmp-minAmp)/(maxCharge-minCharge);
683 if (slope > 10. * averSlope) {
684 ATH_MSG_WARNING( "AVERAGE slope and fit slope do not match "
685 << ros << "/" << drawer
686 << "/" << chan << "/" << gain
687 << " charge " << charge
688 << " amp " << mean
689 << " AVslope " << averSlope
690 << " fit " << slope
691 << " => Put zero slope");
692 slope = 0.;
693 }
694 }
695 }
696
697 m_nDigitalErrors[ros][drawer][chan][gain] = ndigerr;
698
699 // Set this bit if there aren't any digital errors
700 if (ndigerr == 0) {
701 setBit(digiErrorBit, m_qflag[ros][drawer][chan][gain]);
702 }
703
704 m_calib[ros][drawer][chan][gain] = slope;
705 if (fslope->GetNDF() == 0)
706 m_chi2[ros][drawer][chan][gain] = 0.0;
707 else
708 m_chi2[ros][drawer][chan][gain] = fslope->GetChisquare() / fslope->GetNDF();
709
710 // Set this bit if there is a good Chi2 probability
711 if (TMath::Prob(fslope->GetChisquare(), fslope->GetNDF()) > 2.e-6) {
712 setBit(probChi2Bit, m_qflag[ros][drawer][chan][gain]);
713 }
714
715 // update quality flag if calibration is successful
716 if (!badPts && fslope->GetNDF() > 0) setBit(calibratedBit, m_qflag[ros][drawer][chan][gain]);
717
718 // update quality flag if calibration is within 5% of nominal
719 // saved for legacy support
720 ratio = (slope / m_defaultCalib[gain_ind]);
721 if (ratio > 0.95 && ratio < 1.05) setBit(rangeBit, m_qflag[ros][drawer][chan][gain]);
722
723 // update quality flag if calibration if the probability of this calibration
724 // constant, given a 1.6% gaussian-sigma of the calibration constants, is greater
725 // than 1/10000 (number of channels)
726 //
727 // Mathemica code: NSolve[Erf[x/(1.6*Sqrt[2])] == 0.9999, x]
728 // x -> 6.22495
729 ratio = (slope / meanCalib[gain_ind]);
730 if (ratio > 0.9378 && ratio < 1.0623) setBit(probBit, m_qflag[ros][drawer][chan][gain]);
731
732 // If the maximum response in the fit range is less than 600 ADC counts, then
733 // all the response in most likely noise
734 if (maxPointInFitRange > 600) {
735 setBit(noiseBit, m_qflag[ros][drawer][chan][gain]);
736 }
737
738 // RMS criteria. If any collection of injections at a fixed-charge has
739 // an RMS less than 5 ADC counts, then set this bit.
740 if (maxRMS < 5.0) {
741 setBit(injRMSBit, m_qflag[ros][drawer][chan][gain]);
742 }
743
744 // set the sample check bits
745
746 // this bit is set if there were no events found in the fit range
747 // with the maximum sample value in the first or last sample
748
749 if (m_edgeSample[ros][drawer][chan][gain] == 0) {
750 setBit(edgeSamp, m_qflag[ros][drawer][chan][gain]);
751 }
752 // this bit is set if there were no events found in the fit range
753 // with the maximum sample value in the second or sixth sample
754 if (m_nextToEdgeSample[ros][drawer][chan][gain] == 0) {
755 setBit(nextToEdgeSamp, m_qflag[ros][drawer][chan][gain]);
756 }
757
758 // Determine failure/passing of StuckBit quality flag
759 // And store information about bits in an array
760 // which will be written to the ntuple
761 int NoStuckBit = 1;
762 for(int i = 0; i < NBITS; i++) {
763 // If a bit is stuck at zero...
764 if(m_sampleBit[ros][drawer][chan][gain][i] == 0 && (m_numSamp[ros][drawer][chan][gain] != 0)) {
765 // write information to m_bitStatus array of shorts
766 // each bit in short corresponds to a bit in an adc
767 // with 6 short bits left over
768 m_bitStatus[ros][drawer][chan][gain][0] += (1<<i);
769 NoStuckBit = 0;
770 ATH_MSG_DEBUG( "\n\nBIT STUCK AT ZERO: "
771 << ros << " " << drawer << " " << chan << " " << gain << " " << i << "\n");
772
773 }
774 // Same for a bit stuck at one
775 else if (m_sampleBit[ros][drawer][chan][gain][i] == m_numSamp[ros][drawer][chan][gain] && (m_numSamp[ros][drawer][chan][gain] != 0)) {
776 m_bitStatus[ros][drawer][chan][gain][1] += (1<<i);
777 NoStuckBit = 0;
778 ATH_MSG_DEBUG( "\n\nBIT STUCK AT ONE: "
779 << ros << " " << drawer << " " << chan << " " << gain << " " << i << "\n");
780 }
781 } //end bit loop
782
783 // If no stuck bits are found, this adc passes StuckBit m_qflag
784 if(NoStuckBit) {
785 setBit(stuckbitBit, m_qflag[ros][drawer][chan][gain]);
786 }
787
788 gr->SetName("scan_" + arrayString(ros, drawer, chan, gain));
789 grrms->SetName("scan_" + arrayString(ros, drawer, chan, gain));
790
791 m_scanMap->Add(new TObjString("scan" + arrayString(ros, drawer, chan, gain)), gr);
792 m_scanMapRMS->Add(new TObjString("scan" + arrayString(ros, drawer, chan, gain)), grrms);
793 }
794 }
795 return StatusCode::SUCCESS;
796}
797
798StatusCode TileCisDefaultCalibTool::writeNtuple(int runNumber, int runType, TFile* rootFile) {
799
800 ATH_MSG_INFO( "writeNtuple(" << runNumber << "," << runType << "," << rootFile << ")" );
801
802 TTree *t = new TTree(m_ntupleID.c_str(), "TileCalib-Ntuple");
803 t->Branch("RunNumber", &runNumber, "runNo/I");
804 t->Branch("calib", *m_calib, "calib[5][64][48][2]/F");
805 t->Branch("qflag", *m_qflag, "qflag[5][64][48][2]/I");
806 t->Branch("nDAC", *m_nDAC, "nDAC[5][64][48][2]/I");
807 t->Branch("nDigitalErrors", *m_nDigitalErrors, "nDigitalErrors[5][64][48][2]/I");
808 t->Branch("chi2", *m_chi2, "chi2[5][64][48][2]/F");
809 t->Branch("BitStatus", *m_bitStatus, "BitStatus[5][64][48][2][4]/s");
810
811 if (!m_stuckBitsProbs.empty()) {
812 if (m_stuckBitsProbs.retrieve().isFailure()) {
813 ATH_MSG_WARNING("Impossible to get ITileStuckBitsProbsTool and stuck bits probabilities!");
814 } else {
815 m_stuckBitsProbs->saveStuckBitsProbabilities(t);
816 }
817 }
818
819 // Fill with current values (i.e. tree will have only one entry for this whole run)
820 t->Fill();
821 t->Write();
822
823 // Save graphs for all calibrated adc channels
824 m_scanMap->Write("cisScans", TObject::kSingleKey);
825 m_scanMapRMS->Write("cisScansRMS", TObject::kSingleKey);
826
827 return StatusCode::SUCCESS;
828}
829
831
832 ATH_MSG_INFO( "finalize()" );
833
834 return StatusCode::SUCCESS;
835}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
double charge(const T &p)
Definition AtlasPID.h:997
Helpers for checking error return status codes and reporting errors.
#define CHECK(...)
Evaluate an expression and check for errors.
#define gr
Handle class for reading from StoreGate.
#define NBITS
#define NBSTATUS
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
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
virtual bool isValid() override final
Can the handle be successfully dereferenced?
static const unsigned int MAX_ROS
Number of ROSs.
static const unsigned int MAX_GAIN
Number of gains per channel.
static const unsigned int MAX_DRAWER
Number of drawers in ROS 1-4.
static const unsigned int MAX_CHAN
Number of channels in drawer.
int(* m_nextToEdgeSample)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
float(* m_calib)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
const TileCablingService * m_cabling
std::map< HWIdentifier, TDACDoubleMap * >::iterator TAdcDoubleMapIter
virtual StatusCode finalizeCalculations() override
std::vector< int > m_fragIDsDemonstrators
int(* m_sampleBit)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN][NBITS]
int(* m_qflag)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
std::map< uint32_t, int > TDACIntMap
SG::ReadHandleKey< TileRawChannelContainer > m_rawChannelContainerKey
virtual StatusCode initialize() override
int(* m_numSamp)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
int(* m_nDAC)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
std::map< uint32_t, double >::iterator TDACDoubleMapIter
float(* m_chi2)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
int(* m_edgeSample)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
ServiceHandle< TileCablingSvc > m_cablingSvc
std::map< uint32_t, double > TDACDoubleMap
SG::ReadHandleKey< TileDQstatus > m_dqStatusKey
ToolHandle< ITileStuckBitsProbsTool > m_stuckBitsProbs
TString arrayString(int ros, int drawer, int chan, int gain)
void setBit(QualityType qb, int &bitflag)
unsigned short(* m_bitStatus)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN][NBSTATUS]
int(* m_nDigitalErrors)[Tile::MAX_DRAWER][Tile::MAX_CHAN][Tile::MAX_GAIN]
virtual StatusCode execute() override
virtual StatusCode finalize() override
virtual StatusCode writeNtuple(int runNumber, int runType, TFile *rootfile) override
TileCisDefaultCalibTool(const std::string &type, const std::string &name, const IInterface *pParent)
virtual StatusCode initNtuple(int runNumber, int runType, TFile *rootfile) override
SG::ReadHandleKey< TileDigitsContainer > m_digitsContainerKey
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
static int isChEmpty(int partition, int drawer, int ch)
True if channel is not fully implemented.
const uint32_t * cispar() const
CIS parameters.
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
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.