Project
Loading...
Searching...
No Matches
CalibdEdx.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
15
16#ifndef ALICEO2_TPC_CALIBDEDX_H_
17#define ALICEO2_TPC_CALIBDEDX_H_
18
19#include <cstddef>
20#include <gsl/span>
21#include <string_view>
22
23// o2 includes
27#include "DataFormatsTPC/Defs.h"
31
32// boost includes
33#include <boost/histogram.hpp>
34#include "THn.h"
35
36class TLinearFitter;
37
38namespace o2::tpc
39{
40
43{
44 public:
54
55 // Interger histogram axis identifying the GEM stacks, without under and overflow bins.
56 using IntAxis = boost::histogram::axis::integer<int, boost::histogram::use_default, boost::histogram::axis::option::none_t>;
57 // Float axis to store data, without under and overflow bins.
58 using FloatAxis = boost::histogram::axis::regular<float, boost::histogram::use_default, boost::histogram::use_default, boost::histogram::axis::option::none_t>;
59
61
62 // Histogram axes types
63 using AxesType = std::tuple<
64 FloatAxis, // dEdx
65 FloatAxis, // Tgl
66 FloatAxis, // Snp
67 IntAxis, // sector
68 IntAxis, // stackType
69 IntAxis // Charge
70 >;
71
72 using Hist = boost::histogram::histogram<AxesType>;
73
76
79 CalibdEdx(int dEdxBins = 70, float mindEdx = 10, float maxdEdx = 90, int angularBins = 36, bool fitSnp = false);
80
81 void setCuts(const TrackCuts& cuts) { mCuts = cuts; }
82 void setApplyCuts(bool apply) { mApplyCuts = apply; }
83 bool getApplyCuts() const { return mApplyCuts; }
84
86 void setSectorFitThreshold(int minEntries) { mSectorThreshold = minEntries; }
88 void set1DFitThreshold(int minEntries) { m1DThreshold = minEntries; }
91 void set2DFitThreshold(int minEntries) { m2DThreshold = minEntries; }
92
93 int getSectorFitThreshold() const { return mSectorThreshold; }
94 int getTglFitThreshold() const { return m1DThreshold; }
95 int getSnpFitThreshold() const { return m2DThreshold; }
96
102 void setElectronCut(float cut, int passes = 3, float lowCutFactor = 1.5)
103 {
104 mFitCut = cut;
105 mFitPasses = passes;
106 mFitLowCutFactor = lowCutFactor;
107 }
108
110 void setMaterialType(o2::base::Propagator::MatCorrType materialType) { mMatType = materialType; }
111
113 void fill(const TrackTPC& tracks);
114 void fill(const gsl::span<const TrackTPC>);
115 void fill(const std::vector<TrackTPC>& tracks) { fill(gsl::span(tracks)); }
116
117 void fill(const TFIDInfo& tfid, const gsl::span<const TrackTPC> tracks)
118 {
119 mTFID = tfid;
120 fill(tracks);
121 }
122 void fill(const TFIDInfo& tfid, const std::vector<TrackTPC>& tracks) { fill(tfid, gsl::span(tracks)); }
123
125 void merge(const CalibdEdx* other);
126
130 void finalize(const bool useGausFits = true);
131
133 const Hist& getHist() const { return mHist; }
135 THnF* getRootHist() const;
136
137 const CalibdEdxCorrection& getCalib() const { return mCalib; }
138
140 void setCalibrationInput(const CalibdEdxCorrection& calib) { mCalibIn = calib; }
141 const CalibdEdxCorrection& getCalibrationInput() const { return mCalibIn; }
142
144 int minStackEntries() const;
145
148 bool hasEnoughData(float minEntries) const;
149
151 static int findTrackSector(const TrackTPC& track, GEMstack, bool& ok);
152
154 void print() const;
155
157 void writeTTree(std::string_view fileName) const;
158
160 void enableDebugOutput(std::string_view fileName);
161
163 void disableDebugOutput();
164
166 void finalizeDebugOutput() const;
167
169 bool hasDebugOutput() const { return static_cast<bool>(mDebugOutputStreamer); }
170
171 constexpr static float MipScale = 1.0 / 50.0;
172
173 constexpr static float scaleTgl(float tgl, GEMstack rocType) { return tgl / conf_dedx_corr::TglScale[rocType]; }
174 constexpr static float recoverTgl(float scaledTgl, GEMstack rocType) { return scaledTgl * conf_dedx_corr::TglScale[rocType]; }
175
177 void dumpToFile(const char* outFile, const char* outName) const;
178
180 static CalibdEdx readFromFile(const char* inFile, const char* inName);
181
185 void setSigmaFitRange(const float lowerSigma, const float upperSigma);
186
187 private:
188 // ATTENTION: Adjust copy constructor
189 bool mFitSnp{};
190 bool mApplyCuts{true};
191 TrackCuts mCuts{0.3, 0.7, 60};
192 int mSectorThreshold = 1000;
193 int m1DThreshold = 500;
194 int m2DThreshold = 5000;
195 float mFitCut = 0.2;
196 float mFitLowCutFactor = 1.5;
197 int mFitPasses = 3;
198 TFIDInfo mTFID{};
199 float mSigmaUpper = 1.;
200 float mSigmaLower = 1.5;
201
202 Hist mHist;
203 CalibdEdxCorrection mCalib{};
204 CalibdEdxCorrection mCalibIn{};
205
207
208 std::unique_ptr<o2::utils::TreeStreamRedirector> mDebugOutputStreamer;
209
210 THnF* getTHnF() const;
211
213 void fitHistGaus(TLinearFitter& fitter, CalibdEdxCorrection& corr, const CalibdEdxCorrection* stackMean = nullptr);
214
216 void setFromRootHist(const THnF* hist);
217
218 ClassDefNV(CalibdEdx, 5);
219};
220
221} // namespace o2::tpc
222
223#endif
Class that creates dE/dx histograms from a sequence of tracks objects.
Definition CalibdEdx.h:43
const CalibdEdxCorrection & getCalibrationInput() const
Definition CalibdEdx.h:141
int getSectorFitThreshold() const
Definition CalibdEdx.h:93
void disableDebugOutput()
Disable debug output to file. Also writes and closes stored time slots.
void setElectronCut(float cut, int passes=3, float lowCutFactor=1.5)
Params used to remove electron points from the fit. The fit to find the MIP peak will be performed pa...
Definition CalibdEdx.h:102
void finalize(const bool useGausFits=true)
static constexpr float recoverTgl(float scaledTgl, GEMstack rocType)
Definition CalibdEdx.h:174
static constexpr float scaleTgl(float tgl, GEMstack rocType)
Definition CalibdEdx.h:173
void setCalibrationInput(const CalibdEdxCorrection &calib)
calibration used during reconstruction
Definition CalibdEdx.h:140
void print() const
Print statistics info.
std::tuple< FloatAxis, FloatAxis, FloatAxis, IntAxis, IntAxis, IntAxis > AxesType
Definition CalibdEdx.h:70
boost::histogram::histogram< AxesType > Hist
Definition CalibdEdx.h:72
boost::histogram::axis::regular< float, boost::histogram::use_default, boost::histogram::use_default, boost::histogram::axis::option::none_t > FloatAxis
Definition CalibdEdx.h:58
void setSectorFitThreshold(int minEntries)
Definition CalibdEdx.h:86
boost::histogram::axis::integer< int, boost::histogram::use_default, boost::histogram::axis::option::none_t > IntAxis
Definition CalibdEdx.h:56
static int findTrackSector(const TrackTPC &track, GEMstack, bool &ok)
Find the sector of a track at a given GEM stack type.
void merge(const CalibdEdx *other)
Add counts from another container.
bool hasDebugOutput() const
Definition CalibdEdx.h:169
const Hist & getHist() const
Return calib data histogram.
Definition CalibdEdx.h:133
void setMaterialType(o2::base::Propagator::MatCorrType materialType)
setting the material type for track propagation
Definition CalibdEdx.h:110
void dumpToFile(const char *outFile, const char *outName) const
dump this object to a file - the boost histogram is converted to a ROOT histogram -
THnF * getRootHist() const
Return calib data as a THn.
int minStackEntries() const
Return the number of hist entries of the gem stack with less statistics.
const CalibdEdxCorrection & getCalib() const
Definition CalibdEdx.h:137
@ Size
Number of axes.
Definition CalibdEdx.h:52
void finalizeDebugOutput() const
Write debug output to file.
void writeTTree(std::string_view fileName) const
Save the histograms to a TTree.
static CalibdEdx readFromFile(const char *inFile, const char *inName)
read the object from a file
void fill(const TFIDInfo &tfid, const gsl::span< const TrackTPC > tracks)
Definition CalibdEdx.h:117
void set2DFitThreshold(int minEntries)
Definition CalibdEdx.h:91
void fill(const TrackTPC &tracks)
Fill histograms using tracks data.
Definition CalibdEdx.cxx:84
void setSigmaFitRange(const float lowerSigma, const float upperSigma)
void setCuts(const TrackCuts &cuts)
Definition CalibdEdx.h:81
int getTglFitThreshold() const
Definition CalibdEdx.h:94
void enableDebugOutput(std::string_view fileName)
Enable debug output to file of the time slots calibrations outputs and dE/dx histograms.
void fill(const TFIDInfo &tfid, const std::vector< TrackTPC > &tracks)
Definition CalibdEdx.h:122
bool getApplyCuts() const
Definition CalibdEdx.h:83
bool hasEnoughData(float minEntries) const
Check if there are enough data to compute the calibration.
void setApplyCuts(bool apply)
Definition CalibdEdx.h:82
o2::dataformats::TFIDInfo TFIDInfo
Definition CalibdEdx.h:60
static constexpr float MipScale
Inverse of target dE/dx value for MIPs.
Definition CalibdEdx.h:171
int getSnpFitThreshold() const
Definition CalibdEdx.h:95
void fill(const std::vector< TrackTPC > &tracks)
Definition CalibdEdx.h:115
void set1DFitThreshold(int minEntries)
Definition CalibdEdx.h:88
track cut class
Definition TrackCuts.h:37
Global TPC definitions and constants.
Definition SimTraits.h:167
GEMstack
TPC GEM stack types.
Definition Defs.h:53
VectorOfTObjectPtrs other