Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
StandardHypoTestInvDemo.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Standard tutorial macro for performing an inverted hypothesis test for computing an interval

This macro will perform a scan of the p-values for computing the interval or limit

Usage:

root> StandardHypoTestInvDemo("fileName","workspace name","S+B modelconfig name","B model name","data set
name",calculator type, test statistic type, use CLS,
number of points, xmin, xmax, number of toys, use number counting)
type = 0 Freq calculator
type = 1 Hybrid calculator
type = 2 Asymptotic calculator
type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
testStatType = 0 LEP
= 1 Tevatron
= 2 Profile Likelihood two sided
= 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
= 4 Profile Likelihood signed ( pll = -pll if mu < mu_hat)
= 5 Max Likelihood Estimate as test statistic
= 6 Number of observed event as test statistic
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
float xmin
float xmax
0x5608553df440 results/example_combined_GaussExample_model.root
Running HypoTestInverter on the workspace combined
RooWorkspace(combined) combined contents
variables
---------
(Lumi,SigXsecOverSM,alpha_syst1,alpha_syst2,alpha_syst3,channelCat,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1,nominalLumi,obs_x_channel1)
p.d.f.s
-------
RooGaussian::alpha_syst1Constraint[ x=alpha_syst1 mean=nom_alpha_syst1 sigma=1 ] = 1
RooGaussian::alpha_syst2Constraint[ x=alpha_syst2 mean=nom_alpha_syst2 sigma=1 ] = 1
RooGaussian::alpha_syst3Constraint[ x=alpha_syst3 mean=nom_alpha_syst3 sigma=1 ] = 1
RooRealSumPdf::channel1_model[ signal_channel1_scaleFactors * signal_channel1_shapes + background1_channel1_scaleFactors * background1_channel1_shapes + background2_channel1_scaleFactors * background2_channel1_shapes ] = 240
RooPoisson::gamma_stat_channel1_bin_0_constraint[ x=nom_gamma_stat_channel1_bin_0 mean=gamma_stat_channel1_bin_0_poisMean ] = 0.019943
RooPoisson::gamma_stat_channel1_bin_1_constraint[ x=nom_gamma_stat_channel1_bin_1 mean=gamma_stat_channel1_bin_1_poisMean ] = 0.039861
RooGaussian::lumiConstraint[ x=Lumi mean=nominalLumi sigma=0.1 ] = 1
RooProdPdf::model_channel1[ lumiConstraint * alpha_syst1Constraint * alpha_syst2Constraint * alpha_syst3Constraint * gamma_stat_channel1_bin_0_constraint * gamma_stat_channel1_bin_1_constraint * channel1_model(obs_x_channel1) ] = 0.190787
RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.190787
functions
--------
RooHistFunc::background1_channel1_Hist_alphanominal[ depList=(obs_x_channel1) depList=(obs_x_channel1) ] = 100
RooStats::HistFactory::FlexibleInterpVar::background1_channel1_epsilon[ paramList=(alpha_syst2) ] = 1
RooProduct::background1_channel1_scaleFactors[ background1_channel1_epsilon * Lumi ] = 1
RooProduct::background1_channel1_shapes[ background1_channel1_Hist_alphanominal * mc_stat_channel1 * channel1_model_binWidth ] = 200
RooHistFunc::background2_channel1_Hist_alphanominal[ depList=(obs_x_channel1) depList=(obs_x_channel1) ] = 0
RooStats::HistFactory::FlexibleInterpVar::background2_channel1_epsilon[ paramList=(alpha_syst3) ] = 1
RooProduct::background2_channel1_scaleFactors[ background2_channel1_epsilon * Lumi ] = 1
RooProduct::background2_channel1_shapes[ background2_channel1_Hist_alphanominal * mc_stat_channel1 * channel1_model_binWidth ] = 0
RooBinWidthFunction::channel1_model_binWidth[ HistFuncForBinWidth=signal_channel1_Hist_alphanominal HistFuncForBinWidth=signal_channel1_Hist_alphanominal ] = 2
RooProduct::gamma_stat_channel1_bin_0_poisMean[ gamma_stat_channel1_bin_0 * gamma_stat_channel1_bin_0_tau ] = 400
RooProduct::gamma_stat_channel1_bin_1_poisMean[ gamma_stat_channel1_bin_1 * gamma_stat_channel1_bin_1_tau ] = 100
ParamHistFunc::mc_stat_channel1[ ] = 1
RooHistFunc::signal_channel1_Hist_alphanominal[ depList=(obs_x_channel1) depList=(obs_x_channel1) ] = 20
RooStats::HistFactory::FlexibleInterpVar::signal_channel1_epsilon[ paramList=(alpha_syst1) ] = 1
RooProduct::signal_channel1_scaleFactors[ signal_channel1_epsilon * SigXsecOverSM * Lumi ] = 1
RooProduct::signal_channel1_shapes[ signal_channel1_Hist_alphanominal * channel1_model_binWidth ] = 40
datasets
--------
RooDataSet::obsData(obs_x_channel1,channelCat)
RooDataSet::asimovData(obs_x_channel1,channelCat)
embedded datasets (in pdfs and functions)
-----------------------------------------
RooDataHist::signal_channel1_Hist_alphanominalDHist(obs_x_channel1)
RooDataHist::background1_channel1_Hist_alphanominalDHist(obs_x_channel1)
RooDataHist::background2_channel1_Hist_alphanominalDHist(obs_x_channel1)
parameter snapshots
-------------------
NominalParamValues = (nominalLumi=1[C],nom_alpha_syst1=0[C],nom_alpha_syst2=0[C],nom_alpha_syst3=0[C],nom_gamma_stat_channel1_bin_0=400[C],nom_gamma_stat_channel1_bin_1=100[C],obs_x_channel1=1.25,Lumi=1 +/- 0.1[C],alpha_syst1=0 +/- 1[C],alpha_syst2=0 +/- 1,alpha_syst3=0 +/- 1,gamma_stat_channel1_bin_0=1 +/- 0.05,gamma_stat_channel1_bin_1=1 +/- 0.1,SigXsecOverSM=1 +/- 0)
named sets
----------
ModelConfig_GlobalObservables:(nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
ModelConfig_NuisParams:(alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
ModelConfig_Observables:(obs_x_channel1,channelCat)
ModelConfig_POI:(SigXsecOverSM)
globalObservables:(nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
observables:(obs_x_channel1,channelCat)
generic objects
---------------
RooStats::ModelConfig::ModelConfig
Using data set obsData
StandardHypoTestInvDemo : POI initial value: SigXsecOverSM = 1
[#1] INFO:InputArguments -- HypoTestInverter ---- Input models:
using as S+B (null) model : ModelConfig
using as B (alternate) model : ModelConfig_with_poi_0
Doing a fixed scan in interval : 0 , 5
[#1] INFO:Eval -- HypoTestInverter::GetInterval - run a fixed scan
[#0] WARNING:InputArguments -- HypoTestInverter::RunFixedScan - xMax > upper bound, using xmax = 3
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 0
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.158989
Snapshot:
1) 0x56085605e8b0 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.158989
Snapshot:
1) 0x56085600d510 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: 0
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 0
CLs = 1 +/- 0
CLb = 1 +/- 0
CLsplusb = 1 +/- 0
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 0.6
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.178068
Snapshot:
1) 0x560856186e30 RooRealVar:: SigXsecOverSM = 0.6 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.178068
Snapshot:
1) 0x56085605ca70 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: -2.35222
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 0.6
CLs = 0.826039 +/- 0.0187037
CLb = 0.914 +/- 0.0125383
CLsplusb = 0.755 +/- 0.0136006
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 1.2
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.197147
Snapshot:
1) 0x560856171c00 RooRealVar:: SigXsecOverSM = 1.2 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.197147
Snapshot:
1) 0x560856266f90 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: -2.9364
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 1.2
CLs = 0.495652 +/- 0.018325
CLb = 0.92 +/- 0.0121326
CLsplusb = 0.456 +/- 0.01575
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 1.8
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.216225
Snapshot:
1) 0x5608561b21d0 RooRealVar:: SigXsecOverSM = 1.8 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.216225
Snapshot:
1) 0x5608561f8770 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: -2.05075
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 1.8
CLs = 0.209052 +/- 0.0137241
CLb = 0.928 +/- 0.0115599
CLsplusb = 0.194 +/- 0.0125046
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 2.4
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.235304
Snapshot:
1) 0x560856191130 RooRealVar:: SigXsecOverSM = 2.4 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.235304
Snapshot:
1) 0x5608561ef3e0 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: 0.0783908
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 2.4
CLs = 0.0503282 +/- 0.00728062
CLb = 0.914 +/- 0.0125383
CLsplusb = 0.046 +/- 0.0066245
[#1] INFO:ObjectHandling -- RooWorkspace::saveSnapshot(combined) replacing previous snapshot with name ModelConfig__snapshot
[#0] PROGRESS:Eval -- Running for SigXsecOverSM = 3
=== Using the following for ModelConfig ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.254383
Snapshot:
1) 0x56085626c2d0 RooRealVar:: SigXsecOverSM = 3 +/- 0 L(0 - 3) "SigXsecOverSM"
=== Using the following for ModelConfig_with_poi_0 ===
Observables: RooArgSet:: = (obs_x_channel1,channelCat)
Parameters of Interest: RooArgSet:: = (SigXsecOverSM)
Nuisance Parameters: RooArgSet:: = (alpha_syst2,alpha_syst3,gamma_stat_channel1_bin_0,gamma_stat_channel1_bin_1)
Global Observables: RooArgSet:: = (nominalLumi,nom_alpha_syst1,nom_alpha_syst2,nom_alpha_syst3,nom_gamma_stat_channel1_bin_0,nom_gamma_stat_channel1_bin_1)
PDF: RooSimultaneous::simPdf[ indexCat=channelCat channel1=model_channel1 ] = 0.254383
Snapshot:
1) 0x5608562fd5d0 RooRealVar:: SigXsecOverSM = 0 +/- 0 L(0 - 3) "SigXsecOverSM"
[#0] PROGRESS:Generation -- Test Statistic on data: 3.27476
[#1] INFO:InputArguments -- Profiling conditional MLEs for Null.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Null.
[#0] PROGRESS:Generation -- generated toys: 500 / 1000
[#1] INFO:InputArguments -- Profiling conditional MLEs for Alt.
[#1] INFO:InputArguments -- Using a ToyMCSampler. Now configuring for Alt.
[#0] PROGRESS:Eval -- P values for SigXsecOverSM = 3
CLs = 0.00534188 +/- 0.0023838
CLb = 0.936 +/- 0.0109457
CLsplusb = 0.005 +/- 0.00223047
Time to perform limit scan
Real time 0:00:10, CP time 10.490
The computed upper limit is: 2.40438 +/- 0.0566658
Expected upper limits, using the B (alternate) model :
expected limit (median) 1.61927
expected limit (-1 sig) 1.09318
expected limit (+1 sig) 2.24198
expected limit (-2 sig) 0.787047
expected limit (+2 sig) 2.86153
[#0] WARNING:Plotting -- Could not determine xmin and xmax of sampling distribution that was given to plot.
[#0] WARNING:Plotting -- Could not determine xmin and xmax of sampling distribution that was given to plot.
#include "TFile.h"
#include "RooWorkspace.h"
#include "RooAbsPdf.h"
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooRandom.h"
#include "TGraphErrors.h"
#include "TCanvas.h"
#include "TLine.h"
#include "TROOT.h"
#include "TSystem.h"
#include <cassert>
using namespace RooFit;
using namespace RooStats;
using std::cout, std::endl;
// structure defining the options
struct HypoTestInvOptions {
bool plotHypoTestResult = true; // plot test statistic result at each point
bool writeResult = true; // write HypoTestInverterResult in a file
TString resultFileName; // file with results (by default is built automatically using the workspace input file name)
bool optimize = true; // optimize evaluation of test statistic
bool useVectorStore = true; // convert data to use new roofit data store
bool generateBinned = false; // generate binned data sets
bool noSystematics = false; // force all systematics to be off (i.e. set all nuisance parameters as constat
// to their nominal values)
double nToysRatio = 2; // ratio Ntoys S+b/ntoysB
double maxPOI = -1; // max value used of POI (in case of auto scan)
bool enableDetailedOutput =
false; // enable detailed output with all fit information for each toys (output will be written in result file)
bool rebuild = false; // re-do extra toys for computing expected limits and rebuild test stat
// distributions (N.B this requires much more CPU (factor is equivalent to nToyToRebuild)
int nToyToRebuild = 100; // number of toys used to rebuild
int rebuildParamValues = 0; // = 0 do a profile of all the parameters on the B (alt snapshot) before performing a
// rebuild operation (default)
// = 1 use initial workspace parameters with B snapshot values
// = 2 use all initial workspace parameters with B
// Otherwise the rebuild will be performed using
int initialFit = -1; // do a first fit to the model (-1 : default, 0 skip fit, 1 do always fit)
int randomSeed = -1; // random seed (if = -1: use default value, if = 0 always random )
int nAsimovBins = 0; // number of bins in observables used for Asimov data sets (0 is the default and it is given by
// workspace, typically is 100)
bool reuseAltToys = false; // reuse same toys for alternate hypothesis (if set one gets more stable bands)
double confLevel = 0.95; // confidence level value
std::string minimizerType =
""; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
std::string massValue = ""; // extra string to tag output file of result
int printLevel = 0; // print level for debugging PL test statistics and calculators
bool useNLLOffset = false; // use NLL offset when fitting (this increase stability of fits)
};
HypoTestInvOptions optHTInv;
// internal class to run the inverter and more
namespace RooStats {
class HypoTestInvTool {
public:
HypoTestInvTool();
~HypoTestInvTool(){};
HypoTestInverterResult *RunInverter(RooWorkspace *w, const char *modelSBName, const char *modelBName,
const char *dataName, int type, int testStatType, bool useCLs, int npoints,
double poimin, double poimax, int ntoys, bool useNumberCounting = false,
const char *nuisPriorName = nullptr);
void AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType, bool useCLs, int npoints,
const char *fileNameBase = nullptr);
void SetParameter(const char *name, const char *value);
void SetParameter(const char *name, bool value);
void SetParameter(const char *name, int value);
void SetParameter(const char *name, double value);
private:
bool mPlotHypoTestResult;
bool mWriteResult;
bool mOptimize;
bool mUseVectorStore;
bool mGenerateBinned;
bool mRebuild;
bool mReuseAltToys;
bool mEnableDetOutput;
int mNToyToRebuild;
int mRebuildParamValues;
int mPrintLevel;
int mInitialFit;
int mRandomSeed;
double mNToysRatio;
double mMaxPoi;
int mAsimovBins;
std::string mMassValue;
std::string
mMinimizerType; // minimizer type (default is what is in ROOT::Math::MinimizerOptions::DefaultMinimizerType()
TString mResultFileName;
};
} // end namespace RooStats
RooStats::HypoTestInvTool::HypoTestInvTool()
: mPlotHypoTestResult(true), mWriteResult(false), mOptimize(true), mUseVectorStore(true), mGenerateBinned(false),
mEnableDetOutput(false), mRebuild(false), mReuseAltToys(false),
mNToyToRebuild(100), mRebuildParamValues(0), mPrintLevel(0), mInitialFit(-1), mRandomSeed(-1), mNToysRatio(2),
mMaxPoi(-1), mAsimovBins(0), mMassValue(""), mMinimizerType(""), mResultFileName()
{
}
void RooStats::HypoTestInvTool::SetParameter(const char *name, bool value)
{
//
// set boolean parameters
//
std::string s_name(name);
if (s_name.find("PlotHypoTestResult") != std::string::npos)
mPlotHypoTestResult = value;
if (s_name.find("WriteResult") != std::string::npos)
mWriteResult = value;
if (s_name.find("Optimize") != std::string::npos)
mOptimize = value;
if (s_name.find("UseVectorStore") != std::string::npos)
mUseVectorStore = value;
if (s_name.find("GenerateBinned") != std::string::npos)
mGenerateBinned = value;
if (s_name.find("EnableDetailedOutput") != std::string::npos)
mEnableDetOutput = value;
if (s_name.find("Rebuild") != std::string::npos)
mRebuild = value;
if (s_name.find("ReuseAltToys") != std::string::npos)
mReuseAltToys = value;
return;
}
void RooStats::HypoTestInvTool::SetParameter(const char *name, int value)
{
//
// set integer parameters
//
std::string s_name(name);
if (s_name.find("NToyToRebuild") != std::string::npos)
mNToyToRebuild = value;
if (s_name.find("RebuildParamValues") != std::string::npos)
mRebuildParamValues = value;
if (s_name.find("PrintLevel") != std::string::npos)
mPrintLevel = value;
if (s_name.find("InitialFit") != std::string::npos)
mInitialFit = value;
if (s_name.find("RandomSeed") != std::string::npos)
mRandomSeed = value;
if (s_name.find("AsimovBins") != std::string::npos)
mAsimovBins = value;
return;
}
void RooStats::HypoTestInvTool::SetParameter(const char *name, double value)
{
//
// set double precision parameters
//
std::string s_name(name);
if (s_name.find("NToysRatio") != std::string::npos)
mNToysRatio = value;
if (s_name.find("MaxPOI") != std::string::npos)
mMaxPoi = value;
return;
}
void RooStats::HypoTestInvTool::SetParameter(const char *name, const char *value)
{
//
// set string parameters
//
std::string s_name(name);
if (s_name.find("MassValue") != std::string::npos)
mMassValue.assign(value);
if (s_name.find("MinimizerType") != std::string::npos)
mMinimizerType.assign(value);
if (s_name.find("ResultFileName") != std::string::npos)
mResultFileName = value;
return;
}
void StandardHypoTestInvDemo(const char *infile = nullptr, const char *wsName = "combined",
const char *modelSBName = "ModelConfig", const char *modelBName = "",
const char *dataName = "obsData", int calculatorType = 0, int testStatType = 0,
bool useCLs = true, int npoints = 6, double poimin = 0, double poimax = 5,
int ntoys = 1000, bool useNumberCounting = false, const char *nuisPriorName = nullptr)
{
/*
Other Parameter to pass in tutorial
apart from standard for filename, ws, modelconfig and data
type = 0 Freq calculator
type = 1 Hybrid calculator
type = 2 Asymptotic calculator
type = 3 Asymptotic calculator using nominal Asimov data sets (not using fitted parameter values but nominal ones)
testStatType = 0 LEP
= 1 Tevatron
= 2 Profile Likelihood
= 3 Profile Likelihood one sided (i.e. = 0 if mu < mu_hat)
= 4 Profiel Likelihood signed ( pll = -pll if mu < mu_hat)
= 5 Max Likelihood Estimate as test statistic
= 6 Number of observed event as test statistic
useCLs scan for CLs (otherwise for CLs+b)
npoints: number of points to scan , for autoscan set npoints = -1
poimin,poimax: min/max value to scan in case of fixed scans
(if min > max, try to find automatically)
ntoys: number of toys to use
useNumberCounting: set to true when using number counting events
nuisPriorName: name of prior for the nuisance. This is often expressed as constraint term in the global model
It is needed only when using the HybridCalculator (type=1)
If not given by default the prior pdf from ModelConfig is used.
extra options are available as global parameters of the macro. They major ones are:
plotHypoTestResult plot result of tests at each point (TS distributions) (default is true)
writeResult write result of scan (default is true)
rebuild rebuild scan for expected limits (require extra toys) (default is false)
generateBinned generate binned data sets for toys (default is false) - be careful not to activate with
a too large (>=3) number of observables
nToyRatio ratio of S+B/B toys (default is 2)
*/
TString filename(infile);
if (filename.IsNull()) {
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) {
// 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
// if input file was specified but not found, quit
if (!file) {
cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
HypoTestInvTool calc;
// set parameters
calc.SetParameter("PlotHypoTestResult", optHTInv.plotHypoTestResult);
calc.SetParameter("WriteResult", optHTInv.writeResult);
calc.SetParameter("Optimize", optHTInv.optimize);
calc.SetParameter("UseVectorStore", optHTInv.useVectorStore);
calc.SetParameter("GenerateBinned", optHTInv.generateBinned);
calc.SetParameter("NToysRatio", optHTInv.nToysRatio);
calc.SetParameter("MaxPOI", optHTInv.maxPOI);
calc.SetParameter("EnableDetailedOutput", optHTInv.enableDetailedOutput);
calc.SetParameter("Rebuild", optHTInv.rebuild);
calc.SetParameter("ReuseAltToys", optHTInv.reuseAltToys);
calc.SetParameter("NToyToRebuild", optHTInv.nToyToRebuild);
calc.SetParameter("RebuildParamValues", optHTInv.rebuildParamValues);
calc.SetParameter("MassValue", optHTInv.massValue.c_str());
calc.SetParameter("MinimizerType", optHTInv.minimizerType.c_str());
calc.SetParameter("PrintLevel", optHTInv.printLevel);
calc.SetParameter("InitialFit", optHTInv.initialFit);
calc.SetParameter("ResultFileName", optHTInv.resultFileName);
calc.SetParameter("RandomSeed", optHTInv.randomSeed);
calc.SetParameter("AsimovBins", optHTInv.nAsimovBins);
// enable offset for all roostats
if (optHTInv.useNLLOffset)
RooWorkspace *w = dynamic_cast<RooWorkspace *>(file->Get(wsName));
std::cout << w << "\t" << filename << std::endl;
if (w != nullptr) {
r = calc.RunInverter(w, modelSBName, modelBName, dataName, calculatorType, testStatType, useCLs, npoints, poimin,
poimax, ntoys, useNumberCounting, nuisPriorName);
if (!r) {
std::cerr << "Error running the HypoTestInverter - Exit " << std::endl;
return;
}
} else {
// case workspace is not present look for the inverter result
std::cout << "Reading an HypoTestInverterResult with name " << wsName << " from file " << filename << std::endl;
r = dynamic_cast<HypoTestInverterResult *>(file->Get(wsName)); //
if (!r) {
std::cerr << "File " << filename << " does not contain a workspace or an HypoTestInverterResult - Exit "
<< std::endl;
file->ls();
return;
}
}
calc.AnalyzeResult(r, calculatorType, testStatType, useCLs, npoints, infile);
return;
}
void RooStats::HypoTestInvTool::AnalyzeResult(HypoTestInverterResult *r, int calculatorType, int testStatType,
bool useCLs, int npoints, const char *fileNameBase)
{
// analyze result produced by the inverter, optionally save it in a file
double lowerLimit = 0;
double llError = 0;
#if defined ROOT_SVN_VERSION && ROOT_SVN_VERSION >= 44126
if (r->IsTwoSided()) {
lowerLimit = r->LowerLimit();
llError = r->LowerLimitEstimatedError();
}
#else
lowerLimit = r->LowerLimit();
llError = r->LowerLimitEstimatedError();
#endif
double upperLimit = r->UpperLimit();
double ulError = r->UpperLimitEstimatedError();
// std::cout << "DEBUG : [ " << lowerLimit << " , " << upperLimit << " ] " << std::endl;
if (lowerLimit < upperLimit * (1. - 1.E-4) && lowerLimit != 0)
std::cout << "The computed lower limit is: " << lowerLimit << " +/- " << llError << std::endl;
std::cout << "The computed upper limit is: " << upperLimit << " +/- " << ulError << std::endl;
// compute expected limit
std::cout << "Expected upper limits, using the B (alternate) model : " << std::endl;
std::cout << " expected limit (median) " << r->GetExpectedUpperLimit(0) << std::endl;
std::cout << " expected limit (-1 sig) " << r->GetExpectedUpperLimit(-1) << std::endl;
std::cout << " expected limit (+1 sig) " << r->GetExpectedUpperLimit(1) << std::endl;
std::cout << " expected limit (-2 sig) " << r->GetExpectedUpperLimit(-2) << std::endl;
std::cout << " expected limit (+2 sig) " << r->GetExpectedUpperLimit(2) << std::endl;
// detailed output
if (mEnableDetOutput) {
mWriteResult = true;
Info("StandardHypoTestInvDemo", "detailed output will be written in output result file");
}
// write result in a file
if (r != nullptr && mWriteResult) {
// write to a file the results
const char *calcType = (calculatorType == 0) ? "Freq" : (calculatorType == 1) ? "Hybr" : "Asym";
const char *limitType = (useCLs) ? "CLs" : "Cls+b";
const char *scanType = (npoints < 0) ? "auto" : "grid";
if (mResultFileName.IsNull()) {
mResultFileName = TString::Format("%s_%s_%s_ts%d_", calcType, limitType, scanType, testStatType);
// strip the / from the filename
if (!mMassValue.empty()) {
mResultFileName += mMassValue.c_str();
mResultFileName += "_";
}
TString name = fileNameBase;
name.Replace(0, name.Last('/') + 1, "");
mResultFileName += name;
}
// get (if existing) rebuilt UL distribution
TString uldistFile = "RULDist.root";
TObject *ulDist = nullptr;
bool existULDist = !gSystem->AccessPathName(uldistFile);
if (existULDist) {
TFile *fileULDist = TFile::Open(uldistFile);
if (fileULDist)
ulDist = fileULDist->Get("RULDist");
}
TFile *fileOut = new TFile(mResultFileName, "RECREATE");
r->Write();
if (ulDist)
ulDist->Write();
Info("StandardHypoTestInvDemo", "HypoTestInverterResult has been written in the file %s", mResultFileName.Data());
fileOut->Close();
}
// plot the result ( p values vs scan points)
std::string typeName = "";
if (calculatorType == 0)
typeName = "Frequentist";
if (calculatorType == 1)
typeName = "Hybrid";
else if (calculatorType == 2 || calculatorType == 3) {
typeName = "Asymptotic";
mPlotHypoTestResult = false;
}
const char *resultName = r->GetName();
TString plotTitle = TString::Format("%s CL Scan for workspace %s", typeName.c_str(), resultName);
HypoTestInverterPlot *plot = new HypoTestInverterPlot("HTI_Result_Plot", plotTitle, r);
// plot in a new canvas with style
TString c1Name = TString::Format("%s_Scan", typeName.c_str());
TCanvas *c1 = new TCanvas(c1Name);
c1->SetLogy(false);
plot->Draw("CLb 2CL"); // plot all and Clb
// if (useCLs)
// plot->Draw("CLb 2CL"); // plot all and Clb
// else
// plot->Draw(""); // plot all and Clb
const int nEntries = r->ArraySize();
// plot test statistics distributions for the two hypothesis
if (mPlotHypoTestResult) {
TCanvas *c2 = new TCanvas("c2");
if (nEntries > 1) {
int ny = TMath::CeilNint(TMath::Sqrt(nEntries));
int nx = TMath::CeilNint(double(nEntries) / ny);
c2->Divide(nx, ny);
}
for (int i = 0; i < nEntries; i++) {
if (nEntries > 1)
c2->cd(i + 1);
SamplingDistPlot *pl = plot->MakeTestStatPlot(i);
pl->SetLogYaxis(true);
pl->Draw();
}
}
gPad = c1;
}
// internal routine to run the inverter
HypoTestInverterResult *RooStats::HypoTestInvTool::RunInverter(RooWorkspace *w, const char *modelSBName,
const char *modelBName, const char *dataName, int type,
int testStatType, bool useCLs, int npoints,
double poimin, double poimax, int ntoys,
bool useNumberCounting, const char *nuisPriorName)
{
std::cout << "Running HypoTestInverter on the workspace " << w->GetName() << std::endl;
w->Print();
RooAbsData *data = w->data(dataName);
if (!data) {
Error("StandardHypoTestDemo", "Not existing data %s", dataName);
return nullptr;
} else
std::cout << "Using data set " << dataName << std::endl;
if (mUseVectorStore) {
data->convertToVectorStore();
}
// get models from WS
// get the modelConfig out of the file
ModelConfig *bModel = (ModelConfig *)w->obj(modelBName);
ModelConfig *sbModel = (ModelConfig *)w->obj(modelSBName);
if (!sbModel) {
Error("StandardHypoTestDemo", "Not existing ModelConfig %s", modelSBName);
return nullptr;
}
// check the model
if (!sbModel->GetPdf()) {
Error("StandardHypoTestDemo", "Model %s has no pdf ", modelSBName);
return nullptr;
}
if (!sbModel->GetParametersOfInterest()) {
Error("StandardHypoTestDemo", "Model %s has no poi ", modelSBName);
return nullptr;
}
if (!sbModel->GetObservables()) {
Error("StandardHypoTestInvDemo", "Model %s has no observables ", modelSBName);
return nullptr;
}
if (!sbModel->GetSnapshot()) {
Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi", modelSBName);
sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
}
// case of no systematics
// remove nuisance parameters from model
if (optHTInv.noSystematics) {
const RooArgSet *nuisPar = sbModel->GetNuisanceParameters();
if (nuisPar && nuisPar->getSize() > 0) {
std::cout << "StandardHypoTestInvDemo"
<< " - Switch off all systematics by setting them constant to their initial values" << std::endl;
}
if (bModel) {
const RooArgSet *bnuisPar = bModel->GetNuisanceParameters();
if (bnuisPar)
}
}
if (!bModel || bModel == sbModel) {
Info("StandardHypoTestInvDemo", "The background model %s does not exist", modelBName);
Info("StandardHypoTestInvDemo", "Copy it from ModelConfig %s and set POI to zero", modelSBName);
bModel = (ModelConfig *)sbModel->Clone();
bModel->SetName(TString(modelSBName) + TString("_with_poi_0"));
RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
if (!var)
return nullptr;
double oldval = var->getVal();
var->setVal(0);
bModel->SetSnapshot(RooArgSet(*var));
var->setVal(oldval);
} else {
if (!bModel->GetSnapshot()) {
Info("StandardHypoTestInvDemo", "Model %s has no snapshot - make one using model poi and 0 values ",
modelBName);
RooRealVar *var = dynamic_cast<RooRealVar *>(bModel->GetParametersOfInterest()->first());
if (var) {
double oldval = var->getVal();
var->setVal(0);
bModel->SetSnapshot(RooArgSet(*var));
var->setVal(oldval);
} else {
Error("StandardHypoTestInvDemo", "Model %s has no valid poi", modelBName);
return nullptr;
}
}
}
// check model has global observables when there are nuisance pdf
// for the hybrid case the globals are not needed
if (type != 1) {
bool hasNuisParam = (sbModel->GetNuisanceParameters() && sbModel->GetNuisanceParameters()->getSize() > 0);
bool hasGlobalObs = (sbModel->GetGlobalObservables() && sbModel->GetGlobalObservables()->getSize() > 0);
if (hasNuisParam && !hasGlobalObs) {
// try to see if model has nuisance parameters first
RooAbsPdf *constrPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisanceConstraintPdf_sbmodel");
if (constrPdf) {
Warning("StandardHypoTestInvDemo", "Model %s has nuisance parameters but no global observables associated",
sbModel->GetName());
Warning("StandardHypoTestInvDemo",
"\tThe effect of the nuisance parameters will not be treated correctly ");
}
}
}
// save all initial parameters of the model including the global observables
RooArgSet initialParameters;
std::unique_ptr<RooArgSet> allParams{sbModel->GetPdf()->getParameters(*data)};
allParams->snapshot(initialParameters);
// run first a data fit
const RooArgSet *poiSet = sbModel->GetParametersOfInterest();
RooRealVar *poi = (RooRealVar *)poiSet->first();
std::cout << "StandardHypoTestInvDemo : POI initial value: " << poi->GetName() << " = " << poi->getVal()
<< std::endl;
// fit the data first (need to use constraint )
bool doFit = mInitialFit;
if (testStatType == 0 && mInitialFit == -1)
doFit = false; // case of LEP test statistic
if (type == 3 && mInitialFit == -1)
doFit = false; // case of Asymptoticcalculator with nominal Asimov
double poihat = 0;
if (mMinimizerType.empty())
else
Info("StandardHypoTestInvDemo", "Using %s as minimizer for computing the test statistic",
if (doFit) {
// do the fit : By doing a fit the POI snapshot (for S+B) is set to the fit value
// and the nuisance parameters nominal values will be set to the fit value.
// This is relevant when using LEP test statistics
Info("StandardHypoTestInvDemo", " Doing a first fit to the observed data ");
RooArgSet constrainParams;
if (sbModel->GetNuisanceParameters())
constrainParams.add(*sbModel->GetNuisanceParameters());
tw.Start();
std::unique_ptr<RooFitResult> fitres{sbModel->GetPdf()->fitTo(
*data, InitialHesse(false), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0),
PrintLevel(mPrintLevel), Constrain(constrainParams), Save(true), Offset(RooStats::IsNLLOffset()))};
if (fitres->status() != 0) {
Warning("StandardHypoTestInvDemo",
"Fit to the model failed - try with strategy 1 and perform first an Hesse computation");
fitres = std::unique_ptr<RooFitResult>{sbModel->GetPdf()->fitTo(
*data, InitialHesse(true), Hesse(false), Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(1),
PrintLevel(mPrintLevel + 1), Constrain(constrainParams), Save(true), Offset(RooStats::IsNLLOffset()))};
}
if (fitres->status() != 0)
Warning("StandardHypoTestInvDemo", " Fit still failed - continue anyway.....");
poihat = poi->getVal();
std::cout << "StandardHypoTestInvDemo - Best Fit value : " << poi->GetName() << " = " << poihat << " +/- "
<< poi->getError() << std::endl;
std::cout << "Time for fitting : ";
tw.Print();
// save best fit value in the poi snapshot
sbModel->SetSnapshot(*sbModel->GetParametersOfInterest());
std::cout << "StandardHypoTestInvo: snapshot of S+B Model " << sbModel->GetName()
<< " is set to the best fit value" << std::endl;
}
// print a message in case of LEP test statistics because it affects result by doing or not doing a fit
if (testStatType == 0) {
if (!doFit)
Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit is not done and the TS will use "
"the nuisances at the model value");
else
Info("StandardHypoTestInvDemo", "Using LEP test statistic - an initial fit has been done and the TS will use "
"the nuisances at the best fit value");
}
// build test statistics and hypotest calculators for running the inverter
SimpleLikelihoodRatioTestStat slrts(*sbModel->GetPdf(), *bModel->GetPdf());
// null parameters must includes snapshot of poi plus the nuisance values
RooArgSet nullParams(*sbModel->GetSnapshot());
if (sbModel->GetNuisanceParameters())
nullParams.add(*sbModel->GetNuisanceParameters());
if (sbModel->GetSnapshot())
slrts.SetNullParameters(nullParams);
RooArgSet altParams(*bModel->GetSnapshot());
if (bModel->GetNuisanceParameters())
altParams.add(*bModel->GetNuisanceParameters());
if (bModel->GetSnapshot())
slrts.SetAltParameters(altParams);
if (mEnableDetOutput)
slrts.EnableDetailedOutput();
// ratio of profile likelihood - need to pass snapshot for the alt
RatioOfProfiledLikelihoodsTestStat ropl(*sbModel->GetPdf(), *bModel->GetPdf(), bModel->GetSnapshot());
ropl.SetSubtractMLE(false);
if (testStatType == 11)
ropl.SetSubtractMLE(true);
ropl.SetPrintLevel(mPrintLevel);
ropl.SetMinimizer(mMinimizerType.c_str());
if (mEnableDetOutput)
ropl.EnableDetailedOutput();
ProfileLikelihoodTestStat profll(*sbModel->GetPdf());
if (testStatType == 3)
profll.SetOneSided(true);
if (testStatType == 4)
profll.SetSigned(true);
profll.SetMinimizer(mMinimizerType.c_str());
profll.SetPrintLevel(mPrintLevel);
if (mEnableDetOutput)
profll.EnableDetailedOutput();
profll.SetReuseNLL(mOptimize);
slrts.SetReuseNLL(mOptimize);
ropl.SetReuseNLL(mOptimize);
if (mOptimize) {
profll.SetStrategy(0);
ropl.SetStrategy(0);
}
if (mMaxPoi > 0)
poi->setMax(mMaxPoi); // increase limit
MaxLikelihoodEstimateTestStat maxll(*sbModel->GetPdf(), *poi);
AsymptoticCalculator::SetPrintLevel(mPrintLevel);
// create the HypoTest calculator class
if (type == 0)
hc = new FrequentistCalculator(*data, *bModel, *sbModel);
else if (type == 1)
hc = new HybridCalculator(*data, *bModel, *sbModel);
// else if (type == 2 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false, mAsimovBins);
// else if (type == 3 ) hc = new AsymptoticCalculator(*data, *bModel, *sbModel, true, mAsimovBins); // for using
// Asimov data generated with nominal values
else if (type == 2)
hc = new AsymptoticCalculator(*data, *bModel, *sbModel, false);
else if (type == 3)
hc = new AsymptoticCalculator(*data, *bModel, *sbModel,
true); // for using Asimov data generated with nominal values
else {
Error("StandardHypoTestInvDemo", "Invalid - calculator type = %d supported values are only :\n\t\t\t 0 "
"(Frequentist) , 1 (Hybrid) , 2 (Asymptotic) ",
type);
return nullptr;
}
// set the test statistic
TestStatistic *testStat = nullptr;
if (testStatType == 0)
testStat = &slrts;
if (testStatType == 1 || testStatType == 11)
testStat = &ropl;
if (testStatType == 2 || testStatType == 3 || testStatType == 4)
testStat = &profll;
if (testStatType == 5)
testStat = &maxll;
if (testStatType == 6)
testStat = &nevtts;
if (testStat == nullptr) {
Error("StandardHypoTestInvDemo", "Invalid - test statistic type = %d supported values are only :\n\t\t\t 0 (SLR) "
", 1 (Tevatron) , 2 (PLR), 3 (PLR1), 4(MLE)",
testStatType);
return nullptr;
}
if (toymcs && (type == 0 || type == 1)) {
// look if pdf is number counting or extended
if (sbModel->GetPdf()->canBeExtended()) {
if (useNumberCounting)
Warning("StandardHypoTestInvDemo", "Pdf is extended: but number counting flag is set: ignore it ");
} else {
// for not extended pdf
if (!useNumberCounting) {
int nEvents = data->numEntries();
Info("StandardHypoTestInvDemo",
"Pdf is not extended: number of events to generate taken from observed data set is %d", nEvents);
toymcs->SetNEventsPerToy(nEvents);
} else {
Info("StandardHypoTestInvDemo", "using a number counting pdf");
toymcs->SetNEventsPerToy(1);
}
}
toymcs->SetTestStatistic(testStat);
if (data->isWeighted() && !mGenerateBinned) {
Info("StandardHypoTestInvDemo", "Data set is weighted, nentries = %d and sum of weights = %8.1f but toy "
"generation is unbinned - it would be faster to set mGenerateBinned to true\n",
data->numEntries(), data->sumEntries());
}
toymcs->SetGenerateBinned(mGenerateBinned);
toymcs->SetUseMultiGen(mOptimize);
if (mGenerateBinned && sbModel->GetObservables()->getSize() > 2) {
Warning("StandardHypoTestInvDemo", "generate binned is activated but the number of observable is %d. Too much "
"memory could be needed for allocating all the bins",
sbModel->GetObservables()->getSize());
}
// set the random seed if needed
if (mRandomSeed >= 0)
}
// specify if need to re-use same toys
if (mReuseAltToys) {
}
if (type == 1) {
HybridCalculator *hhc = dynamic_cast<HybridCalculator *>(hc);
assert(hhc);
hhc->SetToys(ntoys, ntoys / mNToysRatio); // can use less ntoys for b hypothesis
// remove global observables from ModelConfig (this is probably not needed anymore in 5.32)
// check for nuisance prior pdf in case of nuisance parameters
if (bModel->GetNuisanceParameters() || sbModel->GetNuisanceParameters()) {
// fix for using multigen (does not work in this case)
toymcs->SetUseMultiGen(false);
ToyMCSampler::SetAlwaysUseMultiGen(false);
RooAbsPdf *nuisPdf = nullptr;
if (nuisPriorName)
nuisPdf = w->pdf(nuisPriorName);
// use prior defined first in bModel (then in SbModel)
if (!nuisPdf) {
Info("StandardHypoTestInvDemo",
"No nuisance pdf given for the HybridCalculator - try to deduce pdf from the model");
if (bModel->GetPdf() && bModel->GetObservables())
nuisPdf = RooStats::MakeNuisancePdf(*bModel, "nuisancePdf_bmodel");
else
nuisPdf = RooStats::MakeNuisancePdf(*sbModel, "nuisancePdf_sbmodel");
}
if (!nuisPdf) {
if (bModel->GetPriorPdf()) {
nuisPdf = bModel->GetPriorPdf();
Info("StandardHypoTestInvDemo",
"No nuisance pdf given - try to use %s that is defined as a prior pdf in the B model",
nuisPdf->GetName());
} else {
Error("StandardHypoTestInvDemo", "Cannot run Hybrid calculator because no prior on the nuisance "
"parameter is specified or can be derived");
return nullptr;
}
}
assert(nuisPdf);
Info("StandardHypoTestInvDemo", "Using as nuisance Pdf ... ");
nuisPdf->Print();
const RooArgSet *nuisParams =
std::unique_ptr<RooArgSet> np{nuisPdf->getObservables(*nuisParams)};
if (np->getSize() == 0) {
Warning("StandardHypoTestInvDemo",
"Prior nuisance does not depend on nuisance parameters. They will be smeared in their full range");
}
hhc->ForcePriorNuisanceAlt(*nuisPdf);
hhc->ForcePriorNuisanceNull(*nuisPdf);
}
} else if (type == 2 || type == 3) {
if (testStatType == 3)
((AsymptoticCalculator *)hc)->SetOneSided(true);
if (testStatType != 2 && testStatType != 3)
Warning("StandardHypoTestInvDemo",
"Only the PL test statistic can be used with AsymptoticCalculator - use by default a two-sided PL");
} else if (type == 0) {
((FrequentistCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
// store also the fit information for each poi point used by calculator based on toys
if (mEnableDetOutput)
((FrequentistCalculator *)hc)->StoreFitInfo(true);
} else if (type == 1) {
((HybridCalculator *)hc)->SetToys(ntoys, ntoys / mNToysRatio);
// store also the fit information for each poi point used by calculator based on toys
// if (mEnableDetOutput) ((HybridCalculator*) hc)->StoreFitInfo(true);
}
// Get the result
HypoTestInverter calc(*hc);
calc.SetConfidenceLevel(optHTInv.confLevel);
calc.UseCLs(useCLs);
calc.SetVerbose(true);
if (npoints > 0) {
if (poimin > poimax) {
// if no min/max given scan between MLE and +4 sigma
poimin = int(poihat);
poimax = int(poihat + 4 * poi->getError());
}
std::cout << "Doing a fixed scan in interval : " << poimin << " , " << poimax << std::endl;
calc.SetFixedScan(npoints, poimin, poimax);
} else {
// poi->setMax(10*int( (poihat+ 10 *poi->getError() )/10 ) );
std::cout << "Doing an automatic scan in interval : " << poi->getMin() << " , " << poi->getMax() << std::endl;
}
tw.Start();
HypoTestInverterResult *r = calc.GetInterval();
std::cout << "Time to perform limit scan \n";
tw.Print();
if (mRebuild) {
std::cout << "\n***************************************************************\n";
std::cout << "Rebuild the upper limit distribution by re-generating new set of pseudo-experiment and re-compute "
"for each of them a new upper limit\n\n";
allParams = std::unique_ptr<RooArgSet>{sbModel->GetPdf()->getParameters(*data)};
// define on which value of nuisance parameters to do the rebuild
// default is best fit value for bmodel snapshot
if (mRebuildParamValues != 0) {
// set all parameters to their initial workspace values
allParams->assign(initialParameters);
}
if (mRebuildParamValues == 0 || mRebuildParamValues == 1) {
RooArgSet constrainParams;
if (sbModel->GetNuisanceParameters())
constrainParams.add(*sbModel->GetNuisanceParameters());
const RooArgSet *poiModel = sbModel->GetParametersOfInterest();
bModel->LoadSnapshot();
// do a profile using the B model snapshot
if (mRebuildParamValues == 0) {
RooStats::SetAllConstant(*poiModel, true);
sbModel->GetPdf()->fitTo(*data, InitialHesse(false), Hesse(false),
Minimizer(mMinimizerType.c_str(), "Migrad"), Strategy(0), PrintLevel(mPrintLevel),
Constrain(constrainParams), Offset(RooStats::IsNLLOffset()));
std::cout << "rebuild using fitted parameter value for B-model snapshot" << std::endl;
constrainParams.Print("v");
RooStats::SetAllConstant(*poiModel, false);
}
}
std::cout << "StandardHypoTestInvDemo: Initial parameters used for rebuilding: ";
RooStats::PrintListContent(*allParams, std::cout);
tw.Start();
SamplingDistribution *limDist = calc.GetUpperLimitDistribution(true, mNToyToRebuild);
std::cout << "Time to rebuild distributions " << std::endl;
tw.Print();
if (limDist) {
std::cout << "Expected limits after rebuild distribution " << std::endl;
std::cout << "expected upper limit (median of limit distribution) " << limDist->InverseCDF(0.5) << std::endl;
std::cout << "expected -1 sig limit (0.16% quantile of limit dist) "
<< limDist->InverseCDF(ROOT::Math::normal_cdf(-1)) << std::endl;
std::cout << "expected +1 sig limit (0.84% quantile of limit dist) "
<< limDist->InverseCDF(ROOT::Math::normal_cdf(1)) << std::endl;
std::cout << "expected -2 sig limit (.025% quantile of limit dist) "
<< limDist->InverseCDF(ROOT::Math::normal_cdf(-2)) << std::endl;
std::cout << "expected +2 sig limit (.975% quantile of limit dist) "
<< limDist->InverseCDF(ROOT::Math::normal_cdf(2)) << std::endl;
// Plot the upper limit distribution
SamplingDistPlot limPlot((mNToyToRebuild < 200) ? 50 : 100);
limPlot.AddSamplingDistribution(limDist);
limPlot.GetTH1F()->SetStats(true); // display statistics
limPlot.SetLineColor(kBlue);
new TCanvas("limPlot", "Upper Limit Distribution");
limPlot.Draw();
/// save result in a file
limDist->SetName("RULDist");
TFile *fileOut = new TFile("RULDist.root", "RECREATE");
limDist->Write();
fileOut->Close();
// update r to a new updated result object containing the rebuilt expected p-values distributions
// (it will not recompute the expected limit)
if (r)
delete r; // need to delete previous object since GetInterval will return a cloned copy
r = calc.GetInterval();
} else
std::cout << "ERROR : failed to re-build distributions " << std::endl;
}
return r;
}
void ReadResult(const char *fileName, const char *resultName = "", bool useCLs = true)
{
// read a previous stored result from a file given the result name
StandardHypoTestInvDemo(fileName, resultName, "", "", "", 0, 0, useCLs);
}
#ifdef USE_AS_MAIN
int main()
{
}
#endif
int main()
Definition Prototype.cxx:12
@ kBlue
Definition Rtypes.h:66
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
winID h TVirtualViewer3D TVirtualGLPainter char TVirtualGLPainter plot
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
#define gROOT
Definition TROOT.h:406
R__EXTERN TSystem * gSystem
Definition TSystem.h:561
#define gPad
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
static void SetDefaultStrategy(int strat)
Set the default strategy.
static const std::string & DefaultMinimizerType()
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:263
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooAbsArg * first() const
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
static void setDefaultStorageType(StorageType s)
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:218
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
Definition RooAbsPdf.h:157
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
static RooMsgService & instance()
Return reference to singleton instance.
StreamConfig & getStream(Int_t id)
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:48
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
double getError() const
Definition RooRealVar.h:58
void setMax(const char *name, double value)
Set maximum of name range to given value.
Hypothesis Test Calculator based on the asymptotic formulae for the profile likelihood ratio.
Does a frequentist hypothesis test.
Same purpose as HybridCalculatorOriginal, but different implementation.
virtual void ForcePriorNuisanceNull(RooAbsPdf &priorNuisance)
Override the distribution used for marginalizing nuisance parameters that is inferred from ModelConfi...
virtual void ForcePriorNuisanceAlt(RooAbsPdf &priorNuisance)
void SetToys(int toysNull, int toysAlt)
set number of toys
Common base class for the Hypothesis Test Calculators.
TestStatSampler * GetTestStatSampler(void) const
Returns instance of TestStatSampler.
void UseSameAltToys()
Set this for re-using always the same toys for alternate hypothesis in case of calls at different nul...
Class to plot a HypoTestInverterResult, the output of the HypoTestInverter calculator.
HypoTestInverterResult class holds the array of hypothesis test results and compute a confidence inte...
A class for performing a hypothesis test inversion by scanning the hypothesis test results of a HypoT...
MaxLikelihoodEstimateTestStat: TestStatistic that returns maximum likelihood estimate of a specified ...
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
virtual void SetSnapshot(const RooArgSet &set)
Set parameter values for a particular hypothesis if using a common PDF by saving a snapshot in the wo...
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
ModelConfig * Clone(const char *name="") const override
clone
Definition ModelConfig.h:57
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
void LoadSnapshot() const
load the snapshot from ws if it exists
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
const RooArgSet * GetSnapshot() const
get RooArgSet for parameters for a particular hypothesis (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
RooAbsPdf * GetPriorPdf() const
get parameters prior pdf (return nullptr if not existing)
virtual void SetGlobalObservables(const RooArgSet &set)
Specify the global observables.
NumEventsTestStat is a simple implementation of the TestStatistic interface used for simple number co...
ProfileLikelihoodTestStat is an implementation of the TestStatistic interface that calculates the pro...
TestStatistic that returns the ratio of profiled likelihoods.
This class provides simple and straightforward utilities to plot SamplingDistribution objects.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
void SetLineColor(Color_t color, const SamplingDistribution *samplDist=nullptr)
Sets line color for given sampling distribution and fill color for the associated shaded TH1F.
void SetLogYaxis(bool ly)
changes plot to log scale on y axis
double AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions="NORMALIZE HIST")
adds the sampling distribution and returns the scale factor
TH1F * GetTH1F(const SamplingDistribution *samplDist=nullptr)
Returns the TH1F associated with the give SamplingDistribution.
This class simply holds a sampling distribution of some test statistic.
double InverseCDF(double pvalue)
get the inverse of the Cumulative distribution function
TestStatistic class that returns -log(L[null] / L[alt]) where L is the likelihood.
TestStatistic is an interface class to provide a facility for construction test statistics distributi...
ToyMCSampler is an implementation of the TestStatSampler interface.
virtual void SetTestStatistic(TestStatistic *testStatistic, unsigned int i)
Set the TestStatistic (want the argument to be a function of the data & parameter points.
void SetGenerateBinned(bool binned=true)
control to use bin data generation (=> see RooFit::AllBinned() option)
virtual void SetNEventsPerToy(const Int_t nevents)
Forces the generation of exactly n events even for extended PDFs.
void SetUseMultiGen(bool flag)
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
TObject * Get(const char *namecycle) override
Return pointer to object identified by namecycle.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
void ls(Option_t *option="") const override
List file contents.
Definition TFile.cxx:1454
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4086
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition TH1.cxx:9023
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
Mother of all ROOT objects.
Definition TObject.h:41
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:898
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:292
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition TRandom.cxx:615
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
Basic string class.
Definition TString.h:139
TString & Replace(Ssiz_t pos, Ssiz_t n, const char *s)
Definition TString.h:694
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition TSystem.cxx:1296
RooCmdArg InitialHesse(bool flag=true)
RooCmdArg Offset(std::string const &mode)
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg Minimizer(const char *type, const char *alg=nullptr)
RooCmdArg Hesse(bool flag=true)
RooCmdArg Strategy(Int_t code)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
double normal_cdf(double x, double sigma=1, double x0=0)
Cumulative distribution function of the normal (Gaussian) distribution (lower tail).
return c1
Definition legend1.C:41
return c2
Definition legend2.C:14
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
@ NumIntegration
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
bool SetAllConstant(const RooAbsCollection &coll, bool constant=true)
utility function to set all variable constant in a collection (from G.
void RemoveConstantParameters(RooArgSet *set)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
void UseNLLOffset(bool on)
function to set a global flag in RooStats to use NLL offset when performing nll computations Note tha...
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
useful function to print in one line the content of a set with their values
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:678
void removeTopic(RooFit::MsgTopic oldTopic)
Author
Lorenzo Moneta

Definition in file StandardHypoTestInvDemo.C.