Project
Loading...
Searching...
No Matches
Painter.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_PAINTER_H_
13#define ALICEO2_TPC_PAINTER_H_
14
19
20#include <vector>
21#include <string>
22#include <string_view>
23#include "DataFormatsTPC/Defs.h"
25
26class TH1;
27class TH2;
28class TH3F;
29class TH2Poly;
30class TCanvas;
31class TMultiGraph;
32class TTree;
33
34namespace o2::tpc
35{
36
37template <class T>
38class CalDet;
39
40template <class T>
41class CalArray;
42
50
51struct painter {
52
53 enum class Type : int {
54 Pad,
55 Stack,
56 FEC
57 };
58
59 static std::array<int, 6> colors;
60 static std::array<int, 10> markers;
61
64 std::vector<double> xVals = std::vector<double>(4);
65 std::vector<double> yVals = std::vector<double>(4);
66
67 void rotate(float angDeg)
68 {
69 const auto ang = 0.017453292519943295 * angDeg;
70 const auto cs = std::cos(ang);
71 const auto sn = std::sin(ang);
72 for (size_t i = 0; i < xVals.size(); ++i) {
73 const auto x = xVals[i] * cs - yVals[i] * sn;
74 const auto y = xVals[i] * sn + yVals[i] * cs;
75 xVals[i] = x;
76 yVals[i] = y;
77 }
78 }
79 };
80
82 static std::vector<PadCoordinates> getPadCoordinatesSector();
83
85 static std::vector<PadCoordinates> getStackCoordinatesSector();
86
88 static std::vector<PadCoordinates> getFECCoordinatesSector();
89
91 static std::vector<o2::tpc::painter::PadCoordinates> getCoordinates(const Type type);
92
95 static std::vector<double> getRowBinningCM(uint32_t roc = 72);
96
98 static std::string getROCTitle(const int rocNumber);
99
100 // using T=float;
104 template <class T>
105 static TCanvas* draw(const CalDet<T>& calDet, int nbins1D = 300, float xMin1D = 0, float xMax1D = 0, TCanvas* outputCanvas = nullptr);
106
110 template <class T>
111 static TCanvas* draw(const CalArray<T>& calArray);
112
117 template <class T>
118 static void fillHistogram2D(TH2& h2D, const CalDet<T>& calDet, Side side);
119
123 template <class T>
124 static void fillHistogram2D(TH2& h2D, const CalArray<T>& calArray);
125
130 template <class T>
131 static TH2* getHistogram2D(const CalDet<T>& calDet, Side side);
132
137 template <class T>
138 static TH2* getHistogram2D(const CalArray<T>& calArray);
139
146 static TH2Poly* makeSectorHist(const std::string_view name = "hSector", const std::string_view title = "Sector;local #it{x} (cm);local #it{y} (cm)", const float xMin = 83.65f, const float xMax = 247.7f, const float yMin = -43.7f, const float yMax = 43.7f, const Type type = Type::Pad);
147
150 static TH2Poly* makeSideHist(Side side, const Type type = Type::Pad);
151
156 template <class T>
157 static void fillPoly2D(TH2Poly& h2D, const CalDet<T>& calDet, Side side);
158
170 template <class T>
171 static std::vector<TCanvas*> makeSummaryCanvases(const CalDet<T>& calDet, int nbins1D = 300, float xMin1D = 0, float xMax1D = 0, bool onlyFilled = true, std::vector<TCanvas*>* outputCanvases = nullptr);
172
185 static std::vector<TCanvas*> makeSummaryCanvases(const std::string_view fileName, const std::string_view calPadNames, int nbins1D = 300, float xMin1D = 0, float xMax1D = 0, bool onlyFilled = true);
186
199 static std::vector<TCanvas*> makeSummaryCanvases(TTree& tree, const std::string_view draw, std::string_view cut = "1", int nbins1D = 300, float xMin1D = 0, float xMax1D = 0);
200
202 static void drawSectorsXY(Side side, int sectorLineColor = 920, int sectorTextColor = 1);
203
207 static void drawSectorLocalPadNumberPoly(short padTextColor = kBlack, float lineScalePS = 1);
208
212 static void drawSectorInformationPoly(short regionLineColor = kRed, short rowTextColor = kRed);
213
222 template <typename DataT>
223 static TH3F convertCalDetToTH3(const std::vector<CalDet<DataT>>& calDet, const bool norm = true, const int nRBins = 150, const float rMin = 83.5, const float rMax = 254.5, const int nPhiBins = 720, const float zMax = 1);
224
226 static std::vector<TCanvas*> makeSummaryCanvases(const LtrCalibData& ltr, std::vector<TCanvas*>* outputCanvases = nullptr);
227
229 static TCanvas* makeJunkDetectionCanvas(const TObjArray* data, TCanvas* outputCanvas = nullptr);
230
235 static void adjustPalette(TH1* h, float x2ndc, float tickLength = 0.015);
236
252 static TMultiGraph* makeMultiGraph(TTree& tree, std::string_view varX, std::string_view varsY, std::string_view errVarsY = "", std::string_view cut = "", bool makeSparse = true);
253
254}; // struct painter
255
256} // namespace o2::tpc
257
258#endif // ALICEO2_TPC_PAINTER_H_
int32_t i
float float float & zMax
float & yMax
calibration data from laser track calibration
uint32_t roc
Definition RawData.h:3
uint32_t side
Definition RawData.h:0
Class for time synchronization of RawReader instances.
GLint GLenum GLint x
Definition glcorearb.h:403
GLuint const GLchar * name
Definition glcorearb.h:781
GLint y
Definition glcorearb.h:270
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLboolean * data
Definition glcorearb.h:298
Global TPC definitions and constants.
Definition SimTraits.h:167
Side
TPC readout sidE.
Definition Defs.h:35
Drawing helper functions.
pad corner coordinates
Definition Painter.h:63
void rotate(float angDeg)
Definition Painter.h:67
std::vector< double > yVals
Definition Painter.h:65
std::vector< double > xVals
Definition Painter.h:64
static std::vector< double > getRowBinningCM(uint32_t roc=72)
static std::array< int, 6 > colors
Definition Painter.h:59
static TH3F convertCalDetToTH3(const std::vector< CalDet< DataT > > &calDet, const bool norm=true, const int nRBins=150, const float rMin=83.5, const float rMax=254.5, const int nPhiBins=720, const float zMax=1)
static std::vector< TCanvas * > makeSummaryCanvases(const LtrCalibData &ltr, std::vector< TCanvas * > *outputCanvases=nullptr)
make summary canvases for laser calibration data
static std::vector< PadCoordinates > getStackCoordinatesSector()
create a vector of stack corner coordinate for one full sector
static TH2Poly * makeSectorHist(const std::string_view name="hSector", const std::string_view title="Sector;local #it{x} (cm);local #it{y} (cm)", const float xMin=83.65f, const float xMax=247.7f, const float yMin=-43.7f, const float yMax=43.7f, const Type type=Type::Pad)
static std::array< int, 10 > markers
Definition Painter.h:60
static void fillHistogram2D(TH2 &h2D, const CalDet< T > &calDet, Side side)
static void fillHistogram2D(TH2 &h2D, const CalArray< T > &calArray)
static std::vector< o2::tpc::painter::PadCoordinates > getCoordinates(const Type type)
static TMultiGraph * makeMultiGraph(TTree &tree, std::string_view varX, std::string_view varsY, std::string_view errVarsY="", std::string_view cut="", bool makeSparse=true)
static std::string getROCTitle(const int rocNumber)
ROC title from ROC number.
static std::vector< PadCoordinates > getFECCoordinatesSector()
create a vector of FEC corner coordinates for one full sector
static void drawSectorLocalPadNumberPoly(short padTextColor=kBlack, float lineScalePS=1)
static TCanvas * draw(const CalArray< T > &calArray)
static void adjustPalette(TH1 *h, float x2ndc, float tickLength=0.015)
static void fillPoly2D(TH2Poly &h2D, const CalDet< T > &calDet, Side side)
@ FEC
drawing of FECs
@ Stack
drawing stacks
static TH2 * getHistogram2D(const CalDet< T > &calDet, Side side)
static TH2Poly * makeSideHist(Side side, const Type type=Type::Pad)
static TH2 * getHistogram2D(const CalArray< T > &calArray)
static void drawSectorInformationPoly(short regionLineColor=kRed, short rowTextColor=kRed)
static void drawSectorsXY(Side side, int sectorLineColor=920, int sectorTextColor=1)
draw sector boundaris, side name and sector numbers
static std::vector< PadCoordinates > getPadCoordinatesSector()
create a vector of pad corner coordinate for one full sector
static std::vector< TCanvas * > makeSummaryCanvases(const std::string_view fileName, const std::string_view calPadNames, int nbins1D=300, float xMin1D=0, float xMax1D=0, bool onlyFilled=true)
static TCanvas * makeJunkDetectionCanvas(const TObjArray *data, TCanvas *outputCanvas=nullptr)
make a canvas for junk detection data
static TCanvas * draw(const CalDet< T > &calDet, int nbins1D=300, float xMin1D=0, float xMax1D=0, TCanvas *outputCanvas=nullptr)
static std::vector< TCanvas * > makeSummaryCanvases(const CalDet< T > &calDet, int nbins1D=300, float xMin1D=0, float xMax1D=0, bool onlyFilled=true, std::vector< TCanvas * > *outputCanvases=nullptr)
std::unique_ptr< TTree > tree((TTree *) flIn.Get(std::string(o2::base::NameConf::CTFTREENAME).c_str()))