40 using namespace RooFit;
41 using namespace RooStats;
44 void AddModel(RooWorkspace *);
45 void AddData(RooWorkspace *);
46 void DoHypothesisTest(RooWorkspace *);
47 void MakePlots(RooWorkspace *);
50 void rs102_hypotestwithshapes()
56 RooWorkspace *wspace =
new RooWorkspace(
"myWS");
68 DoHypothesisTest(wspace);
78 void AddModel(RooWorkspace *wks)
86 Double_t lowRange = 60, highRange = 200;
89 RooRealVar invMass(
"invMass",
"M_{inv}", lowRange, highRange,
"GeV");
93 RooRealVar mH(
"mH",
"Higgs Mass", 130, 90, 160);
94 RooRealVar sigma1(
"sigma1",
"Width of Gaussian", 12., 2, 100);
95 RooGaussian sigModel(
"sigModel",
"Signal Model", invMass, mH, sigma1);
103 RooRealVar mZ(
"mZ",
"Z Mass", 91.2, 0, 100);
104 RooRealVar sigma1_z(
"sigma1_z",
"Width of Gaussian", 10., 6, 100);
105 RooGaussian zjjModel(
"zjjModel",
"Z+jets Model", invMass, mZ, sigma1_z);
109 sigma1_z.setConstant();
113 RooRealVar a0(
"a0",
"a0", 0.26, -1, 1);
114 RooRealVar a1(
"a1",
"a1", -0.17596, -1, 1);
115 RooRealVar a2(
"a2",
"a2", 0.018437, -1, 1);
116 RooRealVar a3(
"a3",
"a3", 0.02, -1, 1);
117 RooChebychev qcdModel(
"qcdModel",
"A Polynomial for QCD", invMass, RooArgList(a0, a1, a2));
128 RooRealVar fzjj(
"fzjj",
"fraction of zjj background events", .4, 0., 1);
131 RooRealVar fsigExpected(
"fsigExpected",
"expected fraction of signal events", .2, 0., 1);
132 fsigExpected.setConstant();
136 RooRealVar mu(
"mu",
"signal strength in units of SM expectation", 1, 0., 2);
140 RooRealVar ratioSigEff(
"ratioSigEff",
"ratio of signal efficiency to nominal signal efficiency", 1., 0., 2);
141 ratioSigEff.setConstant(kTRUE);
144 RooProduct fsig(
"fsig",
"fraction of signal events", RooArgSet(mu, ratioSigEff, fsigExpected));
147 RooAddPdf model(
"model",
"sig+zjj+qcd background shapes", RooArgList(sigModel, zjjModel, qcdModel),
148 RooArgList(fsig, fzjj));
158 void AddData(RooWorkspace *wks)
163 RooAbsPdf *model = wks->pdf(
"model");
164 RooRealVar *invMass = wks->var(
"invMass");
166 RooDataSet *data = model->generate(*invMass, nEvents);
168 wks->import(*data, Rename(
"data"));
172 void DoHypothesisTest(RooWorkspace *wks)
177 model.SetWorkspace(*wks);
178 model.SetPdf(
"model");
182 ProfileLikelihoodCalculator plc;
183 plc.SetData(*(wks->data(
"data")));
187 RooRealVar *mu = wks->var(
"mu");
191 RooArgSet *nullParams = (RooArgSet *)poi.snapshot();
192 nullParams->setRealValue(
"mu", 0);
201 plc.SetNullParameters(*nullParams);
204 HypoTestResult *htr = plc.GetHypoTest();
205 cout <<
"-------------------------------------------------" << endl;
206 cout <<
"The p-value for the null is " << htr->NullPValue() << endl;
207 cout <<
"Corresponding to a significance of " << htr->Significance() << endl;
208 cout <<
"-------------------------------------------------\n\n" << endl;
212 void MakePlots(RooWorkspace *wks)
220 RooAbsPdf *model = wks->pdf(
"model");
221 RooAbsPdf *sigModel = wks->pdf(
"sigModel");
222 RooAbsPdf *zjjModel = wks->pdf(
"zjjModel");
223 RooAbsPdf *qcdModel = wks->pdf(
"qcdModel");
225 RooRealVar *mu = wks->var(
"mu");
226 RooRealVar *invMass = wks->var(
"invMass");
227 RooAbsData *data = wks->data(
"data");
232 mu->setConstant(kFALSE);
234 model->fitTo(*data, Save(kTRUE), Minos(kFALSE), Hesse(kFALSE), PrintLevel(-1));
238 RooPlot *frame = invMass->frame();
240 model->plotOn(frame);
241 model->plotOn(frame, Components(*sigModel), LineStyle(kDashed), LineColor(kRed));
242 model->plotOn(frame, Components(*zjjModel), LineStyle(kDashed), LineColor(kBlack));
243 model->plotOn(frame, Components(*qcdModel), LineStyle(kDashed), LineColor(kGreen));
245 frame->SetTitle(
"An example fit to the signal + background model");
253 mu->setConstant(kTRUE);
255 model->fitTo(*data, Save(kTRUE), Minos(kFALSE), Hesse(kFALSE), PrintLevel(-1));
259 RooPlot *xframe2 = invMass->frame();
260 data->plotOn(xframe2, DataError(RooAbsData::SumW2));
261 model->plotOn(xframe2);
262 model->plotOn(xframe2, Components(*zjjModel), LineStyle(kDashed), LineColor(kBlack));
263 model->plotOn(xframe2, Components(*qcdModel), LineStyle(kDashed), LineColor(kGreen));
265 xframe2->SetTitle(
"An example fit to the background-only model");