9#include <cdc/calibration/SpaceResolutionCalibrationAlgorithm.h>
10#include <cdc/geometry/CDCGeometryPar.h>
11#include <cdc/dbobjects/CDCSpaceResols.h>
13#include <framework/database/DBObjPtr.h>
14#include <framework/logging/Logger.h>
22#include <TStopwatch.h>
28typedef std::array<float, 3> array3;
33 " -------------------------- Position Resolution Calibration Algorithm -------------------------\n"
38 B2INFO(
"Creating histograms");
43 for (
int i = 0; i < 50; ++i) {
44 yb.push_back(-0.07 + i * (0.14 / 50));
46 for (
int i = 0; i < 50; ++i) {
47 yu.push_back(-0.08 + i * (0.16 / 50));
53 for (
int i = 1; i < np; ++i) {
57 for (
int il = 0; il < 56; ++il) {
58 for (
int lr = 0; lr < 2; ++lr) {
61 m_hBiased[il][lr][al][th] =
new TH2F(Form(
"hb_%d_%d_%d_%d", il, lr, al, th),
62 Form(
"lay_%d_lr%d_al_%3.0f_th_%3.0f;Drift Length [cm];#DeltaX", il, lr,
m_iAlpha[al],
m_iTheta[th]),
63 xbin.size() - 1, &xbin.at(0), yb.size() - 1, &yb.at(0));
64 m_hUnbiased[il][lr][al][th] =
new TH2F(Form(
"hu_%d_%d_%d_%d", il, lr, al, th),
65 Form(
"lay_%d_lr%d_al_%3.0f_th_%3.0f;Drift Length [cm];#DeltaX", il, lr,
m_iAlpha[al],
m_iTheta[th]),
66 xbin.size() - 1, &xbin.at(0), yu.size() - 1, &yu.at(0));
73 auto tree = getObjectPtr<TTree>(
"tree");
86 tree->SetBranchAddress(
"lay", &lay);
87 tree->SetBranchAddress(
"ndf", &ndf);
88 tree->SetBranchAddress(
"Pval", &Pval);
89 tree->SetBranchAddress(
"x_u", &x_u);
90 tree->SetBranchAddress(
"x_b", &x_b);
91 tree->SetBranchAddress(
"x_mea", &x_mea);
92 tree->SetBranchAddress(
"weight", &w);
93 tree->SetBranchAddress(
"alpha", &alpha);
94 tree->SetBranchAddress(
"theta", &theta);
97 std::vector<TString> list_vars = {
"lay",
"ndf",
"Pval",
"x_u",
"x_b",
"x_mea",
"weight",
"alpha",
"theta"};
98 tree->SetBranchStatus(
"*", 0);
100 for (TString brname : list_vars) {
101 tree->SetBranchStatus(brname, 1);
105 const Long64_t nEntries = tree->GetEntries();
106 B2INFO(
"Number of entries: " << nEntries);
111 for (Long64_t i = 0; i < nEntries; ++i) {
113 if (std::fabs(x_b) < 0.02 || std::fabs(x_u) < 0.02)
continue;
130 int ilr = x_u > 0 ? 1 : 0;
132 if (ial == -99 || ith == -99) {
133 TString command = Form(
"Error in alpha=%3.2f and theta = %3.2f>> error", alpha, theta);
134 B2FATAL(
"ERROR" << command);
137 absRes_u = fabs(x_mea) - fabs(x_u);
138 absRes_b = fabs(x_mea) - fabs(x_b);
140 int ilay =
static_cast<int>(lay);
141 m_hUnbiased[ilay][ilr][ial][ith]->Fill(fabs(x_u), absRes_u, w);
142 m_hBiased[ilay][ilr][ial][ith]->Fill(fabs(x_b), absRes_b, w);
146 B2INFO(
"Time to fill histograms: " << timer.RealTime() <<
"s");
148 B2INFO(
"Start to obtain the biased and unbiased sigmas...");
149 TF1* gb =
new TF1(
"gb",
"gaus", -0.05, 0.05);
150 TF1* gu =
new TF1(
"gu",
"gaus", -0.06, 0.06);
151 TF1* g0b =
new TF1(
"g0b",
"gaus", -0.015, 0.07);
152 TF1* g0u =
new TF1(
"g0u",
"gaus", -0.015, 0.08);
154 std::vector<double> sigma;
155 std::vector<double> dsigma;
156 std::vector<double> s2;
157 std::vector<double> ds2;
158 std::vector<double> xl;
159 std::vector<double> dxl;
160 std::vector<double> dxl0;
162 ofstream ofss(
"IntReso.dat");
166 for (
int il = 0; il < 56; ++il) {
167 for (
int lr = 0; lr < 2; ++lr) {
171 B2DEBUG(21,
"layer-lr-al-th " << il <<
" - " << lr <<
" - " << al <<
" - " << th);
172 if (
m_hBiased[il][lr][al][th]->GetEntries() < 5000) {
177 auto* proYb =
m_hBiased[il][lr][al][th]->ProjectionY();
178 auto* proYu =
m_hUnbiased[il][lr][al][th]->ProjectionY();
180 g0b->SetParLimits(0, 0,
m_hBiased[il][lr][al][th]->GetEntries() * 5);
181 g0u->SetParLimits(0, 0,
m_hUnbiased[il][lr][al][th]->GetEntries() * 5);
182 g0b->SetParLimits(1, -0.01, 0.004);
183 g0u->SetParLimits(1, -0.01, 0.004);
184 g0b->SetParLimits(2, 0.0, proYb->GetRMS() * 5);
185 g0u->SetParLimits(2, 0.0, proYu->GetRMS() * 5);
187 g0b->SetParameter(0,
m_hBiased[il][lr][al][th]->GetEntries());
188 g0u->SetParameter(0,
m_hUnbiased[il][lr][al][th]->GetEntries());
189 g0b->SetParameter(1, 0);
190 g0u->SetParameter(1, 0);
191 g0b->SetParameter(2, proYb->GetRMS());
192 g0u->SetParameter(2, proYu->GetRMS());
194 B2DEBUG(21,
"Nentries: " <<
m_hBiased[il][lr][al][th]->GetEntries());
195 m_hBiased[il][lr][al][th]->SetDirectory(0);
201 m_hBiased[il][lr][al][th]->FitSlicesY(g0b, firstbin, ib1, minEntry);
204 m_hMeanBiased[il][lr][al][th] = (TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_1", il, lr, al, th))->Clone(Form(
"hb_%d_%d_%d_%d_m", il,
208 m_hSigmaBiased[il][lr][al][th] = (TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_2", il, lr, al, th))->Clone(Form(
"hb_%d_%d_%d_%d_s",
215 m_hBiased[il][lr][al][th]->FitSlicesY(gb, ib1 + 1, np, minEntry);
217 m_hMeanBiased[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_1", il, lr, al, th)));
219 m_hSigmaBiased[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_2", il, lr, al, th)));
220 B2DEBUG(21,
"entries (2nd): " <<
m_hSigmaBiased[il][lr][al][th]->GetEntries());
225 m_hUnbiased[il][lr][al][th]->FitSlicesY(g0u, firstbin, ib1, minEntry);
227 m_hMeanUnbiased[il][lr][al][th] = (TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_1", il, lr, al, th))->Clone(Form(
"hu_%d_%d_%d_%d_m",
231 m_hSigmaUnbiased[il][lr][al][th] = (TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_2", il, lr, al, th))->Clone(Form(
"hu_%d_%d_%d_%d_s",
239 m_hUnbiased[il][lr][al][th]->FitSlicesY(gu, ib1 + 1, np, minEntry);
241 m_hMeanUnbiased[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_1", il, lr, al, th)));
243 m_hSigmaUnbiased[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_2", il, lr, al, th)));
245 B2WARNING(
"sliced histo not found");
261 if (
m_hSigmaBiased[il][lr][al][th]->GetBinContent(j) == 0)
continue;
267 double XL =
m_hSigmaBiased[il][lr][al][th]->GetXaxis()->GetBinCenter(j);
268 double dXL = (
m_hSigmaBiased[il][lr][al][th]->GetXaxis()->GetBinWidth(j)) / 2;
269 double s_int = std::sqrt(sb * su);
270 double ds_int = 0.5 * s_int * (dsb / sb + dsu / su);
271 if (ds_int > 0.02)
continue;
275 sigma.push_back(s_int);
276 dsigma.push_back(ds_int);
277 s2.push_back(s_int * s_int);
278 ds2.push_back(2 * s_int * ds_int);
279 ofss << il <<
" " << lr <<
" " << al <<
" " << th <<
" " << j <<
" " << XL <<
" " << dXL <<
" " << s_int <<
" " <<
283 if (xl.size() < 7 || xl.size() >
Max_np) {
285 B2WARNING(
"number of element might out of range");
continue;
289 B2DEBUG(21,
"Create Histo for layer-lr: " << il <<
" " << lr);
290 m_graph[il][lr][al][th] =
new TGraphErrors(xl.size(), &xl.at(0), &sigma.at(0), &dxl.at(0), &dsigma.at(0));
291 m_graph[il][lr][al][th]->SetMarkerSize(0.5);
292 m_graph[il][lr][al][th]->SetMarkerStyle(8);
293 m_graph[il][lr][al][th]->SetTitle(Form(
"Layer_%d lr%d #alpha = %3.0f #theta = %3.0f", il, lr,
m_iAlpha[al],
m_iTheta[th]));
294 m_graph[il][lr][al][th]->SetName(Form(
"lay%d_lr%d_al%d_th%d", il, lr, al, th));
297 m_gFit[il][lr][al][th] =
new TGraphErrors(xl.size(), &xl.at(0), &s2.at(0), &dxl0.at(0), &ds2.at(0));
298 m_gFit[il][lr][al][th]->SetMarkerSize(0.5);
299 m_gFit[il][lr][al][th]->SetMarkerStyle(8);
300 m_gFit[il][lr][al][th]->SetTitle(Form(
"L%d lr%d #alpha = %3.0f #theta = %3.0f ", il, lr,
m_iAlpha[al],
m_iTheta[th]));
301 m_gFit[il][lr][al][th]->SetName(Form(
"sigma2_lay%d_lr%d_al%d_th%d", il, lr, al, th));
303 gDirectory->Delete(
"hu_%d_%d_%d_%d_0");
315 B2INFO(
"Start calibration");
316 gPrintViaErrorHandler =
true;
318 gErrorIgnoreLevel = 3001;
321 B2INFO(
"ExpRun used for DB Geometry : " << exprun.first <<
" " << exprun.second);
329 TF1* func =
new TF1(
"func",
"[0]/(x*x + [1])+[2]* x+[3]+[4]*exp([5]*(x-[6])*(x-[6]))", 0, 1.);
330 TH1F* hprob =
new TH1F(
"h1",
"", 20, 0, 1);
334 for (
int i = 0; i < 56; ++i) {
335 for (
int lr = 0; lr < 2; ++lr) {
338 if (!
m_gFit[i][lr][al][th])
continue;
342 B2DEBUG(199,
"xmax for fitting: " << upFit);
344 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-5, 0.00008, -30, intp6);
345 func->SetParLimits(0, 1
E-7, 1
E-4);
346 func->SetParLimits(1, 0.0045, 0.02);
347 func->SetParLimits(2, 1
E-6, 0.0005);
348 func->SetParLimits(3, 1
E-8, 0.0005);
349 func->SetParLimits(4, 0., 0.001);
350 func->SetParLimits(5, -40, 0.);
351 func->SetParLimits(6, intp6 - 0.5, intp6 + 0.2);
353 B2DEBUG(21,
"Fitting for layer: " << i <<
"lr: " << lr <<
" ial" << al <<
" ith:" << th);
354 B2DEBUG(21,
"Fit status before fit:" <<
m_fitStatus[i][lr][al][th]);
356 for (
int j = 0; j < 10; j++) {
358 B2DEBUG(21,
"loop: " << j);
359 B2DEBUG(21,
"Int p6: " << intp6);
360 B2DEBUG(21,
"Number of Point: " <<
m_gFit[i][lr][al][th]->GetN());
361 Int_t stat =
m_gFit[i][lr][al][th]->Fit(
"func",
"MQE",
"", 0.05, upFit);
362 B2DEBUG(21,
"stat of fit" << stat);
363 std::string Fit_status = gMinuit->fCstatu.Data();
364 B2DEBUG(21,
"FIT STATUS: " << Fit_status);
365 if (Fit_status ==
"OK" || Fit_status ==
"SUCCESSFUL" || Fit_status ==
"CALL LIMIT"
366 || Fit_status ==
"PROBLEMS") {
367 if (fabs(func->Eval(0.3)) > 0.00035 || func->Eval(0.3) < 0) {
368 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-7, 0.0007, -30, intp6 + 0.05 * j);
369 func->SetParLimits(6, intp6 + 0.05 * j - 0.5, intp6 + 0.05 * j + 0.2);
373 B2DEBUG(21,
"Prob of fit: " << func->GetProb());
379 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-7, 0.0007, -30, intp6 + 0.05 * j);
380 func->SetParLimits(6, intp6 + 0.05 * j - 0.5, intp6 + 0.05 * j + 0.2);
391 B2DEBUG(21,
"ProbFit: Lay_lr_al_th: " << i <<
" " << lr <<
" " << al <<
" " << th << func->GetProb());
392 hprob->Fill(func->GetProb());
393 func->GetParameters(
m_sigma[i][lr][al][th]);
405 int nFitCompleted = 0;
406 for (
int l = 0; l < 56; ++l) {
407 for (
int lr = 0; lr < 2; ++lr) {
418 if (
static_cast<double>(nFitCompleted) / nTotal <
m_threshold) {
419 B2WARNING(
"Less than " <<
m_threshold * 100 <<
" % of Sigmas were fitted.");
428 B2INFO(
"saving histograms");
430 TFile* ff =
new TFile(
m_histName.c_str(),
"RECREATE");
431 TDirectory* top = gDirectory;
432 TDirectory* Direct[56];
434 auto hNDF = getObjectPtr<TH1F>(
"hNDF");
435 auto hPval = getObjectPtr<TH1F>(
"hPval");
436 auto hEvtT0 = getObjectPtr<TH1F>(
"hEventT0");
438 if (hNDF && hPval && hEvtT0) {
445 for (
int il = 0; il < 56; ++il) {
447 Direct[il] = gDirectory->mkdir(Form(
"lay_%d", il));
450 for (
int lr = 0; lr < 2; ++lr) {
453 if (!
m_graph[il][lr][al][th])
continue;
454 if (!
m_gFit[il][lr][al][th])
continue;
462 m_graph[il][lr][al][th]->Write();
463 m_gFit[il][lr][al][th]->Write();
470 B2INFO(
"Finish store histogram");
475 B2INFO(
"Writing calibrated sigma's");
481 const float deg2rad = M_PI / 180.0;
484 std::array<float, 3> alpha3 = {
m_lowerAlpha[i]* deg2rad,
493 std::array<float, 3> theta3 = {
m_lowerTheta[i]* deg2rad,
503 for (
int iCL = 0; iCL < 56; ++iCL) {
504 for (
int iLR = 1; iLR >= 0; --iLR) {
505 std::vector<float> sgbuff;
508 for (
int i = 0; i < 7; ++i) {
509 sgbuff.push_back(
m_sigma[iCL][iLR][ialpha][itheta][i]);
515 for (
int i = 0; i < 7; ++i) {
516 sgbuff.push_back(
m_sigmaPost[iCL][iLR][ialpha][itheta][i]);
532 B2RESULT(
"Histos succesfully fitted: " << nfitted);
533 B2RESULT(
"Histos fit failure: " << nfailure);
540 B2INFO(
"Prepare calibration of space resolution");
542 const double rad2deg = 180 / M_PI;
551 array3 alpha = dbSigma->getAlphaBin(i);
559 array3 theta = dbSigma->getThetaBin(i);
566 for (
unsigned short iCL = 0; iCL < 56; ++iCL) {
567 for (
unsigned short iLR = 0; iLR < 2; ++iLR) {
570 const std::vector<float> params = dbSigma->getSigmaParams(iCL, iLR, iA, iT);
571 unsigned short np = params.size();
573 for (
unsigned short i = 0; i < np; ++i) {
Database object for space resolutions.
void setSigmaParams(const SigmaID sigmaID, const std::vector< float > ¶ms)
Set sigma parameters for the specified id.
void setThetaBin(const array3 &theta)
Set theta-angle bin (rad)
void outputToFile(std::string fileName) const
Output the contents in text file format.
void setAlphaBin(const array3 &alpha)
Set alpha-angle bin (rad)
void setSigmaParamMode(unsigned short mode)
Set sigma parameterization mode.
void storeHisto()
store histogram
TH1F * m_hSigmaBiased[56][2][Max_nalpha][Max_ntheta]
sigma histogram of biased residual
void prepare()
Prepare the calibration of space resolution.
float m_lowerTheta[7]
Lower boundays of theta bins.
std::string m_histName
root file name
double m_threshold
minimal requirement for the fraction of fitted results
double m_minNdf
Minimum NDF
unsigned short m_sigmaParamMode
sigma mode for this calibration.
void createHisto()
create histogram
SpaceResolutionCalibrationAlgorithm()
Constructor.
unsigned short m_sigmaParamModePost
sigma mode before this calibration.
TH1F * m_hMeanBiased[56][2][Max_nalpha][Max_ntheta]
mean histogram biased residual
double m_sigmaPost[56][2][18][7][8]
sigma prameters before calibration
int m_nAlphaBins
number of alpha bins
int m_nThetaBins
number of theta bins
double m_minPval
Minimum Prob(chi2) of track.
double m_binWidth
width of each bin, unit cm
float m_iAlpha[18]
represented alphas of alpha bins.
double getUpperBoundaryForFit(TGraphErrors *graph)
search max point at boundary region
TGraphErrors * m_graph[56][2][18][7]
sigma graph.
float m_lowerAlpha[18]
Lower boundays of alpha bins.
TH2F * m_hBiased[56][2][Max_nalpha][Max_ntheta]
2D histogram of biased residual
static const unsigned short Max_np
Maximum number of point =1/binwidth.
double m_sigma[56][2][18][7][8]
new sigma prameters.
TH1F * m_hSigmaUnbiased[56][2][Max_nalpha][Max_ntheta]
sigma histogram of ubiased residual
void write()
save calibration, in text file or db
float m_upperAlpha[18]
Upper boundays of alpha bins.
TH1F * m_hMeanUnbiased[56][2][Max_nalpha][Max_ntheta]
mean histogram of unbiased residual
bool m_textOutput
output text file if true
EResult calibrate() override
Run algo on data.
float m_upperTheta[7]
Upper boundays of theta bins.
TGraphErrors * m_gFit[56][2][18][7]
sigma*sigma graph for fit
int m_fitStatus[56][2][Max_nalpha][Max_ntheta]
Fit flag; 1:OK ; 0:error.
std::string m_outputFileName
Output sigma filename.
float m_iTheta[7]
represented alphas of theta bins.
TH2F * m_hUnbiased[56][2][Max_nalpha][Max_ntheta]
2D histogram of unbiased residual
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_NotEnoughData
Needs more data =2 in Python.
Class for accessing objects in the database.
Abstract base class for different kinds of events.