35void replaceAll(std::string &str,
const std::string &from,
const std::string &to)
40 while ((start_pos = str.find(from, start_pos)) != std::string::npos) {
41 str.replace(start_pos, from.length(), to);
42 start_pos += to.length();
50namespace Experimental {
54 :
RooAbsReal{
name, title}, _params{
"!params",
"List of parameters", this}, _useEvaluator{useEvaluator}
57 _absReal = std::make_unique<RooEvaluatorWrapper>(obj,
const_cast<RooAbsData *
>(
data),
false,
"", simPdf,
false);
80 _params(
"!params", this, other._params),
81 _funcName(other._funcName),
84 _hasGradient(other._hasGradient),
85 _gradientVarBuffer(other._gradientVarBuffer),
86 _observables(other._observables)
94 std::stack<std::vector<double>> vectorBuffers;
95 std::map<RooFit::Detail::DataKey, std::span<const double>> spans;
98 spans = RooFit::BatchModeDataHelpers::getDataSpans(*
data,
"", simPdf,
true,
false, vectorBuffers);
102 for (
auto const &item : spans) {
103 std::size_t
n = item.second.size();
106 for (std::size_t i = 0; i <
n; ++i) {
113 for (
auto *param : paramSet) {
115 std::stringstream errorMsg;
116 errorMsg <<
"In creation of function " << GetName()
117 <<
" wrapper: input param expected to be of type RooAbsReal.";
119 throw std::runtime_error(errorMsg.str().c_str());
121 if (spans.find(param) == spans.end()) {
125 _gradientVarBuffer.resize(_params.size());
128 _nodeOutputSizes = RooFit::BatchModeDataHelpers::determineOutputSizes(
130 auto found = spans.find(key);
131 return found != spans.end() ? found->second.size() : -1;
136std::string RooFuncWrapper::declareFunction(std::string
const &funcBody)
138 static int iFuncWrapper = 0;
139 auto funcName =
"roo_func_wrapper_" + std::to_string(iFuncWrapper++);
142 std::stringstream bodyWithSigStrm;
143 bodyWithSigStrm <<
"double " << funcName <<
"(double* params, double const* obs, double const* xlArr) {\n"
144 << funcBody <<
"\n}";
145 _collectedFunctions.emplace_back(funcName);
146 if (!
gInterpreter->Declare(bodyWithSigStrm.str().c_str())) {
147 std::stringstream errorMsg;
148 errorMsg <<
"Function " << funcName <<
" could not be compiled. See above for details.";
150 throw std::runtime_error(errorMsg.str().c_str());
155void RooFuncWrapper::createGradient()
157 std::string gradName = _funcName +
"_grad_0";
158 std::string requestName = _funcName +
"_req";
161 gInterpreter->Declare(
"#include <Math/CladDerivator.h>\n");
164 std::stringstream requestFuncStrm;
165 requestFuncStrm <<
"#pragma clad ON\n"
166 "void " << requestName <<
"() {\n"
167 " clad::gradient(" << _funcName <<
", \"params\");\n"
171 if (!
gInterpreter->Declare(requestFuncStrm.str().c_str())) {
172 std::stringstream errorMsg;
173 errorMsg <<
"Function " << GetName() <<
" could not be differentiated. See above for details.";
175 throw std::runtime_error(errorMsg.str().c_str());
178 _grad =
reinterpret_cast<Grad>(
gInterpreter->ProcessLine((gradName +
";").c_str()));
182void RooFuncWrapper::gradient(
double *out)
const
184 updateGradientVarBuffer();
185 std::fill(out, out + _params.size(), 0.0);
187 _grad(_gradientVarBuffer.data(), _observables.data(), _xlArr.data(), out);
190void RooFuncWrapper::updateGradientVarBuffer()
const
192 std::transform(_params.begin(), _params.end(), _gradientVarBuffer.begin(),
193 [](
RooAbsArg *obj) { return static_cast<RooAbsReal *>(obj)->getVal(); });
196double RooFuncWrapper::evaluate()
const
199 return _absReal->getVal();
200 updateGradientVarBuffer();
202 return _func(_gradientVarBuffer.data(), _observables.data(), _xlArr.data());
205void RooFuncWrapper::gradient(
const double *
x,
double *
g)
const
207 std::fill(
g,
g + _params.size(), 0.0);
209 _grad(
const_cast<double *
>(
x), _observables.data(), _xlArr.data(),
g);
219 ctx.
addResult(param,
"params[" + std::to_string(idx) +
"]");
223 for (
auto const &item : _obsInfos) {
224 const char *
name = item.first->GetName();
228 if (item.second.size == 1) {
229 ctx.
addResult(
name,
"obs[" + std::to_string(item.second.idx) +
"]");
240void RooFuncWrapper::writeDebugMacro(std::string
const &
filename)
const
242 std::stringstream allCode;
243 std::set<std::string> seenFunctions;
246 for (std::string
const &
name : _collectedFunctions) {
247 if (seenFunctions.count(
name) > 0) {
250 seenFunctions.insert(
name);
251 std::unique_ptr<TInterpreterValue>
v =
gInterpreter->MakeInterpreterValue();
253 std::string s =
v->ToString();
254 for (
int i = 0; i < 2; ++i) {
255 s = s.erase(0, s.find(
"\n") + 1);
257 allCode << s << std::endl;
260 std::ofstream outFile;
262 outFile << R
"(//auto-generated test macro
263#include <RooFit/Detail/MathFuncs.h>
264#include <Math/CladDerivator.h>
266#pragma cling optimize(2)
270void gradient_request() {
272 << _funcName << R"(, "params");
277 updateGradientVarBuffer();
279 auto writeVector = [&](std::string
const &
name, std::span<const double>
vec) {
280 std::stringstream decl;
281 decl <<
"std::vector<double> " <<
name <<
" = {";
282 for (std::size_t i = 0; i <
vec.size(); ++i) {
286 if (i <
vec.size() - 1)
291 std::string declStr = decl.str();
293 replaceAll(declStr,
"inf",
"std::numeric_limits<double>::infinity()");
294 replaceAll(declStr,
"nan",
"NAN");
299 outFile <<
"// clang-format off\n" << std::endl;
300 writeVector(
"parametersVec", _gradientVarBuffer);
301 outFile << std::endl;
302 writeVector(
"observablesVec", _observables);
303 outFile << std::endl;
304 writeVector(
"auxConstantsVec", _xlArr);
305 outFile << std::endl;
306 outFile <<
"// clang-format on\n" << std::endl;
309// To run as a ROOT macro
313 std::vector<double> gradientVec(parametersVec.size());
315 auto func = [&](std::span<double> params) {
317 << _funcName << R"((params.data(), observablesVec.data(), auxConstantsVec.data());
319 auto grad = [&](std::span<double> params, std::span<double> out) {
321 << _funcName << R"(_grad_0(parametersVec.data(), observablesVec.data(), auxConstantsVec.data(),
325 grad(parametersVec, gradientVec);
327 auto numDiff = [&](int i) {
328 const double eps = 1e-6;
329 std::vector<double> p{parametersVec};
330 p[i] = parametersVec[i] - eps;
331 double funcValDown = func(p);
332 p[i] = parametersVec[i] + eps;
333 double funcValUp = func(p);
334 return (funcValUp - funcValDown) / (2 * eps);
337 for (std::size_t i = 0; i < parametersVec.size(); ++i) {
338 std::cout << i << ":" << std::endl;
339 std::cout << " numr : " << numDiff(i) << std::endl;
340 std::cout << " clad : " << gradientVec[i] << std::endl;
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 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
Common abstract base class for objects that represent a value and a "shape" in RooFit.
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...
Abstract base class for binned and unbinned datasets.
Abstract base class for objects that represent a real value and implements functionality common to al...
RooArgSet is a container object that can hold multiple RooAbsArg objects.
A class to maintain the context for squashing of RooFit models into code.
std::string assembleCode(std::string const &returnExpr)
Assemble and return the final code with the return expression and global statements.
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.
void addVecObs(const char *key, int idx)
Since the squashed code represents all observables as a single flattened array, it is important to ke...
std::string const & getResult(RooAbsArg const &arg)
Gets the result for the given node using the node name.
A wrapper class to store a C++ function of type 'double (*)(double*, double*)'.
double(*)(double *, double const *, double const *) Func
std::unique_ptr< RooAbsReal > _absReal
std::string buildCode(RooAbsReal const &head)
std::vector< double > _observables
void loadParamsAndData(RooAbsArg const *head, RooArgSet const ¶mSet, const RooAbsData *data, RooSimultaneous const *simPdf)
std::map< RooFit::Detail::DataKey, ObsInfo > _obsInfos
void(*)(double *, double const *, double const *, double *) Grad
std::string declareFunction(std::string const &funcBody)
RooFuncWrapper(const char *name, const char *title, RooAbsReal &obj, const RooAbsData *data=nullptr, RooSimultaneous const *simPdf=nullptr, bool useEvaluator=false)
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...