Belle II Software development
CDCDedx1DCellAlgorithm.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <TH1D.h>
12#include <TH2D.h>
13
14#include <cdc/dbobjects/CDCDedx1DCell.h>
15#include <calibration/CalibrationAlgorithm.h>
16#include <framework/database/DBObjPtr.h>
17
18namespace Belle2 {
23
29
30 public:
31
32 static constexpr int m_kNGroups = 3;
33
38
43
47 void setEntaRange(double min = -1.0, double max = 1.0) {m_eaMin = min; m_eaMax = max;}
48
52 void setEntaBins(unsigned int value = 316) {m_eaB = value;}
53
57 void setHistBins(int value = 250) {m_dedxBin = value;}
58
62 void setHistRange(double min = 0.0, double max = 5.0) {m_dedxMin = min; m_dedxMax = max;}
63
67 void setPtLimit(double value) {m_ptMax = value;}
68
72 void setCosLimit(double value) {m_cosMax = value;}
73
77 void setTrucationBins(double lowedge, double upedge)
78 {
79 m_truncMin = lowedge; m_truncMax = upedge ;
80 }
81
85 void setSplitFactor(int value) {m_binSplit = value;}
86
90 void setChargeType(int value) {m_chargeType = value;}
91
95 void setAdjustmentFactor(int value) {m_adjustFac = value;}
96
100 void setLayerTrunc(bool value = false) {isFixTrunc = value;}
101
105 void setVariableBins(bool value) {isVarBins = value;}
106
110 void setRotSymmetry(bool value) {isRotSymm = value;}
111
115 void enableExtraPlots(bool value = false) {isMakePlots = value;}
116
120 void setPrintLog(bool value = false) {isPrintLog = value;}
121
126 void setMergePayload(bool value) { isMerge = value;}
127
131 void setSuffix(const std::string& value) {m_suffix = value;}
132
137 unsigned int getRepresentativeLayer(unsigned int il) const
138 {
139 static const std::array<unsigned int, m_kNGroups> repLayer = {1, 9, 17};
140 return repLayer.at(il);
141 }
142
147 void setBaselineFactor(double charge, double factor)
148 {
149
150 m_adjustFac = factor;
151 if (charge < 0)m_chargeType = -1.0;
152 else if (charge > 0)m_chargeType = 1.0;
153 else if (charge == 0)m_chargeType = 0.0;
154 else
155 B2FATAL("Choose charge value either +/-1 or 0");
156 }
157
161 int rotationalBin(int nbin, int ibin)
162 {
163 if (nbin % 4 != 0)return ibin;
164 int jbin = ibin;
165 if (ibin < nbin / 4) jbin = ibin + nbin / 2 ;
166 else if (ibin >= 3 * nbin / 4) jbin = ibin - nbin / 2 ;
167 return jbin;
168 }
169
173 void getExpRunInfo();
174
178 void CreateBinMapping();
179
183 void defineHisto(std::array<std::vector<TH1D*>, 3>& hdedxhit, std::array<TH1D*, 3>& hdedxlay, std::array<TH1D*, 3>& hentalay);
184
188 void getTruncatedBins(TH1D* hist, int& binlow, int& binhigh);
189
193 double getTruncationMean(TH1D* hist, int binlow, int binhigh);
194
198 void createPayload();
199
203 void plotMergeFactor(std::map<int, std::vector<double>> bounds, const std::array<int, 2> nDev,
204 std::map<int, std::vector<int>> steps);
205
209 void plotdedxHist(std::array<std::vector<TH1D*>, m_kNGroups>& hdedxhit);
210
214 void plotLayerDist(std::array<TH1D*, m_kNGroups>& hdedxlay);
215
219 void plotQaPars(std::array<TH1D*, m_kNGroups>& hentalay, TH2D* hptcosth);
220
224 void plotRelConst(std::vector<double>tempconst, std::vector<double>layerconst, int il);
225
229 void plotConstants();
230
234 void plotEventStats();
235
236 protected:
237
241 virtual EResult calibrate() override;
242
243
244 private:
245
246 double m_eaMin;
247 double m_eaMax;
248 double m_eaBW;
250 int m_eaB;
251
252 double m_dedxMin;
253 double m_dedxMax;
255
256 double m_ptMax;
257 double m_cosMax;
258
259 double m_truncMin;
260 double m_truncMax;
261
263
266
272 bool isMerge;
273
274 std::string m_suffix;
275 std::string m_runExp;
276 std::string m_label[m_kNGroups] = {"SL0", "SL1", "SL2-8"};
277
278 std::vector<int> m_eaBinLocal;
279 std::array<std::vector<int>, m_kNGroups> m_binIndex;
280 std::array<std::vector<double>, m_kNGroups>m_binValue;
281
282 std::vector<std::vector<double>> m_onedcors;
283
285
286 };
287
288} // namespace Belle2
void plotLayerDist(std::array< TH1D *, m_kNGroups > &hdedxlay)
function to draw dedx dist.
void setSplitFactor(int value)
set bin split factor for all range
CDCDedx1DCellAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
void setChargeType(int value)
set charge type
double m_eaMax
upper edge of entrance angle
void plotMergeFactor(std::map< int, std::vector< double > > bounds, const std::array< int, 2 > nDev, std::map< int, std::vector< int > > steps)
function to plot merging factor
double m_truncMax
upper threshold on truncation
int m_binSplit
multiply nbins by this factor in full range
double m_truncMin
lower threshold on truncation
double m_adjustFac
factor with that one what to adjust baseline
void setBaselineFactor(double charge, double factor)
adjust baseline based on charge or global overall works for only single charge or both
void setSuffix(const std::string &value)
adding suffix to control plots
void enableExtraPlots(bool value=false)
funtion to set flag active for plotting
void getTruncatedBins(TH1D *hist, int &binlow, int &binhigh)
function to get bins of truncation from histogram
void CreateBinMapping()
class function to create vectors for bin mapping (Var->symm)
double m_chargeType
charge type for baseline adj
void getExpRunInfo()
function to extract calibration run/exp
void setHistRange(double min=0.0, double max=5.0)
function to set min/max range of dedx dist for calibration
unsigned int getRepresentativeLayer(unsigned int il) const
Representative CDC layer for each SL group (used to access group-wise constants): SL0 => 1,...
void setEntaBins(unsigned int value=316)
function to set number of entrance angle bins for calibration
static constexpr int m_kNGroups
SL grouping: inner (SL0), middle (SL1), outer (SL2–8)
DBObjPtr< CDCDedx1DCell > m_DBOneDCell
One cell correction DB object.
double m_cosMax
a limit on cos theta
void setHistBins(int value=250)
function to set nbins of dedx dist for calibration
bool isPrintLog
print more debug information
virtual ~CDCDedx1DCellAlgorithm()
Destructor.
std::string m_suffix
add suffix to all plot name
int m_eaB
reset # of bins for entrance angle for each experiment
double getTruncationMean(TH1D *hist, int binlow, int binhigh)
function to get truncated mean
double m_ptMax
a limit on transverse momentum
void setLayerTrunc(bool value=false)
function to set truncation method (local vs global)
void plotConstants()
function to draw the old/new final constants
void setPrintLog(bool value=false)
funtion to set flag to print log
std::array< std::vector< int >, m_kNGroups > m_binIndex
symm/Var bin numbers
void setMergePayload(bool value)
set false if generating absolute (not relative) payload
bool isFixTrunc
true = fix window for all out/inner layers
bool isVarBins
true: if variable bin size is requested
void plotQaPars(std::array< TH1D *, m_kNGroups > &hentalay, TH2D *hptcosth)
function to draw pt vs costh and entrance angle distribution for Inner/Outer layer
void setPtLimit(double value)
function to set pt limit
double m_eaBW
binwdith of entrance angle bin
bool isRotSymm
if rotation symmetry requested
std::string m_runExp
add run and exp to title of plot
void setAdjustmentFactor(int value)
set adjustment factor
void defineHisto(std::array< std::vector< TH1D * >, 3 > &hdedxhit, std::array< TH1D *, 3 > &hdedxlay, std::array< TH1D *, 3 > &hentalay)
function to define histograms
void setVariableBins(bool value)
Set Var bins flag to on or off.
void setEntaRange(double min=-1.0, double max=1.0)
function to set min/max range of entrance angle for calibration
void plotEventStats()
function to draw the stats plots
int rotationalBin(int nbin, int ibin)
class function to set rotation symmetry
virtual EResult calibrate() override
1D cell algorithm
void setRotSymmetry(bool value)
set rotation sys to copy constants from one region to other
std::array< std::vector< double >, m_kNGroups > m_binValue
enta Var bin values
double m_dedxMax
upper edge of dedxhit
void createPayload()
function to generate final constants
void setTrucationBins(double lowedge, double upedge)
function to set bins of truncation from histogram
bool isMakePlots
produce plots for status
void plotRelConst(std::vector< double >tempconst, std::vector< double >layerconst, int il)
function to draw symm/Var layer constant
std::vector< std::vector< double > > m_onedcors
final vectors of calibration
bool isMerge
print more debug information
void setCosLimit(double value)
function to set costheta limit
double m_dedxMin
lower edge of dedxhit
std::string m_label[m_kNGroups]
add inner/outer superlayer label
void plotdedxHist(std::array< std::vector< TH1D * >, m_kNGroups > &hdedxhit)
function to draw the dE/dx histogram in enta bins
double m_eaMin
lower edge of entrance angle
EResult
The result of calibration.
CalibrationAlgorithm(const std::string &collectorModuleName)
Constructor - sets the prefix for collected objects (won't be accesses until execute(....
Class for accessing objects in the database.
Definition DBObjPtr.h:21
Abstract base class for different kinds of events.