00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #ifndef __CINT__
00021 #include "RooGlobalFunc.h"
00022 #endif
00023 #include "RooDataSet.h"
00024 #include "RooRealVar.h"
00025 #include "RooGaussian.h"
00026 #include "RooAddPdf.h"
00027 #include "RooProdPdf.h"
00028 #include "RooAddition.h"
00029 #include "RooProduct.h"
00030 #include "TCanvas.h"
00031 #include "RooChebychev.h"
00032 #include "RooAbsPdf.h"
00033 #include "RooFit.h"
00034 #include "RooFitResult.h"
00035 #include "RooPlot.h"
00036 #include "RooAbsArg.h"
00037 #include "RooWorkspace.h"
00038 #include "RooStats/ProfileLikelihoodCalculator.h"
00039 #include "RooStats/HypoTestResult.h"
00040 #include <string>
00041
00042
00043
00044 using namespace RooFit;
00045 using namespace RooStats;
00046
00047
00048 void AddModel(RooWorkspace*);
00049 void AddData(RooWorkspace*);
00050 void DoHypothesisTest(RooWorkspace*);
00051 void MakePlots(RooWorkspace*);
00052
00053
00054 void rs102_hypotestwithshapes() {
00055
00056
00057
00058
00059 RooWorkspace* wspace = new RooWorkspace("myWS");
00060
00061
00062 AddModel(wspace);
00063
00064
00065 AddData(wspace);
00066
00067
00068
00069
00070
00071 DoHypothesisTest(wspace);
00072
00073
00074 MakePlots(wspace);
00075
00076
00077 delete wspace;
00078 }
00079
00080
00081 void AddModel(RooWorkspace* wks){
00082
00083
00084
00085
00086
00087
00088 Double_t lowRange = 60, highRange = 200;
00089
00090
00091 RooRealVar invMass("invMass", "M_{inv}", lowRange, highRange,"GeV");
00092
00093
00094
00095
00096 RooRealVar mH("mH","Higgs Mass",130,90,160) ;
00097 RooRealVar sigma1("sigma1","Width of Gaussian",12.,2,100) ;
00098 RooGaussian sigModel("sigModel", "Signal Model", invMass, mH, sigma1);
00099
00100 mH.setConstant();
00101
00102 sigma1.setConstant();
00103
00104
00105
00106 RooRealVar mZ("mZ", "Z Mass", 91.2, 0, 100);
00107 RooRealVar sigma1_z("sigma1_z","Width of Gaussian",10.,6,100) ;
00108 RooGaussian zjjModel("zjjModel", "Z+jets Model", invMass, mZ, sigma1_z);
00109
00110 mZ.setConstant();
00111
00112 sigma1_z.setConstant();
00113
00114
00115
00116
00117 RooRealVar a0("a0","a0",0.26,-1,1) ;
00118 RooRealVar a1("a1","a1",-0.17596,-1,1) ;
00119 RooRealVar a2("a2","a2",0.018437,-1,1) ;
00120 RooRealVar a3("a3","a3",0.02,-1,1) ;
00121 RooChebychev qcdModel("qcdModel","A Polynomail for QCD",invMass,RooArgList(a0,a1,a2)) ;
00122
00123
00124 a0.setConstant();
00125 a1.setConstant();
00126 a2.setConstant();
00127
00128
00129
00130
00131
00132 RooRealVar fzjj("fzjj","fraction of zjj background events",.4,0.,1) ;
00133
00134
00135 RooRealVar fsigExpected("fsigExpected","expected fraction of signal events",.2,0.,1) ;
00136 fsigExpected.setConstant();
00137
00138
00139
00140 RooRealVar mu("mu","signal strength in units of SM expectation",1,0.,2) ;
00141
00142
00143
00144 RooRealVar ratioSigEff("ratioSigEff","ratio of signal efficiency to nominal signal efficiency",1. ,0.,2) ;
00145 ratioSigEff.setConstant(kTRUE);
00146
00147
00148 RooProduct fsig("fsig","fraction of signal events",RooArgSet(mu,ratioSigEff,fsigExpected)) ;
00149
00150
00151 RooAddPdf model("model","sig+zjj+qcd background shapes",RooArgList(sigModel,zjjModel, qcdModel),RooArgList(fsig,fzjj)) ;
00152
00153
00154
00155
00156
00157 wks->import(model);
00158 }
00159
00160
00161 void AddData(RooWorkspace* wks){
00162
00163
00164 Int_t nEvents = 150;
00165 RooAbsPdf* model = wks->pdf("model");
00166 RooRealVar* invMass = wks->var("invMass");
00167
00168 RooDataSet* data = model->generate(*invMass,nEvents);
00169
00170 wks->import(*data, Rename("data"));
00171
00172 }
00173
00174
00175 void DoHypothesisTest(RooWorkspace* wks){
00176
00177
00178
00179 ModelConfig model;
00180 model.SetWorkspace(*wks);
00181 model.SetPdf("model");
00182
00183
00184
00185 ProfileLikelihoodCalculator plc;
00186 plc.SetData( *(wks->data("data") ));
00187
00188
00189
00190 RooRealVar* mu = wks->var("mu");
00191
00192
00193 RooArgSet poi(*mu);
00194 RooArgSet * nullParams = (RooArgSet*) poi.snapshot();
00195 nullParams->setRealValue("mu",0);
00196
00197
00198
00199 plc.SetModel(model);
00200
00201
00202
00203
00204
00205 plc.SetNullParameters( *nullParams);
00206
00207
00208
00209
00210 HypoTestResult* htr = plc.GetHypoTest();
00211 cout << "-------------------------------------------------" << endl;
00212 cout << "The p-value for the null is " << htr->NullPValue() << endl;
00213 cout << "Corresponding to a signifcance of " << htr->Significance() << endl;
00214 cout << "-------------------------------------------------\n\n" << endl;
00215
00216
00217 }
00218
00219
00220 void MakePlots(RooWorkspace* wks) {
00221
00222
00223
00224
00225
00226
00227 RooAbsPdf* model = wks->pdf("model");
00228 RooAbsPdf* sigModel = wks->pdf("sigModel");
00229 RooAbsPdf* zjjModel = wks->pdf("zjjModel");
00230 RooAbsPdf* qcdModel = wks->pdf("qcdModel");
00231
00232 RooRealVar* mu = wks->var("mu");
00233 RooRealVar* invMass = wks->var("invMass");
00234 RooAbsData* data = wks->data("data");
00235
00236
00237
00238
00239
00240 mu->setConstant(kFALSE);
00241
00242 model->fitTo(*data,Save(kTRUE),Minos(kFALSE), Hesse(kFALSE),PrintLevel(-1));
00243
00244
00245 new TCanvas();
00246 RooPlot* frame = invMass->frame() ;
00247 data->plotOn(frame ) ;
00248 model->plotOn(frame) ;
00249 model->plotOn(frame,Components(*sigModel),LineStyle(kDashed), LineColor(kRed)) ;
00250 model->plotOn(frame,Components(*zjjModel),LineStyle(kDashed),LineColor(kBlack)) ;
00251 model->plotOn(frame,Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
00252
00253 frame->SetTitle("An example fit to the signal + background model");
00254 frame->Draw() ;
00255
00256
00257
00258
00259
00260 mu->setVal(0);
00261 mu->setConstant(kTRUE);
00262
00263 model->fitTo(*data, Save(kTRUE), Minos(kFALSE), Hesse(kFALSE),PrintLevel(-1));
00264
00265
00266 new TCanvas();
00267 RooPlot* xframe2 = invMass->frame() ;
00268 data->plotOn(xframe2, DataError(RooAbsData::SumW2)) ;
00269 model->plotOn(xframe2) ;
00270 model->plotOn(xframe2, Components(*zjjModel),LineStyle(kDashed),LineColor(kBlack)) ;
00271 model->plotOn(xframe2, Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
00272
00273 xframe2->SetTitle("An example fit to the background-only model");
00274 xframe2->Draw() ;
00275
00276
00277 }
00278