Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooEvaluatorWrapper.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*
4 * Project: RooFit
5 * Authors:
6 * Jonas Rembser, CERN 2023
7 *
8 * Copyright (c) 2023, CERN
9 *
10 * Redistribution and use in source and binary forms,
11 * with or without modification, are permitted according to the terms
12 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
13 */
14
15/**
16\internal
17\file RooEvaluatorWrapper.cxx
18\class RooEvaluatorWrapper
19\ingroup Roofitcore
20
21Wraps a RooFit::Evaluator that evaluates a RooAbsReal back into a RooAbsReal.
22**/
23
24#include "RooEvaluatorWrapper.h"
25
26#include <RooAbsData.h>
27#include <RooAbsPdf.h>
28#include <RooMsgService.h>
29#include <RooRealVar.h>
30#include <RooSimultaneous.h>
31
33#include "RooFitImplHelpers.h"
34
35#include <TInterpreter.h>
36
37#include <fstream>
38
39namespace RooFit::Experimental {
40
41RooEvaluatorWrapper::RooEvaluatorWrapper(RooAbsReal &topNode, RooAbsData *data, bool useGPU,
42 std::string const &rangeName, RooAbsPdf const *pdf,
44 : RooAbsReal{"RooEvaluatorWrapper", "RooEvaluatorWrapper"},
45 _evaluator{std::make_unique<RooFit::Evaluator>(topNode, useGPU)},
46 _topNode("topNode", "top node", this, topNode, false, false),
47 _data{data},
48 _paramSet("paramSet", "Set of parameters", this),
49 _rangeName{rangeName},
50 _pdf{pdf},
51 _takeGlobalObservablesFromData{takeGlobalObservablesFromData}
52{
53 if (data) {
54 setData(*data, false);
55 }
56 _paramSet.add(_evaluator->getParameters());
57 for (auto const &item : _dataSpans) {
58 _paramSet.remove(*_paramSet.find(item.first->GetName()));
59 }
60}
61
62RooEvaluatorWrapper::RooEvaluatorWrapper(const RooEvaluatorWrapper &other, const char *name)
64 _evaluator{other._evaluator},
65 _topNode("topNode", this, other._topNode),
66 _data{other._data},
67 _paramSet("paramSet", "Set of parameters", this),
68 _rangeName{other._rangeName},
69 _pdf{other._pdf},
70 _takeGlobalObservablesFromData{other._takeGlobalObservablesFromData},
72{
73 _paramSet.add(other._paramSet);
74}
75
76RooEvaluatorWrapper::~RooEvaluatorWrapper() = default;
77
78bool RooEvaluatorWrapper::getParameters(const RooArgSet *observables, RooArgSet &outputSet,
79 bool stripDisconnected) const
80{
81 outputSet.add(_evaluator->getParameters());
82 if (observables) {
83 outputSet.remove(*observables, /*silent*/ false, /*matchByNameOnly*/ true);
84 }
85 // Exclude the data variables from the parameters which are not global observables
86 for (auto const &item : _dataSpans) {
87 if (_data->getGlobalObservables() && _data->getGlobalObservables()->find(item.first->GetName())) {
88 continue;
89 }
90 RooAbsArg *found = outputSet.find(item.first->GetName());
91 if (found) {
92 outputSet.remove(*found);
93 }
94 }
95 // If we take the global observables as data, we have to return these as
96 // parameters instead of the parameters in the model. Otherwise, the
97 // constant parameters in the fit result that are global observables will
98 // not have the right values.
99 if (_takeGlobalObservablesFromData && _data->getGlobalObservables()) {
100 outputSet.replace(*_data->getGlobalObservables());
101 }
102
103 // The disconnected parameters are stripped away in
104 // RooAbsArg::getParametersHook(), that is only called in the original
105 // RooAbsArg::getParameters() implementation. So he have to call it to
106 // identify disconnected parameters to remove.
107 if (stripDisconnected) {
109 _topNode->getParameters(observables, paramsStripped, true);
111 for (RooAbsArg *param : outputSet) {
112 if (!paramsStripped.find(param->GetName())) {
113 toRemove.add(*param);
114 }
115 }
116 outputSet.remove(toRemove, /*silent*/ false, /*matchByNameOnly*/ true);
117 }
118
119 return false;
120}
121
122/// @brief A wrapper class to store a C++ function of type 'double (*)(double*, double*)'.
123/// The parameters can be accessed as params[<relative position of param in paramSet>] in the function body.
124/// The observables can be accessed as obs[i + j], where i represents the observable position and j
125/// represents the data entry.
126class RooFuncWrapper {
127public:
129
130 bool hasGradient() const { return _hasGradient; }
131 void gradient(double *out) const
132 {
134 std::fill(out, out + _params.size(), 0.0);
135
136 _grad(_gradientVarBuffer.data(), _observables.data(), _xlArr.data(), out);
137 }
138
139 void createGradient();
140
141 void writeDebugMacro(std::string const &) const;
142
143 std::vector<std::string> const &collectedFunctions() { return _collectedFunctions; }
144
145 double evaluate() const
146 {
148 return _func(_gradientVarBuffer.data(), _observables.data(), _xlArr.data());
149 }
150
151 void loadData(RooAbsData const &data, RooSimultaneous const *simPdf);
152
153private:
154 void updateGradientVarBuffer() const;
155
157
158 using Func = double (*)(double *, double const *, double const *);
159 using Grad = void (*)(double *, double const *, double const *, double *);
160
161 RooArgList _params;
162 std::string _funcName;
163 Func _func;
164 Grad _grad;
165 bool _hasGradient = false;
166 mutable std::vector<double> _gradientVarBuffer;
167 std::vector<double> _observables;
168 std::unordered_map<RooFit::Detail::DataKey, std::size_t> _obsInfos;
169 std::vector<double> _xlArr;
170 std::vector<std::string> _collectedFunctions;
171};
172
173namespace {
174
175void replaceAll(std::string &str, const std::string &from, const std::string &to)
176{
177 if (from.empty())
178 return;
179 size_t start_pos = 0;
180 while ((start_pos = str.find(from, start_pos)) != std::string::npos) {
181 str.replace(start_pos, from.length(), to);
182 start_pos += to.length(); // In case 'to' contains 'from', like replacing 'x' with 'yx'
183 }
184}
185
187{
190
191 std::unordered_set<RooFit::Detail::DataKey> dependsOnData;
192 for (RooAbsArg *arg : dataObs) {
193 dependsOnData.insert(arg);
194 }
195
196 for (RooAbsArg *arg : serverSet) {
197 if (arg->getAttribute("__obs__")) {
198 dependsOnData.insert(arg);
199 }
200 for (RooAbsArg *server : arg->servers()) {
201 if (server->isValueServer(*arg)) {
202 if (dependsOnData.find(server) != dependsOnData.end() && !arg->isReducerNode()) {
203 dependsOnData.insert(arg);
204 break;
205 }
206 }
207 }
208 }
209
210 return dependsOnData;
211}
212
213} // namespace
214
215RooFuncWrapper::RooFuncWrapper(RooAbsReal &obj, const RooAbsData *data, RooSimultaneous const *simPdf,
216 RooArgSet const &paramSet)
217{
218 // Load the observables from the dataset
219 if (data) {
221 }
222
223 // Define the parameters
224 for (auto *param : paramSet) {
225 if (_obsInfos.find(param) == _obsInfos.end()) {
226 _params.add(*param);
227 }
228 }
229 _gradientVarBuffer.resize(_params.size());
230
231 // Figure out which part of the computation graph depends on data
232 std::unordered_set<RooFit::Detail::DataKey> dependsOnData;
233 if (data) {
234 dependsOnData = getDependsOnData(obj, *data->get());
235 }
236
237 // Set up the code generation context
239
240 // First update the result variable of params in the compute graph to in[<position>].
241 int idx = 0;
242 for (RooAbsArg *param : _params) {
243 ctx.addResult(param, "params[" + std::to_string(idx) + "]");
244 idx++;
245 }
246
247 for (auto const &item : _obsInfos) {
248 const char *obsName = item.first->GetName();
249 ctx.addResult(obsName, "obs");
250 ctx.addVecObs(obsName, item.second);
251 }
252
253 // Declare the function and create its derivative.
254 auto print = [](std::string const &msg) { oocoutI(nullptr, Fitting) << msg << std::endl; };
255 ROOT::Math::Util::TimingScope timingScope(print, "Function JIT time:");
256 _funcName = ctx.buildFunction(obj, dependsOnData);
257
258 // Make sure the codegen implementations are known to the interpreter
259 gInterpreter->Declare("#include <RooFit/CodegenImpl.h>\n");
260
261 if (!gInterpreter->Declare(ctx.collectedCode().c_str())) {
262 std::stringstream errorMsg;
263 std::string debugFileName = "_codegen_" + _funcName + ".cxx";
264 errorMsg << "Function " << _funcName << " could not be compiled. See above for details. Full code dumped to file "
265 << debugFileName << " for debugging";
266 {
267 std::ofstream outFile;
268 outFile.open(debugFileName.c_str());
269 outFile << ctx.collectedCode();
270 }
271 oocoutE(nullptr, InputArguments) << errorMsg.str() << std::endl;
272 throw std::runtime_error(errorMsg.str().c_str());
273 }
274
275 _func = reinterpret_cast<Func>(gInterpreter->ProcessLine((_funcName + ";").c_str()));
276
277 _xlArr = ctx.xlArr();
278 _collectedFunctions = ctx.collectedFunctions();
279}
280
281void RooFuncWrapper::loadData(RooAbsData const &data, RooSimultaneous const *simPdf)
282{
283 // Extract observables
284 std::stack<std::vector<double>> vectorBuffers; // for data loading
285 auto spans = RooFit::BatchModeDataHelpers::getDataSpans(data, "", simPdf, true, false, vectorBuffers);
286
287 _observables.clear();
288 // The first elements contain the sizes of the packed observable arrays
289 std::size_t total = 0;
290 _observables.reserve(2 * spans.size());
291 std::size_t idx = 0;
292 for (auto const &item : spans) {
293 _obsInfos.emplace(item.first, idx);
294 _observables.push_back(total + 2 * spans.size());
295 _observables.push_back(item.second.size());
296 total += item.second.size();
297 idx += 1;
298 }
299 idx = 0;
300 for (auto const &item : spans) {
301 std::size_t n = item.second.size();
302 _observables.reserve(_observables.size() + n);
303 for (std::size_t i = 0; i < n; ++i) {
304 _observables.push_back(item.second[i]);
305 }
306 idx += n;
307 }
308}
309
310void RooFuncWrapper::createGradient()
311{
312#ifdef ROOFIT_CLAD
313 std::string gradName = _funcName + "_grad_0";
314 std::string requestName = _funcName + "_req";
315
316 // Calculate gradient
317 gInterpreter->Declare("#include <Math/CladDerivator.h>\n");
318 // disable clang-format for making the following code unreadable.
319 // clang-format off
320 std::stringstream requestFuncStrm;
321 requestFuncStrm << "#pragma clad ON\n"
322 "void " << requestName << "() {\n"
323 " clad::gradient(" << _funcName << ", \"params\");\n"
324 "}\n"
325 "#pragma clad OFF";
326 // clang-format on
327 auto print = [](std::string const &msg) { oocoutI(nullptr, Fitting) << msg << std::endl; };
328
329 bool cladSuccess = false;
330 {
331 ROOT::Math::Util::TimingScope timingScope(print, "Gradient generation time:");
332 cladSuccess = !gInterpreter->Declare(requestFuncStrm.str().c_str());
333 }
334 if (cladSuccess) {
335 std::stringstream errorMsg;
336 errorMsg << "Function could not be differentiated. See above for details.";
337 oocoutE(nullptr, InputArguments) << errorMsg.str() << std::endl;
338 throw std::runtime_error(errorMsg.str().c_str());
339 }
340
341 // Clad provides different overloads for the gradient, and we need to
342 // resolve to the one that we want. Without the static_cast, getting the
343 // function pointer would be ambiguous.
344 std::stringstream ss;
345 ROOT::Math::Util::TimingScope timingScope(print, "Gradient IR to machine code time:");
346 ss << "static_cast<void (*)(double *, double const *, double const *, double *)>(" << gradName << ");";
347 _grad = reinterpret_cast<Grad>(gInterpreter->ProcessLine(ss.str().c_str()));
348 _hasGradient = true;
349#else
350 _hasGradient = false;
351 std::stringstream errorMsg;
352 errorMsg << "Function could not be differentiated since ROOT was built without Clad support.";
353 oocoutE(nullptr, InputArguments) << errorMsg.str() << std::endl;
354 throw std::runtime_error(errorMsg.str().c_str());
355#endif
356}
357
358void RooFuncWrapper::updateGradientVarBuffer() const
359{
360 std::transform(_params.begin(), _params.end(), _gradientVarBuffer.begin(), [](RooAbsArg *obj) {
361 return obj->isCategory() ? static_cast<RooAbsCategory *>(obj)->getCurrentIndex()
362 : static_cast<RooAbsReal *>(obj)->getVal();
363 });
364}
365
366/// @brief Dumps a macro "filename.C" that can be used to test and debug the generated code and gradient.
367void RooFuncWrapper::writeDebugMacro(std::string const &filename) const
368{
369 std::stringstream allCode;
370 std::set<std::string> seenFunctions;
371
372 // Remove duplicated declared functions
373 for (std::string const &name : _collectedFunctions) {
374 if (seenFunctions.count(name) > 0) {
375 continue;
376 }
377 seenFunctions.insert(name);
378 std::unique_ptr<TInterpreterValue> v = gInterpreter->MakeInterpreterValue();
379 gInterpreter->Evaluate(name.c_str(), *v);
380 std::string s = v->ToString();
381 for (int i = 0; i < 2; ++i) {
382 s = s.erase(0, s.find("\n") + 1);
383 }
384 allCode << s << std::endl;
385 }
386
387 std::ofstream outFile;
388 outFile.open(filename + ".C");
389 outFile << R"(//auto-generated test macro
390#include <RooFit/Detail/MathFuncs.h>
391#include <Math/CladDerivator.h>
392
393)" << allCode.str()
394 << R"(
395#pragma clad ON
396void gradient_request() {
397 clad::gradient()"
398 << _funcName << R"(, "params");
399}
400#pragma clad OFF
401)";
402
404
405 auto writeVector = [&](std::string const &name, std::span<const double> vec) {
406 std::stringstream decl;
407 decl << "std::vector<double> " << name << " = {";
408 for (std::size_t i = 0; i < vec.size(); ++i) {
409 if (i % 10 == 0)
410 decl << "\n ";
411 decl << vec[i];
412 if (i < vec.size() - 1)
413 decl << ", ";
414 }
415 decl << "\n};\n";
416
417 std::string declStr = decl.str();
418
419 replaceAll(declStr, "inf", "std::numeric_limits<double>::infinity()");
420 replaceAll(declStr, "nan", "NAN");
421
422 outFile << declStr;
423 };
424
425 outFile << "// clang-format off\n" << std::endl;
426 writeVector("parametersVec", _gradientVarBuffer);
427 outFile << std::endl;
428 writeVector("observablesVec", _observables);
429 outFile << std::endl;
430 writeVector("auxConstantsVec", _xlArr);
431 outFile << std::endl;
432 outFile << "// clang-format on\n" << std::endl;
433
434 outFile << R"(
435// To run as a ROOT macro
436void )" << filename
437 << R"(()
438{
439 std::vector<double> gradientVec(parametersVec.size());
440
441 auto func = [&](std::span<double> params) {
442 return )"
443 << _funcName << R"((params.data(), observablesVec.data(), auxConstantsVec.data());
444 };
445 auto grad = [&](std::span<double> params, std::span<double> out) {
446 return )"
447 << _funcName << R"(_grad_0(parametersVec.data(), observablesVec.data(), auxConstantsVec.data(),
448 out.data());
449 };
450
451 grad(parametersVec, gradientVec);
452
453 auto numDiff = [&](int i) {
454 const double eps = 1e-6;
455 std::vector<double> p{parametersVec};
456 p[i] = parametersVec[i] - eps;
457 double funcValDown = func(p);
458 p[i] = parametersVec[i] + eps;
459 double funcValUp = func(p);
460 return (funcValUp - funcValDown) / (2 * eps);
461 };
462
463 for (std::size_t i = 0; i < parametersVec.size(); ++i) {
464 std::cout << i << ":" << std::endl;
465 std::cout << " numr : " << numDiff(i) << std::endl;
466 std::cout << " clad : " << gradientVec[i] << std::endl;
467 }
468}
469)";
470}
471
472double RooEvaluatorWrapper::evaluate() const
473{
475 return _funcWrapper->evaluate();
476
477 if (!_evaluator)
478 return 0.0;
479
480 _evaluator->setOffsetMode(hideOffset() ? RooFit::EvalContext::OffsetMode::WithoutOffset
481 : RooFit::EvalContext::OffsetMode::WithOffset);
482
483 return _evaluator->run()[0];
484}
485
486bool RooEvaluatorWrapper::setData(RooAbsData &data, bool /*cloneData*/)
487{
488 // To make things easier for RooFit, we only support resetting with
489 // datasets that have the same structure, e.g. the same columns and global
490 // observables. This is anyway the usecase: resetting same-structured data
491 // when iterating over toys.
492 constexpr auto errMsg = "Error in RooAbsReal::setData(): only resetting with same-structured data is supported.";
493
494 _data = &data;
495 bool isInitializing = _paramSet.empty();
496 const std::size_t oldSize = _dataSpans.size();
497
498 std::stack<std::vector<double>>{}.swap(_vectorBuffers);
499 bool skipZeroWeights = !_pdf || !_pdf->getAttribute("BinnedLikelihoodActive");
500 auto simPdf = dynamic_cast<RooSimultaneous const *>(_pdf);
501 _dataSpans = RooFit::BatchModeDataHelpers::getDataSpans(*_data, _rangeName, simPdf, skipZeroWeights,
502 _takeGlobalObservablesFromData, _vectorBuffers);
503 if (!isInitializing && _dataSpans.size() != oldSize) {
504 coutE(DataHandling) << errMsg << std::endl;
505 throw std::runtime_error(errMsg);
506 }
507 for (auto const &item : _dataSpans) {
508 const char *name = item.first->GetName();
509 _evaluator->setInput(name, item.second, false);
510 if (_paramSet.find(name)) {
511 coutE(DataHandling) << errMsg << std::endl;
512 throw std::runtime_error(errMsg);
513 }
514 }
515 if (_funcWrapper) {
516 _funcWrapper->loadData(*_data, simPdf);
517 }
518 return true;
519}
520
521void RooEvaluatorWrapper::createFuncWrapper()
522{
523 // Get the parameters.
525 this->getParameters(_data ? _data->get() : nullptr, paramSet, /*sripDisconnectedParams=*/false);
526
528 std::make_unique<RooFuncWrapper>(*_topNode, _data, dynamic_cast<RooSimultaneous const *>(_pdf), paramSet);
529}
530
531void RooEvaluatorWrapper::generateGradient()
532{
533 if (!_funcWrapper)
535 _funcWrapper->createGradient();
536}
537
538void RooEvaluatorWrapper::setUseGeneratedFunctionCode(bool flag)
539{
543}
544
545void RooEvaluatorWrapper::gradient(double *out) const
546{
547 _funcWrapper->gradient(out);
548}
549
550bool RooEvaluatorWrapper::hasGradient() const
551{
552 if (!_funcWrapper)
553 return false;
554 return _funcWrapper->hasGradient();
555}
556
557void RooEvaluatorWrapper::writeDebugMacro(std::string const &filename) const
558{
559 if (_funcWrapper)
560 return _funcWrapper->writeDebugMacro(filename);
561}
562
563} // namespace RooFit::Experimental
564
565/// \endcond
#define oocoutE(o, a)
#define oocoutI(o, a)
#define coutE(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
static unsigned int total
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
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
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.
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 buildFunction(RooAbsArg const &arg, std::unordered_set< RooFit::Detail::DataKey > const &dependsOnData={})
Assemble and return the final code with the return expression and global statements.
std::vector< std::string > const & collectedFunctions()
std::vector< double > const & xlArr()
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
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:67
void getSortedComputationGraph(RooAbsArg const &func, RooArgSet &out)
void evaluate(typename Architecture_t::Tensor_t &A, EActivationFunction f)
Apply the given activation function to each value in the given tensor A.
Definition Functions.h:98