Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooFuncWrapper.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 * Authors:
4 * Garima Singh, CERN 2022
5 *
6 * Copyright (c) 2022, CERN
7 *
8 * Redistribution and use in source and binary forms,
9 * with or without modification, are permitted according to the terms
10 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
11 */
12
13#include <RooFuncWrapper.h>
14
15#include <RooAbsData.h>
17#include <RooFit/Evaluator.h>
18#include <RooGlobalFunc.h>
19#include <RooHelpers.h>
20#include <RooMsgService.h>
21#include <RooRealVar.h>
22#include <RooSimultaneous.h>
23
24#include "RooEvaluatorWrapper.h"
26
27#include <TROOT.h>
28#include <TSystem.h>
29
30#include <fstream>
31#include <set>
32
33namespace {
34
35void replaceAll(std::string &str, const std::string &from, const std::string &to)
36{
37 if (from.empty())
38 return;
39 size_t start_pos = 0;
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(); // In case 'to' contains 'from', like replacing 'x' with 'yx'
43 }
44}
45
46} // namespace
47
48namespace RooFit {
49
50namespace Experimental {
51
52RooFuncWrapper::RooFuncWrapper(const char *name, const char *title, RooAbsReal &obj, const RooAbsData *data,
53 RooSimultaneous const *simPdf, bool useEvaluator)
54 : RooAbsReal{name, title}, _params{"!params", "List of parameters", this}, _useEvaluator{useEvaluator}
55{
56 if (_useEvaluator) {
57 _absReal = std::make_unique<RooEvaluatorWrapper>(obj, const_cast<RooAbsData *>(data), false, "", simPdf, false);
58 }
59
60 // Get the parameters.
61 RooArgSet paramSet;
62 obj.getParameters(data ? data->get() : nullptr, paramSet);
63
64 // Load the parameters and observables.
65 auto spans = loadParamsAndData(paramSet, data, simPdf);
66
67 // Set up the code generation context
68 std::map<RooFit::Detail::DataKey, std::size_t> nodeOutputSizes =
69 RooFit::BatchModeDataHelpers::determineOutputSizes(obj, [&spans](RooFit::Detail::DataKey key) -> int {
70 auto found = spans.find(key);
71 return found != spans.end() ? found->second.size() : -1;
72 });
73
75
76 // First update the result variable of params in the compute graph to in[<position>].
77 int idx = 0;
78 for (RooAbsArg *param : _params) {
79 ctx.addResult(param, "params[" + std::to_string(idx) + "]");
80 idx++;
81 }
82
83 for (auto const &item : _obsInfos) {
84 const char *obsName = item.first->GetName();
85 // If the observable is scalar, set name to the start idx. else, store
86 // the start idx and later set the the name to obs[start_idx + curr_idx],
87 // here curr_idx is defined by a loop producing parent node.
88 if (item.second.size == 1) {
89 ctx.addResult(obsName, "obs[" + std::to_string(item.second.idx) + "]");
90 } else {
91 ctx.addResult(obsName, "obs");
92 ctx.addVecObs(obsName, item.second.idx);
93 }
94 }
95
96 gInterpreter->Declare("#pragma cling optimize(2)");
97
98 // Declare the function and create its derivative.
99 _funcName = ctx.buildFunction(obj, nodeOutputSizes);
100 _func = reinterpret_cast<Func>(gInterpreter->ProcessLine((_funcName + ";").c_str()));
101
102 _xlArr = ctx.xlArr();
103 _collectedFunctions = ctx.collectedFunctions();
104}
105
107 : RooAbsReal(other, name),
108 _params("!params", this, other._params),
109 _funcName(other._funcName),
110 _func(other._func),
111 _grad(other._grad),
112 _hasGradient(other._hasGradient),
113 _gradientVarBuffer(other._gradientVarBuffer),
114 _observables(other._observables)
115{
116}
117
118std::map<RooFit::Detail::DataKey, std::span<const double>>
120{
121 // Extract observables
122 std::stack<std::vector<double>> vectorBuffers; // for data loading
123 std::map<RooFit::Detail::DataKey, std::span<const double>> spans;
124
125 if (data) {
126 spans = RooFit::BatchModeDataHelpers::getDataSpans(*data, "", simPdf, true, false, vectorBuffers);
127 }
128
129 std::size_t idx = 0;
130 for (auto const &item : spans) {
131 std::size_t n = item.second.size();
132 _obsInfos.emplace(item.first, ObsInfo{idx, n});
133 _observables.reserve(_observables.size() + n);
134 for (std::size_t i = 0; i < n; ++i) {
135 _observables.push_back(item.second[i]);
136 }
137 idx += n;
138 }
139
140 // Extract parameters
141 for (auto *param : paramSet) {
142 if (!dynamic_cast<RooAbsReal *>(param)) {
143 std::stringstream errorMsg;
144 errorMsg << "In creation of function " << GetName()
145 << " wrapper: input param expected to be of type RooAbsReal.";
146 coutE(InputArguments) << errorMsg.str() << std::endl;
147 throw std::runtime_error(errorMsg.str().c_str());
148 }
149 if (spans.find(param) == spans.end()) {
150 _params.add(*param);
151 }
152 }
153 _gradientVarBuffer.resize(_params.size());
154
155 return spans;
156}
157
158void RooFuncWrapper::createGradient()
159{
160 std::string gradName = _funcName + "_grad_0";
161 std::string requestName = _funcName + "_req";
162
163 // Calculate gradient
164 gInterpreter->Declare("#include <Math/CladDerivator.h>\n");
165 // disable clang-format for making the following code unreadable.
166 // clang-format off
167 std::stringstream requestFuncStrm;
168 requestFuncStrm << "#pragma clad ON\n"
169 "void " << requestName << "() {\n"
170 " clad::gradient(" << _funcName << ", \"params\");\n"
171 "}\n"
172 "#pragma clad OFF";
173 // clang-format on
174 if (!gInterpreter->Declare(requestFuncStrm.str().c_str())) {
175 std::stringstream errorMsg;
176 errorMsg << "Function " << GetName() << " could not be differentiated. See above for details.";
177 oocoutE(nullptr, InputArguments) << errorMsg.str() << std::endl;
178 throw std::runtime_error(errorMsg.str().c_str());
179 }
180
181 _grad = reinterpret_cast<Grad>(gInterpreter->ProcessLine((gradName + ";").c_str()));
182 _hasGradient = true;
183}
184
185void RooFuncWrapper::gradient(double *out) const
186{
187 updateGradientVarBuffer();
188 std::fill(out, out + _params.size(), 0.0);
189
190 _grad(_gradientVarBuffer.data(), _observables.data(), _xlArr.data(), out);
191}
192
193void RooFuncWrapper::updateGradientVarBuffer() const
194{
195 std::transform(_params.begin(), _params.end(), _gradientVarBuffer.begin(),
196 [](RooAbsArg *obj) { return static_cast<RooAbsReal *>(obj)->getVal(); });
197}
198
199double RooFuncWrapper::evaluate() const
200{
201 if (_useEvaluator)
202 return _absReal->getVal();
203 updateGradientVarBuffer();
204
205 return _func(_gradientVarBuffer.data(), _observables.data(), _xlArr.data());
206}
207
208void RooFuncWrapper::gradient(const double *x, double *g) const
209{
210 std::fill(g, g + _params.size(), 0.0);
211
212 _grad(const_cast<double *>(x), _observables.data(), _xlArr.data(), g);
213}
214
215/// @brief Dumps a macro "filename.C" that can be used to test and debug the generated code and gradient.
216void RooFuncWrapper::writeDebugMacro(std::string const &filename) const
217{
218 std::stringstream allCode;
219 std::set<std::string> seenFunctions;
220
221 // Remove duplicated declared functions
222 for (std::string const &name : _collectedFunctions) {
223 if (seenFunctions.count(name) > 0) {
224 continue;
225 }
226 seenFunctions.insert(name);
227 std::unique_ptr<TInterpreterValue> v = gInterpreter->MakeInterpreterValue();
228 gInterpreter->Evaluate(name.c_str(), *v);
229 std::string s = v->ToString();
230 for (int i = 0; i < 2; ++i) {
231 s = s.erase(0, s.find("\n") + 1);
232 }
233 allCode << s << std::endl;
234 }
235
236 std::ofstream outFile;
237 outFile.open(filename + ".C");
238 outFile << R"(//auto-generated test macro
239#include <RooFit/Detail/MathFuncs.h>
240#include <Math/CladDerivator.h>
241
242#pragma cling optimize(2)
243)" << allCode.str()
244 << R"(
245#pragma clad ON
246void gradient_request() {
247 clad::gradient()"
248 << _funcName << R"(, "params");
249}
250#pragma clad OFF
251)";
252
253 updateGradientVarBuffer();
254
255 auto writeVector = [&](std::string const &name, std::span<const double> vec) {
256 std::stringstream decl;
257 decl << "std::vector<double> " << name << " = {";
258 for (std::size_t i = 0; i < vec.size(); ++i) {
259 if (i % 10 == 0)
260 decl << "\n ";
261 decl << vec[i];
262 if (i < vec.size() - 1)
263 decl << ", ";
264 }
265 decl << "\n};\n";
266
267 std::string declStr = decl.str();
268
269 replaceAll(declStr, "inf", "std::numeric_limits<double>::infinity()");
270 replaceAll(declStr, "nan", "NAN");
271
272 outFile << declStr;
273 };
274
275 outFile << "// clang-format off\n" << std::endl;
276 writeVector("parametersVec", _gradientVarBuffer);
277 outFile << std::endl;
278 writeVector("observablesVec", _observables);
279 outFile << std::endl;
280 writeVector("auxConstantsVec", _xlArr);
281 outFile << std::endl;
282 outFile << "// clang-format on\n" << std::endl;
283
284 outFile << R"(
285// To run as a ROOT macro
286void )" << filename
287 << R"(()
288{
289 std::vector<double> gradientVec(parametersVec.size());
290
291 auto func = [&](std::span<double> params) {
292 return )"
293 << _funcName << R"((params.data(), observablesVec.data(), auxConstantsVec.data());
294 };
295 auto grad = [&](std::span<double> params, std::span<double> out) {
296 return )"
297 << _funcName << R"(_grad_0(parametersVec.data(), observablesVec.data(), auxConstantsVec.data(),
298 out.data());
299 };
300
301 grad(parametersVec, gradientVec);
302
303 auto numDiff = [&](int i) {
304 const double eps = 1e-6;
305 std::vector<double> p{parametersVec};
306 p[i] = parametersVec[i] - eps;
307 double funcValDown = func(p);
308 p[i] = parametersVec[i] + eps;
309 double funcValUp = func(p);
310 return (funcValUp - funcValDown) / (2 * eps);
311 };
312
313 for (std::size_t i = 0; i < parametersVec.size(); ++i) {
314 std::cout << i << ":" << std::endl;
315 std::cout << " numr : " << numDiff(i) << std::endl;
316 std::cout << " clad : " << gradientVec[i] << std::endl;
317 }
318}
319)";
320}
321
322} // namespace Experimental
323
324} // namespace RooFit
#define g(i)
Definition RSha256.hxx:105
#define oocoutE(o, a)
#define coutE(a)
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
char name[80]
Definition TGX11.cxx:110
#define gInterpreter
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
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.
Definition RooAbsData.h:57
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
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.
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::vector< std::string > _collectedFunctions
std::map< RooFit::Detail::DataKey, std::span< const double > > loadParamsAndData(RooArgSet const &paramSet, const RooAbsData *data, RooSimultaneous const *simPdf)
std::map< RooFit::Detail::DataKey, ObsInfo > _obsInfos
void(*)(double *, double const *, double const *, double *) Grad
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.
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
void replaceAll(std::string &inOut, std::string_view what, std::string_view with)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
@ InputArguments