Project
Loading...
Searching...
No Matches
CalibPedestal.h
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#ifndef ALICEO2_TPC_CALIBPEDESTAL_H_
13#define ALICEO2_TPC_CALIBPEDESTAL_H_
14
17
18#include <string>
19#include <vector>
20#include <memory>
21#include <unordered_map>
22
23#include "Rtypes.h"
24
25#include "DataFormatsTPC/Defs.h"
26#include "TPCBase/CalDet.h"
27#include "TPCBase/CRU.h"
30
31class TH2;
32
33namespace o2
34{
35namespace tpc
36{
37
44
46{
47 public:
48 using vectorType = std::vector<float>;
49
50 // enum class StatisticsType {
51 // GausFit, ///< Use Gaus fit for pedestal and noise
52 // MeanStdDev ///< Use mean and standard deviation
53 // };
54
57
59 ~CalibPedestal() override = default;
60
62 void init();
63
71 Int_t updateROC(const Int_t roc, const Int_t row, const Int_t pad,
72 const Int_t timeBin, const Float_t signal) final;
73
75 Int_t updateCRU(const CRU& cru, const Int_t row, const Int_t pad,
76 const Int_t timeBin, const Float_t signal) final { return 0; }
77
79 void resetData();
80
82 void setADCRange(int minADC, int maxADC)
83 {
84 mADCMin = minADC;
85 mADCMax = maxADC;
86 mNumberOfADCs = mADCMax - mADCMin + 1;
87 }
88
90 void setStatisticsType(StatisticsType statisticsType) { mStatisticsType = statisticsType; }
91
93 void setTimeBinRange(int first, int last)
94 {
95 mFirstTimeBin = first;
96 mLastTimeBin = last;
97 }
99 void analyse();
100
104 const CalPad& getPedestal() const { return mCalDets.at("Pedestals"); }
105
109 const CalPad& getNoise() const { return mCalDets.at("Noise"); }
110
111 const auto& getCalDets() const { return mCalDets; }
112
114 // const std::vector<const o2::tpc::CalDet<float>*> getCalDets() const { return std::vector<const o2::tpc::CalDet<float>*>{&mCalDets.at("Pedestals"), &mCalDets.at("Noise")}; }
115
117 StatisticsType getStatisticsType() const { return mStatisticsType; }
118
120 void dumpToFile(const std::string filename, uint32_t type = 0) final;
121
123 void endEvent() final{};
124
127
128 private:
129 int mFirstTimeBin;
130 int mLastTimeBin;
131 int mADCMin;
132 int mADCMax;
133 int mNumberOfADCs;
134 StatisticsType mStatisticsType;
135 std::unordered_map<std::string, CalPad> mCalDets;
136
137 std::vector<std::unique_ptr<vectorType>> mADCdata;
138
143 vectorType* getVector(ROC roc, bool create = kFALSE);
144
146 void resetEvent() final {}
147};
148
149} // namespace tpc
150
151} // namespace o2
152#endif
Implementation of the parameter class for the pedestal calibration.
uint32_t roc
Definition RawData.h:3
Pedestal calibration class.
TH2 * createControlHistogram(ROC roc)
generate a control histogram
void setTimeBinRange(int first, int last)
set the time bin range to analyse
void resetData()
Reset pedestal data.
Int_t updateROC(const Int_t roc, const Int_t row, const Int_t pad, const Int_t timeBin, const Float_t signal) final
const CalPad & getNoise() const
void init()
initialize the clusterer from CalibPedestalParam
~CalibPedestal() override=default
default destructor
void dumpToFile(const std::string filename, uint32_t type=0) final
Dump the relevant data to file.
StatisticsType getStatisticsType() const
return all pad clibrations as vector
const CalPad & getPedestal() const
void setStatisticsType(StatisticsType statisticsType)
set the statistics type
const auto & getCalDets() const
void endEvent() final
Dummy end event.
std::vector< float > vectorType
Int_t updateCRU(const CRU &cru, const Int_t row, const Int_t pad, const Int_t timeBin, const Float_t signal) final
not used
void setADCRange(int minADC, int maxADC)
set the adc range
void analyse()
Analyse the buffered adc values and calculate noise and pedestal.
Base class for raw data calibrations.
GLint first
Definition glcorearb.h:399
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
StatisticsType
Statistics type.
Definition Defs.h:94
PadSubset
Definition of the different pad subsets.
Definition Defs.h:63
@ ROC
ROCs (up to 72)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
std::string filename()
std::vector< int > row