Logo ROOT   6.07/09
Reference Guide
StandardProfileInspectorDemo.C File Reference

Detailed Description

View in nbviewer Open in SWAN Standard demo of the ProfileInspector class StandardProfileInspectorDemo

This is a standard demo that can be used with any ROOT file prepared in the standard way. You specify:

With default parameters the macro will attempt to run the standard hist2workspace example and read the ROOT file that it produces.

The actual heart of the demo is only about 10 lines long.

The ProfileInspector plots the conditional maximum likelihood estimate of each nuisance parameter in the model vs. the parameter of interest. (aka. profiled value of nuisance parameter vs. parameter of interest) (aka. best fit nuisance parameter with p.o.i fixed vs. parameter of interest)

pict1_StandardProfileInspectorDemo.C.png
Processing /mnt/vdb/lsf/workspace/root-makedoc-v608/rootspi/rdoc/src/v6-08-00-patches/tutorials/roostats/StandardProfileInspectorDemo.C...
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:Minization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minization -- Including the following contraint terms in minimization: (lumiConstraint,alpha_syst1Constraint,alpha_syst2Constraint,alpha_syst3Constraint,gamma_stat_channel1_bin_0_constraint,gamma_stat_channel1_bin_1_constraint)
[#1] INFO:Minization -- RooProfileLL::evaluate(nll_simPdf_obsData_with_constr_Profile[SigXsecOverSM]) Creating instance of MINUIT
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_simPdf_obsData_with_constr) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minization -- RooProfileLL::evaluate(nll_simPdf_obsData_with_constr_Profile[SigXsecOverSM]) determining minimum likelihood for current configurations w.r.t all observable
RooAbsTestStatistic::initSimMode: creating slave calculator #0 for state channel1 (2 dataset entries)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:Fitting -- RooAbsTestStatistic::initSimMode: created 1 slave calculators.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gamma_stat_channel1_bin_0_constraint_Int[gamma_stat_channel1_bin_0]) using numeric integrator RooIntegrator1D to calculate Int(gamma_stat_channel1_bin_0)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gamma_stat_channel1_bin_1_constraint_Int[gamma_stat_channel1_bin_1]) using numeric integrator RooIntegrator1D to calculate Int(gamma_stat_channel1_bin_1)
[#1] INFO:Minization -- RooProfileLL::evaluate(nll_simPdf_obsData_with_constr_Profile[SigXsecOverSM]) minimum found at (SigXsecOverSM=1.11601)
....................................................................................................
#include "TFile.h"
#include "TROOT.h"
#include "TCanvas.h"
#include "TList.h"
#include "TMath.h"
#include "TSystem.h"
#include "RooWorkspace.h"
#include "RooAbsData.h"
using namespace RooFit;
using namespace RooStats;
void StandardProfileInspectorDemo(const char* infile = "",
const char* workspaceName = "combined",
const char* modelConfigName = "ModelConfig",
const char* dataName = "obsData") {
// -------------------------------------------------------
// First part is just to access a user-defined file
// or create the standard example file if it doesn't exist
const char* filename = "";
if (!strcmp(infile,"")) {
filename = "results/example_combined_GaussExample_model.root";
bool fileExist = !gSystem->AccessPathName(filename); // note opposite return code
// if file does not exists generate with histfactory
if (!fileExist) {
#ifdef _WIN32
cout << "HistFactory file cannot be generated on Windows - exit" << endl;
return;
#endif
// Normally this would be run on the command line
cout <<"will run standard hist2workspace example"<<endl;
gROOT->ProcessLine(".! prepareHistFactory .");
gROOT->ProcessLine(".! hist2workspace config/example.xml");
cout <<"\n\n---------------------"<<endl;
cout <<"Done creating example input"<<endl;
cout <<"---------------------\n\n"<<endl;
}
}
else
filename = infile;
// Try to open the file
TFile *file = TFile::Open(filename);
// if input file was specified byt not found, quit
if(!file ){
cout <<"StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
// -------------------------------------------------------
// Tutorial starts here
// -------------------------------------------------------
// get the workspace out of the file
RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
if(!w){
cout <<"workspace not found" << endl;
return;
}
// get the modelConfig out of the file
ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
// get the modelConfig out of the file
RooAbsData* data = w->data(dataName);
// make sure ingredients are found
if(!data || !mc){
w->Print();
cout << "data or ModelConfig was not found" <<endl;
return;
}
// -----------------------------
// now use the profile inspector
TList* list = p.GetListOfProfilePlots(*data,mc);
// now make plots
TCanvas* c1 = new TCanvas("c1","ProfileInspectorDemo",800,200);
if(list->GetSize()>4){
double n = list->GetSize();
int nx = (int)sqrt(n) ;
int ny = TMath::CeilNint(n/nx);
nx = TMath::CeilNint( sqrt(n) );
c1->Divide(ny,nx);
} else
c1->Divide(list->GetSize());
for(int i=0; i<list->GetSize(); ++i){
c1->cd(i+1);
list->At(i)->Draw("al");
}
cout << endl;
}
Author
Kyle Cranmer

Definition in file StandardProfileInspectorDemo.C.