Lambdas used to check and fit the result of the H1 analysis.
auto checkH1 = [](
TList *out) {
if (!out) {
std::cout << "checkH1 >>> Test failure: output list not found\n";
return -1;
}
auto hdmd =
dynamic_cast<TH1F *
>(out->FindObject(
"hdmd"));
if (!hdmd) {
std::cout << "checkH1 >>> Test failure: 'hdmd' histo not found\n";
return -1;
}
if ((
Int_t)(hdmd->GetEntries()) != 7525) {
std::cout << "checkH1 >>> Test failure: 'hdmd' histo: wrong number"
" of entries ("
<< (
Int_t)(hdmd->GetEntries()) <<
": expected 7525) \n";
return -1;
}
if (
TMath::Abs((hdmd->GetMean() - 0.15512023) / 0.15512023) > 0.001) {
std::cout << "checkH1 >>> Test failure: 'hdmd' histo: wrong mean (" << hdmd->GetMean()
<< ": expected 0.15512023) \n";
return -1;
}
auto h2 =
dynamic_cast<TH2F *
>(out->FindObject(
"h2"));
if (!h2) {
std::cout << "checkH1 >>> Test failure: 'h2' histo not found\n";
return -1;
}
if ((
Int_t)(h2->GetEntries()) != 7525) {
std::cout << "checkH1 >>> Test failure: 'h2' histo: wrong number"
" of entries ("
<< (
Int_t)(h2->GetEntries()) <<
": expected 7525) \n";
return -1;
}
if (
TMath::Abs((h2->GetMean() - 0.15245688) / 0.15245688) > 0.001) {
std::cout << "checkH1 >>> Test failure: 'h2' histo: wrong mean (" << h2->GetMean() << ": expected 0.15245688) \n";
return -1;
}
return 0;
};
auto doFit = [](
TList *out,
const char *lfn = 0) ->
Int_t {
if (lfn)
gSystem->RedirectOutput(lfn,
"a", &redH);
if (hdmd == 0 || h2 == 0) {
std::cout << "doFit: hdmd = " << hdmd << " , h2 = " << h2 << "\n";
return -1;
if (lfn)
gSystem->RedirectOutput(0, 0, &redH);
}
c1->SetBottomMargin(0.15);
hdmd->GetXaxis()->SetTitle("m_{K#pi#pi} - m_{K#pi}[GeV/c^{2}]");
hdmd->GetXaxis()->SetTitleOffset(1.4);
if (
gROOT->GetListOfFunctions()->FindObject(
"f5"))
delete gROOT->GetFunction(
"f5");
return 0;
par[2] / 2.5066 / par[4] *
TMath::Exp(-xp3 / 2 / par[4] / par[4]));
return res;
};
auto f5 =
new TF1(
"f5",
fdm5, 0.139, 0.17, 5);
f5->SetParameters(1000000, .25, 2000, .1454, .001);
hdmd->Fit("f5", "lr");
Double_t ref_f5[4] = {959915.0, 0.351114, 1185.03, 0.145569};
for (int i : {0, 1, 2, 3}) {
if ((
TMath::Abs((f5->GetParameters())[i] - ref_f5[i]) / ref_f5[i]) > 0.001) {
std::cout << "\n >>> Test failure: fit to 'f5': parameter '" << f5->GetParName(i) << "' has wrong value ("
<< (f5->GetParameters())[i] << ": expected" << ref_f5[i] << ") \n";
if (lfn)
gSystem->RedirectOutput(0, 0, &redH);
return -1;
}
}
auto c2 =
new TCanvas(
"c2",
"tauD0", 100, 100, 800, 600);
c2->SetBottomMargin(0.15);
if (
gROOT->GetListOfFunctions()->FindObject(
"f2"))
delete gROOT->GetFunction(
"f2");
const auto dxbin = (0.17 - 0.13) / 40;
const auto sigma = 0.0012;
return 0;
return res;
};
auto f2 =
new TF1(
"f2",
fdm2, 0.139, 0.17, 2);
f2->SetParameters(10000, 10);
std::cout << "doFit: restricting fit to two bins only in this example...\n";
h2->FitSlicesX(f2, 10, 20, 10, "g5 l");
Double_t ref_f2[2] = {52432.2, 105.481};
for (int i : {0, 1}) {
if ((
TMath::Abs((f2->GetParameters())[i] - ref_f2[i]) / ref_f2[i]) > 0.001) {
std::cout << "\n >>> Test failure: fit to 'f2': parameter '" << f2->GetParName(i) << "' has wrong value ("
<< (f2->GetParameters())[i] << ": expected" << ref_f2[i] << ") \n";
if (lfn)
gSystem->RedirectOutput(0, 0, &redH);
return -1;
}
}
h2_1->GetXaxis()->SetTitle("#tau[ps]");
h2_1->SetMarkerStyle(21);
h2_1->Draw();
auto psdmd = (
TPaveStats *)hdmd->GetListOfFunctions()->FindObject(
"stats");
psdmd->SetOptStat(1110);
if (lfn)
gSystem->RedirectOutput(0, 0, &redH);
return 0;
};
auto hdmd =
new TH1F(
"hdmd",
"Dm_d", 40, 0.13, 0.17);
auto h2 =
new TH2F(
"h2",
"ptD0 vs Dm_d", 30, 0.135, 0.165, 30, -3, 6);
while (reader.Next()) {
continue;
if (*fPtds_d <= 2.5)
continue;
continue;
(*fIk)--;
(*fIpi)--;
if (fNhitrp.At(*fIk) * fNhitrp.At(*fIpi) <= 1)
continue;
if (fRend.At(*fIk) - fRstart.At(*fIk) <= 22)
continue;
if (fRend.At(*fIpi) - fRstart.At(*fIpi) <= 22)
continue;
if (fNlhk.At(*fIk) <= 0.1)
continue;
if (fNlhpi.At(*fIpi) <= 0.1)
continue;
(*fIpis)--;
if (fNlhpi.At(*fIpis) <= 0.1)
continue;
if (*fNjets < 1)
continue;
hdmd->Fill(*fDm_d);
h2->Fill(*fDm_d, *fRpd0_t / 0.029979 * 1.8646 / *fPtd0_d);
}
};
auto elist =
new TEntryList(
"elist",
"H1 selection from Cut");
while (reader.Next()) {
continue;
if (*fPtds_d <= 2.5)
continue;
continue;
(*fIk)--;
(*fIpi)--;
if (fNhitrp.At(*fIk) * fNhitrp.At(*fIpi) <= 1)
continue;
if (fRend.At(*fIk) - fRstart.At(*fIk) <= 22)
continue;
if (fRend.At(*fIpi) - fRstart.At(*fIpi) <= 22)
continue;
if (fNlhk.At(*fIk) <= 0.1)
continue;
if (fNlhpi.At(*fIpi) <= 0.1)
continue;
(*fIpis)--;
if (fNlhpi.At(*fIpis) <= 0.1)
continue;
if (*fNjets < 1)
continue;
elist->Enter(reader.GetCurrentEntry(), reader.GetTree());
}
return elist;
};
auto hdmd =
new TH1F(
"hdmd",
"Dm_d", 40, 0.13, 0.17);
auto h2 =
new TH2F(
"h2",
"ptD0 vs Dm_d", 30, 0.135, 0.165, 30, -3, 6);
while (reader.Next()) {
hdmd->Fill(*fDm_d);
h2->Fill(*fDm_d, *fRpd0_t / 0.029979 * 1.8646 / *fPtd0_d);
}
};
int Int_t
Signed integer 4 bytes (int).
double Double_t
Double 8 bytes.
1-D histogram with a double per channel (see TH1 documentation)
1-D histogram with a float per channel (see TH1 documentation)
2-D histogram with a float per channel (see TH1 documentation)
Use the TLine constructor to create a simple line.
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
An interface for reading collections stored in ROOT columnar datasets.
An interface for reading values stored in ROOT columnar datasets.
A simple, robust and fast interface to read values from ROOT columnar datasets such as TTree,...
Double_t fdm5(Double_t *xx, Double_t *par)
Double_t fdm2(Double_t *xx, Double_t *par)
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.