Belle II Software development
T0Correction Class Reference

Class for T0 Correction . More...

#include <T0Correction.h>

Public Member Functions

 T0Correction ()
 Constructor.
 
virtual ~T0Correction ()
 Destructor.
 
virtual void setDebug (bool debug=false)
 turn on/off debug.
 
virtual void setUseDB (bool useDB=false)
 use DB or text mode.
 
virtual void storeHisto (bool storeHist=false)
 store Hisotgram or not.
 
void setMinimumNDF (double minndf)
 minimum ndf require for track.
 
void setMinimumPval (double minPval)
 minimum pvalue requirement.
 
void inputFileNames (std::string inputname)
 input root file name.
 
void outputFileName (std::string outputname)
 output xt T0 file name (for text mode)
 
void execute ()
 run t0 correction.
 

Protected Member Functions

virtual bool calibrate ()
 Run algo on data.
 
virtual void CreateHisto ()
 create histo for each channel
 
virtual void Write ()
 write outut or store db
 

Private Attributes

TH1F * m_hTotal
 1D histogram of delta T whole channel
 
TH1F * m_h1 [56][385]
 1D histogram for each channel
 
TH1F * m_hT0b [300]
 1D histogram for each board
 
double m_xmin = 0.07
 minimum drift length
 
double m_ndfmin = 5
 minimum ndf required
 
double m_Pvalmin = 0.
 minimum pvalue required
 
double dt [56][385] = {{0.}}
 dt of each channel
 
double err_dt [56][385] = {{0.}}
 error of dt of each channel
 
double dtb [300] = {0.}
 dt of each board
 
double err_dtb [300] = {0.}
 error of dt of board
 
bool m_debug
 debug.
 
bool m_storeHisto
 store histo or not
 
bool m_useDB
 use DB or text mode
 
std::string m_outputT0FileName = "t0_new.dat"
 output t0 file name for text file
 
std::string m_inputRootFileName = "rootfile/output*"
 input file names
 
int m_firstExperiment
 First experiment.
 
int m_firstRun
 First run.
 
int m_lastExperiment
 Last experiment.
 
int m_lastRun
 Last run.
 

Detailed Description

Class for T0 Correction .

Definition at line 21 of file T0Correction.h.

Constructor & Destructor Documentation

◆ T0Correction()

Constructor.

Definition at line 28 of file T0Correction.cc.

28 :
31{
32 /*
33 setDescription(
34 " -------------------------- Test Calibration Algoritm -------------------------\n"
35 " \n"
36 " Testing algorithm which just gets mean of a test histogram collected by \n"
37 " CaTest module and provides a DB object with another histogram with one \n"
38 " entry at calibrated value. \n"
39 " ------------------------------------------------------------------------------\n"
40 );
41 */
42}
int m_lastExperiment
Last experiment.
Definition: T0Correction.h:74
int m_firstExperiment
First experiment.
Definition: T0Correction.h:72
int m_firstRun
First run.
Definition: T0Correction.h:73

◆ ~T0Correction()

virtual ~T0Correction ( )
inlinevirtual

Destructor.

Definition at line 26 of file T0Correction.h.

26{}

Member Function Documentation

◆ calibrate()

bool calibrate ( )
protectedvirtual

Run algo on data.

Definition at line 116 of file T0Correction.cc.

117{
118 B2INFO("Start calibration");
119
120 gROOT->SetBatch(1);
121 gErrorIgnoreLevel = 3001;
122
123 CreateHisto();
124 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
125
126 TF1* g1 = new TF1("g1", "gaus", -100, 100);
127 vector<double> b, db, Sb, dSb;
128 vector<double> c[56];
129 vector<double> dc[56];
130 vector<double> s[56];
131 vector<double> ds[56];
132
133 B2INFO("Gaus fitting for whole channel");
134 double par[3];
135 m_hTotal->SetDirectory(0);
136 double mean = m_hTotal->GetMean();
137 m_hTotal->Fit("g1", "Q", "", mean - 15, mean + 15);
138 g1->GetParameters(par);
139
140 B2INFO("Gaus fitting for each board");
141 for (int ib = 1; ib < 300; ++ib) {
142 if (m_hT0b[ib]->GetEntries() < 10) continue;
143 mean = m_hT0b[ib]->GetMean();
144 m_hT0b[ib]->SetDirectory(0);
145 m_hT0b[ib]->Fit("g1", "Q", "", mean - 15, mean + 15);
146 g1->GetParameters(par);
147 b.push_back(ib);
148 db.push_back(0.0);
149 Sb.push_back(par[1]);
150 dSb.push_back(g1->GetParError(1));
151 dtb[ib] = par[1];
152 err_dtb[ib] = g1->GetParError(1);
153 }
154 B2INFO("Gaus fitting for each cell");
155 for (int ilay = 0; ilay < 56; ++ilay) {
156 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
157 const int n = m_h1[ilay][iwire]->GetEntries();
158 B2DEBUG(99, "layer " << ilay << " wire " << iwire << " entries " << n);
159 if (n < 10) continue;
160 mean = m_h1[ilay][iwire]->GetMean();
161 m_h1[ilay][iwire]->SetDirectory(0);
162 m_h1[ilay][iwire]->Fit("g1", "Q", "", mean - 15, mean + 15);
163 g1->GetParameters(par);
164 c[ilay].push_back(iwire);
165 dc[ilay].push_back(0.0);
166 s[ilay].push_back(par[1]);
167 ds[ilay].push_back(g1->GetParError(1));
168
169 dt[ilay][iwire] = par[1];
170 err_dt[ilay][iwire] = g1->GetParError(1);
171 }
172 }
173
174
175 if (m_storeHisto) {
176 B2INFO("Store histo");
177 TFile* fout = new TFile("Correct_T0.root", "RECREATE");
178 fout->cd();
179 TGraphErrors* gr[56];
180 TDirectory* top = gDirectory;
181 m_hTotal->Write();
182 TDirectory* subDir[56];
183 for (int ilay = 0; ilay < 56; ++ilay) {
184 subDir[ilay] = top ->mkdir(Form("lay_%d", ilay));
185 subDir[ilay]->cd();
186 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
187 m_h1[ilay][iwire]->Write();
188 }
189 }
190 top->cd();
191 TDirectory* corrT0 = top->mkdir("DeltaT0");
192 corrT0->cd();
193
194 if (b.size() > 2) {
195 TGraphErrors* grb = new TGraphErrors(b.size(), &b.at(0), &Sb.at(0), &db.at(0), &dSb.at(0));
196 grb->SetMarkerColor(2);
197 grb->SetMarkerSize(1.0);
198 grb->SetTitle("#DeltaT0;BoardID;#DeltaT0[ns]");
199 grb->SetMaximum(10);
200 grb->SetMinimum(-10);
201 grb->SetName("Board");
202 grb->Write();
203 }
204 for (int sl = 0; sl < 56; ++sl) {
205 if (c[sl].size() < 2) continue;
206 gr[sl] = new TGraphErrors(c[sl].size(), &c[sl].at(0), &s[sl].at(0), &dc[sl].at(0), &ds[sl].at(0));
207 gr[sl]->SetMarkerColor(2);
208 gr[sl]->SetMarkerSize(1.0);
209 gr[sl]->SetTitle(Form("Layer_%d;IWire;#LT t_{mea}-t_{fit} #GT [ns]", sl));
210 gr[sl]->SetMaximum(10);
211 gr[sl]->SetMinimum(-10);
212 gr[sl]->SetName(Form("lay%d", sl));
213 gr[sl]->Write();
214 }
215 fout->Close();
216 }
217 B2INFO("Write constants");
218 Write();
219 return true;
220}
The Class for CDC Geometry Parameters.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
double dtb[300]
dt of each board
Definition: T0Correction.h:64
TH1F * m_hTotal
1D histogram of delta T whole channel
Definition: T0Correction.h:55
double dt[56][385]
dt of each channel
Definition: T0Correction.h:62
bool m_storeHisto
store histo or not
Definition: T0Correction.h:68
virtual void Write()
write outut or store db
virtual void CreateHisto()
create histo for each channel
Definition: T0Correction.cc:44
TH1F * m_h1[56][385]
1D histogram for each channel
Definition: T0Correction.h:56
double err_dt[56][385]
error of dt of each channel
Definition: T0Correction.h:63
TH1F * m_hT0b[300]
1D histogram for each board
Definition: T0Correction.h:57
double err_dtb[300]
error of dt of board
Definition: T0Correction.h:65

◆ CreateHisto()

void CreateHisto ( )
protectedvirtual

create histo for each channel

Definition at line 44 of file T0Correction.cc.

45{
46
47 B2INFO("CreateHisto");
48 double x;
49 double t_mea;
50 double w;
51 double t_fit;
52 double ndf;
53 double Pval;
54 int IWire;
55 int lay;
56
57 // auto& tree = getObject<TTree>("cdcCalib");
58 TChain* tree = new TChain("tree");
59 tree->Add(m_inputRootFileName.c_str());
60 tree->SetBranchAddress("lay", &lay);
61 tree->SetBranchAddress("IWire", &IWire);
62 tree->SetBranchAddress("x_u", &x);
63 tree->SetBranchAddress("t", &t_mea);
64 tree->SetBranchAddress("t_fit", &t_fit);
65 tree->SetBranchAddress("weight", &w);
66 tree->SetBranchAddress("ndf", &ndf);
67 tree->SetBranchAddress("Pval", &Pval);
68
69 /* Disable unused branch */
70 std::vector<TString> list_vars = {"lay", "IWire", "x_u", "t", "t_fit", "weight", "Pval", "ndf"};
71 tree->SetBranchStatus("*", 0);
72
73 for (TString brname : list_vars) {
74 tree->SetBranchStatus(brname, 1);
75 }
76
77
78 double halfCSize[56];
80 for (int i = 0; i < 56; ++i) {
81 double R = cdcgeo.senseWireR(i);
82 double nW = cdcgeo.nWiresInLayer(i);
83 halfCSize[i] = M_PI * R / nW;
84 }
85
86 m_hTotal = new TH1F("hTotal", "hTotal", 30, -15, 15);
87 //for each channel
88 for (int il = 0; il < 56; il++) {
89 const int nW = cdcgeo.nWiresInLayer(il);
90 for (int ic = 0; ic < nW; ++ic) {
91 m_h1[il][ic] = new TH1F(Form("h1%d@%d", il, ic), Form("L%d_cell%d", il, ic), 30, -15, 15);
92 }
93 }
94 //for each board
95 for (int ib = 0; ib < 300; ++ib) {
96 m_hT0b[ib] = new TH1F(Form("hT0b%d", ib), Form("boardID_%d", ib), 100, -20, 20);
97 }
98 //read data
99 const int nEntries = tree->GetEntries();
100 B2INFO("Number of entries: " << nEntries);
101 for (int i = 0; i < nEntries; ++i) {
102 tree->GetEntry(i);
103 double xmax = halfCSize[lay] - 0.1;
104 if ((fabs(x) < m_xmin) || (fabs(x) > xmax)
105 || (ndf < m_ndfmin)
106 || (Pval < m_Pvalmin)) continue; /*select good region*/
107 //each channel
108 m_hTotal->Fill(t_mea - t_fit);
109 m_h1[lay][IWire]->Fill(t_mea - t_fit);
110 //each board
111 m_hT0b[cdcgeo.getBoardID(WireID(lay, IWire))]->Fill(t_mea - t_fit);
112 }
113 B2INFO("Finish making histogram for all channels");
114}
double R
typedef autogenerated by FFTW
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
double m_xmin
minimum drift length
Definition: T0Correction.h:58
std::string m_inputRootFileName
input file names
Definition: T0Correction.h:71
double m_ndfmin
minimum ndf required
Definition: T0Correction.h:59
double m_Pvalmin
minimum pvalue required
Definition: T0Correction.h:60
Class to identify a wire inside the CDC.
Definition: WireID.h:34

◆ execute()

void execute ( )
inline

run t0 correction.

Definition at line 42 of file T0Correction.h.

43 {
44 calibrate();
45 }
virtual bool calibrate()
Run algo on data.

◆ inputFileNames()

void inputFileNames ( std::string  inputname)
inline

input root file name.

Definition at line 38 of file T0Correction.h.

38{m_inputRootFileName.assign(inputname);}

◆ outputFileName()

void outputFileName ( std::string  outputname)
inline

output xt T0 file name (for text mode)

Definition at line 40 of file T0Correction.h.

40{m_outputT0FileName.assign(outputname);}
std::string m_outputT0FileName
output t0 file name for text file
Definition: T0Correction.h:70

◆ setDebug()

virtual void setDebug ( bool  debug = false)
inlinevirtual

turn on/off debug.

Definition at line 28 of file T0Correction.h.

28{m_debug = debug; }

◆ setMinimumNDF()

void setMinimumNDF ( double  minndf)
inline

minimum ndf require for track.

Definition at line 34 of file T0Correction.h.

34{m_ndfmin = minndf;}

◆ setMinimumPval()

void setMinimumPval ( double  minPval)
inline

minimum pvalue requirement.

Definition at line 36 of file T0Correction.h.

36{m_Pvalmin = minPval;}

◆ setUseDB()

virtual void setUseDB ( bool  useDB = false)
inlinevirtual

use DB or text mode.

Definition at line 30 of file T0Correction.h.

30{m_useDB = useDB; }
bool m_useDB
use DB or text mode
Definition: T0Correction.h:69

◆ storeHisto()

virtual void storeHisto ( bool  storeHist = false)
inlinevirtual

store Hisotgram or not.

Definition at line 32 of file T0Correction.h.

32{m_storeHisto = storeHist;}

◆ Write()

void Write ( )
protectedvirtual

write outut or store db

Definition at line 221 of file T0Correction.cc.

222{
223 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
224 ofstream ofs(m_outputT0FileName.c_str());
226 tz.construct();
227
228 double T0;
229 TH1F* T0B[300];
230 for (int ib = 0; ib < 300; ++ib) {
231 T0B[ib] = new TH1F(Form("T0B%d", ib), Form("boardID_%d", ib), 9000, 0, 9000);
232 }
233 //for calculate T0 mean of each board
234 for (int ilay = 0; ilay < 56; ++ilay) {
235 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
236 WireID wireid(ilay, iwire);
237 T0 = cdcgeo.getT0(wireid);
238 T0B[cdcgeo.getBoardID(wireid)]->Fill(T0);
239 }
240 }
241
242 //correct T0 and write
243 for (int ilay = 0; ilay < 56; ++ilay) {
244 for (unsigned int iwire = 0; iwire < cdcgeo.nWiresInLayer(ilay); ++iwire) {
245 WireID wireid(ilay, iwire);
246 int bID = cdcgeo.getBoardID(wireid);
247 if (abs(err_dt[ilay][iwire]) > 2 || abs(dt[ilay][iwire]) > 1.e3) {
248 // if (abs(err_dt[ilay][iwire]) > 2) {
249 T0 = T0B[bID]->GetMean();
250 dt[ilay][iwire] = dtb[bID];
251 } else {
252 T0 = cdcgeo.getT0(wireid);
253 }
254 ofs << ilay << "\t" << iwire << "\t" << T0 - dt[ilay][iwire] << std::endl;
255 if (m_useDB) {
256 tz->setT0(wireid, T0 - dt[ilay][iwire]);
257 }
258 }
259 }
260 ofs.close();
261 if (m_useDB) {
264 tz.import(iov);
265 B2RESULT("T0 was calibrated and imported to database.");
266 }
267
268}
float getT0(const WireID &wireID) const
Returns t0 parameter of the specified sense wire.
bool import(const IntervalOfValidity &iov)
Import the object to database.
Definition: DBImportBase.cc:36
Class for importing a single object to the database.
void construct(Args &&... params)
Construct an object of type T in this DBImportObjPtr using the provided constructor arguments.
A class that describes the interval of experiments/runs for which an object in the database is valid.

Member Data Documentation

◆ dt

double dt[56][385] = {{0.}}
private

dt of each channel

Definition at line 62 of file T0Correction.h.

◆ dtb

double dtb[300] = {0.}
private

dt of each board

Definition at line 64 of file T0Correction.h.

◆ err_dt

double err_dt[56][385] = {{0.}}
private

error of dt of each channel

Definition at line 63 of file T0Correction.h.

◆ err_dtb

double err_dtb[300] = {0.}
private

error of dt of board

Definition at line 65 of file T0Correction.h.

◆ m_debug

bool m_debug
private

debug.

Definition at line 67 of file T0Correction.h.

◆ m_firstExperiment

int m_firstExperiment
private

First experiment.

Definition at line 72 of file T0Correction.h.

◆ m_firstRun

int m_firstRun
private

First run.

Definition at line 73 of file T0Correction.h.

◆ m_h1

TH1F* m_h1[56][385]
private

1D histogram for each channel

Definition at line 56 of file T0Correction.h.

◆ m_hT0b

TH1F* m_hT0b[300]
private

1D histogram for each board

Definition at line 57 of file T0Correction.h.

◆ m_hTotal

TH1F* m_hTotal
private

1D histogram of delta T whole channel

Definition at line 55 of file T0Correction.h.

◆ m_inputRootFileName

std::string m_inputRootFileName = "rootfile/output*"
private

input file names

Definition at line 71 of file T0Correction.h.

◆ m_lastExperiment

int m_lastExperiment
private

Last experiment.

Definition at line 74 of file T0Correction.h.

◆ m_lastRun

int m_lastRun
private

Last run.

Definition at line 75 of file T0Correction.h.

◆ m_ndfmin

double m_ndfmin = 5
private

minimum ndf required

Definition at line 59 of file T0Correction.h.

◆ m_outputT0FileName

std::string m_outputT0FileName = "t0_new.dat"
private

output t0 file name for text file

Definition at line 70 of file T0Correction.h.

◆ m_Pvalmin

double m_Pvalmin = 0.
private

minimum pvalue required

Definition at line 60 of file T0Correction.h.

◆ m_storeHisto

bool m_storeHisto
private

store histo or not

Definition at line 68 of file T0Correction.h.

◆ m_useDB

bool m_useDB
private

use DB or text mode

Definition at line 69 of file T0Correction.h.

◆ m_xmin

double m_xmin = 0.07
private

minimum drift length

Definition at line 58 of file T0Correction.h.


The documentation for this class was generated from the following files: