14#if !defined(R__ALPHA) && !defined(R__SOLARIS) && !defined(R__ACC) && !defined(R__FBSD)
31using std::endl, std::vector;
35 template<
class listT,
class stringT>
void getParameterNames(
const listT*
l,std::vector<stringT>& names){
38 for (
auto const *obj : *
l) {
39 names.push_back(obj->GetName());
43 for(
const auto&
p:names){
63 double sb2 = sigma_b*sigma_b;
68 double bpsb2 =
b + sb2;
71 double za2 = 2.*( (spb)* std::log( ( spb)*(bpsb2)/(b2+ spb*sb2) ) -
72 (b2/sb2) * std::log(1. + ( sb2 * s)/(
b * bpsb2) ) );
77 double za2 = 2.*( (s+
b) * std::log(1. + s/
b) -s );
78 return std::sqrt(za2);
94 if (
auto prod =
dynamic_cast<RooProdPdf *
>(&pdf)) {
96 for (
int i = 0,
n = list.
size(); i <
n; ++i) {
102 auto iter = pdf.
servers().begin();
103 assert(iter != pdf.
servers().end());
107 assert(sim !=
nullptr);
109 for (
int ic = 0, nc = cat->
numBins((
const char *)
nullptr); ic < nc; ++ic) {
113 if (catPdf !=
nullptr)
FactorizePdf(observables, *catPdf, obsTerms, constraints);
119 if (!constraints.
contains(pdf)) constraints.
add(pdf);
128 oocoutE(
nullptr,InputArguments) <<
"RooStatsUtils::FactorizePdf - invalid input model: missing observables" << endl;
140 if(constraints.
empty()) {
141 oocoutW(
nullptr, Eval) <<
"RooStatsUtils::MakeNuisancePdf - no constraints found on nuisance parameters in the input model" << endl;
150 oocoutE(
nullptr, InputArguments) <<
"RooStatsUtils::MakeNuisancePdf - invalid input model: missing pdf and/or observables" << endl;
158 if (
auto prod =
dynamic_cast<RooProdPdf *
>(&pdf)) {
162 for (
int i = 0,
n = list.
size(); i <
n; ++i) {
165 if(newPdfi !=
nullptr) newList.
add(*newPdfi);
168 if (newList.
empty()) {
171 }
else if (newList.
size() == 1) {
181 auto iter = pdf.
servers().begin();
183 auto uPdf =
dynamic_cast<RooAbsPdf *
>(*(iter++));
184 auto extended_term =
dynamic_cast<RooAbsReal *
>(*(iter++));
185 assert(uPdf !=
nullptr);
186 assert(extended_term !=
nullptr);
187 assert(iter == pdf.
servers().end());
190 if(newUPdf ==
nullptr)
return nullptr;
196 assert(sim !=
nullptr);
200 for (
int ic = 0, nc = cat->
numBins((
const char *)
nullptr); ic < nc; ++ic) {
206 if (newPdf ==
nullptr) {
delete cat;
return nullptr; }
207 pdfList.
add(*newPdf);
223 if(!unconstrainedPdf) {
224 oocoutE(
nullptr, InputArguments) <<
"RooStats::MakeUnconstrainedPdf - invalid observable list passed (observables not found in original pdf) or invalid pdf passed (without observables)" << endl;
228 return unconstrainedPdf;
234 oocoutE(
nullptr, InputArguments) <<
"RooStatsUtils::MakeUnconstrainedPdf - invalid input model: missing pdf and/or observables" << endl;
247 BranchStore(
const vector<TString> ¶ms = vector<TString>(),
double _inval = -999.) :
fInval(_inval)
249 for (
unsigned int i = 0; i < params.size(); i++)
255 for(std::map<TString, double>::iterator it =
fVarVals.begin();it!=
fVarVals.end();++it) {
264 for(std::map<TString, double>::iterator it =
fVarVals.begin();it!=
fVarVals.end();++it) {
270 for(std::map<TString, double>::iterator it =
fVarVals.begin();it!=
fVarVals.end();++it) {
278 if (
data.numEntries() == 0) {
282 for (
auto *rvar : dynamic_range_cast<RooRealVar *>(*
data.get(0))) {
285 V.push_back(rvar->GetName());
286 if (rvar->hasAsymError()) {
290 else if (rvar->hasError()) {
301 for(
int entry = 0;entry<
data.numEntries();entry++) {
303 for (
auto const *rvar : dynamic_range_cast<RooRealVar *>(*
data.get(entry))) {
306 bs->
fVarVals[rvar->GetName()] = rvar->getValV();
307 if (rvar->hasAsymError()) {
311 else if (rvar->hasError()) {
332 for (std::size_t i = 0; i<
l.size(); ++i) {
348 bool copySnapshots,
const char *mcname,
349 const char *newmcname) {
354 for (
auto it : objects) {
360 throw std::runtime_error(
"unable to retrieve ModelConfig");
365 std::vector<TString> poilist;
366 std::vector<TString> nplist;
367 std::vector<TString> obslist;
368 std::vector<TString> globobslist;
392 for (
auto* snap : static_range_cast<RooArgSet const*>(oldWS->
getSnapshots())) {
401 for (
auto d :
data) {
406 ::getArgs(newWS, poilist, poiset);
408 ::getArgs(newWS, nplist, npset);
410 ::getArgs(newWS, obslist, obsset);
412 ::getArgs(newWS, globobslist, globobsset);
#define NamespaceImp(name)
winID h TVirtualViewer3D TVirtualGLPainter p
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 Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void on
TIterator Use servers() and begin()
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
virtual TObject * clone(const char *newname=nullptr) const =0
void SetName(const char *name) override
Set the name of the TNamed.
Abstract base class for objects that represent a discrete value that can be set from the outside,...
Int_t numBins(const char *rangeName=nullptr) const override
Return the number of fit bins ( = number of types )
void setBin(Int_t ibin, const char *rangeName=nullptr) override
Set category to i-th fit bin, which is the i-th registered state.
virtual const char * getCurrentLabel() const
Return label string of current state.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
Abstract interface for all probability density functions.
Abstract base class for objects that represent a real value and implements functionality common to al...
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Container class to hold unbinned data.
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Efficient implementation of a product of PDFs of the form.
Variable that can be changed from the outside.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
std::map< TString, double > fVarVals
void AssignToTTree(TTree &myTree)
BranchStore(const vector< TString > ¶ms=vector< TString >(), double _inval=-999.)
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
virtual void SetObservables(const RooArgSet &set)
Specify the observables.
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
virtual void SetParametersOfInterest(const RooArgSet &set)
Specify parameters of interest.
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)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
virtual void SetNuisanceParameters(const RooArgSet &set)
Specify the nuisance parameters (parameters that are not POI).
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the workspace if not already there.
virtual void SetGlobalObservables(const RooArgSet &set)
Specify the global observables.
Persistable container for RooFit projects.
TObject * obj(RooStringView name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
bool saveSnapshot(RooStringView, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of given parameters.
std::list< RooAbsData * > allData() const
Return list of all dataset in the workspace.
RooLinkedList const & getSnapshots() const
std::list< TObject * > allGenericObjects() const
Return list of all generic objects in the workspace.
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.
RooRealVar * var(RooStringView name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
static void autoImportClassCode(bool flag)
If flag is true, source code of classes not the ROOT distribution is automatically imported if on obj...
A TTree is a list of TBranches.
virtual void ResetAddress()
Reset the address of the branch.
const char * GetName() const override
Returns name of object.
const char * GetTitle() const override
Returns title of object.
const char * Data() const
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
A TTree represents a columnar dataset.
virtual Int_t Fill()
Fill all branches.
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
RooCmdArg RecycleConflictNodes(bool flag=true)
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...
void FillTree(TTree &myTree, const RooDataSet &data)
BranchStore * CreateBranchStore(const RooDataSet &data)
RooAbsPdf * StripConstraints(RooAbsPdf &pdf, const RooArgSet &observables)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
extract constraint terms from pdf
void FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints)
RooStatsConfig & GetGlobalRooStatsConfig()
Retrieve the config object which can be used to set flags for things like offsetting the likelihood o...
double AsimovSignificance(double s, double b, double sigma_b=0.0)
Compute the Asimov Median significance for a Poisson process with s = expected number of signal event...
void UseNLLOffset(bool on)
function to set a global flag in RooStats to use NLL offset when performing nll computations Note tha...
RooAbsPdf * MakeUnconstrainedPdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name=nullptr)
remove constraints from pdf and return the unconstrained pdf
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set
RooWorkspace * MakeCleanWorkspace(RooWorkspace *oldWS, const char *newName, bool copySnapshots, const char *mcname, const char *newmcname)
function that clones a workspace, copying all needed components and discarding all others
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
bool useLikelihoodOffset
Offset the likelihood by passing RooFit::Offset to fitTo().