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 <algorithm>
12#include <iostream>
13#include <fstream>
14
15#include <TMath.h>
16#include <TH1D.h>
17#include <TCanvas.h>
18#include <TH2D.h>
19#include <TTree.h>
20#include <TStyle.h>
21#include <TLegend.h>
22#include <TRandom.h>
23
24
25#include <cdc/dbobjects/CDCDedx1DCell.h>
26#include <calibration/CalibrationAlgorithm.h>
27#include <framework/database/DBObjPtr.h>
28
29namespace Belle2 {
40
41 public:
42
47
52
56 void setSuffix(const std::string& value) {m_suffix = value;}
57
61 void setVariableBins(bool value) {isVarBins = value;}
62
66 void setLayerTrunc(bool value = false) {isFixTrunc = value;}
67
72 void setMergePayload(bool value) { isMerge = value;}
73
77 void setSplitFactor(int value) {m_binSplit = value;}
78
82 void setRotSymmetry(bool value) {isRotSymm = value;}
83
87 void setTrucationBins(double lowedge, double upedge)
88 {
89 m_truncMin = lowedge; m_truncMax = upedge ;
90 }
91
95 void enableExtraPlots(bool value = false) {isMakePlots = value;}
96
100 void setPtLimit(double value) {m_ptMax = value;}
101
105 void setCosLimit(double value) {m_cosMax = value;}
106
111 void setBaselineFactor(double charge, double factor)
112 {
113
114 m_adjustFac = factor;
115 if (charge < 0)m_chargeType = -1.0;
116 else if (charge > 0)m_chargeType = 1.0;
117 else if (charge == 0)m_chargeType = 0.0;
118 else
119 B2FATAL("Choose charge value either +/-1 or 0");
120 }
121
125 int rotationalBin(int nbin, int ibin)
126 {
127 if (nbin % 4 != 0)return ibin;
128 int jbin = ibin;
129 if (ibin < nbin / 4) jbin = ibin + nbin / 2 ;
130 else if (ibin >= 3 * nbin / 4) jbin = ibin - nbin / 2 ;
131 return jbin;
132 }
133
137 void getExpRunInfo();
138
142 void CreateBinMapping();
143
147 void defineHisto(std::vector<TH1D*> hdedxhit[2], TH1D* hdedxlay[2], TH1D* hentalay[2]);
148
152 void getTruncatedBins(TH1D* hist, int& binlow, int& binhigh);
153
157 double getTruncationMean(TH1D* hist, int binlow, int binhigh);
158
162 void createPayload();
163
167 void plotMergeFactor(std::map<int, std::vector<double>> bounds, const std::array<int, 2> nDev,
168 std::map<int, std::vector<int>> steps);
169
173 void plotdedxHist(std::vector<TH1D*> hdedxhit[2]);
174
178 void plotLayerDist(TH1D* hdedxL[2]);
179
183 void plotQaPars(TH1D* hentalay[2], TH2D* hptcosth);
184
188 void plotRelConst(std::vector<double>tempconst, std::vector<double>layerconst, int il);
189
193 void plotConstants();
194
198 void plotEventStats();
199
200 protected:
204 virtual EResult calibrate() override;
205
206
207 private:
208
209 double m_eaMin;
210 double m_eaMax;
211 double m_eaBW;
213 int m_eaB;
215 double m_dedxMin;
216 double m_dedxMax;
219 double m_ptMax;
220 double m_cosMax;
222 double m_truncMin;
223 double m_truncMax;
235 bool isMerge;
237 std::string m_suffix;
238 std::string m_runExp;
239 std::string m_label[2] = {"IL", "OL"};
241 std::vector<int> m_eaBinLocal;
242 std::array<std::vector<int>, 2> m_binIndex;
243 std::array<std::vector<double>, 2>m_binValue;
245 std::vector<std::vector<double>> m_onedcors;
250 };
252} // namespace Belle2
A calibration algorithm for CDC dE/dx electron: 1D enta cleanup correction.
void setSplitFactor(int value)
set bin split factor for all range
CDCDedx1DCellAlgorithm()
Constructor: Sets the description, the properties and the parameters of the algorithm.
std::string m_label[2]
add inner/outer layer label
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
std::array< std::vector< int >, 2 > m_binIndex
symm/Var bin numbers
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)
function 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 get extract calibration run/exp
DBObjPtr< CDCDedx1DCell > m_DBOneDCell
One cell correction DB object.
double m_cosMax
a limit on cos theta
bool isPrintLog
print more debug information
std::array< std::vector< double >, 2 > m_binValue
enta Var bin values
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 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 plotdedxHist(std::vector< TH1D * > hdedxhit[2])
function to draw the dE/dx histogram in enta bins
void setPtLimit(double value)
function to set pt limit
void defineHisto(std::vector< TH1D * > hdedxhit[2], TH1D *hdedxlay[2], TH1D *hentalay[2])
function to define histograms
double m_eaBW
binwdith of entrance angle bin
bool isRotSymm
if rotation symmetry requested
std::string m_runExp
add suffix to all plot name
void setVariableBins(bool value)
Set Var bins flag to on or off.
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
void plotQaPars(TH1D *hentalay[2], TH2D *hptcosth)
function to draw pt vs costh and entrance angle distribution for Inner/Outer layer
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 cos limit
double m_dedxMin
lower edge of dedxhit
void plotLayerDist(TH1D *hdedxL[2])
function to draw dedx dist.
double m_eaMin
lower edge of entrance angle
Base class for calibration algorithms.
EResult
The result of calibration.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Abstract base class for different kinds of events.