rs102_hypotestwithshapes.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////
00002 //
00003 // rs102_hypotestwithshapes for RooStats project
00004 // Author: Kyle Cranmer <cranmer@cern.ch>
00005 // 
00006 // Modified from version of February 29, 2008
00007 //
00008 // This tutorial macro shows a typical search for a new particle 
00009 // by studying an invariant mass distribution.  
00010 // The macro creates a simple signal model and two background models, 
00011 // which are added to a RooWorkspace.
00012 // The macro creates a toy dataset, and then uses a RooStats 
00013 // ProfileLikleihoodCalculator to do a hypothesis test of the 
00014 // background-only and signal+background hypotheses.
00015 // In this example, shape uncertainties are not taken into account, but
00016 // normalization uncertainties are.  
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 // use this order for safety on library loading
00044 using namespace RooFit;
00045 using namespace RooStats;
00046 
00047 // see below for implementation
00048 void AddModel(RooWorkspace*);
00049 void AddData(RooWorkspace*);
00050 void DoHypothesisTest(RooWorkspace*);
00051 void MakePlots(RooWorkspace*);
00052 
00053 //____________________________________
00054 void rs102_hypotestwithshapes() {
00055 
00056   // The main macro.
00057 
00058   // Create a workspace to manage the project.
00059   RooWorkspace* wspace = new RooWorkspace("myWS");
00060 
00061   // add the signal and background models to the workspace
00062   AddModel(wspace);
00063 
00064   // add some toy data to the workspace
00065   AddData(wspace);
00066 
00067   // inspect the workspace if you wish
00068   //  wspace->Print();
00069 
00070   // do the hypothesis test
00071   DoHypothesisTest(wspace);
00072 
00073   // make some plots
00074   MakePlots(wspace);
00075 
00076   // cleanup
00077   delete wspace;
00078 }
00079  
00080 //____________________________________
00081 void AddModel(RooWorkspace* wks){
00082 
00083   // Make models for signal (Higgs) and background (Z+jets and QCD)
00084   // In real life, this part requires an intellegent modeling 
00085   // of signal and background -- this is only an example.  
00086 
00087   // set range of observable
00088   Double_t lowRange = 60, highRange = 200;
00089 
00090   // make a RooRealVar for the observable
00091   RooRealVar invMass("invMass", "M_{inv}", lowRange, highRange,"GeV");
00092  
00093 
00094   /////////////////////////////////////////////
00095   // make a simple signal model. 
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   // we will test this specific mass point for the signal
00100   mH.setConstant();
00101   // and we assume we know the mass resolution
00102   sigma1.setConstant();
00103 
00104   /////////////////////////////////////////////
00105   // make zjj model.  Just like signal model
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   // we know Z mass
00110   mZ.setConstant();
00111   // assume we know resolution too
00112   sigma1_z.setConstant();
00113   
00114 
00115   //////////////////////////////////////////////
00116   // make QCD model
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   // let's assume this shape is known, but the normalization is not
00124   a0.setConstant();
00125   a1.setConstant();
00126   a2.setConstant();
00127 
00128   //////////////////////////////////////////////
00129   // combined model
00130 
00131   // Setting the fraction of Zjj to be 40% for initial guess.
00132   RooRealVar fzjj("fzjj","fraction of zjj background events",.4,0.,1) ; 
00133 
00134   // Set the expected fraction of signal to 20%.
00135   RooRealVar fsigExpected("fsigExpected","expected fraction of signal events",.2,0.,1) ; 
00136   fsigExpected.setConstant(); // use mu as main parameter, so fix this.
00137 
00138   // Introduce mu: the signal strength in units of the expectation.  
00139   // eg. mu = 1 is the SM, mu = 0 is no signal, mu=2 is 2x the SM
00140   RooRealVar mu("mu","signal strength in units of SM expectation",1,0.,2) ; 
00141 
00142   // Introduce ratio of signal efficiency to nominal signal efficiency. 
00143   // This is useful if you want to do limits on cross section.
00144   RooRealVar ratioSigEff("ratioSigEff","ratio of signal efficiency to nominal signal efficiency",1. ,0.,2) ; 
00145   ratioSigEff.setConstant(kTRUE);  
00146 
00147   // finally the signal fraction is the product of the terms above.
00148   RooProduct fsig("fsig","fraction of signal events",RooArgSet(mu,ratioSigEff,fsigExpected)) ; 
00149 
00150   // full model
00151   RooAddPdf model("model","sig+zjj+qcd background shapes",RooArgList(sigModel,zjjModel, qcdModel),RooArgList(fsig,fzjj)) ; 
00152 
00153   // interesting for debugging and visualizing the model
00154   //  model.printCompactTree("","fullModel.txt");
00155   //  model.graphVizTree("fullModel.dot");
00156 
00157   wks->import(model);
00158 }
00159 
00160 //____________________________________
00161 void AddData(RooWorkspace* wks){
00162   // Add a toy dataset
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   // Use a RooStats ProfileLikleihoodCalculator to do the hypothesis test.
00179   ModelConfig model; 
00180   model.SetWorkspace(*wks);
00181   model.SetPdf("model");
00182 
00183   //plc.SetData("data");
00184  
00185   ProfileLikelihoodCalculator plc; 
00186   plc.SetData( *(wks->data("data") )); 
00187 
00188   // here we explicitly set the value of the parameters for the null.  
00189   // We want no signal contribution, eg. mu = 0
00190   RooRealVar* mu = wks->var("mu");
00191 //   RooArgSet* nullParams = new RooArgSet("nullParams");
00192 //   nullParams->addClone(*mu);
00193   RooArgSet poi(*mu);
00194   RooArgSet * nullParams = (RooArgSet*) poi.snapshot(); 
00195   nullParams->setRealValue("mu",0); 
00196 
00197   
00198   //plc.SetNullParameters(*nullParams);
00199   plc.SetModel(model);
00200   // NOTE: using snapshot will import nullparams 
00201   // in the WS and merge with existing "mu" 
00202   // model.SetSnapshot(*nullParams);
00203   
00204   //use instead setNuisanceParameters
00205   plc.SetNullParameters( *nullParams);
00206 
00207  
00208 
00209   // We get a HypoTestResult out of the calculator, and we can query it.
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   // Make plots of the data and the best fit model in two cases:
00223   // first the signal+background case
00224   // second the background-only case.
00225 
00226   // get some things out of workspace
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   // Make plots for the Alternate hypothesis, eg. let mu float
00239 
00240   mu->setConstant(kFALSE);
00241 
00242   model->fitTo(*data,Save(kTRUE),Minos(kFALSE), Hesse(kFALSE),PrintLevel(-1));
00243   
00244   //plot sig candidates, full model, and individual componenets
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   //  cdata->SaveAs("alternateFit.gif");
00256 
00257   //////////////////////////////////////////////////////////
00258   // Do Fit to the Null hypothesis.  Eg. fix mu=0
00259 
00260   mu->setVal(0); // set signal fraction to 0
00261   mu->setConstant(kTRUE); // set constant 
00262 
00263   model->fitTo(*data, Save(kTRUE), Minos(kFALSE), Hesse(kFALSE),PrintLevel(-1));
00264 
00265   // plot signal candidates with background model and components
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   //  cbkgonly->SaveAs("nullFit.gif");
00276 
00277 }
00278 

Generated on Tue Jul 5 15:45:10 2011 for ROOT_528-00b_version by  doxygen 1.5.1