REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDataSetGainMap.h
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see https://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see https://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
23#ifndef REST_TRestDataSetGainMap
24#define REST_TRestDataSetGainMap
25
26#include <TCanvas.h>
27#include <TFile.h>
28#include <TGraph.h>
29#include <TH1.h>
30#include <TH2.h>
31#include <TRestStringOutput.h>
32#include <TSpectrum.h>
33#include <TTree.h>
34#include <TVector2.h>
35
36#include "TRestCut.h"
37#include "TRestDataSet.h"
38#include "TRestMetadata.h"
39
42 public:
43 class Module;
44
45 private:
47 std::string fObservable = ""; //<
48
50 std::string fSpatialObservableX = ""; //<
51
53 std::string fSpatialObservableY = ""; //<
54
56 std::vector<Module> fModulesCal = {}; //<
57
59 std::string fCalibFileName = ""; //<
60
62 std::string fOutputFileName = ""; //<
63
65 TRestCut* fCut = nullptr; //<
66
67 void Initialize() override;
68 void InitFromConfigFile() override;
69
70 public:
71 std::set<int> GetPlaneIDs() const;
72 std::set<int> GetModuleIDs(const int planeId) const;
73 std::map<int, std::set<int>> GetModuleIDs() const;
74 Int_t GetNumberOfPlanes() const { return GetPlaneIDs().size(); }
75 Int_t GetNumberOfModules() const {
76 int sum = 0;
77 for (auto pID : GetPlaneIDs()) sum += GetModuleIDs(pID).size();
78 return sum;
79 }
80
81 std::string GetCalibrationFileName() const { return fCalibFileName; }
82 std::string GetOutputFileName() const { return fOutputFileName; }
83 std::string GetObservable() const { return fObservable; }
84 std::string GetSpatialObservableX() const { return fSpatialObservableX; }
85 std::string GetSpatialObservableY() const { return fSpatialObservableY; }
86 TRestCut* GetCut() const { return fCut; }
87
88 Module* GetModule(const size_t index = 0);
89 Module* GetModule(const int planeID, const int moduleID);
90 double GetSlopeParameter(const int planeID, const int moduleID, const double x, const double y);
91 double GetInterceptParameter(const int planeID, const int moduleID, const double x, const double y);
92 double GetSlopeParameterFullSpc(const int planeID, const int moduleID);
93 double GetInterceptParameterFullSpc(const int planeID, const int moduleID);
94
95 void SetCalibrationFileName(const std::string& fileName) { fCalibFileName = fileName; }
96 void SetOutputFileName(const std::string& fileName) { fOutputFileName = fileName; }
97 void SetModule(const Module& moduleCal);
98 void SetObservable(const std::string& observable) { fObservable = observable; }
99 void SetSpatialObservableX(const std::string& spatialObservableX) {
100 fSpatialObservableX = spatialObservableX;
101 }
102 void SetSpatialObservableY(const std::string& spatialObservableY) {
103 fSpatialObservableY = spatialObservableY;
104 }
105 void SetCut(TRestCut* cut) { fCut = cut; }
106
107 void Import(const std::string& fileName);
108 void Export(const std::string& fileName = "");
109
111
112 void PrintMetadata() override;
113
114 void GenerateGainMap();
115 void CalibrateDataSet(const std::string& dataSetFileName, std::string outputFileName = "",
116 std::vector<std::string> excludeColumns = {});
117
119 TRestDataSetGainMap(const char* configFilename, std::string name = "");
121
122 ClassDefOverride(TRestDataSetGainMap, 2);
123
124 class Module {
125 private:
127 const TRestDataSetGainMap* p = nullptr; //<!
128
129 std::pair<double, double> FitPeaks(TH1F* hSeg, TGraph* gr);
130 std::pair<double, double> UpdateCalibrationFits(TH1F* hSeg, TGraph* gr);
131
132 public:
135 Int_t fPlaneId = -1; //<
136
137 // Module ID (unique identifier). Although it is not linked to any TRestDetectorReadout it is
138 // recommended to use the same.
139 Int_t fModuleId = -1; //<
140
142 std::vector<double> fEnergyPeaks = {}; //<
143
145 std::vector<TVector2> fRangePeaks = {}; //<
146
148 TVector2 fCalibRange = TVector2(0, 0); //<
149
151 Int_t fNBins = 100; //<
152
154 std::string fDefinitionCut = ""; //<
155
157 Int_t fNumberOfSegmentsX = 1; //<
158
160 Int_t fNumberOfSegmentsY = 1; //<
161
163 TVector2 fReadoutRange = TVector2(-1, 246.24); //<
164
166 std::set<double> fSplitX = {}; //<
167
169 std::set<double> fSplitY = {}; //<
170
173 std::string fDataSetFileName = ""; //<
174
176 std::vector<std::vector<double>> fSlope = {}; //<
177
179 double fFullSlope = 0; //<
180
182 double fFullIntercept = 0; //<
183
185 std::vector<std::vector<double>> fIntercept = {}; //<
186
189 bool fZeroPoint = false; //<
190
192 bool fAutoRangePeaks = true; //<
193
195 std::vector<std::vector<TH1F*>> fSegSpectra = {}; //<
196
198 std::vector<std::vector<TGraph*>> fSegLinearFit = {}; //<
199
201 TH1F* fFullSpectrum = nullptr; //<
202
204 TGraph* fFullLinearFit = nullptr; //<
205
206 public:
207 void AddPeak(const double& energyPeak, const TVector2& rangePeak = TVector2(0, 0)) {
208 fEnergyPeaks.push_back(energyPeak);
209 fRangePeaks.push_back(rangePeak);
210 }
211
212 void LoadConfigFromTiXmlElement(const TiXmlElement* module);
213
214 TRestDataSetGainMap* GetParent() const { return const_cast<TRestDataSetGainMap*>(p); }
215 std::pair<int, int> GetIndexMatrix(const double x, const double y) const;
216 double GetSlope(const double x, const double y) const;
217 double GetIntercept(const double x, const double y) const;
218 double GetSlopeFullSpc() const { return fFullSlope; };
219 double GetInterceptFullSpc() const { return fFullIntercept; };
220
221 Int_t GetPlaneId() const { return fPlaneId; }
222 Int_t GetModuleId() const { return fModuleId; }
223 std::string GetObservable() const { return p->fObservable; }
224 std::string GetSpatialObservableX() const { return p->fSpatialObservableX; }
225 std::string GetSpatialObservableY() const { return p->fSpatialObservableY; }
226 inline std::string GetModuleDefinitionCut() const { return fDefinitionCut; }
227 Int_t GetNumberOfSegmentsX() const { return fNumberOfSegmentsX; }
228 Int_t GetNumberOfSegmentsY() const { return fNumberOfSegmentsY; }
229
230 std::set<double> GetSplitX() const { return fSplitX; }
231 std::set<double> GetSplitY() const { return fSplitY; }
232 std::string GetDataSetFileName() const { return fDataSetFileName; }
233 TVector2 GetReadoutRange() const { return fReadoutRange; }
234
235 void DrawSpectrum(const bool drawFits = true, const int color = -1, TCanvas* c = nullptr);
236 void DrawSpectrum(const TVector2& position, bool drawFits = true, int color = -1,
237 TCanvas* c = nullptr);
238 void DrawSpectrum(const int index_x, const int index_y, bool drawFits = true, int color = -1,
239 TCanvas* c = nullptr);
240 void DrawFullSpectrum(const bool drawFits = true, const int color = -1, TCanvas* c = nullptr);
241
242 void DrawLinearFit(TCanvas* c = nullptr);
243 void DrawLinearFit(const TVector2& position, TCanvas* c = nullptr);
244 void DrawLinearFit(const int index_x, const int index_y, TCanvas* c = nullptr);
245
246 void DrawGainMap(const int peakNumber = 0, const bool fullModuleAsRef = true);
247
248 void Refit(const TVector2& position, const double energy, const TVector2& range);
249 void Refit(const size_t x, const size_t y, const size_t peakNumber, const TVector2& range);
250 void RefitFullSpc(const double energy, const TVector2& range);
251 void RefitFullSpc(const size_t peakNumber, const TVector2& range);
252 void UpdateCalibrationFits(const size_t x, const size_t y);
254
255 void SetPlaneId(const Int_t& planeId) { fPlaneId = planeId; }
256 void SetModuleId(const Int_t& moduleId) { fModuleId = moduleId; }
257 void SetModuleDefinitionCut(const std::string& moduleDefinitionCut) {
258 fDefinitionCut = moduleDefinitionCut;
259 }
260 void SetCalibrationRange(const TVector2& calibrationRange) { fCalibRange = calibrationRange; }
261 void SetNBins(const Int_t& nBins) { fNBins = nBins; }
262 void SetSplitX();
263 void SetSplitY();
264 void SetSplitX(const std::set<double>& splitX);
265 void SetSplitY(const std::set<double>& splitY);
266 void SetSplits();
267 void SetSplits(const std::set<double>& splitXandY) {
268 SetSplitX(splitXandY);
269 SetSplitY(splitXandY);
270 }
271
272 void SetNumberOfSegmentsX(const Int_t& numberOfSegmentsX) { fNumberOfSegmentsX = numberOfSegmentsX; }
273 void SetNumberOfSegmentsY(const Int_t& numberOfSegmentsY) { fNumberOfSegmentsY = numberOfSegmentsY; }
274
275 void SetDataSetFileName(const std::string& dataSetFileName) { fDataSetFileName = dataSetFileName; }
276 void SetReadoutRange(const TVector2& readoutRange) { fReadoutRange = readoutRange; }
277 void SetZeroPoint(const bool& ZeroPoint) { fZeroPoint = ZeroPoint; }
278 void SetAutoRangePeaks(const bool& autoRangePeaks) { fAutoRangePeaks = autoRangePeaks; }
279
280 void Print() const;
281
282 void GenerateGainMap();
283 void Initialize();
284
285 Module() {}
286 Module(const TRestDataSetGainMap& parent) : p(&parent){};
287 Module(const TRestDataSetGainMap& parent, const Int_t planeId, const Int_t moduleId) : p(&parent) {
288 SetPlaneId(planeId);
289 SetModuleId(moduleId);
290 };
291 ~Module(){};
292 };
293};
294#endif
A class to help on cuts definitions. To be used with TRestAnalysisTree.
Definition: TRestCut.h:31
double fFullIntercept
Intercept of the calibration linear fit of whole module.
std::vector< std::vector< double > > fIntercept
Array containing the intercept of the linear fit for each segment.
Int_t fNumberOfSegmentsY
Number of segments in the y direction.
void SetSplitY()
Function to set the class members for segmentation of the detector plane along the Y axis.
std::vector< std::vector< TGraph * > > fSegLinearFit
Array containing the calibration linear fit for each segment.
void SetSplits()
Function to set the class members for segmentation of the detector plane along the X and Y axis.
void SetSplitX()
Function to set the class members for segmentation of the detector plane along the X axis.
std::vector< std::vector< TH1F * > > fSegSpectra
Array containing the observable spectrum for each segment.
void DrawSpectrum(const bool drawFits=true, const int color=-1, TCanvas *c=nullptr)
Function to draw the spectrum for each segment of the module on the same canvas. The canvas is divide...
void LoadConfigFromTiXmlElement(const TiXmlElement *module)
Function to read the parameters from the RML element (TiXmlElement) and set those class members.
void GenerateGainMap()
Function that calculates the calibration parameters for each segment defined at fSplitX and fSplitY a...
double fFullSlope
Slope of the calibration linear fit of whole module.
double GetSlope(const double x, const double y) const
Function to get the calibration parameter slope for a given x and y position on the detector plane.
std::set< double > fSplitX
Split points in the x direction.
void Print() const
Prints on screen the information about the members of Module.
void Refit(const TVector2 &position, const double energy, const TVector2 &range)
Function to fit again manually a peak for a given segment of the module.
std::set< double > fSplitY
Split points in the y direction.
void RefitFullSpc(const double energy, const TVector2 &range)
Function to fit again manually a peak for the whole module spectrum. The calibration curve is updated...
std::vector< TVector2 > fRangePeaks
Range of the peaks to be used for the calibration. If empty it will be automatically calculated.
const TRestDataSetGainMap * p
Pointer to the parent class.
void UpdateCalibrationFitsFullSpc()
Function to update the calibration curve for the whole module. The calibration curve is cleared and t...
std::vector< std::vector< double > > fSlope
Array containing the slope of the linear fit for each segment.
Int_t fNumberOfSegmentsX
Number of segments in the x direction.
void DrawGainMap(const int peakNumber=0, const bool fullModuleAsRef=true)
Function to draw the relative gain map for a given energy peak of the module.
Int_t fNBins
Number of bins for the spectrum histograms.
TVector2 fCalibRange
Calibration range. If fCalibRange.X()>=fCalibRange.Y() the range will be automatically calculated.
std::string fDefinitionCut
Cut that defines which events are from this module.
TGraph * fFullLinearFit
Calibration linear fit for the whole module.
TH1F * fFullSpectrum
Spectrum of the observable for the whole module.
bool fAutoRangePeaks
Automatic range for the peaks fitting. See GenerateGainMap() for more information of the logic.
double GetIntercept(const double x, const double y) const
Function to get the calibration parameter intercept for a given x and y position on the detector plan...
TVector2 fReadoutRange
Readout dimensions.
std::pair< int, int > GetIndexMatrix(const double x, const double y) const
Function to get the index of the matrix of calibration parameters for a given x and y position on the...
std::vector< double > fEnergyPeaks
Energy of the peaks to be used for the calibration.
Metadata class to calculate,store and apply the gain corrected calibration of a group of detectors.
double GetInterceptParameter(const int planeID, const int moduleID, const double x, const double y)
Function to get the intercept parameter of the module with planeID and moduleID at physical position ...
double GetSlopeParameter(const int planeID, const int moduleID, const double x, const double y)
Function to get the slope parameter of the module with planeID and moduleID at physical position (x,...
std::string fObservable
Observable that will be used to calculate the gain map.
double GetInterceptParameterFullSpc(const int planeID, const int moduleID)
Function to get the intercept parameter of the whole module with planeID and moduleID.
void CalibrateDataSet(const std::string &dataSetFileName, std::string outputFileName="", std::vector< std::string > excludeColumns={})
Function to calibrate a dataset with this gain map.
~TRestDataSetGainMap()
Default destructor.
void GenerateGainMap()
Function to calculate the calibration parameters of all modules.
TRestDataSetGainMap()
Default constructor.
void Initialize() override
Making default settings.
std::vector< Module > fModulesCal
List of modules.
void SetModule(const Module &moduleCal)
Function to set a module calibration. If the module calibration already exists (same planeId and modu...
void InitFromConfigFile() override
Initialization of TRestDataSetGainMap members through a RML file.
double GetSlopeParameterFullSpc(const int planeID, const int moduleID)
Function to get the slope parameter of the whole module with planeID and moduleID.
void Export(const std::string &fileName="")
Function to export the calibration to the file fileName.
void PrintMetadata() override
Prints on screen the information about the metadata members.
std::map< int, std::set< int > > GetModuleIDs() const
Function to get the map of the module IDs for each plane ID.
std::string fSpatialObservableY
Observable that will be used to segmentize the gain map in the y direction.
std::string fCalibFileName
Name of the file that contains the calibration data.
TRestCut * fCut
Cut to be applied to the calibration data.
std::string fOutputFileName
Name of the file where the gain map was (or will be) exported.
std::set< int > GetPlaneIDs() const
Function to get a list (set) of the plane IDs.
std::string fSpatialObservableX
Observable that will be used to segmentize the gain map in the x direction.
Module * GetModule(const size_t index=0)
Function to retrieve the module calibration by index. Default is 0.
void Import(const std::string &fileName)
Function to import the calibration parameters from the root file fileName.
A base class for any REST metadata class.
Definition: TRestMetadata.h:70