00001
00002
00003
00004 #include "RooRandom.h"
00005 #include "RooRealVar.h"
00006 #include "RooGaussian.h"
00007 #include "RooPolynomial.h"
00008 #include "RooArgSet.h"
00009 #include "RooAddPdf.h"
00010 #include "RooDataSet.h"
00011 #include "RooExtendPdf.h"
00012 #include "RooConstVar.h"
00013
00014 #ifndef __CINT__ // problem including this file with CINT
00015 #include "RooGlobalFunc.h"
00016 #endif
00017
00018 #include "RooStats/HybridCalculatorOriginal.h"
00019 #include "RooStats/HybridResult.h"
00020 #include "RooStats/HybridPlot.h"
00021
00022 void HybridOriginalDemo(int ntoys = 1000)
00023 {
00024
00025
00026
00027
00028
00029
00030
00031
00032 using namespace RooFit;
00033 using namespace RooStats;
00034
00035
00036 RooRandom::randomGenerator()->SetSeed(3007);
00037
00038
00039 RooRealVar x("x","",-3,3);
00040 RooArgList observables(x);
00041
00042
00043
00044
00045
00046 RooGaussian sig_pdf("sig_pdf","",x, RooConst(0.0),RooConst(0.8));
00047 RooRealVar sig_yield("sig_yield","",20,0,300);
00048
00049
00050
00051
00052 RooPolynomial bkg_pdf("bkg_pdf","", x, RooConst(0));
00053 RooRealVar bkg_yield("bkg_yield","",40,0,300);
00054 RooExtendPdf bkg_ext_pdf("bkg_ext_pdf","",bkg_pdf,bkg_yield);
00055
00056
00057 sig_yield.setConstant(kTRUE);
00058
00059
00060 RooAddPdf tot_pdf("tot_pdf","",RooArgList(sig_pdf,bkg_pdf),RooArgList(sig_yield,bkg_yield));
00061
00062
00063
00064 RooGaussian bkg_yield_prior("bkg_yield_prior","",bkg_yield,RooConst(bkg_yield.getVal()),RooConst(10.));
00065
00066 RooArgSet nuisance_parameters(bkg_yield);
00067
00068
00069 RooDataSet* data = tot_pdf.generate(observables,RooFit::Extended());
00070
00071
00072
00073
00074
00075
00076
00077 HybridCalculatorOriginal myHybridCalc(*data, tot_pdf , bkg_ext_pdf ,
00078 &nuisance_parameters, &bkg_yield_prior);
00079
00080
00081 myHybridCalc.SetTestStatistic(1);
00082
00083
00084 myHybridCalc.SetNumberOfToys(ntoys);
00085 myHybridCalc.UseNuisance(true);
00086
00087
00088 myHybridCalc.SetGenerateBinned(false);
00089
00090
00091 HybridResult* myHybridResult = myHybridCalc.GetHypoTest();
00092
00093 if (! myHybridResult) {
00094 std::cerr << "\nError returned from Hypothesis test" << std::endl;
00095 return;
00096 }
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119 HybridPlot* myHybridPlot = myHybridResult->GetPlot("myHybridPlot","Plot of results with HybridCalculatorOriginal",100);
00120 myHybridPlot->Draw();
00121
00122
00123 double clsb_data = myHybridResult->CLsplusb();
00124 double clb_data = myHybridResult->CLb();
00125 double cls_data = myHybridResult->CLs();
00126 double data_significance = myHybridResult->Significance();
00127 double min2lnQ_data = myHybridResult->GetTestStat_data();
00128
00129
00130 double mean_sb_toys_test_stat = myHybridPlot->GetSBmean();
00131 myHybridResult->SetDataTestStatistics(mean_sb_toys_test_stat);
00132 double toys_significance = myHybridResult->Significance();
00133
00134 std::cout << "Completed HybridCalculatorOriginal example:\n";
00135 std::cout << " - -2lnQ = " << min2lnQ_data << endl;
00136 std::cout << " - CL_sb = " << clsb_data << std::endl;
00137 std::cout << " - CL_b = " << clb_data << std::endl;
00138 std::cout << " - CL_s = " << cls_data << std::endl;
00139 std::cout << " - significance of data = " << data_significance << std::endl;
00140 std::cout << " - mean significance of toys = " << toys_significance << std::endl;
00141 }
00142
00143