10#include <cdc/calibration/SpaceResolutionCalibration.h>
11#include <cdc/geometry/CDCGeometryPar.h>
12#include <cdc/calibration/CDCDatabaseImporter.h>
14#include <framework/database/DBObjPtr.h>
15#include <framework/logging/Logger.h>
16#include <framework/utilities/FileSystem.h>
40 B2INFO(
"createHisto");
45 TChain* tree =
new TChain(
"tree");
48 if (!tree->GetBranch(
"ndf")) {
49 B2FATAL(
"input data do not exits, please check!");
50 gSystem->Exec(
"echo rootfile do not exits or something wrong >> error");
65 tree->SetBranchAddress(
"lay", &lay);
66 tree->SetBranchAddress(
"ndf", &ndf);
67 tree->SetBranchAddress(
"Pval", &Pval);
68 tree->SetBranchAddress(
"x_u", &x_u);
69 tree->SetBranchAddress(
"x_b", &x_b);
70 tree->SetBranchAddress(
"x_mea", &x_mea);
71 tree->SetBranchAddress(
"weight", &w);
72 tree->SetBranchAddress(
"alpha", &alpha);
73 tree->SetBranchAddress(
"theta", &theta);
76 std::vector<TString> list_vars = {
"lay",
"ndf",
"Pval",
"x_u",
"x_b",
"x_mea",
"weight",
"alpha",
"theta"};
77 tree->SetBranchStatus(
"*", 0);
79 for (TString brname : list_vars) {
80 tree->SetBranchStatus(brname, 1);
86 for (
int i = 0; i < 50; ++i) {
87 yb.push_back(-0.07 + i * (0.14 / 50));
89 for (
int i = 0; i < 50; ++i) {
90 yu.push_back(-0.08 + i * (0.16 / 50));
96 for (
int i = 1; i < m_np; ++i) {
100 for (
int il = 0; il < 56; ++il) {
101 for (
int lr = 0; lr < 2; ++lr) {
102 for (
int al = 0; al <
m_nalpha; ++al) {
103 for (
int th = 0; th <
m_ntheta; ++th) {
104 hist_b[il][lr][al][th] =
new TH2F(Form(
"hb_%d_%d_%d_%d", il, lr, al, th),
105 Form(
"lay_%d_lr%d_al_%3.0f_th_%3.0f;Drift Length [cm];#DeltaX", il, lr,
ialpha[al],
itheta[th]),
106 xbin.size() - 1, &xbin.at(0), yb.size() - 1, &yb.at(0));
107 hist_u[il][lr][al][th] =
new TH2F(Form(
"hu_%d_%d_%d_%d", il, lr, al, th),
108 Form(
"lay_%d_lr%d_al_%3.0f_th_%3.0f;Drift Length [cm];#DeltaX", il, lr,
ialpha[al],
itheta[th]),
109 xbin.size() - 1, &xbin.at(0), yu.size() - 1, &yu.at(0));
116 const int nEntries = tree->GetEntries();
117 B2INFO(
"Number of entries: " << nEntries);
121 for (
int i = 0; i < nEntries; ++i) {
124 if (std::fabs(x_b) < 0.02 || std::fabs(x_u) < 0.02)
continue;
127 for (
int k = 0; k <
m_nalpha; ++k) {
134 for (
int j = 0; j <
m_ntheta; ++j) {
141 ilr = x_u > 0 ? 1 : 0;
143 if (ial == -99 || ith == -99 || ilr == -99) {
144 TString command = Form(
"Error in alpha=%3.2f and theta = %3.2f>> error", alpha, theta);
146 B2FATAL(
"ERROR" << command);
148 absRes_u = fabs(x_mea) - fabs(x_u);
149 absRes_b = fabs(x_mea) - fabs(x_b);
151 hist_u[lay][ilr][ial][ith]->Fill(fabs(x_u), absRes_u, w);
152 hist_b[lay][ilr][ial][ith]->Fill(fabs(x_b), absRes_b, w);
156 B2INFO(
"Finish reading data");
158 TF1* gb =
new TF1(
"gb",
"gaus", -0.05, 0.05);
159 TF1* gu =
new TF1(
"gu",
"gaus", -0.06, 0.06);
160 TF1* g0b =
new TF1(
"g0b",
"gaus", -0.015, 0.07);
161 TF1* g0u =
new TF1(
"g0u",
"gaus", -0.015, 0.08);
162 g0b->SetParLimits(1, -0.01, 0.004);
163 g0u->SetParLimits(1, -0.01, 0.004);
165 std::vector<double> sigma;
166 std::vector<double> dsigma;
167 std::vector<double> s2;
168 std::vector<double> ds2;
169 std::vector<double> xl;
170 std::vector<double> dxl;
171 std::vector<double> dxl0;
176 for (
int il = 0; il < 56; ++il) {
177 for (
int lr = 0; lr < 2; ++lr) {
178 for (
int al = 0; al <
m_nalpha; ++al) {
179 for (
int th = 0; th <
m_ntheta; ++th) {
181 B2DEBUG(199,
"layer-lr-al-th " << il <<
" - " << lr <<
" - " << al <<
" - " << th);
182 if (
hist_b[il][lr][al][th]->GetEntries() < 20000) {
186 B2DEBUG(199,
"Nentries: " <<
hist_b[il][lr][al][th]->GetEntries());
187 hist_b[il][lr][al][th]->SetDirectory(0);
188 hist_u[il][lr][al][th]->SetDirectory(0);
190 hist_b[il][lr][al][th]->FitSlicesY(g0b, firstbin, ib1, minEntry);
192 TH1F* hm1 = (TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_1", il, lr, al, th));
193 TH1F* hs1 = (TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_2", il, lr, al, th));
199 hb_m[il][lr][al][th] = (TH1F*)hm1->Clone(Form(
"hb_%d_%d_%d_%d_m", il, lr, al, th));
200 hb_s[il][lr][al][th] = (TH1F*)hs1->Clone(Form(
"hb_%d_%d_%d_%d_s", il, lr, al, th));
202 hb_m[il][lr][al][th]->SetDirectory(0);
203 hb_s[il][lr][al][th]->SetDirectory(0);
207 hist_b[il][lr][al][th]->FitSlicesY(gb, ib1 + 1, m_np, minEntry);
208 hb_m[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_1", il, lr, al, th)));
209 hb_s[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hb_%d_%d_%d_%d_2", il, lr, al, th)));
210 B2DEBUG(199,
"entries (2nd): " <<
hb_s[il][lr][al][th]->GetEntries());
213 hist_u[il][lr][al][th]->FitSlicesY(g0u, firstbin, ib1, minEntry);
214 hu_m[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", il, lr, al,
216 hu_s[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", il, lr, al,
218 hu_m[il][lr][al][th]->SetDirectory(0);
219 hu_s[il][lr][al][th]->SetDirectory(0);
222 hist_u[il][lr][al][th]->FitSlicesY(gu, ib1 + 1, m_np, minEntry);
223 hu_m[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_1", il, lr, al, th)));
224 hu_s[il][lr][al][th]->Add((TH1F*)gDirectory->Get(Form(
"hu_%d_%d_%d_%d_2", il, lr, al, th)));
225 if (!
hu_s[il][lr][al][th] || !
hb_s[il][lr][al][th]) {
226 B2WARNING(
"slice histo do not found");
231 xl.clear(); dxl.clear(); dxl0.clear(); sigma.clear(); dsigma.clear(); s2.clear(); ds2.clear();
232 for (Int_t j = 1; j <
hu_s[il][lr][al][th]->GetNbinsX(); j++) {
233 if (
hu_s[il][lr][al][th]->GetBinContent(j) == 0)
continue;
234 if (
hb_s[il][lr][al][th]->GetBinContent(j) == 0)
continue;
235 double sb =
hb_s[il][lr][al][th]->GetBinContent(j);
236 double su =
hu_s[il][lr][al][th]->GetBinContent(j);
238 double dsb =
hb_s[il][lr][al][th]->GetBinError(j);
239 double dsu =
hu_s[il][lr][al][th]->GetBinError(j);
240 double XL =
hb_s[il][lr][al][th]->GetXaxis()->GetBinCenter(j);
241 double dXL = (
hb_s[il][lr][al][th]->GetXaxis()->GetBinWidth(j)) / 2;
242 double s_int = std::sqrt(sb * su);
243 double ds_int = 0.5 * s_int * (dsb / sb + dsu / su);
244 if (ds_int > 0.02)
continue;
248 sigma.push_back(s_int);
249 dsigma.push_back(ds_int);
250 s2.push_back(s_int * s_int);
251 ds2.push_back(2 * s_int * ds_int);
254 if (xl.size() < 7 || xl.size() >
Max_np) {
256 B2WARNING(
"number of element might out of range");
continue;
260 B2DEBUG(199,
"Create Histo for layer-lr: " << il <<
" " << lr);
261 gr[il][lr][al][th] =
new TGraphErrors(xl.size(), &xl.at(0), &sigma.at(0), &dxl.at(0), &dsigma.at(0));
262 gr[il][lr][al][th]->SetMarkerSize(0.5);
263 gr[il][lr][al][th]->SetMarkerStyle(8);
266 gr[il][lr][al][th]->SetTitle(Form(
"Layer_%d_lr%d | #alpha = %3.0f | #theta = %3.0f", il, lr,
ialpha[al],
itheta[th]));
267 gr[il][lr][al][th]->SetName(Form(
"lay%d_lr%d_al%d_th%d", il, lr, al, th));
270 gfit[il][lr][al][th] =
new TGraphErrors(xl.size(), &xl.at(0), &s2.at(0), &dxl0.at(0), &ds2.at(0));
271 gfit[il][lr][al][th]->SetMarkerSize(0.5);
272 gfit[il][lr][al][th]->SetMarkerStyle(8);
275 gfit[il][lr][al][th]->SetTitle(Form(
"L%d-lr%d | #alpha = %3.0f | #theta = %3.0f ", il, lr,
ialpha[al],
itheta[th]));
276 gfit[il][lr][al][th]->SetName(Form(
"sigma2_lay%d_lr%d_al%d_th%d", il, lr, al, th));
279 gDirectory->Delete(
"hu_%d_%d_%d_%d_0");
289 gErrorIgnoreLevel = 3001;
291 B2INFO(
"Start to calibrate");
292 gSystem->Exec(
"mkdir -p Sigma_Fit_Err");
294 TF1* func =
new TF1(
"func",
"[0]/(x*x + [1])+[2]* x+[3]+[4]*exp([5]*(x-[6])*(x-[6]))", 0, 1.);
295 TH1F* hprob =
new TH1F(
"h1",
"", 20, 0, 1);
299 for (
int i = 0; i < 56; ++i) {
300 for (
int lr = 0; lr < 2; ++lr) {
301 for (
int al = 0; al <
m_nalpha; ++al) {
302 for (
int th = 0; th <
m_ntheta; ++th) {
308 B2DEBUG(199,
"xmax for fitting: " << upFit);
310 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-5, 0.00008, -30, intp6);
311 func->SetParLimits(0, 1
E-7, 1
E-4);
312 func->SetParLimits(1, 0.00001, 0.02);
313 func->SetParLimits(2, 1
E-6, 0.0005);
314 func->SetParLimits(3, 1
E-8, 0.0005);
315 func->SetParLimits(4, 0., 0.001);
316 func->SetParLimits(5, -40, 0.);
317 func->SetParLimits(6, intp6 - 0.5, intp6 + 0.3);
318 B2DEBUG(199,
"FITTING for layer: " << i <<
"lr: " << lr <<
" ial" << al <<
" ith:" << th);
319 B2DEBUG(199,
"Fit flag before fit:" <<
m_fitflag[i][lr][al][th]);
322 for (
int j = 0; j < 10; j++) {
324 B2DEBUG(199,
"loop: " << j);
325 B2DEBUG(199,
"Int p6: " << intp6);
326 B2DEBUG(199,
"Number of Point: " <<
gfit[i][lr][al][th]->GetN());
327 Int_t stat =
gfit[i][lr][al][th]->Fit(
"func",
"MQE",
"", 0.05, upFit);
328 B2DEBUG(199,
"stat of fit" << stat);
329 std::string Fit_status = gMinuit->fCstatu.Data();
330 B2DEBUG(199,
"FIT STATUS: " << Fit_status);
332 if (Fit_status ==
"OK" ||
333 Fit_status ==
"SUCCESSFUL" ||
334 Fit_status ==
"CALL LIMIT" ||
335 Fit_status ==
"PROBLEMS") {
336 if (fabs(func->Eval(0.3)) > 0.00035 || func->Eval(0.3) < 0) {
337 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-7, 0.0007, -30, intp6 + 0.05 * j);
341 B2DEBUG(199,
"Prob of fit: " << func->GetProb());
347 func->SetParameters(5
E-6, 0.007, 1
E-4, 1
E-7, 0.0007, -30, intp6 + 0.05 * j);
350 TCanvas* c1 =
new TCanvas(
"c1",
"", 600, 600);
351 gfit[i][lr][al][th]->Draw();
352 c1->SaveAs(Form(
"Sigma_Fit_Err/%d_%d_al%d_th%d_%s.png", i, lr, al, th, Fit_status.c_str()));
353 B2WARNING(
"Fit error: " << i <<
" " << lr <<
" " << al <<
" " << th);
358 B2DEBUG(199,
"ProbFit: Lay_lr_al_th: " << i <<
" " << lr <<
" " << al <<
" " << th << func->GetProb());
359 hprob->Fill(func->GetProb());
360 func->GetParameters(
sigma_new[i][lr][al][th]);
374 B2INFO(
"storeHisto");
375 TFile* ff =
new TFile(
"sigma_histo.root",
"RECREATE");
376 TDirectory* top = gDirectory;
377 TDirectory* Direct[56];
379 for (
int il = 0; il < 56; ++il) {
381 Direct[il] = gDirectory->mkdir(Form(
"lay_%d", il));
383 for (
int lr = 0; lr < 2; ++lr) {
384 for (
int al = 0; al <
m_nalpha; ++al) {
385 for (
int th = 0; th <
m_ntheta; ++th) {
386 if (!
gr[il][lr][al][th])
continue;
387 if (!
gfit[il][lr][al][th])
continue;
389 hist_b[il][lr][al][th]->Write();
390 hist_u[il][lr][al][th]->Write();
391 hb_m[il][lr][al][th]->Write();
392 hb_s[il][lr][al][th]->Write();
393 hu_m[il][lr][al][th]->Write();
394 hu_s[il][lr][al][th]->Write();
395 gr[il][lr][al][th]->Write();
396 gfit[il][lr][al][th]->Write();
403 B2INFO(
"Finish store histogram");
408 B2INFO(
"Exporting parameters...");
415 for (
int i = 0; i <
m_nalpha; ++i) {
416 ofs << std::setprecision(4) <<
l_alpha[i] <<
" " << std::setprecision(4) <<
u_alpha[i] <<
" " << std::setprecision(
421 for (
int i = 0; i <
m_ntheta; ++i) {
422 ofs << std::setprecision(4) <<
l_theta[i] <<
" " << std::setprecision(4) <<
u_theta[i] <<
" " << std::setprecision(
426 ofs << 0 <<
" " << 7 << endl;
427 for (
int al = 0; al <
m_nalpha; ++al) {
428 for (
int th = 0; th <
m_ntheta; ++th) {
429 for (
int i = 0; i < 56; ++i) {
430 for (
int lr = 1; lr >= 0; --lr) {
432 ofs << i << std::setw(4) <<
itheta[th] << std::setw(4) <<
ialpha[al] << std::setw(4) << lr << std::setw(15);
435 for (
int p = 0; p < 7; ++p) {
436 if (p != 6) { ofs << std::setprecision(7) <<
sigma_new[i][lr][al][th][p] << std::setw(15);}
437 if (p == 6) { ofs << std::setprecision(7) <<
sigma_new[i][lr][al][th][p] << std::endl;}
440 B2WARNING(
"Fitting error and old sigma will be used. (Layer " << i <<
") (lr = " << lr <<
") (al = " << al <<
") (th = " << th <<
451 for (
int p = 0; p < 7; ++p) {
452 if (p != 6) { ofs << std::setprecision(7) <<
sigma_old[i][lr][ial_old][ith_old][p] << std::setw(15);}
453 if (p == 6) { ofs << std::setprecision(7) <<
sigma_old[i][lr][ial_old][ith_old][p] << std::endl;}
462 B2RESULT(
"Histos succesfully fitted: " << nfitted);
463 B2RESULT(
"Histos fit failure: " << nfailure);
475 B2INFO(
"reading sigma from DB");
478 B2INFO(
"Number of theta bins from input sigma: " <<
ntheta_old);
480 B2INFO(
"Read sigma from text");
482 B2INFO(
"number of alpha bins from input sigma: " <<
nalpha_old);
491 if (fileName ==
"") {
494 if (fileName ==
"") {
495 cout <<
"CDCGeometryPar: " << fileName <<
" not exist!" << endl;
497 cout <<
"CDCGeometryPar: " << fileName <<
" exists." << endl;
499 if (!ifs) cout <<
"CDCGeometryPar: cannot open " << fileName <<
" !" << endl;
506 cout <<
"Fail to read alpha bins !" << endl;
return;
508 double alpha0, alpha1, alpha2;
509 for (
unsigned short i = 0; i <
nalpha_old; ++i) {
510 ifs >> alpha0 >> alpha1 >> alpha2;
518 if (
ntheta_old == 0 ||
ntheta_old > 7) cout <<
"CDCGeometryPar: fail to read theta bins !" << endl;
520 cout <<
"CDCGeometryPar: fail to read theta bins !" << endl;
522 double theta0, theta1, theta2;
523 for (
unsigned short i = 0; i <
ntheta_old; ++i) {
524 ifs >> theta0 >> theta1 >> theta2;
530 unsigned short np = 0;
531 unsigned short iCL, iLR;
539 const double epsi = 0.1;
541 ifs >> theta >> alpha >> iLR;
542 for (
int i = 0; i < np; ++i) {
547 for (
unsigned short i = 0; i <
ntheta_old; ++i) {
553 if (ith < 0) cout <<
"CDCGeometryPar: thetas in sigma.dat are inconsistent !" << endl;
556 for (
unsigned short i = 0; i <
nalpha_old; ++i) {
562 if (ial < 0) cout <<
"CDCGeometryPar: alphas in sigma.dat are inconsistent !" << endl;
564 for (
int i = 0; i < np; ++i) {
565 sigma_old[iCL][iLR][ial][ith][i] = sigma[i];
572 typedef std::array<float, 3> array3;
574 nalpha_old = (*m_sResolFromDB)->getNoOfAlphaBins();
575 double rad2deg = 180 / M_PI;
576 for (
unsigned short i = 0; i <
nalpha_old; ++i) {
578 array3 alpha = (*m_sResolFromDB)->getAlphaBin(i);
585 ntheta_old = (*m_sResolFromDB)->getNoOfThetaBins();
587 for (
unsigned short i = 0; i <
ntheta_old; ++i) {
589 array3 theta = (*m_sResolFromDB)->getThetaBin(i);
596 for (
unsigned short iCL = 0; iCL < c_maxNSenseLayers; ++iCL) {
597 for (
unsigned short iLR = 0; iLR < 2; ++iLR) {
598 for (
unsigned short iA = 0; iA <
nalpha_old; ++iA) {
599 for (
unsigned short iT = 0; iT <
ntheta_old; ++iT) {
600 const std::vector<float> params = (*m_sResolFromDB)->getSigmaParams(iCL, iLR, iA, iT);
601 unsigned short np = params.size();
603 for (
unsigned short i = 0; i < np; ++i) {
604 sigma_old[iCL][iLR][iA][iT][i] = params[i];
614 B2INFO(
"readProfile");
617 B2INFO(
"use Sigma bining from input Sigma");
620 B2INFO(
"Number of alpha bins: " <<
m_nalpha);
621 for (
int i = 0; i <
m_nalpha; ++i) {
625 B2INFO(
"Number of theta bins: " <<
m_ntheta);
626 for (
int i = 0; i <
m_ntheta; ++i) {
631 B2INFO(
"use Sigma bining from profile file");
636 double dumy1, dumy2, dumy3;
638 B2INFO(
"Number of alpha bins: " <<
m_nalpha);
640 for (
int i = 0; i <
m_nalpha; ++i) {
641 proxt >> dumy1 >> dumy2 >> dumy3;
646 B2INFO(
"Number of theta bins: " <<
m_ntheta);
648 for (
int i = 0; i <
m_ntheta; ++i) {
649 proxt >> dumy1 >> dumy2 >> dumy3;
654 B2INFO(
"Finish asssign sigma bining");
void importSigma(std::string fileName)
Import sigma table to the database.
virtual void storeHisto()
store histogram
double ialpha[18]
represented alphas of alpha bins.
std::string m_inputRootFileNames
Input root file names.
TH1F * hb_s[56][2][Max_nalpha][Max_ntheta]
sigma histogram of ubiased residual
int m_nalpha
number of alpha bins
TGraphErrors * gr[56][2][18][7]
sigma graph.
int m_fitflag[56][2][Max_nalpha][Max_ntheta]
Fit flag; 1:OK ; 0:error.
static const int Max_ntheta
maximum theta bin
TH2F * hist_u[56][2][Max_nalpha][Max_ntheta]
2D histogram of unbiased residual
bool m_useProfileFromInputSigma
Use binning from old sigma or new one form input.
TH1F * hu_m[56][2][Max_nalpha][Max_ntheta]
mean histogram biased residual
double u_alpha[18]
Upper boundays of alpha bins.
virtual void readProfile()
read sigma bining (alpha, theta bining)
double m_ndfmin
Minimum NDF
virtual void createHisto()
create histogram
double sigma_old[56][2][18][7][8]
old sigma prameters.
TGraphErrors * gfit[56][2][18][7]
sigma*sigma graph for fit
DBObjPtr< CDCSpaceResols > * m_sResolFromDB
Database for sigma.
std::string m_ProfileFileName
Profile file name.
double l_theta_old[7]
Lower boundays of theta bins from input.
double l_alpha[18]
Lower boundays of alpha bins.
virtual void readSigmaFromDB()
read sigma from DB
double m_binWidth
width of each bin, unit cm
double getUpperBoundaryForFit(TGraphErrors *graph)
search max point at boundary region
int m_ntheta
number of theta bins
TH2F * hist_b[56][2][Max_nalpha][Max_ntheta]
2D histogram of biased residual
int nalpha_old
number of alpha bins from input
std::string m_outputSigmaFileName
Output sigma file name.
TH1F * hu_s[56][2][Max_nalpha][Max_ntheta]
sigma histogram of biased residual
TH1F * hb_m[56][2][Max_nalpha][Max_ntheta]
mean histogram of unbiased residual
double itheta_old[7]
represented alphas of theta bins from input.
static const unsigned short Max_np
Maximum number of point =1/binwidth.
std::string m_sigmafile
Sigma file name, for text mode.
double l_alpha_old[18]
Lower boundays of alpha bins from input.
double u_theta_old[7]
Upper boundays of theta bins from input.
virtual bool calibrate()
Run algo on data.
virtual void write()
save calibration, in text file or db
double itheta[7]
represented alphas of theta bins.
double ialpha_old[18]
represented alphas of alpha bins from input.
SpaceResolutionCalibration()
Constructor.
double l_theta[7]
Lower boundays of theta bins.
double sigma_new[56][2][18][7][8]
new sigma prameters.
virtual void readSigmaFromText()
read sigma from text file
static const int Max_nalpha
Maximum alpha bin.
double m_Pvalmin
Minimum Prob(chi2) of track.
double u_alpha_old[18]
Upper boundays of alpha bins from input.
virtual void readSigma()
read sigma from previous calibration, (input sigma)
double u_theta[7]
Upper boundays of theta bins.
unsigned short m_sigmaParamMode_old
sigma mode from input.
int ntheta_old
number of theta bins from input
bool m_useDB
use db or text mode
Class for accessing objects in the database.
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Abstract base class for different kinds of events.