PrepareWorkspace_GaussOverFlat_withSystematics.C File Reference

Go to the source code of this file.

Functions

void PrepareWorkspace_GaussOverFlat_withSystematics (TString fileName="WS_GaussOverFlat_withSystematics.root")

Function Documentation

void PrepareWorkspace_GaussOverFlat_withSystematics ( TString  fileName = "WS_GaussOverFlat_withSystematics.root"  ) 

Definition at line 3 of file PrepareWorkspace_GaussOverFlat_withSystematics.C.

00004 {
00005   // In this macro a PDF model is built assuming signal has a Gaussian
00006   // PDF and the background a flat PDF.  The parameter of interest is
00007   // the signal yield and we assume for it a flat prior.  It is shown
00008   // how two types of systematics uncertainties can be expressed;
00009   // those are a sytematic uncertainty on the background yield and
00010   // another on one of the parameters (sigma) of the signal shape.
00011   // All needed objects are stored in a ROOT file (within a
00012   // RooWorkspace container); this ROOT file can then be fed as input
00013   // to various statistical methods.
00014 
00015   using namespace RooFit;
00016   using namespace RooStats;
00017 
00018   // use an observable for this shape-based analysis
00019   RooRealVar* mass = new RooRealVar("mass","mass",0,500,"GeV/c^{2}");
00020   mass->setBins(100);
00021   RooArgSet* observables = new RooArgSet(*mass,"observables");
00022 
00023   // signal (Gaussian) and background (flat) PDFs
00024   RooRealVar* sigSigma = new RooRealVar("sigSigma","sigma in signal PDF",0,100);
00025   RooAbsPdf* sigPdf = new RooGaussian("sigPdf","signal PDF",*mass,RooConst(200),*sigSigma);
00026   RooAbsPdf* bkgPdf = new RooPolynomial("bkgPdf","background PDF",*mass,RooFit::RooConst(0));
00027   
00028   // S+B model: the sum of both shapes weighted with the yields
00029   RooRealVar* S = new RooRealVar("S","signal yield",1000,0,100000);
00030   RooRealVar* B = new RooRealVar("B","background yield",1000,0,100000);
00031   RooAbsPdf* model = new RooAddPdf("model","S+B PDF",RooArgList(*sigPdf,*bkgPdf),RooArgList(*S,*B));
00032   
00033   // B-only model: the same as with a signal yield fixed to 0
00034   RooAbsPdf* modelBkg = new RooExtendPdf("modelBkg","B-only PDF",*bkgPdf,*B);
00035 
00036   // assume a Gaussian uncertainty on the detector resolution affecting the signal width (of 10%)
00037   // another nuisance parameter is the background yield (apply an uncertainty of 20%)
00038   RooAbsPdf* prior_sigSigma = new RooGaussian("prior_sigSigma","prior probability on sigSigma",*sigSigma,RooConst(sigSigma->getVal()),RooConst(sigSigma->getVal()*0.10));
00039   RooAbsPdf* prior_B = new RooGaussian("prior_B","prior probability on B",*B,RooConst(B->getVal()),RooConst(B->getVal()*0.20));
00040   RooAbsPdf* priorNuisance = new RooProdPdf("priorNuisance","prior on the nuisance parameters",*prior_B,*prior_sigSigma);
00041   RooArgSet* parameters = new RooArgSet(*B,*sigSigma,"parameters");
00042 
00043   // assume a flat prior on our parameter of interest (POI) which is the signal yield
00044   RooAbsPdf* priorPOI = new RooPolynomial("priorPOI","flat prior on the POI",*S,RooFit::RooConst(0));
00045   RooArgSet* POI = new RooArgSet(*S,"POI");
00046   
00047   // different options are shown for the data generation from the model
00048   
00049   // unbinned data with Poisson fluctuations
00050 //   RooAbsData* data = (RooDataSet*) model->generate(*observables,RooFit::Extended(),Name("data"));
00051 
00052   // binned data with Poisson fluctuations
00053 //   RooAbsData* data = (RooDataHist*) model->generateBinned(*observables,Extended(),Name("data"));
00054   
00055   // binned without any fluctuations (average case)
00056   RooAbsData* data = (RooDataHist*) model->generateBinned(*observables,Name("data"),ExpectedData());
00057 
00058   // control plot of the generated data
00059 //   RooPlot* plot = mass->frame();
00060 //   data->plotOn(plot);
00061 //   plot->Draw();
00062 
00063   // use a RooWorkspace to store the pdf models, prior informations, list of parameters,...
00064   RooWorkspace myWS("myWS");
00065   myWS.import(*data,Rename("data"));
00066   myWS.import(*model,RecycleConflictNodes());
00067   myWS.import(*modelBkg,RecycleConflictNodes());
00068   myWS.import(*priorPOI,RecycleConflictNodes());
00069   myWS.import(*priorNuisance,RecycleConflictNodes());  
00070   myWS.defineSet("observables",*observables,kTRUE);
00071   myWS.defineSet("parameters",*parameters,kTRUE);
00072   myWS.defineSet("POI",*POI,kTRUE);
00073 
00074   // store the workspace in a ROOT file  
00075   TFile file(fileName,"RECREATE");
00076   file.cd();
00077   myWS.Write();
00078   file.Write();
00079   file.Close();
00080   
00081   std::cout << "\nRooFit model initialized and stored in " << fileName << std::endl;
00082 }


Generated on Tue Oct 6 09:48:21 2009 for Bayesian Analysis Toolkit by  doxygen 1.6.1