Project
Loading...
Searching...
No Matches
TimeCalibrationParams.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
14
15#include <fairlogger/Logger.h>
16
17#include <TH1.h>
18
19#include <iostream>
20
21using namespace o2::emcal;
22
23void TimeCalibrationParams::addTimeCalibParam(unsigned short cellID, short time, bool isLowGain)
24{
25 if (cellID >= mTimeCalibParamsHG.size()) {
26 throw CalibContainerIndexException(cellID);
27 }
28 if (!isLowGain) {
29 mTimeCalibParamsHG[cellID] = time;
30 } else {
31 mTimeCalibParamsLG[cellID] = time;
32 }
33}
34
35short TimeCalibrationParams::getTimeCalibParam(unsigned short cellID, bool isLowGain) const
36{
37 if (cellID >= mTimeCalibParamsHG.size()) {
38 throw CalibContainerIndexException(cellID);
39 }
40 if (isLowGain) {
41 return mTimeCalibParamsLG[cellID];
42 } else {
43 return mTimeCalibParamsHG[cellID];
44 }
45}
46
48{
49
50 if (!isLowGain) {
51 auto hist = new TH1S("TimeCalibrationParams", "Time Calibration Params HG", 17664, 0, 17664);
52 hist->SetDirectory(nullptr);
53
54 for (std::size_t icell{0}; icell < mTimeCalibParamsHG.size(); ++icell) {
55 hist->SetBinContent(icell + 1, mTimeCalibParamsHG[icell]);
56 }
57
58 return hist;
59 } else {
60 auto hist = new TH1S("TimeCalibrationParams", "Time Calibration Params LG", 17664, 0, 17664);
61 hist->SetDirectory(nullptr);
62
63 for (std::size_t icell{0}; icell < mTimeCalibParamsLG.size(); ++icell) {
64 hist->SetBinContent(icell + 1, mTimeCalibParamsLG[icell]);
65 }
66
67 return hist;
68 }
69}
70
72{
73 return mTimeCalibParamsHG == other.mTimeCalibParamsHG && mTimeCalibParamsLG == other.mTimeCalibParamsLG;
74}
int16_t time
Definition RawEventData.h:4
Error handling for invalid index in calibration request.
short getTimeCalibParam(unsigned short cellID, bool isLowGain) const
Get the time calibration coefficient for a certain cell.
bool operator==(const TimeCalibrationParams &other) const
Comparison of two time calibration coefficients.
TH1 * getHistogramRepresentation(bool isLowGain) const
Convert the time calibration coefficient array to a histogram.
VectorOfTObjectPtrs other