8#include <cdc/calibration/TimeWalkCalibrationAlgorithm.h>
9#include <cdc/dbobjects/CDCTimeWalks.h>
10#include <cdc/geometry/CDCGeometryPar.h>
11#include <cdc/dataobjects/WireID.h>
15#include <TDirectory.h>
18#include <TStopwatch.h>
20#include <framework/database/DBObjPtr.h>
21#include <framework/database/IntervalOfValidity.h>
22#include <framework/logging/Logger.h>
31 " -------------------------- Time Walk Calibration Algorithm -------------------------\n"
37 B2INFO(
"Creating and filling histograms");
42 for (
int i = 0; i < 56; ++i) {
45 halfCSize[i] = M_PI *
R / nW;
49 for (
int i = 0; i < 300; ++i) {
50 m_h2[i] =
new TH2D(Form(
"board_%d", i), Form(
"board_%d", i), 50, 0., 500, 60, -30, 30);
55 auto tree = getObjectPtr<TTree>(
"tree");
66 tree->SetBranchAddress(
"lay", &lay);
67 tree->SetBranchAddress(
"IWire", &IWire);
68 tree->SetBranchAddress(
"x_u", &x);
69 tree->SetBranchAddress(
"t", &t_mea);
70 tree->SetBranchAddress(
"t_fit", &t_fit);
71 tree->SetBranchAddress(
"ndf", &ndf);
72 tree->SetBranchAddress(
"Pval", &Pval);
73 tree->SetBranchAddress(
"adc", &adc);
76 std::vector<TString> list_vars = {
"lay",
"IWire",
"x_u",
"t",
"t_fit",
"Pval",
"ndf",
"adc"};
77 tree->SetBranchStatus(
"*", 0);
79 for (TString brname : list_vars) {
80 tree->SetBranchStatus(brname, 1);
83 const Long64_t nEntries = tree->GetEntries();
84 B2INFO(
"Number of entries: " << nEntries);
87 for (Long64_t i = 0; i < nEntries; ++i) {
89 const double xmax = halfCSize[lay] - 0.12;
90 if ((fabs(x) <
m_xmin) || (fabs(x) > xmax)
97 B2INFO(
"Time to fill histograms: " << time.RealTime() <<
"s");
103 B2INFO(
"Start calibration");
107 B2INFO(
"ExpRun used for DB Geometry : " << exprun.first <<
" " << exprun.second);
109 B2INFO(
"Creating CDCGeometryPar object");
117 fold =
new TF1(
"fold",
"[0]/sqrt(x)", 0, 500);
119 fold =
new TF1(
"fold",
"[0]*exp(-1*[1]*x)", 0, 500);
121 if (fold ==
nullptr) {
122 B2FATAL(
"Old fitting function is not defined.");
126 B2INFO(
"time walk formular: ");
127 [](TF1 * f) {
auto title = f->GetTitle(); B2INFO(
"Expression : " << title);}(fold);
130 for (
int ib = 1; ib < 300; ++ib) {
132 B2DEBUG(21,
"Board ID:" << ib);
133 m_h2[ib]->SetDirectory(0);
136 if (
m_h2[ib]->GetEntries() < 500) {
140 m_h2[ib]->FitSlicesY(0, 1, -1, 10);
141 TString name =
m_h2[ib]->GetName();
142 TString hm_name = name +
"_1";
143 m_h1[ib] = (TH1D*)gDirectory->Get(hm_name);
145 m_h1[ib]->SetDirectory(0);
146 if (
m_h1[ib]->GetEntries() < 5) {
148 B2WARNING(
"Low statistic, number of points after slice fit: " <<
m_h1[ib]->GetEntries());
154 fold->SetParameter(p,
m_tw_old[ib][p]);
159 TF1* func =
new TF1(
"func",
"[0]+[1]/sqrt(x)", 0, 500);
160 func->SetParameters(-4, 28);
161 m_h1[ib]->Fit(
"func",
"MQ",
"", 20, 150);
166 B2FATAL(
"Mode " <<
m_twParamMode <<
" is not available, please check again");
170 TF1* f1 =
m_h1[ib]->GetFunction(
"func");
171 if (!f1) {
m_flag[ib] = 0;
continue;}
174 m_tw_new[ib][i - 1] = f1->GetParameter(i);
178 B2DEBUG(21,
"Prob of fitting:" << f1->GetProb());
179 B2DEBUG(21,
"Fitting Param 0-1:" << f1->GetParameter(0) <<
" - " << f1->GetParameter(1));
196 B2INFO(
"Storing histogram");
197 B2DEBUG(21,
"Store 1D histogram");
198 TFile* fhist =
new TFile(
m_histName.c_str(),
"RECREATE");
199 auto hNDF = getObjectPtr<TH1F>(
"hNDF");
200 auto hPval = getObjectPtr<TH1F>(
"hPval");
201 auto hEvtT0 = getObjectPtr<TH1F>(
"hEventT0");
203 if (hNDF && hPval && hEvtT0) {
209 TDirectory* old = gDirectory;
210 TDirectory* h1D = old->mkdir(
"h1D");
211 TDirectory* h2D = old->mkdir(
"h2D");
213 for (
int ib = 1; ib < 300; ++ib) {
215 if (
m_h1[ib]->GetEntries() < 5)
continue;
216 m_h1[ib]->SetMinimum(-5);
217 m_h1[ib]->SetMaximum(15);
221 B2DEBUG(21,
"Store 2D histogram");
223 for (
int ib = 1; ib < 300; ++ib) {
224 if (
m_h2[ib] ==
nullptr)
continue;
225 if (
m_h2[ib]->GetEntries() < 5)
continue;
231 B2INFO(
"Hitograms were stored");
236 TH1F* hDtw =
new TH1F(
"hDtw",
"", 100, -1, 1);
237 for (
int ib = 0; ib < 300; ++ib) {
239 if (std::isnan(dtw) == 0) {
244 B2INFO(hDtw->GetRMS());
246 if (hDtw->GetRMS() < 0.02) {
256 B2INFO(
"Save to the local DB");
261 for (
int ib = 0; ib < 300; ++ib) {
266 const int num =
static_cast<int>(
m_tw_old[ib].size());
267 for (
int i = 0; i < num; ++i) {
277 B2RESULT(
"Time-walk coeff. table has been written to " <<
m_outputTWFileName.c_str());
281 B2RESULT(
"Failure to calibrate time-walk for " << nfailure <<
" board");
287 B2INFO(
"Prepare calibration");
293 B2FATAL(
"Mode " <<
m_twParamMode <<
" haven't implemented yet.");
297 const int nEntries = dbTw->getEntries();
298 for (
int ib = 0; ib < nEntries; ++ib) {
299 m_tw_old[ib] = dbTw->getTimeWalkParams(ib);
306 B2INFO(
"Add constant term into T0 database");
310 for (
int ilay = 0; ilay < 56; ++ilay) {
312 for (
unsigned int iwire = 0; iwire < nW; ++iwire) {
313 WireID wireid(ilay, iwire);
315 T0 = cdcgeo.
getT0(wireid);
319 tz->
setT0(wireid, T0);
332 int max = h1->GetMaximumBin();
333 double maxX = h1->GetBinCenter(max);
334 double maxY = h1->GetBinContent(max);
335 B2DEBUG(21,
"Max: id - x - y : " << max <<
" " << maxX <<
" " << maxY);
339 h1->Fit(
"pol0",
"MQE",
"", maxX + 125, 500);
340 if (h1->GetFunction(
"pol0")) {
341 p0 = h1->GetFunction(
"pol0")->GetParameter(0);
345 TH1D* hshift =
new TH1D(
"hshift",
"shift", h1->GetNbinsX(), 0, 500);
346 hshift->SetDirectory(0);
347 for (
int i = 0; i <= h1->GetNbinsX(); ++i) {
348 hshift->SetBinContent(i, h1->GetBinContent(i) - p0);
350 hshift->Fit(
"expo",
"MQE",
"", 0, 400);
351 double p1 = maxY + p0;
353 if (hshift->GetFunction(
"expo")) {
354 p1 = exp(hshift->GetFunction(
"expo")->GetParameter(0));
355 p2 = hshift->GetFunction(
"expo")->GetParameter(1);
359 TF1* func =
new TF1(
"func",
"[0]+ [1]*exp(-1*[2]*x)", 0, 500);
360 func->SetParameters(p0, p1, -1 * p2);
361 func->SetParLimits(0, -5, 5);
362 func->SetParLimits(1, -5, 500);
363 func->SetParLimits(2, -0.01, 0.1);
364 h1->Fit(
"func",
"MQE",
"", 10, 400);
Database object for time-walk.
void outputToFile(std::string fileName) const
Output the contents in text file format.
void setTwParamMode(unsigned short mode)
Set tw parameterization mode mode=0: tw (in ns) = p0/sqrt(FADCsum); mode=1: tw( in ns) = p0*exp(-p1*F...
void setTimeWalkParams(unsigned short boardID, const std::vector< float > ¶ms)
Set the time-walk coefficients in the list.
Database object for timing offset (t0).
void outputToFile(std::string fileName) const
Output the contents in text file format.
void setT0(unsigned short iCLayer, unsigned short iWire, float t0)
Set t0 in the list.
The Class for CDC Geometry Parameters.
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
float getT0(const WireID &wireID) const
Returns t0 parameter of the specified sense wire.
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 senseWireR(int layerId) const
Returns radius of sense wire in each layer.
double m_xmin
minimum value cut of drift length.
void prepare()
prepare calibration.
std::string m_outputT0FileName
t0 file name after calibration.
TH2D * m_h2[300]
2D histogram of residual vs ADC for each board.
unsigned short m_twParamMode
=0 for P0/Sqrt(ADC); =1 for P0*Exp(-P1*ADC).
void storeHist()
Store histograms.
std::string m_histName
root file name
double m_constTerm[300]
const term in fitting, it will be added to T0 instead tw
bool m_storeHisto
Store all Histogram or not.
TimeWalkCalibrationAlgorithm()
Constructor.
double m_minNdf
minimum number of degree of freedom required for track.
unsigned short m_flag[300]
flag for fit status
void createHisto()
Create histo for calibrate.
double m_minPval
minimum number of Prob(chi2) of fitted track.
void updateT0()
update constant term to t0 database.
std::vector< float > m_tw_old[300]
tw list before calibration.
DBObjPtr< CDCGeometry > m_cdcGeo
Geometry of CDC.
std::string m_outputTWFileName
tw file name after calibration.
std::vector< float > m_tw_new[300]
tw list after calibration.
void write()
save calibration.
unsigned short m_nTwParams
No.
void fitToExponentialFunc(TH1D *h1)
fit tw histogram
bool m_textOutput
output text file if true
EResult calibrate() override
Run algo on data.
EResult checkConvergence()
check convergence.
TH1D * m_h1[300]
Mean of residual as function of ADC of each board.
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
@ c_Iterate
Needs iteration =1 in Python.
Class for accessing objects in the database.
Class to identify a wire inside the CDC.
Abstract base class for different kinds of events.