Project
Loading...
Searching...
No Matches
ThresholdCalibratorSpec.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
13
17
18#ifdef WITH_OPENMP
19#include <omp.h>
20#endif
21
22namespace o2
23{
24namespace its
25{
26
28// Define error function for ROOT fitting
29double erf(double* xx, double* par)
30{
31 return (nInjScaled / 2) * TMath::Erf((xx[0] - par[0]) / (sqrt(2) * par[1])) + (nInjScaled / 2);
32}
33
34// ITHR erf is reversed
35double erf_ithr(double* xx, double* par)
36{
37 return (nInjScaled / 2) * (1 - TMath::Erf((xx[0] - par[0]) / (sqrt(2) * par[1])));
38}
39
41// Default constructor
43 : mChipModSel(inpConf.chipModSel), mChipModBase(inpConf.chipModBase)
44{
45 mSelfName = o2::utils::Str::concat_string(ChipMappingITS::getName(), "ITSThresholdCalibrator");
46}
47
49// Default deconstructor
51{
52 // Clear dynamic memory
53
54 delete[] this->mX;
55 this->mX = nullptr;
56
57 if (this->mFitType == FIT) {
58 delete this->mFitHist;
59 this->mFitHist = nullptr;
60 delete this->mFitFunction;
61 this->mFitFunction = nullptr;
62 }
63}
64
67{
68 LOGF(info, "ITSThresholdCalibrator init...", mSelfName);
69
70 mPercentageCut = ic.options().get<short int>("percentage-cut");
71
72 mColStep = ic.options().get<short int>("s-curve-col-step");
73 if (mColStep >= N_COL) {
74 LOG(warning) << "mColStep = " << mColStep << ": saving s-curves of only 1 pixel (pix 0) per row";
75 }
76
77 std::string fittype = ic.options().get<std::string>("fittype");
78 if (fittype == "derivative") {
79 this->mFitType = DERIVATIVE;
80
81 } else if (fittype == "fit") {
82 this->mFitType = FIT;
83
84 } else if (fittype == "hitcounting") {
85 this->mFitType = HITCOUNTING;
86
87 } else {
88 LOG(error) << "fittype " << fittype
89 << " not recognized, please use 'derivative', 'fit', or 'hitcounting'";
90 throw fittype;
91 }
92
93 // Get metafile directory from input
94 try {
95 this->mMetafileDir = ic.options().get<std::string>("meta-output-dir");
96 } catch (std::exception const& e) {
97 LOG(warning) << "Input parameter meta-output-dir not found"
98 << "\n*** Setting metafile output directory to /dev/null";
99 }
100 if (this->mMetafileDir != "/dev/null") {
101 this->mMetafileDir = o2::utils::Str::rectifyDirectory(this->mMetafileDir);
102 }
103
104 // Get ROOT output directory from input
105 try {
106 this->mOutputDir = ic.options().get<std::string>("output-dir");
107 } catch (std::exception const& e) {
108 LOG(warning) << "Input parameter output-dir not found"
109 << "\n*** Setting ROOT output directory to ./";
110 }
111 this->mOutputDir = o2::utils::Str::rectifyDirectory(this->mOutputDir);
112
113 // Get metadata data type from input
114 try {
115 this->mMetaType = ic.options().get<std::string>("meta-type");
116 } catch (std::exception const& e) {
117 LOG(warning) << "Input parameter meta-type not found"
118 << "\n*** Disabling 'type' in metadata output files";
119 }
120
121 this->mVerboseOutput = ic.options().get<bool>("verbose");
122
123 // Get number of threads
124 this->mNThreads = ic.options().get<int>("nthreads");
125
126 // Check fit type vs nthreads (fit option is not thread safe!)
127 if (mFitType == FIT && mNThreads > 1) {
128 throw std::runtime_error("Multiple threads are requested with fit method which is not thread safe");
129 }
130
131 // Machine hostname
132 this->mHostname = boost::asio::ip::host_name();
133
134 // check flag to tag single noisy pix in digital and analog scans
135 this->mTagSinglePix = ic.options().get<bool>("enable-single-pix-tag");
136
137 // get min and max ithr and vcasn (default if not specified)
138 inMinVcasn = ic.options().get<short int>("min-vcasn");
139 inMaxVcasn = ic.options().get<short int>("max-vcasn");
140 inMinIthr = ic.options().get<short int>("min-ithr");
141 inMaxIthr = ic.options().get<short int>("max-ithr");
142 if (inMinVcasn > inMaxVcasn || inMinIthr > inMaxIthr) {
143 throw std::runtime_error("Min VCASN/ITHR is larger than Max VCASN/ITHR: check the settings, analysis not possible");
144 }
145
146 // Get flag to enable most-probable value calculation
147 isMpv = ic.options().get<bool>("enable-mpv");
148
149 // Parameters to operate in manual mode (when run type is not recognized automatically)
150 isManualMode = ic.options().get<bool>("manual-mode");
151 if (isManualMode) {
152 try {
153 manualMin = ic.options().get<short int>("manual-min");
154 } catch (std::exception const& e) {
155 throw std::runtime_error("Min value of the scan parameter not found, mandatory in manual mode");
156 }
157
158 try {
159 manualMax = ic.options().get<short int>("manual-max");
160 } catch (std::exception const& e) {
161 throw std::runtime_error("Max value of the scan parameter not found, mandatory in manual mode");
162 }
163
164 try {
165 manualScanType = ic.options().get<std::string>("manual-scantype");
166 } catch (std::exception const& e) {
167 throw std::runtime_error("Scan type not found, mandatory in manual mode");
168 }
169
170 try {
171 saveTree = ic.options().get<bool>("save-tree");
172 } catch (std::exception const& e) {
173 throw std::runtime_error("Please specify if you want to save the ROOT trees, mandatory in manual mode");
174 }
175
176 // this is not mandatory since it's 1 by default
177 manualStep = ic.options().get<short int>("manual-step");
178
179 // this is not mandatory since it's 0 by default
180 manualMin2 = ic.options().get<short int>("manual-min2");
181
182 // this is not mandatory since it's 0 by default
183 manualMax2 = ic.options().get<short int>("manual-max2");
184
185 // this is not mandatory since it's 1 by default
186 manualStep2 = ic.options().get<short int>("manual-step2");
187
188 // this is not mandatory since it's 5 by default
189 manualStrobeWindow = ic.options().get<short int>("manual-strobewindow");
190
191 // Flag to scale the number of injections by 3 in case --meb-select is used
192 scaleNinj = ic.options().get<bool>("scale-ninj");
193 }
194
195 // Flag to enable the analysis of CRU_ITS data
196 isCRUITS = ic.options().get<bool>("enable-cru-its");
197
198 // Number of injections
199 nInj = ic.options().get<int>("ninj");
201
202 // flag to set the url ccdb mgr
203 this->mCcdbMgrUrl = ic.options().get<std::string>("ccdb-mgr-url");
204 // FIXME: Temporary solution to retrieve ConfDBmap
205 long int ts = o2::ccdb::getCurrentTimestamp();
206 LOG(info) << "Getting confDB map from ccdb - timestamp: " << ts;
208 mgr.setURL(mCcdbMgrUrl);
209 mgr.setTimestamp(ts);
210 mConfDBmap = mgr.get<std::vector<int>>("ITS/Calib/Confdbmap");
211
212 // Parameters to dump s-curves on disk
213 isDumpS = ic.options().get<bool>("dump-scurves");
214 maxDumpS = ic.options().get<int>("max-dump");
215 chipDumpS = ic.options().get<std::string>("chip-dump"); // comma-separated list of chips
216 chipDumpList = getIntegerVect(chipDumpS);
217 if (isDumpS && mFitType != FIT) {
218 LOG(error) << "S-curve dump enabled but `fittype` is not fit. Please check";
219 }
220 if (isDumpS) {
221 fileDumpS = TFile::Open(Form("s-curves_%d.root", mChipModSel), "RECREATE"); // in case of multiple processes, every process will have it's own file
222 if (maxDumpS < 0) {
223 LOG(info) << "`max-dump` " << maxDumpS << ". Dumping all s-curves";
224 } else {
225 LOG(info) << "`max-dump` " << maxDumpS << ". Dumping " << maxDumpS << " s-curves";
226 }
227 if (!chipDumpList.size()) {
228 LOG(info) << "Dumping s-curves for all chips";
229 } else {
230 LOG(info) << "Dumping s-curves for chips: " << chipDumpS;
231 }
232 }
233
234 // flag to enable the calculation of the slope in 2d pulse shape scans
235 doSlopeCalculation = ic.options().get<bool>("calculate-slope");
236 if (doSlopeCalculation) {
237 try {
238 chargeA = ic.options().get<int>("charge-a");
239 } catch (std::exception const& e) {
240 throw std::runtime_error("You want to do the slop calculation but you did not specify charge-a");
241 }
242
243 try {
244 chargeB = ic.options().get<int>("charge-b");
245 } catch (std::exception const& e) {
246 throw std::runtime_error("You want to do the slop calculation but you did not specify charge-b");
247 }
248 }
249
250 // Variable to select from which multi-event buffer select the hits
251 mMeb = ic.options().get<int>("meb-select");
252 if (mMeb > 2) {
253 LOG(error) << "MEB cannot be greater than 2. Please check your command line.";
254 }
255
256 return;
257}
258
260// Get number of active links for a given RU
261short int ITSThresholdCalibrator::getNumberOfActiveLinks(bool* links)
262{
263 int nL = 0;
264 for (int i = 0; i < 3; i++) {
265 if (links[i]) {
266 nL++;
267 }
268 }
269 return nL;
270}
271
273// Get link ID: 0,1,2 for IB RUs / 0,1 for OB RUs
274short int ITSThresholdCalibrator::getLinkID(short int chipID, short int ruID)
275{
276 if (chipID < 432) {
277 return (chipID - ruID * 9) / 3;
278 } else if (chipID >= 432 && chipID < 6480) {
279 return (chipID - 48 * 9 - (ruID - 48) * 112) / 56;
280 } else {
281 return (chipID - 48 * 9 - 54 * 112 - (ruID - 102) * 196) / 98;
282 }
283}
284
286// Get list of chipID (from 0 to 24119) attached to a RU based on the links which are active
287std::vector<short int> ITSThresholdCalibrator::getChipListFromRu(short int ruID, bool* links)
288{
289 std::vector<short int> cList;
290 int a, b;
291 if (ruID < 48) {
292 a = ruID * 9;
293 b = a + 9 - 1;
294 } else if (ruID >= 48 && ruID < 102) {
295 a = 48 * 9 + (ruID - 48) * 112;
296 b = a + 112 - 1;
297 } else {
298 a = 48 * 9 + 54 * 112 + (ruID - 102) * 196;
299 b = a + 196 - 1;
300 }
301
302 for (int c = a; c <= b; c++) {
303 short int lid = getLinkID(c, ruID);
304 if (links[lid]) {
305 cList.push_back(c);
306 }
307 }
308
309 return cList;
310}
311
313// Get RU ID (from 0 to 191) from a given O2ChipID (from 0 to 24119)
314short int ITSThresholdCalibrator::getRUID(short int chipID)
315{
316 // below there are the inverse of the formulas in getChipListFromRu(...)
317 if (chipID < 432) { // IB
318 return chipID / 9;
319 } else if (chipID >= 432 && chipID < 6480) { // ML
320 return (chipID - 48 * 9 + 112 * 48) / 112;
321 } else { // OL
322 return (chipID - 48 * 9 - 54 * 112 + 102 * 196) / 196;
323 }
324}
325
327// Convert comma-separated list of integers to a vector of int
328std::vector<short int> ITSThresholdCalibrator::getIntegerVect(std::string& s)
329{
330 std::stringstream ss(s);
331 std::vector<short int> result;
332 char ch;
333 short int tmp;
334 while (ss >> tmp) {
335 result.push_back(tmp);
336 ss >> ch;
337 }
338 return result;
339}
340
342// Open a new ROOT file and threshold TTree for that file
343void ITSThresholdCalibrator::initThresholdTree(bool recreate /*=true*/)
344{
345
346 // Create output directory to store output
347 std::string dir = this->mOutputDir + fmt::format("{}_{}/", mDataTakingContext.envId, mDataTakingContext.runNumber);
349 LOG(info) << "Created " << dir << " directory for ROOT trees output";
350
351 std::string filename = dir + mDataTakingContext.runNumber + '_' +
352 std::to_string(this->mFileNumber) + '_' + this->mHostname + "_modSel" + std::to_string(mChipModSel) + ".root.part";
353
354 // Check if file already exists
355 struct stat buffer;
356 if (recreate && stat(filename.c_str(), &buffer) == 0) {
357 LOG(warning) << "File " << filename << " already exists, recreating";
358 }
359
360 // Initialize ROOT output file
361 // to prevent premature external usage, use temporary name
362 const char* option = recreate ? "RECREATE" : "UPDATE";
363 mRootOutfile = new TFile(filename.c_str(), option);
364
365 // Tree containing the s-curves points
366 mScTree = new TTree("s-curve-points", "s-curve-points");
367 mScTree->Branch("chipid", &vChipid, "vChipID[1024]/S");
368 mScTree->Branch("row", &vRow, "vRow[1024]/S");
369
370 // Initialize output TTree branches
371 mThresholdTree = new TTree("ITS_calib_tree", "ITS_calib_tree");
372 mThresholdTree->Branch("chipid", &vChipid, "vChipID[1024]/S");
373 mThresholdTree->Branch("row", &vRow, "vRow[1024]/S");
374 if (mScanType == 'T' || mScanType == 'V' || mScanType == 'I') {
375 std::string bName = mScanType == 'T' ? "thr" : mScanType == 'V' ? "vcasn"
376 : "ithr";
377 mThresholdTree->Branch(bName.c_str(), &vThreshold, "vThreshold[1024]/S");
378 mThresholdTree->Branch("noise", &vNoise, "vNoise[1024]/F");
379 mThresholdTree->Branch("spoints", &vPoints, "vPoints[1024]/b");
380 mThresholdTree->Branch("success", &vSuccess, "vSuccess[1024]/O");
381
382 mScTree->Branch("chg", &vCharge, "vCharge[1024]/b");
383 mScTree->Branch("hits", &vHits, "vHits[1024]/b");
384 } else if (mScanType == 'D' || mScanType == 'A') { // this->mScanType == 'D' and this->mScanType == 'A'
385 mThresholdTree->Branch("n_hits", &vThreshold, "vThreshold[1024]/S");
386 } else if (mScanType == 'P') {
387 mThresholdTree->Branch("n_hits", &vThreshold, "vThreshold[1024]/S");
388 mThresholdTree->Branch("strobedel", &vMixData, "vMixData[1024]/S");
389 } else if (mScanType == 'p' || mScanType == 't') {
390 mThresholdTree->Branch("n_hits", &vThreshold, "vThreshold[1024]/S");
391 mThresholdTree->Branch("strobedel", &vMixData, "vMixData[1024]/S");
392 mThresholdTree->Branch("charge", &vCharge, "vCharge[1024]/b");
393 if (doSlopeCalculation) {
394 mSlopeTree = new TTree("line_tree", "line_tree");
395 mSlopeTree->Branch("chipid", &vChipid, "vChipID[1024]/S");
396 mSlopeTree->Branch("row", &vRow, "vRow[1024]/S");
397 mSlopeTree->Branch("slope", &vSlope, "vSlope[1024]/F");
398 mSlopeTree->Branch("intercept", &vIntercept, "vIntercept[1024]/F");
399 }
400 } else if (mScanType == 'R') {
401 mThresholdTree->Branch("n_hits", &vThreshold, "vThreshold[1024]/S");
402 mThresholdTree->Branch("vresetd", &vMixData, "vMixData[1024]/S");
403 } else if (mScanType == 'r') {
404 mThresholdTree->Branch("thr", &vThreshold, "vThreshold[1024]/S");
405 mThresholdTree->Branch("noise", &vNoise, "vNoise[1024]/F");
406 mThresholdTree->Branch("success", &vSuccess, "vSuccess[1024]/O");
407 mThresholdTree->Branch("vresetd", &vMixData, "vMixData[1024]/S");
408 }
409
410 return;
411}
412
414// Returns upper / lower limits for threshold determination.
415// data is the number of trigger counts per charge injected;
416// x is the array of charge injected values;
417// NPoints is the length of both arrays.
418bool ITSThresholdCalibrator::findUpperLower(
419 std::vector<std::vector<unsigned short int>> data, const short int& NPoints,
420 short int& lower, short int& upper, bool flip, int iloop2)
421{
422 // Initialize (or re-initialize) upper and lower
423 upper = -1;
424 lower = -1;
425
426 if (flip) { // ITHR case. lower is at large mX[i], upper is at small mX[i]
427
428 for (int i = 0; i < NPoints; i++) {
429 int comp = mScanType != 'r' ? data[iloop2][i] : data[i][iloop2];
430 if (comp == 0) {
431 upper = i;
432 break;
433 }
434 }
435
436 if (upper == -1) {
437 return false;
438 }
439 for (int i = upper; i > 0; i--) {
440 int comp = mScanType != 'r' ? data[iloop2][i] : data[i][iloop2];
441 if (comp >= nInjScaled) {
442 lower = i;
443 break;
444 }
445 }
446
447 } else { // not flipped
448
449 for (int i = 0; i < NPoints; i++) {
450 int comp = mScanType != 'r' ? data[iloop2][i] : data[i][iloop2];
451 if (comp >= nInjScaled) {
452 upper = i;
453 break;
454 }
455 }
456
457 if (upper == -1) {
458 return false;
459 }
460 for (int i = upper; i > 0; i--) {
461 int comp = mScanType != 'r' ? data[iloop2][i] : data[i][iloop2];
462 if (comp == 0) {
463 lower = i;
464 break;
465 }
466 }
467 }
468
469 // If search was successful, return central x value
470 if ((lower == -1) || (upper < lower)) {
471 return false;
472 }
473 return true;
474}
475
477// Main findThreshold function which calls one of the three methods
478bool ITSThresholdCalibrator::findThreshold(
479 const short int& chipID, std::vector<std::vector<unsigned short int>> data, const float* x, short int& NPoints,
480 float& thresh, float& noise, int& spoints, int iloop2)
481{
482 bool success = false;
483
484 switch (this->mFitType) {
485 case DERIVATIVE: // Derivative method
486 success = this->findThresholdDerivative(data, x, NPoints, thresh, noise, spoints, iloop2);
487 break;
488
489 case FIT: // Fit method
490 success = this->findThresholdFit(chipID, data, x, NPoints, thresh, noise, spoints, iloop2);
491 break;
492
493 case HITCOUNTING: // Hit-counting method
494 success = this->findThresholdHitcounting(data, x, NPoints, thresh, iloop2);
495 // noise = 0;
496 break;
497 }
498
499 return success;
500}
501
503// Use ROOT to find the threshold and noise via S-curve fit
504// data is the number of trigger counts per charge injected;
505// x is the array of charge injected values;
506// NPoints is the length of both arrays.
507// thresh, noise, chi2 pointers are updated with results from the fit
508// spoints: number of points in the S of the S-curve (with n_hits between 0 and 50, excluding first and last point)
509// iloop2 is 0 for thr scan but is equal to vresetd index in 2D vresetd scan
510bool ITSThresholdCalibrator::findThresholdFit(
511 const short int& chipID, std::vector<std::vector<unsigned short int>> data, const float* x, const short int& NPoints,
512 float& thresh, float& noise, int& spoints, int iloop2)
513{
514 // Find lower & upper values of the S-curve region
515 short int lower, upper;
516 bool flip = (this->mScanType == 'I');
517 auto fndVal = std::find(chipDumpList.begin(), chipDumpList.end(), chipID);
518
519 if (!this->findUpperLower(data, NPoints, lower, upper, flip, iloop2) || lower == upper) {
520 if (this->mVerboseOutput) {
521 LOG(warning) << "Start-finding unsuccessful: (lower, upper) = ("
522 << lower << ", " << upper << ")";
523 }
524
525 if (isDumpS && (dumpCounterS[chipID] < maxDumpS || maxDumpS < 0) && (fndVal != chipDumpList.end() || !chipDumpList.size())) { // save bad s-curves
526 for (int i = 0; i < NPoints; i++) {
527 this->mFitHist->SetBinContent(i + 1, mScanType != 'r' ? data[iloop2][i] : data[i][iloop2]);
528 }
529 fileDumpS->cd();
530 mFitHist->Write();
531 }
532 if (isDumpS) {
533 dumpCounterS[chipID]++;
534 }
535
536 return false;
537 }
538 float start = (this->mX[upper] + this->mX[lower]) / 2;
539
540 if (start < 0) {
541 if (this->mVerboseOutput) {
542 LOG(warning) << "Start-finding unsuccessful: Start = " << start;
543 }
544 return false;
545 }
546
547 for (int i = 0; i < NPoints; i++) {
548 this->mFitHist->SetBinContent(i + 1, mScanType != 'r' ? data[iloop2][i] : data[i][iloop2]);
549 }
550
551 // Initialize starting parameters
552 this->mFitFunction->SetParameter(0, start);
553 this->mFitFunction->SetParameter(1, 8);
554
555 this->mFitHist->Fit("mFitFunction", "RQL");
556 if (isDumpS && (dumpCounterS[chipID] < maxDumpS || maxDumpS < 0) && (fndVal != chipDumpList.end() || !chipDumpList.size())) { // save good s-curves
557 fileDumpS->cd();
558 mFitHist->Write();
559 }
560 if (isDumpS) {
561 dumpCounterS[chipID]++;
562 }
563
564 noise = this->mFitFunction->GetParameter(1);
565 thresh = this->mFitFunction->GetParameter(0);
566 float chi2 = this->mFitFunction->GetChisquare() / this->mFitFunction->GetNDF();
567 spoints = upper - lower - 1;
568
569 // Clean up histogram for next time it is used
570 this->mFitHist->Reset();
571
572 return (chi2 < 5);
573}
574
576// Use ROOT to find the threshold and noise via derivative method
577// data is the number of trigger counts per charge injected;
578// x is the array of charge injected values;
579// NPoints is the length of both arrays.
580// spoints: number of points in the S of the S-curve (with n_hits between 0 and 50, excluding first and last point)
581// iloop2 is 0 for thr scan but is equal to vresetd index in 2D vresetd scan
582bool ITSThresholdCalibrator::findThresholdDerivative(std::vector<std::vector<unsigned short int>> data, const float* x, const short int& NPoints,
583 float& thresh, float& noise, int& spoints, int iloop2)
584{
585 // Find lower & upper values of the S-curve region
586 short int lower, upper;
587 bool flip = (this->mScanType == 'I');
588 if (!this->findUpperLower(data, NPoints, lower, upper, flip, iloop2) || lower == upper) {
589 if (this->mVerboseOutput) {
590 LOG(warning) << "Start-finding unsuccessful: (lower, upper) = (" << lower << ", " << upper << ")";
591 }
592 return false;
593 }
594
595 int deriv_size = upper - lower;
596 float deriv[deriv_size];
597 float xfx = 0, fx = 0;
598
599 // Fill array with derivatives
600 for (int i = lower; i < upper; i++) {
601 deriv[i - lower] = std::abs(mScanType != 'r' ? (data[iloop2][i + 1] - data[iloop2][i]) : (data[i + 1][iloop2] - data[i][iloop2])) / (this->mX[i + 1] - mX[i]);
602 xfx += this->mX[i] * deriv[i - lower];
603 fx += deriv[i - lower];
604 }
605
606 if (fx > 0.) {
607 thresh = xfx / fx;
608 }
609 float stddev = 0;
610 for (int i = lower; i < upper; i++) {
611 stddev += std::pow(this->mX[i] - thresh, 2) * deriv[i - lower];
612 }
613
614 stddev /= fx;
615 noise = std::sqrt(stddev);
616 spoints = upper - lower - 1;
617
618 return fx > 0.;
619}
620
622// Use ROOT to find the threshold and noise via derivative method
623// data is the number of trigger counts per charge injected;
624// x is the array of charge injected values;
625// NPoints is the length of both arrays.
626// iloop2 is 0 for thr scan but is equal to vresetd index in 2D vresetd scan
627bool ITSThresholdCalibrator::findThresholdHitcounting(
628 std::vector<std::vector<unsigned short int>> data, const float* x, const short int& NPoints, float& thresh, int iloop2)
629{
630 unsigned short int numberOfHits = 0;
631 bool is50 = false;
632 for (unsigned short int i = 0; i < NPoints; i++) {
633 numberOfHits += (mScanType != 'r') ? data[iloop2][i] : data[i][iloop2];
634 int comp = (mScanType != 'r') ? data[iloop2][i] : data[i][iloop2];
635 if (!is50 && comp == nInjScaled) {
636 is50 = true;
637 }
638 }
639
640 // If not enough counts return a failure
641 if (!is50) {
642 if (this->mVerboseOutput) {
643 LOG(warning) << "Calculation unsuccessful: too few hits. Skipping this pixel";
644 }
645 return false;
646 }
647
648 if (this->mScanType == 'T') {
649 thresh = this->mX[N_RANGE - 1] - numberOfHits / float(nInjScaled);
650 } else if (this->mScanType == 'V') {
651 thresh = (this->mX[N_RANGE - 1] * nInjScaled - numberOfHits) / float(nInjScaled);
652 } else if (this->mScanType == 'I') {
653 thresh = (numberOfHits + nInjScaled * this->mX[0]) / float(nInjScaled);
654 } else {
655 LOG(error) << "Unexpected runtype encountered in findThresholdHitcounting()";
656 return false;
657 }
658
659 return true;
660}
661
663// Run threshold extraction on completed row and update memory
664void ITSThresholdCalibrator::extractThresholdRow(const short int& chipID, const short int& row)
665{
666 if (this->mScanType == 'D' || this->mScanType == 'A') {
667 // Loop over all columns (pixels) in the row
668 for (short int col_i = 0; col_i < this->N_COL; col_i++) {
669 vChipid[col_i] = chipID;
670 vRow[col_i] = row;
671 vThreshold[col_i] = this->mPixelHits[chipID][row][col_i][0][0];
672 if (vThreshold[col_i] > nInj) {
673 this->mNoisyPixID[chipID].push_back(col_i * 1000 + row);
674 } else if (vThreshold[col_i] > 0 && vThreshold[col_i] < nInj) {
675 this->mIneffPixID[chipID].push_back(col_i * 1000 + row);
676 } else if (vThreshold[col_i] == 0) {
677 this->mDeadPixID[chipID].push_back(col_i * 1000 + row);
678 }
679 }
680 } else if (this->mScanType == 'P' || this->mScanType == 'p' || mScanType == 'R' || mScanType == 't') {
681 // Loop over all columns (pixels) in the row
682 for (short int var1_i = 0; var1_i < this->N_RANGE; var1_i++) {
683 for (short int chg_i = 0; chg_i < this->N_RANGE2; chg_i++) {
684 for (short int col_i = 0; col_i < this->N_COL; col_i++) {
685 vChipid[col_i] = chipID;
686 vRow[col_i] = row;
687 vThreshold[col_i] = this->mPixelHits[chipID][row][col_i][chg_i][var1_i];
688 vMixData[col_i] = (var1_i * this->mStep) + mMin;
689 if (mScanType != 'R') {
690 vMixData[col_i]++; // +1 because a delay of n correspond to a real delay of n+1 (from ALPIDE manual)
691 }
692 vCharge[col_i] = (unsigned char)(chg_i * this->mStep2 + mMin2);
693 }
694 this->saveThreshold();
695 }
696 }
697
698 if (doSlopeCalculation) {
699 int delA = -1, delB = -1;
700 for (short int col_i = 0; col_i < this->N_COL; col_i++) {
701 for (short int chg_i = 0; chg_i < 2; chg_i++) {
702 bool isFound = false;
703 int checkchg = !chg_i ? chargeA / mStep2 : chargeB / mStep2;
704 for (short int sdel_i = N_RANGE - 1; sdel_i >= 0; sdel_i--) {
705 if (mPixelHits[chipID][row][col_i][checkchg - 1][sdel_i] == nInj) {
706 if (!chg_i) {
707 delA = mMin + sdel_i * mStep + mStep / 2;
708 isFound = true;
709 } else {
710 delB = mMin + sdel_i * mStep + mStep / 2;
711 isFound = true;
712 }
713 break;
714 }
715 } // end loop on strobe delays
716
717 if (!isFound) { // if not found, take the first point with hits starting from the left (i.e. in principle the closest point to the one with MAX n_hits)
718 for (short int sdel_i = 0; sdel_i < N_RANGE; sdel_i++) {
719 if (mPixelHits[chipID][row][col_i][checkchg - 1][sdel_i] > 0) {
720 if (!chg_i) {
721 delA = mMin + sdel_i * mStep + mStep / 2;
722 } else {
723 delB = mMin + sdel_i * mStep + mStep / 2;
724 }
725 break;
726 }
727 } // end loop on strobe delays
728 } // end if on isFound
729 } // end loop on the two charges
730
731 if (delA > 0 && delB > 0 && delA != delB) {
732 vSlope[col_i] = ((float)(chargeA - chargeB) / (float)(delA - delB));
733 vIntercept[col_i] = (float)chargeA - (float)(vSlope[col_i] * delA);
734 if (vSlope[col_i] < 0) { // protection for non expected slope
735 vSlope[col_i] = 0.;
736 vIntercept[col_i] = 0.;
737 }
738 } else {
739 vSlope[col_i] = 0.;
740 vIntercept[col_i] = 0.;
741 }
742 } // end loop on pix
743
744 mSlopeTree->Fill();
745 }
746
747 } else { // threshold, vcasn, ithr, vresetd_2d
748
749 short int iRU = getRUID(chipID);
750#ifdef WITH_OPENMP
751 omp_set_num_threads(mNThreads);
752#pragma omp parallel for schedule(dynamic)
753#endif
754 // Loop over all columns (pixels) in the row
755 for (short int col_i = 0; col_i < this->N_COL; col_i++) {
756 // Do the threshold fit
757 float thresh = 0., noise = 0.;
758 bool success = false;
759 int spoints = 0;
760 int scan_i = mScanType == 'r' ? (mLoopVal[iRU][row] - mMin) / mStep : 0;
761 if (isDumpS) { // already protected for multi-thread in the init
762 mFitHist->SetName(Form("scurve_chip%d_row%d_col%d_scani%d", chipID, row, col_i, scan_i));
763 }
764
765 success = this->findThreshold(chipID, mPixelHits[chipID][row][col_i],
766 this->mX, mScanType == 'r' ? N_RANGE2 : N_RANGE, thresh, noise, spoints, scan_i);
767
768 vChipid[col_i] = chipID;
769 vRow[col_i] = row;
770 vThreshold[col_i] = (mScanType == 'T' || mScanType == 'r') ? (short int)(thresh * 10.) : (short int)(thresh);
771 vNoise[col_i] = (float)(noise * 10.); // always factor 10 also for ITHR/VCASN to not have all zeros
772 vSuccess[col_i] = success;
773 vPoints[col_i] = spoints > 0 ? (unsigned char)(spoints) : 0;
774
775 if (mScanType == 'r') {
776 vMixData[col_i] = mLoopVal[iRU][row];
777 }
778 }
779 if (mScanType == 'r') {
780 this->saveThreshold(); // save before moving to the next vresetd
781 }
782
783 // Fill the ScTree tree
784 if (mScanType == 'T' || mScanType == 'V' || mScanType == 'I') { // TODO: store also for other scans?
785 for (int ichg = mMin; ichg <= mMax; ichg += mStep) {
786 for (short int col_i = 0; col_i < this->N_COL; col_i += mColStep) {
787 vCharge[col_i] = ichg;
788 vHits[col_i] = mPixelHits[chipID][row][col_i][0][(ichg - mMin) / mStep];
789 }
790 mScTree->Fill();
791 }
792 }
793 } // end of the else
794
795 // Saves threshold information to internal memory
796 if (mScanType != 'P' && mScanType != 'p' && mScanType != 't' && mScanType != 'R' && mScanType != 'r') {
797 this->saveThreshold();
798 }
799}
800
802void ITSThresholdCalibrator::saveThreshold()
803{
804 // write to TTree
805 this->mThresholdTree->Fill();
806
807 if (this->mScanType == 'V' || this->mScanType == 'I' || this->mScanType == 'T') {
808 // Save info in a map for later averaging
809 int sumT = 0, sumSqT = 0, sumN = 0, sumSqN = 0;
810 int countSuccess = 0, countUnsuccess = 0;
811 for (int i = 0; i < this->N_COL; i++) {
812 if (vSuccess[i]) {
813 sumT += vThreshold[i];
814 sumN += (int)vNoise[i];
815 sumSqT += (vThreshold[i]) * (vThreshold[i]);
816 sumSqN += ((int)vNoise[i]) * ((int)vNoise[i]);
817 countSuccess++;
818 if (vThreshold[i] >= mMin && vThreshold[i] <= mMax && (mScanType == 'I' || mScanType == 'V')) {
819 mpvCounter[vChipid[0]][vThreshold[i] - mMin]++;
820 }
821 } else {
822 countUnsuccess++;
823 }
824 }
825 short int chipID = vChipid[0];
826 std::array<long int, 6> dataSum{{sumT, sumSqT, sumN, sumSqN, countSuccess, countUnsuccess}};
827 if (!(this->mThresholds.count(chipID))) {
828 this->mThresholds[chipID] = dataSum;
829 } else {
830 std::array<long int, 6> dataAll{{this->mThresholds[chipID][0] + dataSum[0], this->mThresholds[chipID][1] + dataSum[1], this->mThresholds[chipID][2] + dataSum[2], this->mThresholds[chipID][3] + dataSum[3], this->mThresholds[chipID][4] + dataSum[4], this->mThresholds[chipID][5] + dataSum[5]}};
831 this->mThresholds[chipID] = dataAll;
832 }
833 }
834
835 return;
836}
837
839// Perform final operations on output objects. In the case of a full threshold
840// scan, rename ROOT file and create metadata file for writing to EOS
841void ITSThresholdCalibrator::finalizeOutput()
842{
843 // Check that objects actually exist in memory
844 if (!(mScTree) || !(this->mRootOutfile) || !(this->mThresholdTree) || (doSlopeCalculation && !(this->mSlopeTree))) {
845 return;
846 }
847
848 // Ensure that everything has been written to the ROOT file
849 this->mRootOutfile->cd();
850 this->mThresholdTree->Write(nullptr, TObject::kOverwrite);
851 this->mScTree->Write(nullptr, TObject::kOverwrite);
852
853 if (doSlopeCalculation) {
854 this->mSlopeTree->Write(nullptr, TObject::kOverwrite);
855 }
856
857 // Clean up the mThresholdTree, mScTree and ROOT output file
858 delete this->mThresholdTree;
859 this->mThresholdTree = nullptr;
860 delete mScTree;
861 mScTree = nullptr;
862 if (doSlopeCalculation) {
863 delete this->mSlopeTree;
864 this->mSlopeTree = nullptr;
865 }
866
867 this->mRootOutfile->Close();
868 delete this->mRootOutfile;
869 this->mRootOutfile = nullptr;
870
871 // Check that expected output directory exists
872 std::string dir = this->mOutputDir + fmt::format("{}_{}/", mDataTakingContext.envId, mDataTakingContext.runNumber);
873 if (!std::filesystem::exists(dir)) {
874 LOG(error) << "Cannot find expected output directory " << dir;
875 return;
876 }
877
878 // Expected ROOT output filename
879 std::string filename = mDataTakingContext.runNumber + '_' +
880 std::to_string(this->mFileNumber) + '_' + this->mHostname + "_modSel" + std::to_string(mChipModSel);
881 std::string filenameFull = dir + filename;
882 try {
883 std::rename((filenameFull + ".root.part").c_str(),
884 (filenameFull + ".root").c_str());
885 } catch (std::exception const& e) {
886 LOG(error) << "Failed to rename ROOT file " << filenameFull
887 << ".root.part, reason: " << e.what();
888 }
889
890 // Create metadata file
892 mdFile->fillFileData(filenameFull + ".root");
893 mdFile->setDataTakingContext(mDataTakingContext);
894 if (!(this->mMetaType.empty())) {
895 mdFile->type = this->mMetaType;
896 }
897 mdFile->priority = "high";
898 mdFile->lurl = filenameFull + ".root";
899 auto metaFileNameTmp = fmt::format("{}{}.tmp", this->mMetafileDir, filename);
900 auto metaFileName = fmt::format("{}{}.done", this->mMetafileDir, filename);
901 try {
902 std::ofstream metaFileOut(metaFileNameTmp);
903 metaFileOut << mdFile->asString() << '\n';
904 metaFileOut.close();
905 std::filesystem::rename(metaFileNameTmp, metaFileName);
906 } catch (std::exception const& e) {
907 LOG(error) << "Failed to create threshold metadata file "
908 << metaFileName << ", reason: " << e.what();
909 }
910 delete mdFile;
911
912 // Next time a file is created, use a larger number
913 this->mFileNumber++;
914
915 return;
916
917} // finalizeOutput
918
920// Set the run_type for this run
921// Initialize the memory needed for this specific type of run
922void ITSThresholdCalibrator::setRunType(const short int& runtype)
923{
924
925 // Save run type info for future evaluation
926 this->mRunType = runtype;
927
928 if (runtype == THR_SCAN) {
929 // full_threshold-scan -- just extract thresholds for each pixel and write to TTree
930 // 512 rows per chip
931 this->mScanType = 'T';
932 this->initThresholdTree();
933 this->mMin = 0;
934 this->mMax = 50;
935 this->N_RANGE = 51;
936 this->mCheckExactRow = true;
937
938 } else if (runtype == THR_SCAN_SHORT || runtype == THR_SCAN_SHORT_100HZ ||
939 runtype == THR_SCAN_SHORT_200HZ || runtype == THR_SCAN_SHORT_33 || runtype == THR_SCAN_SHORT_2_10HZ || runtype == THR_SCAN_SHORT_150INJ) {
940 // threshold_scan_short -- just extract thresholds for each pixel and write to TTree
941 // 10 rows per chip
942 this->mScanType = 'T';
943 this->initThresholdTree();
944 this->mMin = 0;
945 this->mMax = 50;
946 this->N_RANGE = 51;
947 this->mCheckExactRow = true;
948 if (runtype == THR_SCAN_SHORT_150INJ) {
949 nInj = 150;
950 if (mMeb >= 0) {
951 nInjScaled = nInj / 3;
952 }
953 }
954 } else if (runtype == VCASN150 || runtype == VCASN100 || runtype == VCASN100_100HZ || runtype == VCASN130 || runtype == VCASNBB) {
955 // VCASN tuning for different target thresholds
956 // Store average VCASN for each chip into CCDB
957 // ATTENTION: with back bias (VCASNBB) put max vcasn to 130 (default is 80)
958 // 4 rows per chip
959 this->mScanType = 'V';
960 this->initThresholdTree();
961 this->mMin = inMinVcasn; // 30 is the default
962 this->mMax = inMaxVcasn; // 80 is the default
963 this->N_RANGE = mMax - mMin + 1;
964 this->mCheckExactRow = true;
965
966 } else if (runtype == ITHR150 || runtype == ITHR100 || runtype == ITHR100_100HZ || runtype == ITHR130) {
967 // ITHR tuning -- average ITHR per chip
968 // S-curve is backwards from VCASN case, otherwise same
969 // 4 rows per chip
970 this->mScanType = 'I';
971 this->initThresholdTree();
972 this->mMin = inMinIthr; // 25 is the default
973 this->mMax = inMaxIthr; // 100 is the default
974 this->N_RANGE = mMax - mMin + 1;
975 this->mCheckExactRow = true;
976
977 } else if (runtype == DIGITAL_SCAN || runtype == DIGITAL_SCAN_100HZ || runtype == DIGITAL_SCAN_NOMASK) {
978 // Digital scan -- only storing one value per chip, no fit needed
979 this->mScanType = 'D';
980 this->initThresholdTree();
981 this->mFitType = NO_FIT;
982 this->mMin = 0;
983 this->mMax = 0;
984 this->N_RANGE = mMax - mMin + 1;
985 this->mCheckExactRow = false;
986
987 } else if (runtype == ANALOGUE_SCAN) {
988 // Analogue scan -- only storing one value per chip, no fit needed
989 this->mScanType = 'A';
990 this->initThresholdTree();
991 this->mFitType = NO_FIT;
992 this->mMin = 0;
993 this->mMax = 0;
994 this->N_RANGE = mMax - mMin + 1;
995 this->mCheckExactRow = false;
996
997 } else if (runtype == PULSELENGTH_SCAN) {
998 // Pulse length scan
999 this->mScanType = 'P';
1000 this->initThresholdTree();
1001 this->mFitType = NO_FIT;
1002 this->mMin = 0;
1003 this->mMax = 400; // strobe delay goes from 0 to 400 (included) in steps of 4
1004 this->mStep = 1;
1005 this->mStrobeWindow = 1; // it's 0 but it corresponds to 0+1 (as from alpide manual)
1006 this->N_RANGE = (mMax - mMin) / mStep + 1;
1007 this->mCheckExactRow = true;
1008 } else if (runtype == TOT_CALIBRATION_1_ROW) {
1009 // Pulse length scan 2D (charge vs strobe delay)
1010 this->mScanType = 'p'; // small p, just to distinguish from capital P
1011 this->initThresholdTree();
1012 this->mFitType = NO_FIT;
1013 this->mMin = 0;
1014 this->mMax = 2000; // strobe delay goes from 0 to 2000 in steps of 10
1015 this->mStep = 10;
1016 this->mStrobeWindow = 10; // it's 9 but it corresponds to 9+1 (as from alpide manual)
1017 this->N_RANGE = (mMax - mMin) / mStep + 1;
1018 this->mMin2 = 0; // charge min
1019 this->mMax2 = 170; // charge max
1020 this->mStep2 = 1; // step for the charge
1021 this->N_RANGE2 = (mMax2 - mMin2) / mStep2 + 1;
1022 this->mCheckExactRow = true;
1023 } else if (runtype == TOT_CALIBRATION) {
1024 // TOT calibration (like pulse shape 2D but with a reduced range in both strobe delay and charge)
1025 this->mScanType = 't';
1026 this->initThresholdTree();
1027 this->mFitType = NO_FIT;
1028 this->mMin = 300;
1029 this->mMax = 1100; // strobe delay goes from 300 to 1100 (included) in steps of 10
1030 this->mStep = 10;
1031 this->mStrobeWindow = 10; // it's 9 but it corresponds to 9+1 (as from alpide manual)
1032 this->N_RANGE = (mMax - mMin) / mStep + 1;
1033 this->mMin2 = 30; // charge min
1034 this->mMax2 = 60; // charge max
1035 this->mStep2 = 30; // step for the charge
1036 this->mCalculate2DParams = false; // do not calculate time over threshold, pulse length, etc..
1037 this->N_RANGE2 = (mMax2 - mMin2) / mStep2 + 1;
1038 this->mCheckExactRow = true;
1039
1040 } else if (runtype == VRESETD_150 || runtype == VRESETD_300 || runtype == VRESETD_2D) {
1041 this->mScanType = 'R'; // capital R is for 1D scan
1042 if (runtype == VRESETD_150 || runtype == VRESETD_300) {
1043 this->mFitType = NO_FIT;
1044 }
1045 this->mMin = 100;
1046 this->mMax = 240; // vresetd goes from 100 to 240 in steps of 5
1047 this->mStep = 5;
1048 this->N_RANGE = (mMax - mMin) / mStep + 1;
1049 if (runtype == VRESETD_2D) {
1050 this->mScanType = 'r'; // small r, just to distinguish from capital R
1051 this->mMin2 = 0; // charge min
1052 this->mMax2 = 50; // charge max
1053 this->mStep2 = 1; // step for the charge
1054 this->N_RANGE2 = (mMax2 - mMin2) / mStep2 + 1;
1055 }
1056 this->mCheckExactRow = true;
1057 this->initThresholdTree();
1058 } else {
1059 // No other run type recognized by this workflow
1060 LOG(error) << "Runtype " << runtype << " not recognized by calibration workflow (ignore if you are in manual mode)";
1061 if (isManualMode) {
1062 LOG(info) << "Entering manual mode: be sure to have set all parameters correctly";
1063 this->mScanType = manualScanType[0];
1064 this->mMin = manualMin;
1065 this->mMax = manualMax;
1066 this->mMin2 = manualMin2;
1067 this->mMax2 = manualMax2;
1068 this->mStep = manualStep; // 1 by default
1069 this->mStep2 = manualStep2; // 1 by default
1070 this->mStrobeWindow = manualStrobeWindow; // 5 = 125 ns by default
1071 this->N_RANGE = (mMax - mMin) / mStep + 1;
1072 this->N_RANGE2 = (mMax2 - mMin2) / mStep2 + 1;
1073 if (saveTree) {
1074 this->initThresholdTree();
1075 }
1076 this->mFitType = (mScanType == 'D' || mScanType == 'A' || mScanType == 'P' || mScanType == 'p' || mScanType == 't') ? NO_FIT : mFitType;
1077 this->mCheckExactRow = (mScanType == 'D' || mScanType == 'A') ? false : true;
1078 if (scaleNinj) {
1079 nInjScaled = nInj / 3;
1080 }
1081 } else {
1082 throw runtype;
1083 }
1084 }
1085
1086 this->mX = new float[mScanType == 'r' ? N_RANGE2 : N_RANGE];
1087 for (short int i = ((mScanType == 'r') ? mMin2 : mMin); i <= ((mScanType == 'r') ? mMax2 / mStep2 : mMax / mStep); i++) {
1088 this->mX[i - (mScanType == 'r' ? mMin2 : mMin)] = (float)i + 0.5;
1089 }
1090
1091 // Initialize objects for doing the threshold fits
1092 if (this->mFitType == FIT) {
1093 // Initialize the histogram used for error function fits
1094 // Will initialize the TF1 in setRunType (it is different for different runs)
1095 this->mFitHist = new TH1F(
1096 "mFitHist", "mFitHist", mScanType == 'r' ? N_RANGE2 : N_RANGE, mX[0] - 1., mX[(mScanType == 'r' ? N_RANGE2 : N_RANGE) - 1]);
1097
1098 // Initialize correct fit function for the scan type
1099 this->mFitFunction = (this->mScanType == 'I')
1100 ? new TF1("mFitFunction", erf_ithr, mMin, mMax, 2)
1101 : new TF1("mFitFunction", erf, (mScanType == 'T' || mScanType == 'r') ? 3 : mMin, mScanType == 'r' ? mMax2 : mMax, 2);
1102 this->mFitFunction->SetParName(0, "Threshold");
1103 this->mFitFunction->SetParName(1, "Noise");
1104 }
1105
1106 return;
1107}
1108
1110// Calculate pulse parameters in 1D scan: time over threshold, rise time, ...
1111std::vector<float> ITSThresholdCalibrator::calculatePulseParams(const short int& chipID)
1112{
1113
1114 int rt_mindel = -1, rt_maxdel = -1, tot_mindel = -1, tot_maxdel = -1;
1115 int sumRt = 0, sumSqRt = 0, countRt = 0, sumTot = 0, sumSqTot = 0, countTot = 0;
1116
1117 for (auto itrow = mPixelHits[chipID].begin(); itrow != mPixelHits[chipID].end(); itrow++) { // loop over the chip rows
1118 short int row = itrow->first;
1119 for (short int col_i = 0; col_i < this->N_COL; col_i++) { // loop over the pixels on the row
1120 for (short int sdel_i = 0; sdel_i < this->N_RANGE; sdel_i++) { // loop over the strobe delays
1121 if (mPixelHits[chipID][row][col_i][0][sdel_i] > 0.1 * nInj && mPixelHits[chipID][row][col_i][0][sdel_i] < nInj && rt_mindel < 0) { // from left, first bin with 10% hits and 90% hits
1122 rt_mindel = (sdel_i * mStep) + 1; // + 1 because if delay = n, we get n+1 in reality (ALPIDE feature)
1123 }
1124 if (mPixelHits[chipID][row][col_i][0][sdel_i] >= 0.9 * nInj) { // for Rt max take the 90% point
1125 rt_maxdel = (sdel_i * mStep) + 1;
1126 break;
1127 }
1128 }
1129 for (short int sdel_i = 0; sdel_i < N_RANGE; sdel_i++) {
1130 if (mPixelHits[chipID][row][col_i][0][sdel_i] >= 0.5 * nInj) { // for ToT take the 50% point
1131 tot_mindel = (sdel_i * mStep) + 1;
1132 break;
1133 }
1134 }
1135
1136 for (short int sdel_i = N_RANGE - 1; sdel_i >= 0; sdel_i--) { // from right, the first bin with 50% nInj hits
1137 if (mPixelHits[chipID][row][col_i][0][sdel_i] >= 0.5 * nInj) {
1138 tot_maxdel = (sdel_i * mStep) + 1;
1139 break;
1140 }
1141 }
1142
1143 if (tot_maxdel > tot_mindel && tot_mindel >= 0 && tot_maxdel >= 0) {
1144 sumTot += tot_maxdel - tot_mindel - mStrobeWindow;
1145 sumSqTot += (tot_maxdel - tot_mindel - mStrobeWindow) * (tot_maxdel - tot_mindel - mStrobeWindow);
1146 countTot++;
1147 }
1148
1149 if (rt_maxdel > rt_mindel && rt_maxdel > 0 && rt_mindel > 0) {
1150 sumRt += rt_maxdel - rt_mindel + mStrobeWindow;
1151 sumSqRt += (rt_maxdel - rt_mindel + mStrobeWindow) * (rt_maxdel - rt_mindel + mStrobeWindow);
1152 countRt++;
1153 }
1154
1155 rt_mindel = -1;
1156 rt_maxdel = -1;
1157 tot_maxdel = -1;
1158 tot_mindel = -1;
1159 } // end loop over col_i
1160 } // end loop over chip rows
1161
1162 std::vector<float> output; // {avgRt, rmsRt, avgTot, rmsTot}
1163 // Avg Rt
1164 output.push_back(!countRt ? 0. : (float)sumRt / (float)countRt);
1165 // Rms Rt
1166 output.push_back(!countRt ? 0. : (std::sqrt((float)sumSqRt / (float)countRt - output[0] * output[0])) * 25.);
1167 output[0] *= 25.;
1168 // Avg ToT
1169 output.push_back(!countTot ? 0. : (float)sumTot / (float)countTot);
1170 // Rms ToT
1171 output.push_back(!countTot ? 0. : (std::sqrt((float)sumSqTot / (float)countTot - output[2] * output[2])) * 25.);
1172 output[2] *= 25.;
1173
1174 return output;
1175}
1176
1178// Calculate pulse parameters in 2D scan
1179std::vector<float> ITSThresholdCalibrator::calculatePulseParams2D(const short int& chipID)
1180{
1181 long int sumTot = 0, sumSqTot = 0, countTot = 0;
1182 long int sumMinThr = 0, sumSqMinThr = 0, countMinThr = 0;
1183 long int sumMinThrDel = 0, sumSqMinThrDel = 0;
1184 long int sumMaxPl = 0, sumSqMaxPl = 0, countMaxPl = 0;
1185 long int sumMaxPlChg = 0, sumSqMaxPlChg = 0;
1186
1187 for (auto itrow = mPixelHits[chipID].begin(); itrow != mPixelHits[chipID].end(); itrow++) { // loop over the chip rows
1188 short int row = itrow->first;
1189 for (short int col_i = 0; col_i < this->N_COL; col_i++) { // loop over the pixels on the row
1190 int minThr = 1e7, minThrDel = 1e7, maxPl = -1, maxPlChg = -1;
1191 int tot_mindel = 1e7;
1192 bool isFound = false;
1193 for (short int chg_i = 0; chg_i < this->N_RANGE2; chg_i++) { // loop over charges
1194 for (short int sdel_i = 0; sdel_i < this->N_RANGE; sdel_i++) { // loop over the strobe delays
1195 if (mPixelHits[chipID][row][col_i][chg_i][sdel_i] == nInj) { // minimum threshold charge and delay
1196 minThr = chg_i * mStep2;
1197 minThrDel = (sdel_i * mStep) + 1; // +1 because n->n+1 (as from alpide manual)
1198 isFound = true;
1199 break;
1200 }
1201 }
1202 if (isFound) {
1203 break;
1204 }
1205 }
1206 isFound = false;
1207 for (short int sdel_i = this->N_RANGE - 1; sdel_i >= 0; sdel_i--) { // loop over the strobe delays
1208 for (short int chg_i = this->N_RANGE2 - 1; chg_i >= 0; chg_i--) { // loop over charges
1209 if (mPixelHits[chipID][row][col_i][chg_i][sdel_i] == nInj) { // max pulse length charge and delay
1210 maxPl = (sdel_i * mStep) + 1;
1211 maxPlChg = chg_i * mStep2;
1212 isFound = true;
1213 break;
1214 }
1215 }
1216 if (isFound) {
1217 break;
1218 }
1219 }
1220 isFound = false;
1221 for (short int sdel_i = 0; sdel_i < this->N_RANGE; sdel_i++) { // loop over the strobe delays
1222 for (short int chg_i = 0; chg_i < this->N_RANGE2; chg_i++) { // loop over charges
1223 if (mPixelHits[chipID][row][col_i][chg_i][sdel_i] == nInj) { // min delay for the ToT calculation
1224 tot_mindel = (sdel_i * mStep) + 1;
1225 isFound = true;
1226 break;
1227 }
1228 }
1229 if (isFound) {
1230 break;
1231 }
1232 }
1233
1234 if (maxPl > tot_mindel && tot_mindel < 1e7 && maxPl >= 0) { // ToT
1235 sumTot += maxPl - tot_mindel - mStrobeWindow;
1236 sumSqTot += (maxPl - tot_mindel - mStrobeWindow) * (maxPl - tot_mindel - mStrobeWindow);
1237 countTot++;
1238 }
1239
1240 if (minThr < 1e7) { // minimum threshold
1241 sumMinThr += minThr;
1242 sumSqMinThr += minThr * minThr;
1243 sumMinThrDel += minThrDel;
1244 sumSqMinThrDel += minThrDel * minThrDel;
1245 countMinThr++;
1246 }
1247
1248 if (maxPl >= 0) { // pulse length
1249 sumMaxPl += maxPl;
1250 sumSqMaxPl += maxPl * maxPl;
1251 sumMaxPlChg += maxPlChg;
1252 sumSqMaxPlChg += maxPlChg * maxPlChg;
1253 countMaxPl++;
1254 }
1255 } // end loop over col_i
1256 } // end loop over chip rows
1257
1258 // Pulse shape 2D output: avgToT, rmsToT, MTC, rmsMTC, avgMTCD, rmsMTCD, avgMPL, rmsMPL, avgMPLC, rmsMPLC
1259 std::vector<long int> values = {sumTot, sumSqTot, countTot, sumMinThr, sumSqMinThr, countMinThr, sumMinThrDel, sumSqMinThrDel, countMinThr, sumMaxPl, sumSqMaxPl, countMaxPl, sumMaxPlChg, sumSqMaxPlChg, countMaxPl};
1260 std::vector<float> output;
1261
1262 for (int i = 0; i < values.size(); i += 3) {
1263 // Avg
1264 output.push_back(!values[i + 2] ? 0. : (float)values[i] / (float)values[i + 2]);
1265 // Rms
1266 output.push_back(!values[i + 2] ? 0. : std::sqrt((float)values[i + 1] / (float)values[i + 2] - output[output.size() - 1] * output[output.size() - 1]));
1267 if (i == 0 || i == 6 || i == 9) {
1268 output[output.size() - 1] *= 25.;
1269 output[output.size() - 2] *= 25.;
1270 }
1271 }
1272
1273 return output;
1274}
1276// Extract thresholds and update memory
1277void ITSThresholdCalibrator::extractAndUpdate(const short int& chipID, const short int& row)
1278{
1279 // In threshold scan case, reset mThresholdTree before writing to a new file
1280 if ((this->mRowCounter)++ == N_ROWS_PER_FILE) {
1281 // Finalize output and create a new TTree and ROOT file
1282 this->finalizeOutput();
1283 this->initThresholdTree();
1284 // Reset data counter for the next output file
1285 this->mRowCounter = 1;
1286 }
1287
1288 // Extract threshold values and save to memory
1289 this->extractThresholdRow(chipID, row);
1290
1291 return;
1292}
1293
1295// Main running function
1296// Get info from previous stf decoder workflow, then loop over readout frames
1297// (ROFs) to count hits and extract thresholds
1299{
1300 if (mRunStopRequested) { // give up when run stop request arrived
1301 return;
1302 }
1303
1304 updateTimeDependentParams(pc);
1305
1306 // Calibration vector
1307 const auto calibs = pc.inputs().get<gsl::span<o2::itsmft::GBTCalibData>>("calib");
1308 const auto digits = pc.inputs().get<gsl::span<o2::itsmft::Digit>>("digits");
1309 const auto ROFs = pc.inputs().get<gsl::span<o2::itsmft::ROFRecord>>("digitsROF");
1310
1311 // Store some lengths for convenient looping
1312 const unsigned int nROF = (unsigned int)ROFs.size();
1313
1314 // Loop over readout frames (usually only 1, sometimes 2)
1315 for (unsigned int iROF = 0; iROF < nROF; iROF++) {
1316
1317 unsigned int rofIndex = ROFs[iROF].getFirstEntry();
1318 unsigned int rofNEntries = ROFs[iROF].getNEntries();
1319
1320 // Find the correct charge, row, cw counter values for this ROF
1321 short int loopval = -1, realcharge = 0;
1322 short int row = -1;
1323 short int cwcnt = -1;
1324 bool isAllZero = true;
1325 short int ruIndex = -1;
1326 for (short int iRU = 0; iRU < this->N_RU; iRU++) {
1327 const auto& calib = calibs[iROF * this->N_RU + iRU];
1328 if (calib.calibUserField != 0) {
1329 mRuSet.insert(iRU);
1330 ruIndex = iRU;
1331 isAllZero = false;
1332
1333 if (loopval >= 0) {
1334 LOG(warning) << "More than one charge detected!";
1335 }
1336
1337 if (this->mRunType == -1) {
1338 mCdwVersion = isCRUITS ? 0 : ((short int)(calib.calibUserField >> 45)) & 0x7;
1339 LOG(info) << "CDW version: " << mCdwVersion;
1340 short int runtype = isCRUITS ? -2 : !mCdwVersion ? ((short int)(calib.calibUserField >> 24)) & 0xff
1341 : ((short int)(calib.calibUserField >> 9)) & 0x7f;
1342 mConfDBv = !mCdwVersion ? ((short int)(calib.calibUserField >> 32)) & 0xffff : ((short int)(calib.calibUserField >> 32)) & 0x1fff; // confDB version
1343 this->setRunType(runtype);
1344 LOG(info) << "Calibrator will ship these run parameters to aggregator:";
1345 LOG(info) << "Run type : " << mRunType;
1346 LOG(info) << "Scan type : " << mScanType;
1347 LOG(info) << "Fit type : " << std::to_string(mFitType);
1348 LOG(info) << "DB version (ignore in TOT_CALIB & VRESET2D): " << mConfDBv;
1349 }
1350 this->mRunTypeUp = isCRUITS ? -1 : !mCdwVersion ? ((short int)(calib.calibUserField >> 24)) & 0xff
1351 : ((short int)(calib.calibUserField >> 9)) & 0x7f;
1352
1353 // count the zeros
1354 if (!mRunTypeUp) {
1355 mRunTypeRU[iRU]++;
1356 mRunTypeRUCopy[iRU]++;
1357 }
1358 // Divide calibration word (24-bit) by 2^16 to get the first 8 bits
1359 if (this->mScanType == 'T') {
1360 // For threshold scan have to subtract from 170 to get charge value
1361 loopval = isCRUITS ? (short int)((calib.calibUserField >> 16) & 0xff) : !mCdwVersion ? (short int)(170 - (calib.calibUserField >> 16) & 0xff)
1362 : (short int)(170 - (calib.calibUserField >> 16) & 0xffff);
1363 } else if (this->mScanType == 'D' || this->mScanType == 'A') { // Digital scan
1364 loopval = 0;
1365 } else { // VCASN / ITHR tuning and Pulse length scan (it's the strobe delay in this case), and vresetd scan
1366 loopval = !mCdwVersion ? (short int)((calib.calibUserField >> 16) & 0xff) : (short int)((calib.calibUserField >> 16) & 0xffff);
1367 }
1368
1369 if (this->mScanType == 'p' || this->mScanType == 't' || this->mScanType == 'r') {
1370 realcharge = 170 - ((short int)(calib.calibUserField >> 32)) & 0x1fff; // not existing with CDW v0
1371 }
1372
1373 // Last 16 bits should be the row (only uses up to 9 bits)
1374 row = !mCdwVersion ? (short int)(calib.calibUserField & 0xffff) : (short int)(calib.calibUserField & 0x1ff);
1375 // cw counter
1376 cwcnt = (short int)(calib.calibCounter);
1377 // count the last N injections
1378 short int checkVal = (mScanType == 'I') ? mMin : mMax;
1379 if ((mScanType != 'r' && mScanType != 'p' && mScanType != 't' && loopval == checkVal) ||
1380 (mScanType == 'r' && realcharge == mMax2) ||
1381 (mScanType == 'p' && realcharge == mMin2) ||
1382 (mScanType == 't' && loopval == checkVal && realcharge == mMax2)) {
1383 mCdwCntRU[iRU][row]++;
1384 mLoopVal[iRU][row] = loopval; // keep loop val (relevant for VRESET2D and TOT_1ROW scan only)
1385 }
1386 if (this->mVerboseOutput) {
1387 LOG(info) << "RU: " << iRU << " CDWcounter: " << cwcnt << " row: " << row << " Loopval: " << loopval << " realcharge: " << realcharge << " confDBv: " << mCdwVersion;
1388 LOG(info) << "NDIGITS: " << digits.size();
1389 }
1390
1391 break;
1392 }
1393 }
1394
1395 if (isCRUITS && isAllZero) {
1396 if (mRunType == -1) {
1397 short int runtype = -2;
1398 mConfDBv = 0;
1399 this->setRunType(runtype);
1400 LOG(info) << "Running with CRU_ITS data - Calibrator will ship these run parameters to aggregator:";
1401 LOG(info) << "Run type (non-sense) : " << mRunType;
1402 LOG(info) << "Scan type : " << mScanType;
1403 LOG(info) << "Fit type : " << std::to_string(mFitType);
1404 LOG(info) << "DB version (non-sense): " << mConfDBv;
1405 }
1406 loopval = 0;
1407 realcharge = 0;
1408 row = 0;
1409 cwcnt = 0;
1410 }
1411
1412 if (loopval > this->mMax || loopval < this->mMin || ((mScanType == 'p' || mScanType == 't' || mScanType == 'r') && (realcharge > this->mMax2 || realcharge < this->mMin2))) {
1413 if (this->mVerboseOutput) {
1414 LOG(warning) << "CW issues - loopval value " << loopval << " out of range for min " << this->mMin
1415 << " and max " << this->mMax << " (range: " << N_RANGE << ")";
1416 if (mScanType == 'p' || mScanType == 'r' || mScanType == 't') {
1417 LOG(warning) << " and/or realcharge value " << realcharge << " out of range from min " << this->mMin2
1418 << " and max " << this->mMax2 << " (range: " << N_RANGE2 << ")";
1419 }
1420 }
1421 } else {
1422 std::vector<short int> mChips;
1423 // loop to retrieve list of chips and start tagging bad dcols if the hits does not come from this row
1424 for (unsigned int idig = rofIndex; idig < rofIndex + rofNEntries; idig++) { // gets chipid
1425 auto& d = digits[idig];
1426 short int chipID = (short int)d.getChipIndex();
1427 if ((chipID % mChipModBase) != mChipModSel) {
1428 continue;
1429 }
1430 if (d.getRow() != row && mVerboseOutput) {
1431 LOG(info) << "iROF: " << iROF << " ChipID " << chipID << ": current row is " << d.getRow() << " (col = " << d.getColumn() << ") but the one in CW is " << row;
1432 }
1433 if (std::find(mChips.begin(), mChips.end(), chipID) == mChips.end()) {
1434 mChips.push_back(chipID);
1435 }
1436 }
1437 // loop to allocate memory only for allowed rows
1438 for (auto& chipID : mChips) {
1439 // mark active RU links
1440 short int ru = getRUID(chipID);
1441 mActiveLinks[ru][getLinkID(chipID, ru)] = true;
1442 // check rows and allocate memory
1443 if (!this->mPixelHits.count(chipID)) {
1444 if (mScanType == 'D' || mScanType == 'A') { // for digital and analog scan initialize the full matrix for each chipID
1445 for (int irow = 0; irow < 512; irow++) {
1446 this->mPixelHits[chipID][irow] = std::vector<std::vector<std::vector<unsigned short int>>>(this->N_COL, std::vector<std::vector<unsigned short int>>(N_RANGE2, std::vector<unsigned short int>(N_RANGE, 0)));
1447 }
1448 } else {
1449 this->mPixelHits[chipID][row] = std::vector<std::vector<std::vector<unsigned short int>>>(this->N_COL, std::vector<std::vector<unsigned short int>>(N_RANGE2, std::vector<unsigned short int>(N_RANGE, 0)));
1450 }
1451 } else if (!this->mPixelHits[chipID].count(row)) { // allocate memory for chip = chipID or for a row of this chipID
1452 this->mPixelHits[chipID][row] = std::vector<std::vector<std::vector<unsigned short int>>>(this->N_COL, std::vector<std::vector<unsigned short int>>(N_RANGE2, std::vector<unsigned short int>(N_RANGE, 0)));
1453 }
1454 }
1455
1456 // loop to count hits from digits
1457 short int loopPoint = (loopval - this->mMin) / mStep;
1458 short int chgPoint = (realcharge - this->mMin2) / mStep2;
1459 for (unsigned int idig = rofIndex; idig < rofIndex + rofNEntries; idig++) {
1460 auto& d = digits[idig];
1461 short int chipID = (short int)d.getChipIndex();
1462 short int col = (short int)d.getColumn();
1463
1464 if ((chipID % mChipModBase) != mChipModSel) {
1465 continue;
1466 }
1467
1468 if ((!mCheckExactRow || d.getRow() == row) && (mMeb < 0 || cwcnt % 3 == mMeb)) { // row has NOT to be forbidden and we ignore hits coming from other rows (potential masking issue on chip)
1469 // Increment the number of counts for this pixel
1470 this->mPixelHits[chipID][d.getRow()][col][chgPoint][loopPoint]++;
1471 }
1472 }
1473 } // if (charge)
1474
1476 // Prepare the ChipDone object for QC + extract data if the row is completed
1477 if (ruIndex < 0) {
1478 continue;
1479 }
1480 short int nL = 0;
1481 for (int iL = 0; iL < 3; iL++) {
1482 if (mActiveLinks[ruIndex][iL]) {
1483 nL++; // count active links
1484 }
1485 }
1486 std::vector<short int> chipEnabled = getChipListFromRu(ruIndex, mActiveLinks[ruIndex]); // chip boundaries
1487 // Fill the chipDone info string
1488 if (mRunTypeRUCopy[ruIndex] == nInjScaled * nL) {
1489 for (short int iChip = 0; iChip < chipEnabled.size(); iChip++) {
1490 if ((chipEnabled[iChip] % mChipModBase) != mChipModSel) {
1491 continue;
1492 }
1493 addDatabaseEntry(chipEnabled[iChip], "", std::vector<float>(), true);
1494 }
1495 mRunTypeRUCopy[ruIndex] = 0; // reset here is safer (the other counter is reset in finalize)
1496 }
1497 // Check if scan of a row is finished: only for specific scans!
1498 bool passCondition = (mCdwCntRU[ruIndex][row] >= nInjScaled * nL);
1499 if (mScanType == 'p' || mScanType == 't') {
1500 passCondition = passCondition && (mLoopVal[ruIndex][row] == mMax);
1501 if (mVerboseOutput) {
1502 LOG(info) << "PassCondition: " << passCondition << " - (mCdwCntRU,mLoopVal) of RU" << ruIndex << " row " << row << " = (" << mCdwCntRU[ruIndex][row] << ", " << mLoopVal[ruIndex][row] << ")";
1503 }
1504 } else if (mVerboseOutput) {
1505 LOG(info) << "PassCondition: " << passCondition << " - mCdwCntRU of RU" << ruIndex << " row " << row << " = " << mCdwCntRU[ruIndex][row];
1506 }
1507
1508 if (mScanType != 'D' && mScanType != 'A' && mScanType != 'P' && mScanType != 'R' && passCondition) {
1509 // extract data from the row
1510 for (short int iChip = 0; iChip < chipEnabled.size(); iChip++) {
1511 short int chipID = chipEnabled[iChip];
1512 if ((chipID % mChipModBase) != mChipModSel) {
1513 continue;
1514 }
1515 if (!isDumpS || (std::find(chipDumpList.begin(), chipDumpList.end(), chipID) != chipDumpList.end() || !chipDumpList.size())) { // to dump s-curves as histograms
1516 if (mPixelHits.count(chipID)) {
1517 if (mPixelHits[chipID].count(row)) { // make sure the row exists
1518 extractAndUpdate(chipID, row);
1519 if (mScanType != 'p' && (mScanType != 'r' || mLoopVal[ruIndex][row] == mMax)) { // do not erase for scantype = p because in finalize() we have calculate2Dparams
1520 mPixelHits[chipID].erase(row);
1521 }
1522 }
1523 }
1524 }
1525 }
1526 mCdwCntRU[ruIndex][row] = 0; // reset
1527 }
1528 } // for (ROFs)
1529
1530 if (!(this->mRunTypeUp)) {
1531 finalize();
1532 LOG(info) << "Shipping all outputs to aggregator (before endOfStream arrival!)";
1533 pc.outputs().snapshot(Output{"ITS", "TSTR", (unsigned int)mChipModSel}, this->mTuning);
1534 pc.outputs().snapshot(Output{"ITS", "PIXTYP", (unsigned int)mChipModSel}, this->mPixStat);
1535 pc.outputs().snapshot(Output{"ITS", "RUNT", (unsigned int)mChipModSel}, this->mRunType);
1536 pc.outputs().snapshot(Output{"ITS", "SCANT", (unsigned int)mChipModSel}, this->mScanType);
1537 pc.outputs().snapshot(Output{"ITS", "FITT", (unsigned int)mChipModSel}, this->mFitType);
1538 pc.outputs().snapshot(Output{"ITS", "CONFDBV", (unsigned int)mChipModSel}, this->mConfDBv);
1539 pc.outputs().snapshot(Output{"ITS", "QCSTR", (unsigned int)mChipModSel}, this->mChipDoneQc);
1540 // reset the DCSconfigObject_t before next ship out
1541 mTuning.clear();
1542 mPixStat.clear();
1543 mChipDoneQc.clear();
1544 } else if (pc.transitionState() == TransitionHandlingState::Requested) {
1545 LOG(info) << "Run stop requested during the scan, sending output to aggregator and then stopping to process new data";
1546 mRunStopRequested = true;
1547 finalize(); // calculating average thresholds based on what's collected up to this moment
1548 pc.outputs().snapshot(Output{"ITS", "TSTR", (unsigned int)mChipModSel}, this->mTuning); // dummy here
1549 pc.outputs().snapshot(Output{"ITS", "PIXTYP", (unsigned int)mChipModSel}, this->mPixStat);
1550 pc.outputs().snapshot(Output{"ITS", "RUNT", (unsigned int)mChipModSel}, this->mRunType);
1551 pc.outputs().snapshot(Output{"ITS", "SCANT", (unsigned int)mChipModSel}, this->mScanType);
1552 pc.outputs().snapshot(Output{"ITS", "FITT", (unsigned int)mChipModSel}, this->mFitType);
1553 pc.outputs().snapshot(Output{"ITS", "CONFDBV", (unsigned int)mChipModSel}, this->mConfDBv);
1554 pc.outputs().snapshot(Output{"ITS", "QCSTR", (unsigned int)mChipModSel}, this->mChipDoneQc);
1555 mChipDoneQc.clear();
1556 mPixStat.clear();
1557 mTuning.clear();
1558 }
1559
1560 return;
1561}
1562
1564// Retrieve conf DB map from production ccdb
1566{
1567 if (matcher == ConcreteDataMatcher("ITS", "CONFDBMAP", 0)) {
1568 LOG(info) << "Conf DB map retrieved from CCDB";
1569 mConfDBmap = (std::vector<int>*)obj;
1570 }
1571}
1572
1574// Calculate the average threshold given a vector of threshold objects
1575void ITSThresholdCalibrator::findAverage(const std::array<long int, 6>& data, float& avgT, float& rmsT, float& avgN, float& rmsN)
1576{
1577 avgT = (!data[4]) ? 0. : (float)data[0] / (float)data[4];
1578 rmsT = (!data[4]) ? 0. : std::sqrt((float)data[1] / (float)data[4] - avgT * avgT);
1579 avgN = (!data[4]) ? 0. : (float)data[2] / (float)data[4];
1580 rmsN = (!data[4]) ? 0. : std::sqrt((float)data[3] / (float)data[4] - avgN * avgN);
1581 return;
1582}
1583
1585void ITSThresholdCalibrator::addDatabaseEntry(
1586 const short int& chipID, const char* name, std::vector<float> data, bool isQC)
1587{
1588
1589 /*--------------------------------------------------------------------
1590 // Format of *data
1591 // - Threshold scan: avgT, rmsT, avgN, rmsN, status
1592 // - ITHR/VCASN scan: avg value, rms value, 0, 0, status (0 are just placeholder since they will not be used)
1593 // - Dig/Ana scans: empty vector
1594 // - Pulse shape 1D: avgRt, rmsRt, avgToT, rmsToT
1595 // - Pulse shape 2D: avgToT, rmsToT, MTC, rmsMTC, avgMTCD, rmsMTCD, avgMPL, rmsMPL, avgMPLC, rmsMPLC
1596 */
1597 // Obtain specific chip information from the chip ID (layer, stave, ...)
1598 int lay, sta, ssta, mod, chipInMod; // layer, stave, sub stave, module, chip
1599 this->mp.expandChipInfoHW(chipID, lay, sta, ssta, mod, chipInMod);
1600
1601 char stave[6];
1602 snprintf(stave, 6, "L%d_%02d", lay, sta);
1603
1604 if (isQC) {
1605 o2::dcs::addConfigItem(this->mChipDoneQc, "O2ChipID", std::to_string(chipID));
1606 return;
1607 }
1608
1609 // Get ConfDB id for the chip chipID
1610 int confDBid = (*mConfDBmap)[chipID];
1611
1612 // Bad pix list and bad dcols for dig and ana scan
1613 if (this->mScanType == 'D' || this->mScanType == 'A') {
1614 short int vPixDcolCounter[512] = {0}; // count #bad_pix per dcol
1615 std::string dcolIDs = "";
1616 std::string pixIDs_Noisy = "";
1617 std::string pixIDs_Dead = "";
1618 std::string pixIDs_Ineff = "";
1619 std::vector<int>& v = PixelType == "Noisy" ? mNoisyPixID[chipID] : PixelType == "Dead" ? mDeadPixID[chipID]
1620 : mIneffPixID[chipID];
1621 // Number of pixel types
1622 int n_pixel = v.size(), nDcols = 0;
1623 std::string ds = "-1"; // dummy string
1624 // find bad dcols and add them one by one
1625 if (PixelType == "Noisy") {
1626 for (int i = 0; i < v.size(); i++) {
1627 short int dcol = ((v[i] - v[i] % 1000) / 1000) / 2;
1628 vPixDcolCounter[dcol]++;
1629 }
1630 for (int i = 0; i < 512; i++) {
1631 if (vPixDcolCounter[i] > N_PIX_DCOL) {
1632 o2::dcs::addConfigItem(this->mTuning, "O2ChipID", std::to_string(chipID));
1633 o2::dcs::addConfigItem(this->mTuning, "ChipDbID", std::to_string(confDBid));
1634 o2::dcs::addConfigItem(this->mTuning, "Dcol", std::to_string(i));
1635 o2::dcs::addConfigItem(this->mTuning, "Row", ds);
1636 o2::dcs::addConfigItem(this->mTuning, "Col", ds);
1637
1638 dcolIDs += std::to_string(i) + '|'; // prepare string for second object for ccdb prod
1639 nDcols++;
1640 }
1641 }
1642 }
1643
1644 if (this->mTagSinglePix) {
1645 if (PixelType == "Noisy") {
1646 for (int i = 0; i < v.size(); i++) {
1647 short int dcol = ((v[i] - v[i] % 1000) / 1000) / 2;
1648 if (vPixDcolCounter[dcol] > N_PIX_DCOL) { // single pixels must not be already in dcolIDs
1649 continue;
1650 }
1651
1652 // Noisy pixel IDs
1653 pixIDs_Noisy += std::to_string(v[i]);
1654 if (i + 1 < v.size()) {
1655 pixIDs_Noisy += '|';
1656 }
1657
1658 o2::dcs::addConfigItem(this->mTuning, "O2ChipID", std::to_string(chipID));
1659 o2::dcs::addConfigItem(this->mTuning, "ChipDbID", std::to_string(confDBid));
1660 o2::dcs::addConfigItem(this->mTuning, "Dcol", ds);
1661 o2::dcs::addConfigItem(this->mTuning, "Row", std::to_string(v[i] % 1000));
1662 o2::dcs::addConfigItem(this->mTuning, "Col", std::to_string(int((v[i] - v[i] % 1000) / 1000)));
1663 }
1664 }
1665
1666 if (PixelType == "Dead") {
1667 for (int i = 0; i < v.size(); i++) {
1668 pixIDs_Dead += std::to_string(v[i]);
1669 if (i + 1 < v.size()) {
1670 pixIDs_Dead += '|';
1671 }
1672 }
1673 }
1674
1675 if (PixelType == "Ineff") {
1676 for (int i = 0; i < v.size(); i++) {
1677 pixIDs_Ineff += std::to_string(v[i]);
1678 if (i + 1 < v.size()) {
1679 pixIDs_Ineff += '|';
1680 }
1681 }
1682 }
1683 }
1684
1685 if (!dcolIDs.empty()) {
1686 dcolIDs.pop_back(); // remove last pipe from the string
1687 } else {
1688 dcolIDs = "-1";
1689 }
1690
1691 if (pixIDs_Noisy.empty()) {
1692 pixIDs_Noisy = "-1";
1693 }
1694
1695 if (pixIDs_Dead.empty()) {
1696 pixIDs_Dead = "-1";
1697 }
1698
1699 if (pixIDs_Ineff.empty()) {
1700 pixIDs_Ineff = "-1";
1701 }
1702
1703 o2::dcs::addConfigItem(this->mPixStat, "O2ChipID", std::to_string(chipID));
1704 if (PixelType == "Dead" || PixelType == "Ineff") {
1705 o2::dcs::addConfigItem(this->mPixStat, "PixelType", PixelType);
1706 o2::dcs::addConfigItem(this->mPixStat, "PixelNos", n_pixel);
1707 o2::dcs::addConfigItem(this->mPixStat, "DcolNos", "-1");
1708 } else {
1709 o2::dcs::addConfigItem(this->mPixStat, "PixelType", PixelType);
1710 o2::dcs::addConfigItem(this->mPixStat, "PixelNos", n_pixel);
1711 o2::dcs::addConfigItem(this->mPixStat, "DcolNos", nDcols);
1712 }
1713 }
1714 if (this->mScanType != 'D' && this->mScanType != 'A' && this->mScanType != 'P' && this->mScanType != 'p') {
1715 o2::dcs::addConfigItem(this->mTuning, "O2ChipID", std::to_string(chipID));
1716 o2::dcs::addConfigItem(this->mTuning, "ChipDbID", std::to_string(confDBid));
1717 o2::dcs::addConfigItem(this->mTuning, name, (strcmp(name, "ITHR") == 0 || strcmp(name, "VCASN") == 0) ? std::to_string((int)data[0]) : std::to_string(data[0]));
1718 o2::dcs::addConfigItem(this->mTuning, "Rms", std::to_string(data[1]));
1719 o2::dcs::addConfigItem(this->mTuning, "Status", std::to_string((int)data[4])); // percentage of unsuccess
1720 }
1721 if (this->mScanType == 'T') {
1722 o2::dcs::addConfigItem(this->mTuning, "Noise", std::to_string(data[2]));
1723 o2::dcs::addConfigItem(this->mTuning, "NoiseRms", std::to_string(data[3]));
1724 }
1725
1726 if (this->mScanType == 'P') {
1727 o2::dcs::addConfigItem(this->mTuning, "O2ChipID", std::to_string(chipID));
1728 o2::dcs::addConfigItem(this->mTuning, "ChipDbID", std::to_string(confDBid));
1729 o2::dcs::addConfigItem(this->mTuning, "Tot", std::to_string(data[2])); // time over threshold
1730 o2::dcs::addConfigItem(this->mTuning, "TotRms", std::to_string(data[3])); // time over threshold rms
1731 o2::dcs::addConfigItem(this->mTuning, "Rt", std::to_string(data[0])); // rise time
1732 o2::dcs::addConfigItem(this->mTuning, "RtRms", std::to_string(data[1])); // rise time rms
1733 }
1734
1735 //- Pulse shape 2D: avgToT, rmsToT, MTC, rmsMTC, avgMTCD, rmsMTCD, avgMPL, rmsMPL, avgMPLC, rmsMPLC
1736 if (this->mScanType == 'p') {
1737 o2::dcs::addConfigItem(this->mTuning, "O2ChipID", std::to_string(chipID));
1738 o2::dcs::addConfigItem(this->mTuning, "ChipDbID", std::to_string(confDBid));
1739 o2::dcs::addConfigItem(this->mTuning, "Tot", std::to_string(data[0])); // time over threshold
1740 o2::dcs::addConfigItem(this->mTuning, "TotRms", std::to_string(data[1])); // time over threshold rms
1741 o2::dcs::addConfigItem(this->mTuning, "MinThrChg", std::to_string(data[2])); // Min threshold charge
1742 o2::dcs::addConfigItem(this->mTuning, "MinThrChgRms", std::to_string(data[3])); // Min threshold charge rms
1743 o2::dcs::addConfigItem(this->mTuning, "MinThrChgDel", std::to_string(data[4])); // Min threshold charge delay
1744 o2::dcs::addConfigItem(this->mTuning, "MinThrChgDelRms", std::to_string(data[5])); // Min threshold charge delay rms
1745 o2::dcs::addConfigItem(this->mTuning, "MaxPulLen", std::to_string(data[6])); // Max pulse length
1746 o2::dcs::addConfigItem(this->mTuning, "MaxPulLenRms", std::to_string(data[7])); // Max pulse length rms
1747 o2::dcs::addConfigItem(this->mTuning, "MaxPulLenChg", std::to_string(data[8])); // Max pulse length charge
1748 o2::dcs::addConfigItem(this->mTuning, "MaxPulLenChgRms", std::to_string(data[9])); // Max pulse length charge rms
1749 }
1750
1751 return;
1752}
1753
1754//___________________________________________________________________
1755void ITSThresholdCalibrator::updateTimeDependentParams(ProcessingContext& pc)
1756{
1757 static bool initOnceDone = false;
1758 if (!initOnceDone) {
1759 initOnceDone = true;
1760 mDataTakingContext = pc.services().get<DataTakingContext>();
1761 }
1762 mTimingInfo = pc.services().get<o2::framework::TimingInfo>();
1763}
1764
1767{
1768 // Add configuration item to output strings for CCDB
1769 const char* name = nullptr;
1770 std::set<int> thisRUs;
1771
1772 if (mScanType == 'V' || mScanType == 'I' || mScanType == 'T') {
1773 // Loop over each chip and calculate avg and rms
1774 name = mScanType == 'V' ? "VCASN" : mScanType == 'I' ? "ITHR"
1775 : "THR";
1776 if (mScanType == 'I') {
1777 // Only ITHR scan: assign default ITHR = 50 if chip has no avg ITHR
1778 for (auto& iRU : mRuSet) {
1779 if (mRunTypeRU[iRU] >= nInjScaled * getNumberOfActiveLinks(mActiveLinks[iRU]) || mRunStopRequested) {
1780 std::vector<short int> chipList = getChipListFromRu(iRU, mActiveLinks[iRU]);
1781 for (size_t i = 0; i < chipList.size(); i++) {
1782 if ((chipList[i] % mChipModBase) != mChipModSel) {
1783 continue;
1784 }
1785 if (!mThresholds.count(chipList[i])) {
1786 if (mVerboseOutput) {
1787 LOG(info) << "Setting ITHR = 50 for chip " << chipList[i];
1788 }
1789 std::vector<float> data = {50, 0, 0, 0, 0};
1790 addDatabaseEntry(chipList[i], name, data, false);
1791 }
1792 }
1793 }
1794 }
1795 }
1796
1797 auto it = this->mThresholds.cbegin();
1798 while (it != this->mThresholds.cend()) {
1799 short int iRU = getRUID(it->first);
1800 if (!isCRUITS && (mRunTypeRU[iRU] < nInjScaled * getNumberOfActiveLinks(mActiveLinks[iRU]) && !mRunStopRequested)) {
1801 ++it;
1802 continue;
1803 }
1804 thisRUs.insert(iRU);
1805 float avgT, rmsT, avgN, rmsN, mpvT, outVal;
1806 this->findAverage(it->second, avgT, rmsT, avgN, rmsN);
1807 outVal = avgT;
1808 if (isMpv) {
1809 mpvT = std::distance(mpvCounter[it->first].begin(), std::max_element(mpvCounter[it->first].begin(), mpvCounter[it->first].end())) + mMin;
1810 outVal = mpvT;
1811 }
1812 if (mVerboseOutput) {
1813 LOG(info) << "Average or mpv " << name << " of chip " << it->first << " = " << outVal << " e-";
1814 }
1815 float status = ((float)it->second[4] / (float)(it->second[4] + it->second[5])) * 100.; // percentage of successful threshold extractions
1816 if (status < mPercentageCut && (mScanType == 'I' || mScanType == 'V')) {
1817 if (mScanType == 'I') { // default ITHR if percentage of success < mPercentageCut
1818 outVal = 50.;
1819 if (mVerboseOutput) {
1820 LOG(info) << "Chip " << it->first << " status is " << status << ". Setting ITHR = 50";
1821 }
1822 } else { // better to not set any VCASN if the percentage of success < mPercentageCut
1823 it = this->mThresholds.erase(it);
1824 if (mVerboseOutput) {
1825 LOG(info) << "Chip " << it->first << " status is " << status << ". Ignoring this chip.";
1826 }
1827 continue;
1828 }
1829 }
1830 std::vector<float> data = {outVal, rmsT, avgN, rmsN, status};
1831 this->addDatabaseEntry(it->first, name, data, false);
1832 it = this->mThresholds.erase(it);
1833 }
1834 } else if (this->mScanType == 'D' || this->mScanType == 'A') {
1835 // Loop over each chip and calculate avg and rms
1836 name = "PixID";
1837 // Extract hits from the full matrix
1838 auto itchip = this->mPixelHits.cbegin();
1839 while (itchip != this->mPixelHits.cend()) { // loop over chips collected
1840 short int iRU = getRUID(itchip->first);
1841 if (!isCRUITS && (mRunTypeRU[iRU] < nInjScaled * getNumberOfActiveLinks(mActiveLinks[iRU]) && !mRunStopRequested)) {
1842 ++itchip;
1843 continue;
1844 }
1845 thisRUs.insert(iRU);
1846 if (mVerboseOutput) {
1847 LOG(info) << "Extracting hits for the full matrix of chip " << itchip->first;
1848 }
1849 for (short int irow = 0; irow < 512; irow++) {
1850 this->extractAndUpdate(itchip->first, irow);
1851 }
1852 if (this->mVerboseOutput) {
1853 LOG(info) << "Chip " << itchip->first << " hits extracted";
1854 }
1855 ++itchip;
1856 }
1857
1858 auto it = this->mNoisyPixID.cbegin();
1859 while (it != this->mNoisyPixID.cend()) {
1860 PixelType = "Noisy";
1861 if (mVerboseOutput) {
1862 LOG(info) << "Extracting noisy pixels in the full matrix of chip " << it->first;
1863 }
1864 this->addDatabaseEntry(it->first, name, std::vector<float>(), false); // all zeros are not used here
1865 if (this->mVerboseOutput) {
1866 LOG(info) << "Chip " << it->first << " done";
1867 }
1868 it = this->mNoisyPixID.erase(it);
1869 }
1870
1871 auto it_d = this->mDeadPixID.cbegin();
1872 while (it_d != this->mDeadPixID.cend()) {
1873 if (mVerboseOutput) {
1874 LOG(info) << "Extracting dead pixels in the full matrix of chip " << it_d->first;
1875 }
1876 PixelType = "Dead";
1877 this->addDatabaseEntry(it_d->first, name, std::vector<float>(), false); // all zeros are not used here
1878 it_d = this->mDeadPixID.erase(it_d);
1879 }
1880
1881 auto it_ineff = this->mIneffPixID.cbegin();
1882 while (it_ineff != this->mIneffPixID.cend()) {
1883 if (mVerboseOutput) {
1884 LOG(info) << "Extracting inefficient pixels in the full matrix of chip " << it_ineff->first;
1885 }
1886 PixelType = "Ineff";
1887 this->addDatabaseEntry(it_ineff->first, name, std::vector<float>(), false);
1888 it_ineff = this->mIneffPixID.erase(it_ineff);
1889 }
1890 } else if (this->mScanType == 'P' || this->mScanType == 'p' || mScanType == 'R') { // pulse length scan 1D and 2D, vresetd scan 1D (2D already extracted in run())
1891 name = "Pulse";
1892 // extract hits for the available row(s)
1893 auto itchip = this->mPixelHits.cbegin();
1894 while (itchip != mPixelHits.cend()) {
1895 int iRU = getRUID(itchip->first);
1896 if (!mRunStopRequested && mRunTypeRU[iRU] < nInjScaled * getNumberOfActiveLinks(mActiveLinks[iRU])) {
1897 ++itchip;
1898 continue;
1899 }
1900 thisRUs.insert(iRU);
1901 if (mVerboseOutput) {
1902 LOG(info) << "Extracting hits from pulse shape scan or vresetd scan, chip " << itchip->first;
1903 }
1904
1905 if (mScanType != 'p') { // done already in run()
1906 auto itrow = this->mPixelHits[itchip->first].cbegin();
1907 while (itrow != mPixelHits[itchip->first].cend()) { // in case there are multiple rows, for now it's 1 row
1908 this->extractAndUpdate(itchip->first, itrow->first); // fill the tree - for mScanType = p, it is done already in run()
1909 ++itrow;
1910 }
1911 }
1912
1913 if (mCalculate2DParams && (mScanType == 'P' || mScanType == 'p')) {
1914 this->addDatabaseEntry(itchip->first, name, mScanType == 'P' ? calculatePulseParams(itchip->first) : calculatePulseParams2D(itchip->first), false);
1915 }
1916 if (this->mVerboseOutput) {
1917 LOG(info) << "Chip " << itchip->first << " hits extracted";
1918 }
1919 ++itchip;
1920 }
1921 // reset RU counters so that the chips which are done will not appear again in the DCSConfigObject
1922 }
1923
1924 for (auto& ru : thisRUs) {
1925 mRunTypeRU[ru] = 0; // reset
1926 }
1927
1928 return;
1929}
1930
1932// O2 functionality allowing to do post-processing when the upstream device
1933// tells that there will be no more input data
1935{
1936 if (!isEnded && !mRunStopRequested) {
1937 LOGF(info, "endOfStream report:", mSelfName);
1938 if (isCRUITS) {
1939 finalize();
1940 }
1941 this->finalizeOutput();
1942 isEnded = true;
1943 }
1944 return;
1945}
1946
1948// DDS stop method: simply close the latest tree
1950{
1951 if (!isEnded) {
1952 LOGF(info, "stop() report:", mSelfName);
1953 this->finalizeOutput();
1954 isEnded = true;
1955 }
1956 return;
1957}
1958
1961{
1963 std::vector<InputSpec> inputs;
1964 inputs.emplace_back("digits", detOrig, "DIGITS", 0, Lifetime::Timeframe);
1965 inputs.emplace_back("digitsROF", detOrig, "DIGITSROF", 0, Lifetime::Timeframe);
1966 inputs.emplace_back("calib", detOrig, "GBTCALIB", 0, Lifetime::Timeframe);
1967 // inputs.emplace_back("confdbmap", detOrig, "CONFDBMAP", 0, Lifetime::Condition,
1968 // o2::framework::ccdbParamSpec("ITS/Calib/Confdbmap"));
1969
1970 std::vector<OutputSpec> outputs;
1971 outputs.emplace_back("ITS", "TSTR", inpConf.chipModSel);
1972 outputs.emplace_back("ITS", "PIXTYP", inpConf.chipModSel);
1973 outputs.emplace_back("ITS", "RUNT", inpConf.chipModSel);
1974 outputs.emplace_back("ITS", "SCANT", inpConf.chipModSel);
1975 outputs.emplace_back("ITS", "FITT", inpConf.chipModSel);
1976 outputs.emplace_back("ITS", "CONFDBV", inpConf.chipModSel);
1977 outputs.emplace_back("ITS", "QCSTR", inpConf.chipModSel);
1978
1979 return DataProcessorSpec{
1980 "its-calibrator_" + std::to_string(inpConf.chipModSel),
1981 inputs,
1982 outputs,
1983 AlgorithmSpec{adaptFromTask<ITSThresholdCalibrator>(inpConf)},
1984 Options{{"fittype", VariantType::String, "derivative", {"Fit type to extract thresholds, with options: fit, derivative (default), hitcounting"}},
1985 {"verbose", VariantType::Bool, false, {"Use verbose output mode"}},
1986 {"output-dir", VariantType::String, "./", {"ROOT trees output directory"}},
1987 {"meta-output-dir", VariantType::String, "/dev/null", {"Metadata output directory"}},
1988 {"meta-type", VariantType::String, "", {"metadata type"}},
1989 {"nthreads", VariantType::Int, 1, {"Number of threads, default is 1"}},
1990 {"enable-cw-cnt-check", VariantType::Bool, false, {"Use to enable the check of the calib word counter row by row in addition to the hits"}},
1991 {"enable-single-pix-tag", VariantType::Bool, false, {"Use to enable tagging of single noisy pix in digital and analogue scan"}},
1992 {"ccdb-mgr-url", VariantType::String, "", {"CCDB url to download confDBmap"}},
1993 {"min-vcasn", VariantType::Int, 30, {"Min value of VCASN in vcasn scan, default is 30"}},
1994 {"max-vcasn", VariantType::Int, 70, {"Max value of VCASN in vcasn scan, default is 70"}},
1995 {"min-ithr", VariantType::Int, 25, {"Min value of ITHR in ithr scan, default is 25"}},
1996 {"max-ithr", VariantType::Int, 100, {"Max value of ITHR in ithr scan, default is 100"}},
1997 {"manual-mode", VariantType::Bool, false, {"Flag to activate the manual mode in case run type is not recognized"}},
1998 {"manual-min", VariantType::Int, 0, {"Min value of the variable used for the scan: use only in manual mode"}},
1999 {"manual-max", VariantType::Int, 50, {"Max value of the variable used for the scan: use only in manual mode"}},
2000 {"manual-min2", VariantType::Int, 0, {"Min2 value of the 2nd variable (if any) used for the scan (ex: charge in tot_calib): use only in manual mode"}},
2001 {"manual-max2", VariantType::Int, 50, {"Max2 value of the 2nd variable (if any) used for the scan (ex: charge in tot_calib): use only in manual mode"}},
2002 {"manual-step", VariantType::Int, 1, {"Step value: defines the steps between manual-min and manual-max. Default is 1. Use only in manual mode"}},
2003 {"manual-step2", VariantType::Int, 1, {"Step2 value: defines the steps between manual-min2 and manual-max2. Default is 1. Use only in manual mode"}},
2004 {"manual-scantype", VariantType::String, "T", {"scan type, can be D, T, I, V, P, p: use only in manual mode"}},
2005 {"manual-strobewindow", VariantType::Int, 5, {"strobe duration in clock cycles, default is 5 = 125 ns: use only in manual mode"}},
2006 {"save-tree", VariantType::Bool, false, {"Flag to save ROOT tree on disk: use only in manual mode"}},
2007 {"scale-ninj", VariantType::Bool, false, {"Flag to activate the scale of the number of injects to be used to count hits from specific MEBs: use only in manual mode and in combination with --meb-select"}},
2008 {"enable-mpv", VariantType::Bool, false, {"Flag to enable calculation of most-probable value in vcasn/ithr scans"}},
2009 {"enable-cru-its", VariantType::Bool, false, {"Flag to enable the analysis of raw data on disk produced by CRU_ITS IB commissioning tools"}},
2010 {"ninj", VariantType::Int, 50, {"Number of injections per change, default is 50"}},
2011 {"dump-scurves", VariantType::Bool, false, {"Dump any s-curve to disk in ROOT file. Works only with fit option."}},
2012 {"max-dump", VariantType::Int, -1, {"Maximum number of s-curves to dump in ROOT file per chip. Works with fit option and dump-scurves flag enabled. Default: dump all"}},
2013 {"chip-dump", VariantType::String, "", {"Dump s-curves only for these Chip IDs (0 to 24119). If multiple IDs, write them separated by comma. Default is empty string: dump all"}},
2014 {"calculate-slope", VariantType::Bool, false, {"For Pulse Shape 2D: if enabled it calculate the slope of the charge vs strobe delay trend for each pixel and fill it in the output tree"}},
2015 {"charge-a", VariantType::Int, 0, {"To use with --calculate-slope, it defines the charge (in DAC) for the 1st point used for the slope calculation"}},
2016 {"charge-b", VariantType::Int, 0, {"To use with --calculate-slope, it defines the charge (in DAC) for the 2nd point used for the slope calculation"}},
2017 {"meb-select", VariantType::Int, -1, {"Select from which multi-event buffer consider the hits: 0,1 or 2"}},
2018 {"s-curve-col-step", VariantType::Int, 8, {"save s-curves points to tree every s-curve-col-step pixels on 1 row"}},
2019 {"percentage-cut", VariantType::Int, 25, {"discard chip in ITHR/VCASN scan if the percentage of success is less than this cut"}}}};
2020}
2021} // namespace its
2022} // namespace o2
int32_t i
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
uint32_t col
Definition RawData.h:4
uint32_t c
Definition RawData.h:2
static BasicCCDBManager & instance()
void snapshot(const Output &spec, T const &object)
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
InputRecord & inputs()
The inputs associated with this processing context.
ServiceRegistryRef services()
The services registry associated with this processing context.
TransitionHandlingState transitionState() const
void run(ProcessingContext &pc) final
void stop() final
This is invoked on stop.
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
ITSThresholdCalibrator(const ITSCalibInpConf &inpConf)
void expandChipInfoHW(int idSW, int &lay, int &sta, int &ssta, int &mod, int &chipInMod) const
convert global SW chip ID to name in HW conventions
static constexpr std::string_view getName()
GLint GLenum GLint x
Definition glcorearb.h:403
GLint GLsizei count
Definition glcorearb.h:399
GLuint64EXT * result
Definition glcorearb.h:5662
GLuint buffer
Definition glcorearb.h:655
GLsizeiptr size
Definition glcorearb.h:659
const GLdouble * v
Definition glcorearb.h:832
GLuint const GLchar * name
Definition glcorearb.h:781
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLenum GLsizei GLsizei GLint * values
Definition glcorearb.h:1576
GLboolean * data
Definition glcorearb.h:298
GLuint start
Definition glcorearb.h:469
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
constexpr o2::header::DataOrigin gDataOriginITS
Definition DataHeader.h:570
long getCurrentTimestamp()
returns the timestamp in long corresponding to "now"
void addConfigItem(DCSconfigObject_t &configVector, std::string key, const T value)
std::vector< ConfigParamSpec > Options
double erf_ithr(double *xx, double *par)
double erf(double *xx, double *par)
o2::framework::DataProcessorSpec getITSThresholdCalibratorSpec(const ITSCalibInpConf &inpConf)
return * this
Enum< T >::Iterator begin(Enum< T >)
Definition Defs.h:173
void createDirectoriesIfAbsent(std::string const &path)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
std::string filename()
std::string asString() const
bool fillFileData(const std::string &fname, bool fillmd5=false, const std::string &tmpEnding="")
void setDataTakingContext(const o2::framework::DataTakingContext &dtc)
std::string envId
The environment ID for the deployment.
std::string runNumber
The current run number.
static std::string rectifyDirectory(const std::string_view p)
static std::string concat_string(Ts const &... ts)
o2::mch::DsIndex ds
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< Digit > digits
std::vector< int > row