Demonstrate Z_Bi = Z_Gamma
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b_X_px]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([[py_X_prior_b]_Norm[b]_X_px_NORM[x]]_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([[py_X_prior_b]_Norm[b]_X_px_cdf_NORM[x_prime]]_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([[py_X_prior_b]_Norm[b]_X_px_cdf_Int[x_prime|CDF]_Norm[x_prime]]_Int[b|CDF]) using numeric integrator RooIntegrator1D to calculate Int(b)
[#1] INFO:NumericIntegration -- RooRealIntegral::init([py_X_prior_b]_Norm[b]_denominator_Int[b]) using numeric integrator RooIntegrator1D to calculate Int(b)
Hybrid p-value = 0.999226
Z_Gamma Significance = 3.1655
Z_Bi significance estimation: 3.10804
{
w1->factory(
"Poisson::px(x[150,0,500],sum::splusb(s[0,0,100],b[100,0,300]))");
w1->factory(
"Poisson::py(y[100,0,500],prod::taub(tau[1.],b))");
w1->factory(
"Uniform::prior_b(b)");
w1->factory(
"PROJ::averagedModel(PROD::foo(px|b,py,prior_b),b)");
w1->pdf(
"averagedModel")->plotOn(frame);
w1->pdf(
"px")->plotOn(frame, LineColor(
kRed));
w1->var(
"y")->setVal(100);
w1->var(
"x")->setVal(150);
std::unique_ptr<RooAbsReal>
cdf{
w1->pdf(
"averagedModel")->createCdf(*
w1->var(
"x"))};
cout <<
"Hybrid p-value = " <<
cdf->getVal() << endl;
cout <<
"Z_Gamma Significance = " << PValueToSignificance(1 -
cdf->getVal()) << endl;
double Z_Bi = NumberCountingUtils::BinomialWithTauObsZ(150, 100, 1);
std::cout <<
"Z_Bi significance estimation: " <<
Z_Bi << std::endl;
}
Plot frame and a container for graphics objects within that frame.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Persistable container for RooFit projects.
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
- Authors
- Kyle Cranmer, Wouter Verkerke
Definition in file Zbi_Zgamma.C.