Project
Loading...
Searching...
No Matches
Clusters.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 <string>
13
14// root includes
15#include "TH1.h"
16#include "TH2.h"
17#include "TFile.h"
18
19// o2 includes
20#include "TPCQC/Clusters.h"
21#include "TPCBase/Painter.h"
22#include "TPCBase/ROC.h"
23#include "TPCBase/CRU.h"
24#include "TPCBase/Mapper.h"
29
31
32using namespace o2::tpc::qc;
33
34//______________________________________________________________________________
35template <class T>
36bool Clusters::processCluster(const T& cluster, const o2::tpc::Sector sector, const int row)
37{
38 if (mIsNormalized) {
40 LOGP(warning, "calling denormalize() before filling");
41 }
42
43 const auto& mapper = Mapper::instance();
44
45 const int nROC = row < 63 ? int(sector) : int(sector) + 36;
46 const int rocRow = row < 63 ? row : row - 63;
47
48 const float pad = cluster.getPad();
49
50 size_t position = mapper.getPadNumber(mNClusters.getCalArray(nROC).getPadSubset(), mNClusters.getCalArray(nROC).getPadSubsetNumber(), rocRow, pad);
51
52 if ((nROC < mapper.getNumberOfIROCs())) { // IROC
53 if (position >= mapper.getPadsInIROC()) {
54 LOG(error) << "IROC out of range: pad position " << position << "\tROC " << nROC << "\tIROC has " << mapper.getPadsInIROC() << " pads";
55 return false;
56 }
57 } else { // OROC
58 if (position >= mapper.getPadsInOROC()) {
59 LOG(error) << "OROC out of range: pad position " << position << "\tROC " << nROC << "\tOROC has " << mapper.getPadsInOROC() << " pads";
60 return false;
61 }
62 }
63
64 const auto qMax = cluster.getQmax();
65 const auto qTot = cluster.getQtot();
66 const auto sigmaPad = cluster.getSigmaPad();
67 const auto sigmaTime = cluster.getSigmaTime();
68 const auto timeBin = cluster.getTime();
69
70 float count = mNClusters.getCalArray(nROC).getValue(rocRow, pad);
71 mNClusters.getCalArray(nROC).setValue(rocRow, pad, count + 1);
72
73 float charge = mQMax.getCalArray(nROC).getValue(rocRow, pad);
74 mQMax.getCalArray(nROC).setValue(rocRow, pad, charge + qMax);
75
76 charge = mQTot.getCalArray(nROC).getValue(rocRow, pad);
77 mQTot.getCalArray(nROC).setValue(rocRow, pad, charge + qTot);
78
79 count = mSigmaTime.getCalArray(nROC).getValue(rocRow, pad);
80 mSigmaTime.getCalArray(nROC).setValue(rocRow, pad, count + sigmaTime);
81
82 count = mSigmaPad.getCalArray(nROC).getValue(rocRow, pad);
83 mSigmaPad.getCalArray(nROC).setValue(rocRow, pad, count + sigmaPad);
84
85 count = mTimeBin.getCalArray(nROC).getValue(rocRow, pad);
86 mTimeBin.getCalArray(nROC).setValue(rocRow, pad, count + timeBin);
87
88 return true;
89}
90
91//______________________________________________________________________________
92void Clusters::fillADCValue(int cru, int rowInSector, int padInRow, int timeBin, float adcValue)
93{
94 if (mIsNormalized) {
96 LOGP(warning, "calling denormalize() before filling");
97 }
98
99 const CRU cruID(cru);
100 float val;
101 val = mNClusters.getValue(cruID.sector(), rowInSector, padInRow);
102 mNClusters.setValue(cruID.sector(), rowInSector, padInRow, val + 1);
103
104 val = mQMax.getValue(cruID.sector(), rowInSector, padInRow);
105 mQMax.setValue(cruID.sector(), rowInSector, padInRow, val + adcValue);
106
107 val = mTimeBin.getValue(cruID.sector(), rowInSector, padInRow);
108 mTimeBin.setValue(cruID.sector(), rowInSector, padInRow, val + timeBin);
109}
110
111//______________________________________________________________________________
113{
114 if (mIsNormalized) {
115 return;
116 }
117
118 mQMax /= mNClusters;
119 mQTot /= mNClusters;
120 mSigmaTime /= mNClusters;
121 mSigmaPad /= mNClusters;
122 mTimeBin /= mNClusters;
123
124 mIsNormalized = true;
125}
126
127//______________________________________________________________________________
129{
130 if (!mIsNormalized) {
131 return;
132 }
133
134 mQMax *= mNClusters;
135 mQTot *= mNClusters;
136 mSigmaTime *= mNClusters;
137 mSigmaPad *= mNClusters;
138 mTimeBin *= mNClusters;
139
140 mIsNormalized = false;
141}
142
143//______________________________________________________________________________
145{
146 mNClusters = 0;
147 mQMax = 0;
148 mQTot = 0;
149 mSigmaTime = 0;
150 mSigmaPad = 0;
151 mTimeBin = 0;
152
153 mIsNormalized = false;
154 mProcessedTFs = 0;
155}
156
157//______________________________________________________________________________
159{
160 o2::tpc::CalPad occupancy = mNClusters;
161 occupancy /= float(mProcessedTFs * (o2::constants::lhc::LHCMaxBunches * nHBFPerTF) / float(o2::tpc::ParameterElectronics::TIMEBININBC));
162 return occupancy;
163}
164//______________________________________________________________________________
166{
167 const bool isThisNormalized = mIsNormalized;
168 const bool isOtherNormalized = clusters.mIsNormalized;
169
170 if (isThisNormalized) {
171 denormalize();
172 }
173 if (isOtherNormalized) {
174 clusters.denormalize();
175 }
176
177 mNClusters += clusters.mNClusters;
178 mQMax += clusters.mQMax;
179 mQTot += clusters.mQTot;
180 mSigmaTime += clusters.mSigmaTime;
181 mSigmaPad += clusters.mSigmaPad;
182 mTimeBin += clusters.mTimeBin;
183
184 if (isThisNormalized) {
185 normalize();
186 }
187 if (isOtherNormalized) {
188 clusters.normalize();
189 }
190}
191
192//______________________________________________________________________________
193void Clusters::dumpToFile(std::string filename, int type)
194{
195 if (filename.find(".root") != std::string::npos) {
196 filename.resize(filename.size() - 5);
197 }
198
199 if (type == 0) {
200 const std::string canvasFile = filename + "_canvas.root";
201 auto f = std::unique_ptr<TFile>(TFile::Open(canvasFile.c_str(), "recreate"));
202 f->WriteObject(o2::tpc::painter::draw(mNClusters), mNClusters.getName().data());
203 f->WriteObject(o2::tpc::painter::draw(mQMax), mQMax.getName().data());
204 f->WriteObject(o2::tpc::painter::draw(mQTot), mQTot.getName().data());
205 f->WriteObject(o2::tpc::painter::draw(mSigmaTime), mSigmaTime.getName().data());
206 f->WriteObject(o2::tpc::painter::draw(mSigmaPad), mSigmaPad.getName().data());
207 f->WriteObject(o2::tpc::painter::draw(mTimeBin), mTimeBin.getName().data());
208 f->Close();
209 }
210
211 if (type == 0 || type == 1) {
212 const std::string calPadFile = filename + ".root";
213 auto f = std::unique_ptr<TFile>(TFile::Open(calPadFile.c_str(), "recreate"));
214 TNamed nTFs("processedTFs", std::to_string(mProcessedTFs).data());
215 f->WriteObject(&mNClusters, mNClusters.getName().data());
216 f->WriteObject(&mQMax, mQMax.getName().data());
217 f->WriteObject(&mQTot, mQTot.getName().data());
218 f->WriteObject(&mSigmaTime, mSigmaTime.getName().data());
219 f->WriteObject(&mSigmaPad, mSigmaPad.getName().data());
220 f->WriteObject(&mTimeBin, mTimeBin.getName().data());
221 nTFs.Write();
222 f->Close();
223 }
224
225 if (type == 2) {
226 const std::string calPadFile = filename + ".root";
227 auto f = std::unique_ptr<TFile>(TFile::Open(calPadFile.c_str(), "recreate"));
228 f->WriteObject(this, "ClusterQC");
229 }
230}
231
232// ===| explicit instantiations |===============================================
233template bool Clusters::processCluster<o2::tpc::ClusterNative>(const o2::tpc::ClusterNative&, const o2::tpc::Sector, const int);
234template bool Clusters::processCluster<o2::tpc::KrCluster>(const o2::tpc::KrCluster&, const o2::tpc::Sector, const int);
Class of a TPC cluster in TPC-native coordinates (row, time)
ClassImp(o2::tpc::qc::Clusters)
int16_t charge
Definition RawEventData.h:5
Header to collect LHC related constants.
Definition of the parameter class for the detector electronics.
Struct for Krypton and X-ray clusters.
const Sector sector() const
Definition CRU.h:65
void setValue(const size_t channel, const T &value)
Definition CalArray.h:95
PadSubset getPadSubset() const
Definition CalArray.h:89
int getPadSubsetNumber() const
Definition CalArray.h:93
const T getValue(const size_t channel) const
Definition CalArray.h:96
const CalType & getCalArray(size_t position) const
Definition CalDet.h:63
void setValue(const int sec, const int globalPadInSector, const T &value)
Definition CalDet.h:257
const T getValue(const int sec, const int globalPadInSector) const
Definition CalDet.h:154
const std::string & getName() const
Definition CalDet.h:85
static Mapper & instance(const std::string mappingDir="")
Definition Mapper.h:44
CalPad getOccupancy(int nHBFPerTF=32)
Definition Clusters.cxx:158
void fillADCValue(int cru, int rowInSector, int padInRow, int timeBin, float adcValue)
Definition Clusters.cxx:92
void merge(Clusters &clusters)
Definition Clusters.cxx:165
void dumpToFile(std::string filename, int type=0)
Definition Clusters.cxx:193
GLint GLsizei count
Definition glcorearb.h:399
GLdouble f
Definition glcorearb.h:310
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLboolean * data
Definition glcorearb.h:298
GLuint GLfloat * val
Definition glcorearb.h:1582
constexpr int LHCMaxBunches
std::string to_string(gsl::span< T, Size > span)
Definition common.h:52
std::string filename()
static TCanvas * draw(const CalDet< T > &calDet, int nbins1D=300, float xMin1D=0, float xMax1D=0, TCanvas *outputCanvas=nullptr)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< Cluster > clusters
std::vector< int > row