50void rs101_limitexample()
63 wspace.
factory(
"Poisson::countingModel(obs[150,0,300], "
64 "sum(s[50,0,120]*ratioSigEff[1.,0,3.],b[100]*ratioBkgEff[1.,0.,3.]))");
67 wspace.
factory(
"Gaussian::sigConstraint(gSigEff[1,0,3],ratioSigEff,0.05)");
68 wspace.
factory(
"Gaussian::bkgConstraint(gSigBkg[1,0,3],ratioBkgEff,0.2)");
69 wspace.
factory(
"PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)");
72 RooAbsPdf *modelWithConstraints = wspace.
pdf(
"modelWithConstraints");
91 RooDataSet dataOrig{
"exampleData",
"exampleData", *obs};
94 RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
103 modelConfig.SetPdf(*modelWithConstraints);
104 modelConfig.SetParametersOfInterest(paramOfInterest);
105 modelConfig.SetNuisanceParameters(constrainedParams);
106 modelConfig.SetObservables(*obs);
107 modelConfig.SetGlobalObservables(
RooArgSet(*gSigEff, *gSigBkg));
108 modelConfig.SetName(
"ModelConfig");
109 wspace.
import(modelConfig);
120 plc.SetTestSize(.05);
121 std::unique_ptr<LikelihoodInterval> lrinterval{
static_cast<LikelihoodInterval*
>(plc.GetInterval())};
124 auto dataCanvas =
new TCanvas(
"dataCanvas");
125 dataCanvas->Divide(2, 1);
129 plotInt.SetTitle(
"Profile Likelihood Ratio and Posterior for S");
134 fc.UseAdaptiveSampling(
true);
135 fc.FluctuateNumDataEntries(
false);
139 std::unique_ptr<PointSetInterval> fcint{
static_cast<PointSetInterval*
>(fc.GetInterval())};
141 std::unique_ptr<RooFitResult> fit{modelWithConstraints->
fitTo(data,
Save(
true),
PrintLevel(-1))};
154 mc.SetNumIters(20000);
156 mc.SetNumBurnInSteps(40);
157 mc.SetProposalFunction(*pdfProp);
158 mc.SetLeftSideTailFraction(0.5);
159 std::unique_ptr<MCMCInterval> mcInt{
static_cast<MCMCInterval *
>(mc.GetInterval())};
162 std::cout <<
"Profile lower limit on s = " << lrinterval->LowerLimit(*s) << std::endl;
163 std::cout <<
"Profile upper limit on s = " << lrinterval->UpperLimit(*s) << std::endl;
167 double fcul = fcint->UpperLimit(*s);
168 double fcll = fcint->LowerLimit(*s);
169 std::cout <<
"FC lower limit on s = " << fcll << std::endl;
170 std::cout <<
"FC upper limit on s = " << fcul << std::endl;
171 auto fcllLine =
new TLine(fcll, 0, fcll, 1);
172 auto fculLine =
new TLine(fcul, 0, fcul, 1);
173 fcllLine->SetLineColor(
kRed);
174 fculLine->SetLineColor(
kRed);
175 fcllLine->Draw(
"same");
176 fculLine->Draw(
"same");
177 dataCanvas->Update();
183 mcPlot.SetLineWidth(2);
186 double mcul = mcInt->UpperLimit(*s);
187 double mcll = mcInt->LowerLimit(*s);
188 std::cout <<
"MCMC lower limit on s = " << mcll << std::endl;
189 std::cout <<
"MCMC upper limit on s = " << mcul << std::endl;
190 std::cout <<
"MCMC Actual confidence level: " << mcInt->GetActualConfidenceLevel() << std::endl;
195 std::unique_ptr<RooDataSet> chainData{mcInt->GetChainAsDataSet()};
198 std::cout <<
"plotting the chain data - nentries = " << chainData->numEntries() << std::endl;
201 chain->SetMarkerStyle(6);
202 chain->SetMarkerColor(
kRed);
204 chain->Draw(
"s:ratioSigEff:ratioBkgEff",
"nll_MarkovChain_local_",
"box");
209 std::cout <<
"plotting the scanned points used in the frequentist construction - npoints = "
216 for (
int ievt = 0; ievt < parScanData->
numEntries(); ++ievt) {
221 gr->SetPoint(ievt,
x,
y, z);
224 gr->SetMarkerStyle(24);
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract interface for all probability density functions.
RooFit::OwningPtr< RooFitResult > fitTo(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Fit PDF to given dataset.
void setConstant(bool value=true)
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Container class to hold unbinned data.
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
Variable that can be changed from the outside.
void setVal(double value) override
Set value of variable to 'value'.
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
This class provides simple and straightforward utilities to plot a LikelihoodInterval object.
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
This class provides simple and straightforward utilities to plot a MCMCInterval object.
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
PointSetInterval is a concrete implementation of the ConfInterval interface.
The ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface cla...
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
virtual void SetCovMatrix(const TMatrixDSym &covMatrix)
set the covariance matrix to use for a multi-variate Gaussian proposal
virtual ProposalFunction * GetProposalFunction()
Get the ProposalFunction that we've been designing.
virtual void SetVariables(RooArgList &vars)
virtual void SetCacheSize(Int_t size)
virtual void SetUpdateProposalParameters(bool updateParams)
Persistable container for RooFit projects.
void Print(Option_t *opts=nullptr) const override
Print contents of the workspace.
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooFactoryWSTool & factory()
Return instance to factory tool.
RooRealVar * var(RooStringView name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
RooAbsData * data(RooStringView name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
bool import(const RooAbsArg &arg, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={})
Import a RooAbsArg object, e.g.
Use the TLine constructor to create a simple line.
const char * GetName() const override
Returns name of object.
virtual void SetName(const char *name)
Set the name of the TNamed.
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Stop()
Stop the stopwatch.
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
A TTree represents a columnar dataset.
RooCmdArg Constrain(const RooArgSet ¶ms)
RooCmdArg Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooStats::ModelConfig ModelConfig
Namespace for the RooStats classes.
TTree * GetAsTTree(TString name, TString desc, const RooDataSet &data)
Create a TTree with the given name and description. All RooRealVars in the RooDataSet are represented...