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