57 _pdfObsList(
"pdfObs",
"List of p.d.f. observables",this),
69 coutE(InputArguments) <<
"RooHistPdf::ctor(" <<
GetName()
70 <<
") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
73 for (
const auto arg : vars) {
74 if (!dvars->
find(arg->GetName())) {
75 coutE(InputArguments) <<
"RooHistPdf::ctor(" <<
GetName()
76 <<
") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
107 _pdfObsList(
"pdfObs",
"List of p.d.f. observables",this),
111 _cdfBoundaries(false)
118 if (histObs.
size()!=dvars->
size()) {
119 coutE(InputArguments) <<
"RooHistPdf::ctor(" <<
GetName()
120 <<
") ERROR histogram variable list and RooDataHist must contain the same variables." << std::endl ;
121 throw(std::string(
"RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
124 for (
const auto arg : histObs) {
125 if (!dvars->
find(arg->GetName())) {
126 coutE(InputArguments) <<
"RooHistPdf::ctor(" <<
GetName()
127 <<
") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
128 throw(std::string(
"RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
130 if (!arg->isFundamental()) {
131 coutE(InputArguments) <<
"RooHistPdf::ctor(" <<
GetName()
132 <<
") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << std::endl ;
133 throw(std::string(
"RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
151 std::unique_ptr<RooDataHist> dhist,
int intOrder)
157 std::unique_ptr<RooDataHist> dhist,
int intOrder)
169 _pdfObsList(
"pdfObs",this,other._pdfObsList),
170 _dataHist(other._dataHist),
171 _codeReg(other._codeReg),
172 _intOrder(other._intOrder),
173 _cdfBoundaries(other._cdfBoundaries),
174 _totVolume(other._totVolume),
175 _unitNorm(other._unitNorm)
235 return std::max(ret, 0.0);
242 ooccoutE(klass, InputArguments) <<
"RooHistPdf::weight(" << klass->
GetName()
243 <<
") ERROR: Code Squashing currently only supports non-interpolation cases."
250 ctx.
addResult(klass, weightName +
"[" + idxName +
"]");
290 if (!_x || !_y)
return false;
291 if (!range || !strlen(range) || !_x->
hasRange(range) ||
305 if(lobs ==
nullptr)
return false;
307 bool isOkayForAnalyticalInt =
false;
311 if(!lobs->isJacobianOK(*var))
return false;
312 isOkayForAnalyticalInt =
true;
316 return isOkayForAnalyticalInt;
324 const char* rangeName,
334 for (
unsigned int n=0;
n < pdfObsList.
size() &&
n < histObsList.
size(); ++
n) {
335 const auto pa = pdfObsList[
n];
336 const auto ha = histObsList[
n];
338 if (okayForAnalytical(*pa, allVars)) {
341 if (fullRange(*pa, *ha, rangeName)) {
347 if (code == frcode) {
356 if (intOrder > 1 && !(code & 1)) {
360 return (code >= 2) ? code : 0;
365 const char* rangeName,
371 if (((2 << histObsList.
size()) - 1) == code) {
378 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
379 for (
unsigned int n=0;
n < pdfObsList.
size() &&
n < histObsList.
size(); ++
n) {
380 const auto pa = pdfObsList[
n];
381 const auto ha = histObsList[
n];
383 if (code & (2 <<
n)) {
393 ha->copyCache(pa,
true);
397 double ret = (code & 1) ?
dataHist.
sum(intSet,histObsList,
true,!histFuncMode) :
398 dataHist.
sum(intSet,histObsList,
true,!histFuncMode, ranges);
406 if (((2 << obs.
size()) - 1) != code) {
408 <<
"RooHistPdf::integral(" << klass->
GetName()
409 <<
") ERROR: AD currently only supports integrating over all histogram observables." << std::endl;
412 return std::to_string(
dataHist->
sum(histFuncMode));
447 bool isOkayForAnalyticalInt =
false;
454 if(!(lvalue && lvalue->isJacobianOK(dep)))
return false;
455 isOkayForAnalyticalInt =
true;
459 return isOkayForAnalyticalInt;
495 for (
unsigned int i=0; i < pdfObsList.
size(); ++i) {
515 std::span<const double> boundaries{binning->
array(),
static_cast<std::size_t
>(binning->
numBoundaries())};
543 double* boundaries = binning->
array() ;
545 auto hint =
new std::list<double> ;
550 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
551 hint->push_back(boundaries[i]) ;
584 if (wgt>max) max=wgt ;
638 coutE(ObjectHandling) <<
" RooHistPdf::importWorkspaceHook(" <<
GetName() <<
") unable to import clone of underlying RooDataHist with unique name " << uniqueName <<
", abort" << std::endl ;
650 coutE(ObjectHandling) <<
" RooHistPdf::importWorkspaceHook(" <<
GetName() <<
") unable to import clone of underlying RooDataHist with unique name " << uniqueName <<
", abort" << std::endl ;
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Common abstract base class for objects that represent a value and a "shape" in RooFit.
virtual void copyCache(const RooAbsArg *source, bool valueOnly=false, bool setValDirty=true)=0
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 void syncCache(const RooArgSet *nset=nullptr)=0
friend void RooRefArray::Streamer(TBuffer &)
virtual bool inRange(const char *) const
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
virtual bool isParameterized() const
Interface function.
virtual Int_t numBoundaries() const =0
virtual double * array() const =0
Int_t numTypes(const char *=nullptr) const
Return number of types defined (in range named rangeName if rangeName!=nullptr)
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
bool selectCommon(const RooAbsCollection &refColl, RooAbsCollection &outColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for binned and unbinned datasets.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for objects that are lvalues, i.e.
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Abstract interface for all probability density functions.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
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.
const RooAbsBinning * getBinningPtr(const char *rangeName) const override
bool hasRange(const char *name) const override
Check if variable has a binning with given name.
virtual void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Object to represent discrete states.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
static std::list< double > * plotSamplingHintForBinBoundaries(std::span< const double > boundaries, double xlo, double xhi)
Returns sampling hints for a histogram with given boundaries.
The RooDataHist is a container class to hold N-dimensional binned data.
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
void weights(double *output, std::span< double const > xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of RooDataHist::weight() for one dimensional histograms with up to one dimension...
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
std::string calculateTreeIndexForCodeSquash(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
double weight(std::size_t i) const
Return weight of i-th bin.
double weightFast(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A faster version of RooDataHist::weight that assumes the passed arguments are aligned with the histog...
std::string declWeightArrayForCodeSquash(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, bool correctForBinSize) const
const RooArgSet * get() const override
Get bin centre of current bin.
double sumEntries() const override
Sum the weights of all bins.
A class to maintain the context for squashing of RooFit models into code.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
RooHistPdf implements a propability density function sampled from a multidimensional histogram.
static void rooHistTranslateImpl(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, int intOrder, RooDataHist const *dataHist, const RooArgSet &obs, bool correctForBinSize)
RooArgSet _histObsList
List of observables defining dimensions of histogram.
Int_t _intOrder
Interpolation order.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
bool areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
RooDataHist * _dataHist
Unowned pointer to underlying histogram.
bool _cdfBoundaries
Use boundary conditions for CDFs.
static std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const *klass, RooDataHist const *dataHist, const RooArgSet &obs, bool histFuncMode)
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
double totVolume() const
Return the total volume spanned by the observables of the RooHistPdf.
bool importWorkspaceHook(RooWorkspace &ws) override
Check if our datahist is already in the workspace.
std::string buildCallToAnalyticIntegral(int code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
RooSetProxy _pdfObsList
List of observables mapped onto histogram observables.
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const override
Base function for computing multiple values of a RooAbsReal.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Return integral identified by 'code'.
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
Int_t getMaxVal(const RooArgSet &vars) const override
Only handle case of maximum in all variables.
double _totVolume
! Total volume of space (product of ranges of observables)
RooDataHist * cloneAndOwnDataHist(const char *newname="")
Replaces underlying RooDataHist with a clone, which is now owned, and returns the clone.
std::unique_ptr< RooDataHist > _ownedDataHist
! Owned pointer to underlying histogram
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
~RooHistPdf() override
Destructor.
double evaluate() const override
Return the current value: The value of the bin enclosing the current coordinates of the observables,...
bool _unitNorm
Assume contents is unit normalized (for use as pdf cache)
RooRealVar represents a variable that can be changed from the outside.
void setRange(const char *name, double min, double max)
Set a fit or plotting range.
Persistable container for RooFit projects.
RooAbsData * embeddedData(RooStringView name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
std::list< RooAbsData * > allData() const
Return list of all dataset 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.
Buffer base class used for serializing objects.
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
const char * GetName() const override
Returns name of object.
RooCmdArg Rename(const char *suffix)
RooCmdArg Embedded(bool flag=true)
std::pair< double, double > getRangeOrBinningInterval(RooAbsArg const *arg, const char *rangeName)