Project
Loading...
Searching...
No Matches
CalibTOFapi.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
12#include "TOFBase/CalibTOFapi.h"
13#include <fairlogger/Logger.h> // for LOG
14#include <TH2F.h>
15
16using namespace o2::tof;
17
19
20o2::tof::Diagnostic CalibTOFapi::doDRMerrCalibFromQCHisto(const TH2F* histo, const char* file_output_name)
21{
22 // this is a method which translate the QC output in qc/TOF/MO/TaskRaw/DRMCounter (TH2F) into a Diagnotic object for DRM (patter(crate, error), frequency)
23 // note that, differently from TRM errors, DRM ones are not stored in CTF by design (since very rare, as expected). Such an info is available only at the level of raw sync QC
25
26 for (int j = 1; j <= Geo::kNCrate; j++) {
27 drmDia.fillDRM(j - 1, histo->GetBinContent(1, j));
28 for (int i = 2; i <= histo->GetXaxis()->GetNbins(); i++) {
29 if (histo->GetBinContent(1, j)) {
30 if (histo->GetBinContent(i, j) > 0) {
31 drmDia.fillDRMerror(j - 1, i - 1, histo->GetBinContent(i, j));
32 }
33 }
34 }
35 }
36
37 TFile* fo = new TFile(file_output_name, "RECREATE");
38 fo->WriteObjectAny(&drmDia, drmDia.Class_Name(), "ccdb_object");
39 fo->Close();
40 LOG(debug) << "DRM error ccdb object created in " << file_output_name << " with this content";
41 drmDia.print(true);
42
43 return drmDia;
44}
45
46//______________________________________________________________________
47
49{
50 memset(mEmptyCrateProb, 0., Geo::kNCrate * 4);
51 mTRMerrorProb.clear();
52 mTRMmask.clear();
53 mNoisy.clear();
54}
55
56//______________________________________________________________________
57
58CalibTOFapi::CalibTOFapi(const std::string url)
59{
60 // setting the URL to the CCDB manager
61
62 setURL(url);
63}
64
65//______________________________________________________________________
67{
68 auto& mgr = CcdbManager::instance();
69 long timems = long(mTimeStamp) * 1000;
70 LOG(debug) << "TOF get active map with timestamp (ms) = " << timems;
71 auto fee = mgr.getForTimeStamp<TOFFEElightInfo>("TOF/Calib/FEELIGHT", timems);
72 loadActiveMap(fee);
73}
74
75//______________________________________________________________________
77{
78 // getting the active TOF map
79 memset(mIsOffCh, false, Geo::NCHANNELS);
80
81 if (fee) {
82 LOG(info) << "read Active map (TOFFEELIGHT) for TOF ";
83 for (int ich = 0; ich < TOFFEElightInfo::NCHANNELS; ich++) {
84 //printf("%d) Enabled= %d\n",ich,fee->mChannelEnabled[ich]);
85 if (!fee->getChannelEnabled(ich)) {
86 mIsOffCh[ich] = true;
87 }
88 }
89 } else {
90 LOG(info) << "Active map (TOFFEELIGHT) not available in ccdb";
91 }
92}
93
94//______________________________________________________________________
95
97{
98
99 // getting the LHCphase calibration
100
101 auto& mgr = CcdbManager::instance();
102 long timems = long(mTimeStamp) * 1000;
103 LOG(info) << "TOF get LHCphase with timestamp (ms) = " << timems;
104 mLHCphase = mgr.getForTimeStamp<LhcPhase>("TOF/Calib/LHCphase", timems);
105 if (mLHCphase) {
106 LOG(info) << "read LHCphase for TOF " << mLHCphase->getLHCphase(mTimeStamp);
107 } else {
108 LOG(info) << "LHC phase not available in ccdb";
109 }
110}
111
112//______________________________________________________________________
113
115{
116
117 // getting the TimeSlewing calibration
118 // it includes also offset and information on problematic
119
120 auto& mgr = CcdbManager::instance();
121 long timems = long(mTimeStamp) * 1000;
122 LOG(info) << "TOF get time calibrations with timestamp (ms) = " << timems;
123 mSlewParam = mgr.getForTimeStamp<SlewParam>("TOF/Calib/ChannelCalib", timems);
124 if (mSlewParam) {
125 LOG(info) << "read TimeSlewingParam for TOF";
126 } else {
127 LOG(info) << "TimeSlewingParam for TOF not available in ccdb";
128 }
129}
130
131//______________________________________________________________________
133{
134 TFile* f = TFile::Open(filename);
135 if (f) {
136 mSlewParam = (SlewParam*)f->Get("ccdb_object");
137 } else {
138 LOG(info) << "File " << filename << " not found";
139 }
140}
141
142//______________________________________________________________________
143
145{
146 auto& mgr = CcdbManager::instance();
147 long timems = long(mTimeStamp) * 1000;
148 LOG(info) << "TOF get TRM Diagnostics with timestamp (ms) = " << timems;
149 mDiaFreq = mgr.getForTimeStamp<Diagnostic>("TOF/Calib/Diagnostic", timems);
150
152}
153
154//______________________________________________________________________
155
157{
158 auto& mgr = CcdbManager::instance();
159 long timems = long(mTimeStamp) * 1000;
160 LOG(info) << "TOF get DRM Diagnostics with timestamp (ms) = " << timems;
161 mDiaFreq = mgr.getForTimeStamp<Diagnostic>("TOF/Calib/TRMerrors", timems);
162
164}
165//______________________________________________________________________
166
168{
169 mDiaFreq->print();
170
171 static const int NCH_PER_CRATE = Geo::NSTRIPXSECTOR * Geo::NPADS;
172 // getting the Diagnostic Frequency calibration
173 // needed for simulation
174
175 memset(mIsNoisy, false, Geo::NCHANNELS);
176
177 resetDia();
178
179 if (!mDiaFreq->getFrequencyROW()) {
180 return;
181 }
182
183 float nrow = (float)mDiaFreq->getFrequencyROW();
184 mEmptyTOF = mDiaFreq->getFrequencyEmptyTOF() / nrow;
185
186 nrow -= mDiaFreq->getFrequencyEmptyTOF();
187
188 if (nrow < 1) {
189 return;
190 }
191
192 // fill empty crates
193 int ncrate[Geo::kNCrate];
194 for (int i = 0; i < Geo::kNCrate; i++) {
195 ncrate[i] = mDiaFreq->getFrequencyEmptyCrate(i) - mDiaFreq->getFrequencyEmptyTOF(); // counts of empty crate for non-empty event
196 mEmptyCrateProb[i] = ncrate[i] / nrow;
197 }
198
199 const auto vectorDia = mDiaFreq->getVector();
200 // fill TRM errors and noisy
201 for (auto pair : vectorDia) {
202 auto key = pair.first;
203 int slot = mDiaFreq->getSlot(key);
204
205 if (slot < 13 && slot > 2) { // is TRM
206 int icrate = mDiaFreq->getCrate(key);
207 int crateslot = icrate * 100 + slot;
208 mTRMerrorProb.push_back(std::make_pair(crateslot, pair.second / (nrow - ncrate[icrate])));
209 mTRMmask.push_back(key - mDiaFreq->getTRMKey(icrate, slot)); // remove crate and slot from the key (28 bit errors remaining)
210 continue;
211 }
212
213 int channel = mDiaFreq->getChannel(key);
214 if (channel > -1) { // noisy
215 if (!mDiaFreq->isNoisyChannel(channel, mNoisyThreshold)) {
216 continue;
217 }
218
219 int crate = channel / NCH_PER_CRATE;
220 float prob = pair.second / (nrow - ncrate[crate]);
221 mNoisy.push_back(std::make_pair(channel, prob));
222 continue;
223 }
224 }
225
226 std::sort(mTRMerrorProb.begin(), mTRMerrorProb.end(), [](const auto& a, const auto& b) {
227 return a.first < b.first;
228 });
229
230 std::sort(mNoisy.begin(), mNoisy.end(), [](const auto& a, const auto& b) {
231 return a.first < b.first;
232 });
233
234 int ich = -1;
235 float prob = 0;
236 for (auto [ch, p] : mNoisy) {
237 if (ch != ich) { // new channel
238 if (ich != -1 && prob > 0.5) {
239 mIsNoisy[ich] = true;
240 }
241 ich = ch;
242 prob = p;
243 } else {
244 prob += p;
245 }
246 }
247 if (ich != -1 && prob > 0.5) {
248 mIsNoisy[ich] = true;
249 }
250}
251
252//______________________________________________________________________
253
255{
256 mDiaDRMFreq->print();
257
258 for (int ic = 0; ic < Geo::kNCrate; ic++) { // loop over crates
259 float DRMcounters = mDiaDRMFreq->getFrequencyDRM(ic);
260
261 if (DRMcounters < 1) {
262 for (int ie = 0; ie < N_DRM_ERRORS; ie++) {
263 mErrorInDRM[ic][ie] = 0.;
264 }
265 continue;
266 }
267
268 for (int ie = 0; ie < N_DRM_ERRORS; ie++) { // loop over error types
269 int ie_shifted = ie + DRM_ERRINDEX_SHIFT;
270
271 float frequency = mDiaDRMFreq->getFrequencyDRMerror(ic, ie_shifted) * 1. / DRMcounters; // error frequency
272 if (frequency > 1) {
273 frequency = 1.;
274 }
275 if (frequency > 1E-6) {
276 LOG(debug) << "DRMmap: Crate = " << ic << " - error = " << ie << " - frequency = " << frequency;
277 }
278 mErrorInDRM[ic][ie] = frequency;
279 }
280 }
281}
282
283//______________________________________________________________________
284
285void CalibTOFapi::writeLHCphase(LhcPhase* phase, std::map<std::string, std::string> metadataLHCphase, uint64_t minTimeStamp, uint64_t maxTimeStamp)
286{
287
288 // write LHCphase object to CCDB
289
290 auto& mgr = CcdbManager::instance();
291 CcdbApi api;
292 api.init(mgr.getURL());
293 api.storeAsTFileAny(phase, "TOF/Calib/LHCphase", metadataLHCphase, minTimeStamp, maxTimeStamp);
294}
295
296//______________________________________________________________________
297
298void CalibTOFapi::writeTimeSlewingParam(SlewParam* param, std::map<std::string, std::string> metadataChannelCalib, uint64_t minTimeStamp, uint64_t maxTimeStamp)
299{
300
301 // write TiemSlewing object to CCDB (it includes offset + problematic)
302
303 auto& mgr = CcdbManager::instance();
304 CcdbApi api;
305 api.init(mgr.getURL());
306 if (maxTimeStamp == 0) {
307 api.storeAsTFileAny(param, "TOF/Calib/ChannelCalib", metadataChannelCalib, minTimeStamp);
308 } else {
309 api.storeAsTFileAny(param, "TOF/Calib/ChannelCalib", metadataChannelCalib, minTimeStamp, maxTimeStamp);
310 }
311}
312
313//______________________________________________________________________
314
316{
317
318 // method to know if the channel was problematic or not
319
320 return (mSlewParam->getFractionUnderPeak(ich) < 0.5 || mSlewParam->getSigmaPeak(ich) > 1000);
321 // return mSlewParam->isProblematic(ich);
322}
323
324//______________________________________________________________________
325
327{
328 return mIsNoisy[ich];
329}
330
331//______________________________________________________________________
332
334{
335 return mIsOffCh[ich];
336}
337
338//______________________________________________________________________
339
340float CalibTOFapi::getTimeCalibration(int ich, float tot) const
341{
342
343 // time calibration to correct measured TOF times
344
345 float corr = 0;
346 if (!mLHCphase || !mSlewParam) {
347 LOG(warning) << "Either LHC phase or slewing object null: mLHCphase = " << mLHCphase << ", mSlewParam = " << mSlewParam;
348 return corr;
349 }
350 // printf("LHC phase apply\n");
351 // LHCphase
352 corr += mLHCphase->getLHCphase(mTimeStamp); // timestamp that we use in LHCPhase is in seconds
353 // time slewing + channel offset
354 //printf("eval time sleweing calibration: ch=%d tot=%f (lhc phase = %f)\n",ich,tot,corr);
355 corr += mSlewParam->evalTimeSlewing(ich, tot);
356 //printf("corr = %f\n",corr);
357 return corr;
358}
359
360//______________________________________________________________________
361
362float CalibTOFapi::getTimeCalibration(int ich, float tot, float phase) const
363{
364
365 // time calibration to correct measured TOF times
366
367 float corr = 0;
368 if (!mSlewParam) {
369 LOG(warning) << "slewing object null: mSlewParam = " << mSlewParam;
370 return corr;
371 }
372 // printf("LHC phase apply\n");
373 // LHCphase
374 corr += phase; // timestamp that we use in LHCPhase is in seconds
375 // time slewing + channel offset
376 //printf("eval time sleweing calibration: ch=%d tot=%f (lhc phase = %f)\n",ich,tot,corr);
377 corr += mSlewParam->evalTimeSlewing(ich, tot);
378 //printf("corr = %f\n",corr);
379 return corr;
380}
381
382//______________________________________________________________________
383
384float CalibTOFapi::getTimeDecalibration(int ich, float tot) const
385{
386
387 // time decalibration for simulation (it is just the opposite of the calibration)
388
389 return -getTimeCalibration(ich, tot);
390}
391
392//______________________________________________________________________
393
395{
396 for (auto index : mFillErrChannel) {
397 mIsErrorCh[index] = false;
398 }
399
400 mFillErrChannel.clear();
401}
402
403//______________________________________________________________________
404
406{
407 for (auto index : mFillErrDRMChannel) {
408 mIsErrorDRMCh[index] = false;
409 }
410
411 mFillErrDRMChannel.clear();
412}
413
414//______________________________________________________________________
415
416void CalibTOFapi::processError(int crate, int trm, int mask)
417{
418 if (checkTRMPolicy(mask)) { // check the policy of TRM -> true=good TRM
419 return;
420 }
421 int ech = (crate << 12) + ((trm - 3) << 8);
422 for (int i = ech; i < ech + 256; i++) {
423 int channel = Geo::getCHFromECH(i);
424 if (channel == -1) {
425 continue;
426 }
427 mIsErrorCh[channel] = true;
428 mFillErrChannel.push_back(channel);
429 }
430}
431
432//______________________________________________________________________
433
434void CalibTOFapi::processErrorDRM(int crate, int codeErr)
435{
436 int mask = 1 << codeErr;
437
438 if (checkDRMPolicy(mask)) {
439 return;
440 }
441
442 LOG(debug) << "DRMmask: crate = " << crate << " - mask = " << mask << " - critical mask = " << mDRMCriticalErrorMask;
443
444 for (int trm = 3; trm < 13; trm++) {
445 int ech = (crate << 12) + ((trm - 3) << 8);
446 for (int i = ech; i < ech + 256; i++) {
447 int channel = Geo::getCHFromECH(i);
448 if (channel == -1 || mIsErrorDRMCh[channel] == true) {
449 continue;
450 }
451
452 mIsErrorDRMCh[channel] = true;
453 mFillErrDRMChannel.push_back(channel);
454 }
455 }
456}
457
458//______________________________________________________________________
459
461{
462 return false;
463}
464
465//______________________________________________________________________
466
468{
469 return !(mDRMCriticalErrorMask & mask);
470}
471
472//______________________________________________________________________
473
474bool CalibTOFapi::isChannelError(int channel) const
475{
476 return mIsErrorCh[channel];
477}
478
479//______________________________________________________________________
480
481bool CalibTOFapi::isChannelDRMError(int channel) const
482{
483 if (mIsErrorDRMCh[channel]) {
484 int detId[5];
485 o2::tof::Geo::getVolumeIndices(channel, detId);
486 }
487 return mIsErrorDRMCh[channel];
488}
std::string url
ClassImp(o2::tof::CalibTOFapi)
Class to use TOF calibration (decalibration, calibration)
std::ostringstream debug
uint64_t phase
Definition RawEventData.h:7
int32_t i
uint32_t j
Definition RawData.h:0
StringRef key
static BasicCCDBManager & instance()
int storeAsTFileAny(const T *obj, std::string const &path, std::map< std::string, std::string > const &metadata, long startValidityTimestamp=-1, long endValidityTimestamp=-1, std::vector< char >::size_type maxSize=0) const
Definition CcdbApi.h:157
void init(std::string const &hosts)
Definition CcdbApi.cxx:160
float getLHCphase(int timestamp) const
float evalTimeSlewing(int channel, float tot) const
float getFractionUnderPeak(int sector, int channel) const
float getSigmaPeak(int sector, int channel) const
void readDiagnosticDRMFrequencies()
bool isChannelError(int channel) const
void writeLHCphase(LhcPhase *phase, std::map< std::string, std::string > metadataLHCphase, uint64_t minTimeSTamp, uint64_t maxTimeStamp)
void readTimeSlewingParamFromFile(const char *filename)
void writeTimeSlewingParam(SlewParam *param, std::map< std::string, std::string > metadataChannelCalib, uint64_t minTimeSTamp, uint64_t maxTimeStamp=0)
void loadDiagnosticDRMFrequencies()
bool checkTRMPolicy(int mask) const
static const int DRM_ERRINDEX_SHIFT
static o2::tof::Diagnostic doDRMerrCalibFromQCHisto(const TH2F *histo, const char *file_output_name)
bool isNoisy(int ich)
void processError(int crate, int trm, int mask)
void loadActiveMap(TOFFEElightInfo *fee)
bool isChannelDRMError(int channel) const
void setURL(const std::string url)
Definition CalibTOFapi.h:70
float getTimeDecalibration(int ich, float tot) const
bool checkDRMPolicy(int mask) const
bool isProblematic(int ich)
float getTimeCalibration(int ich, float tot) const
void processErrorDRM(int crate, int codeErr)
Diagnostic class for TOF.
Definition Diagnostic.h:32
uint32_t getFrequencyROW() const
Definition Diagnostic.h:38
void print(bool longFormat=false) const
uint32_t getFrequencyEmptyCrate(int crate) const
Definition Diagnostic.h:39
uint32_t fillDRMerror(int crate, int error, uint32_t frequency)
Definition Diagnostic.h:52
static int getCrate(ULong64_t pattern)
Definition Diagnostic.h:73
uint32_t getFrequencyDRM(int crate) const
Definition Diagnostic.h:49
static ULong64_t getTRMKey(int crate, int trm)
static int getChannel(ULong64_t pattern)
Definition Diagnostic.h:74
uint32_t getFrequencyDRMerror(int crate, int error) const
Definition Diagnostic.h:50
uint32_t getFrequencyEmptyTOF() const
Definition Diagnostic.h:40
const std::map< ULong64_t, uint32_t > & getVector() const
Definition Diagnostic.h:95
bool isNoisyChannel(int channel, int thr=0) const
static int getSlot(ULong64_t pattern)
Definition Diagnostic.h:72
uint32_t fillDRM(int crate, uint32_t frequency)
Definition Diagnostic.h:51
@ kNCrate
Definition Geo.h:97
static constexpr Int_t NSTRIPXSECTOR
Definition Geo.h:118
static constexpr Int_t NPADS
Definition Geo.h:112
static Int_t getCHFromECH(int echan)
Definition Geo.h:354
static void getVolumeIndices(Int_t index, Int_t *detId)
Definition Geo.cxx:543
static constexpr int NCHANNELS
Definition Geo.h:126
GLuint index
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLenum GLfloat param
Definition glcorearb.h:271
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLint GLuint mask
Definition glcorearb.h:291
std::string filename()
bool getChannelEnabled(int idx) const
static constexpr int NCHANNELS
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"