00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #ifndef __CINT__
00023 #include "RooGlobalFunc.h"
00024 #endif
00025 #include "RooRealVar.h"
00026 #include "RooStats/SPlot.h"
00027 #include "RooDataSet.h"
00028 #include "RooRealVar.h"
00029 #include "RooGaussian.h"
00030 #include "RooExponential.h"
00031 #include "RooChebychev.h"
00032 #include "RooAddPdf.h"
00033 #include "RooProdPdf.h"
00034 #include "RooAddition.h"
00035 #include "RooProduct.h"
00036 #include "TCanvas.h"
00037 #include "RooAbsPdf.h"
00038 #include "RooFit.h"
00039 #include "RooFitResult.h"
00040 #include "RooWorkspace.h"
00041 #include "RooConstVar.h"
00042
00043
00044 using namespace RooFit ;
00045 using namespace RooStats ;
00046
00047
00048
00049 void AddModel(RooWorkspace*);
00050 void AddData(RooWorkspace*);
00051 void DoSPlot(RooWorkspace*);
00052 void MakePlots(RooWorkspace*);
00053
00054 void rs301_splot()
00055 {
00056
00057
00058 RooWorkspace* wspace = new RooWorkspace("myWS");
00059
00060
00061
00062 AddModel(wspace);
00063
00064
00065 AddData(wspace);
00066
00067
00068
00069
00070
00071
00072 DoSPlot(wspace);
00073
00074
00075
00076 MakePlots(wspace);
00077
00078
00079 delete wspace;
00080
00081 }
00082
00083
00084
00085 void AddModel(RooWorkspace* ws){
00086
00087
00088
00089
00090
00091
00092 Double_t lowRange = 00, highRange = 200;
00093
00094
00095 RooRealVar invMass("invMass", "M_{inv}", lowRange, highRange,"GeV");
00096 RooRealVar isolation("isolation", "isolation", 0., 20., "GeV");
00097
00098
00099
00100
00101
00102
00103 std::cout << "make z model" << std::endl;
00104
00105 RooRealVar mZ("mZ", "Z Mass", 91.2, lowRange, highRange);
00106 RooRealVar sigmaZ("sigmaZ", "Width of Gaussian", 2,0,10,"GeV");
00107 RooGaussian mZModel("mZModel", "Z+jets Model", invMass, mZ, sigmaZ);
00108
00109 mZ.setConstant();
00110
00111
00112
00113
00114
00115
00116 RooConstVar zIsolDecayConst("zIsolDecayConst",
00117 "z isolation decay constant", -1);
00118 RooExponential zIsolationModel("zIsolationModel", "z isolation model",
00119 isolation, zIsolDecayConst);
00120
00121
00122 RooProdPdf zModel("zModel", "4-d model for Z",
00123 RooArgSet(mZModel, zIsolationModel));
00124
00125
00126
00127
00128 std::cout << "make qcd model" << std::endl;
00129
00130
00131
00132
00133
00134 RooRealVar qcdMassDecayConst("qcdMassDecayConst",
00135 "Decay const for QCD mass spectrum",
00136 -0.01, -100, 100,"1/GeV");
00137 RooExponential qcdMassModel("qcdMassModel", "qcd Mass Model",
00138 invMass, qcdMassDecayConst);
00139
00140
00141
00142
00143
00144
00145 RooConstVar qcdIsolDecayConst("qcdIsolDecayConst",
00146 "Et resolution constant", -.1);
00147 RooExponential qcdIsolationModel("qcdIsolationModel", "QCD isolation model",
00148 isolation, qcdIsolDecayConst);
00149
00150
00151 RooProdPdf qcdModel("qcdModel", "2-d model for QCD",
00152 RooArgSet(qcdMassModel, qcdIsolationModel));
00153
00154
00155
00156
00157
00158
00159 RooRealVar zYield("zYield","fitted yield for Z",50 ,0.,1000) ;
00160 RooRealVar qcdYield("qcdYield","fitted yield for QCD", 100 ,0.,1000) ;
00161
00162
00163 std::cout << "make full model" << std::endl;
00164 RooAddPdf model("model","z+qcd background models",
00165 RooArgList(zModel, qcdModel),
00166 RooArgList(zYield,qcdYield));
00167
00168
00169
00170 model.graphVizTree("fullModel.dot");
00171
00172 std::cout << "import model" << std::endl;
00173
00174 ws->import(model);
00175 }
00176
00177
00178 void AddData(RooWorkspace* ws){
00179
00180
00181
00182 Int_t nEvents = 1000;
00183
00184
00185 RooAbsPdf* model = ws->pdf("model");
00186 RooRealVar* invMass = ws->var("invMass");
00187 RooRealVar* isolation = ws->var("isolation");
00188
00189
00190 std::cout << "make data set and import to workspace" << std::endl;
00191 RooDataSet* data = model->generate(RooArgSet(*invMass, *isolation),nEvents);
00192
00193
00194 ws->import(*data, Rename("data"));
00195
00196 }
00197
00198
00199 void DoSPlot(RooWorkspace* ws){
00200 std::cout << "Calculate sWeights" << std::endl;
00201
00202
00203 RooAbsPdf* model = ws->pdf("model");
00204 RooRealVar* zYield = ws->var("zYield");
00205 RooRealVar* qcdYield = ws->var("qcdYield");
00206 RooDataSet* data = (RooDataSet*) ws->data("data");
00207
00208
00209 model->fitTo(*data, Extended() );
00210
00211
00212
00213 RooRealVar* sigmaZ = ws->var("sigmaZ");
00214 RooRealVar* qcdMassDecayConst = ws->var("qcdMassDecayConst");
00215 sigmaZ->setConstant();
00216 qcdMassDecayConst->setConstant();
00217
00218
00219 RooMsgService::instance().setSilentMode(true);
00220
00221
00222
00223
00224 RooStats::SPlot* sData = new RooStats::SPlot("sData","An SPlot",
00225 *data, model, RooArgList(*zYield,*qcdYield) );
00226
00227
00228
00229
00230 std::cout << "Check SWeights:" << std::endl;
00231
00232
00233 std::cout << std::endl << "Yield of Z is "
00234 << zYield->getVal() << ". From sWeights it is "
00235 << sData->GetYieldFromSWeight("zYield") << std::endl;
00236
00237
00238 std::cout << "Yield of QCD is "
00239 << qcdYield->getVal() << ". From sWeights it is "
00240 << sData->GetYieldFromSWeight("qcdYield") << std::endl
00241 << std::endl;
00242
00243 for(Int_t i=0; i < 10; i++)
00244 {
00245 std::cout << "z Weight " << sData->GetSWeight(i,"zYield")
00246 << " qcd Weight " << sData->GetSWeight(i,"qcdYield")
00247 << " Total Weight " << sData->GetSumOfEventSWeight(i)
00248 << std::endl;
00249 }
00250
00251 std::cout << std::endl;
00252
00253
00254 std::cout << "import new dataset with sWeights" << std::endl;
00255 ws->import(*data, Rename("dataWithSWeights"));
00256
00257
00258 }
00259
00260 void MakePlots(RooWorkspace* ws){
00261
00262
00263
00264 std::cout << "make plots" << std::endl;
00265
00266
00267 TCanvas* cdata = new TCanvas("sPlot","sPlot demo", 400, 600);
00268 cdata->Divide(1,3);
00269
00270
00271 RooAbsPdf* model = ws->pdf("model");
00272 RooAbsPdf* zModel = ws->pdf("zModel");
00273 RooAbsPdf* qcdModel = ws->pdf("qcdModel");
00274
00275 RooRealVar* isolation = ws->var("isolation");
00276 RooRealVar* invMass = ws->var("invMass");
00277
00278
00279 RooDataSet* data = (RooDataSet*) ws->data("dataWithSWeights");
00280
00281
00282
00283 model->fitTo(*data, Extended() );
00284
00285
00286
00287 cdata->cd(1);
00288 RooPlot* frame = invMass->frame() ;
00289 data->plotOn(frame ) ;
00290 model->plotOn(frame) ;
00291 model->plotOn(frame,Components(*zModel),LineStyle(kDashed), LineColor(kRed)) ;
00292 model->plotOn(frame,Components(*qcdModel),LineStyle(kDashed),LineColor(kGreen)) ;
00293
00294 frame->SetTitle("Fit of model to discriminating variable");
00295 frame->Draw() ;
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307 cdata->cd(2);
00308
00309
00310 RooDataSet * dataw_z = new RooDataSet(data->GetName(),data->GetTitle(),data,*data->get(),0,"zYield_sw") ;
00311
00312 RooPlot* frame2 = isolation->frame() ;
00313 dataw_z->plotOn(frame2, DataError(RooAbsData::SumW2) ) ;
00314
00315 frame2->SetTitle("isolation distribution for Z");
00316 frame2->Draw() ;
00317
00318
00319
00320
00321
00322 cdata->cd(3);
00323 RooDataSet * dataw_qcd = new RooDataSet(data->GetName(),data->GetTitle(),data,*data->get(),0,"qcdYield_sw") ;
00324 RooPlot* frame3 = isolation->frame() ;
00325 dataw_qcd->plotOn(frame3,DataError(RooAbsData::SumW2) ) ;
00326
00327 frame3->SetTitle("isolation distribution for QCD");
00328 frame3->Draw() ;
00329
00330
00331
00332 }