Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooNLLVarNew.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*
4 * Project: RooFit
5 * Authors:
6 * Jonas Rembser, CERN 2021
7 * Emmanouil Michalainas, CERN 2021
8 *
9 * Copyright (c) 2021, CERN
10 *
11 * Redistribution and use in source and binary forms,
12 * with or without modification, are permitted according to the terms
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
14 */
15
16/**
17\file RooNLLVarNew.cxx
18\class RooNLLVarNew
19\ingroup Roofitcore
20
21This is a simple class designed to produce the nll values needed by the fitter.
22This class calls functions from `RooBatchCompute` library to provide faster
23computation times.
24**/
25
27
28#include <RooHistPdf.h>
29#include <RooBatchCompute.h>
30#include <RooDataHist.h>
31#include <RooNaNPacker.h>
32#include <RooConstVar.h>
33#include <RooRealVar.h>
34#include <RooSetProxy.h>
36
37#include "RooFitImplHelpers.h"
38
39#include <ROOT/StringUtils.hxx>
40
41#include <TClass.h>
42#include <TMath.h>
43#include <Math/Util.h>
44
45#include <numeric>
46#include <stdexcept>
47#include <vector>
48
49
50namespace RooFit {
51namespace Detail {
52
53// Declare constexpr static members to make them available if odr-used in C++14.
54constexpr const char *RooNLLVarNew::weightVarName;
55constexpr const char *RooNLLVarNew::weightVarNameSumW2;
56
57namespace {
58
59// Use RooConstVar for dummies such that they don't get included in getParameters().
60std::unique_ptr<RooConstVar> dummyVar(const char *name)
61{
62 return std::make_unique<RooConstVar>(name, name, 1.0);
63}
64
65// Helper class to represent a template pdf based on the fit dataset.
66class RooOffsetPdf : public RooAbsPdf {
67public:
68 RooOffsetPdf(const char *name, const char *title, RooArgSet const &observables, RooAbsReal &weightVar)
69 : RooAbsPdf(name, title),
70 _observables("!observables", "List of observables", this),
71 _weightVar{"!weightVar", "weightVar", this, weightVar, true, false}
72 {
73 for (RooAbsArg *obs : observables) {
74 _observables.add(*obs);
75 }
76 }
77 RooOffsetPdf(const RooOffsetPdf &other, const char *name = nullptr)
79 _observables("!servers", this, other._observables),
80 _weightVar{"!weightVar", this, other._weightVar}
81 {
82 }
83 TObject *clone(const char *newname) const override { return new RooOffsetPdf(*this, newname); }
84
85 void doEval(RooFit::EvalContext &ctx) const override
86 {
87 std::span<double> output = ctx.output();
88 std::size_t nEvents = output.size();
89
90 std::span<const double> weights = ctx.at(_weightVar);
91
92 // Create the template histogram from the data. This operation is very
93 // expensive, but since the offset only depends on the observables it
94 // only has to be done once.
95
96 RooDataHist dataHist{"data", "data", _observables};
97 // Loop over events to fill the histogram
98 for (std::size_t i = 0; i < nEvents; ++i) {
99 for (auto *var : static_range_cast<RooRealVar *>(_observables)) {
100 var->setVal(ctx.at(var)[i]);
101 }
102 dataHist.add(_observables, weights[weights.size() == 1 ? 0 : i]);
103 }
104
105 // Lookup bin weights via RooHistPdf
106 RooHistPdf pdf{"offsetPdf", "offsetPdf", _observables, dataHist};
107 for (std::size_t i = 0; i < nEvents; ++i) {
108 for (auto *var : static_range_cast<RooRealVar *>(_observables)) {
109 var->setVal(ctx.at(var)[i]);
110 }
111 output[i] = pdf.getVal(_observables);
112 }
113 }
114
115private:
116 double evaluate() const override { return 0.0; } // should never be called
117
118 RooSetProxy _observables;
120};
121
122} // namespace
123
124/** Construct a RooNLLVarNew
125\param name the name
126\param title the title
127\param pdf The pdf for which the nll is computed for
128\param observables The observabes of the pdf
129\param isExtended Set to true if this is an extended fit
130**/
131RooNLLVarNew::RooNLLVarNew(const char *name, const char *title, RooAbsPdf &pdf, RooArgSet const &observables,
132 bool isExtended, RooFit::OffsetMode offsetMode)
133 : RooAbsReal(name, title),
134 _pdf{"pdf", "pdf", this, pdf},
135 _weightVar{"weightVar", "weightVar", this, dummyVar(weightVarName)},
137 _binnedL{pdf.getAttribute("BinnedLikelihoodActive")}
138{
139 RooArgSet obs;
140 pdf.getObservables(&observables, obs);
141
142 // In the "BinnedLikelihoodActiveYields" mode, the pdf values can directly
143 // be interpreted as yields and don't need to be multiplied by the bin
144 // widths. That's why we don't need to even fill them in this case.
145 if (_binnedL && !pdf.getAttribute("BinnedLikelihoodActiveYields")) {
147 }
148
149 if (isExtended && !_binnedL) {
150 std::unique_ptr<RooAbsReal> expectedEvents = pdf.createExpectedEventsFunc(&obs);
151 if (expectedEvents) {
152 _expectedEvents =
153 std::make_unique<RooTemplateProxy<RooAbsReal>>("expectedEvents", "expectedEvents", this, *expectedEvents);
154 addOwnedComponents(std::move(expectedEvents));
155 }
156 }
157
159 enableOffsetting(offsetMode == RooFit::OffsetMode::Initial);
161
162 // In the binned likelihood code path, we directly use that data weights for
163 // the offsetting.
164 if (!_binnedL && _doBinOffset) {
165 auto offsetPdf = std::make_unique<RooOffsetPdf>("_offset_pdf", "_offset_pdf", obs, *_weightVar);
166 _offsetPdf = std::make_unique<RooTemplateProxy<RooAbsPdf>>("offsetPdf", "offsetPdf", this, *offsetPdf);
167 addOwnedComponents(std::move(offsetPdf));
168 }
169}
170
171RooNLLVarNew::RooNLLVarNew(const RooNLLVarNew &other, const char *name)
173 _pdf{"pdf", this, other._pdf},
174 _weightVar{"weightVar", this, other._weightVar},
175 _weightSquaredVar{"weightSquaredVar", this, other._weightSquaredVar},
180 _prefix{other._prefix},
181 _binw{other._binw}
182{
183 if (other._expectedEvents) {
184 _expectedEvents = std::make_unique<RooTemplateProxy<RooAbsReal>>("expectedEvents", this, *other._expectedEvents);
185 }
186}
187
188void RooNLLVarNew::fillBinWidthsFromPdfBoundaries(RooAbsReal const &pdf, RooArgSet const &observables)
189{
190 // Check if the bin widths were already filled
191 if (!_binw.empty()) {
192 return;
193 }
194
195 if (observables.size() != 1) {
196 throw std::runtime_error("BinnedPdf optimization only works with a 1D pdf.");
197 } else {
198 auto *var = static_cast<RooRealVar *>(observables.first());
199 std::list<double> *boundaries = pdf.binBoundaries(*var, var->getMin(), var->getMax());
200 std::list<double>::iterator biter = boundaries->begin();
201 _binw.resize(boundaries->size() - 1);
202 double lastBound = (*biter);
203 ++biter;
204 int ibin = 0;
205 while (biter != boundaries->end()) {
206 _binw[ibin] = (*biter) - lastBound;
207 lastBound = (*biter);
208 ibin++;
209 ++biter;
210 }
211 }
212}
213
214void RooNLLVarNew::doEvalBinnedL(RooFit::EvalContext &ctx, std::span<const double> preds,
215 std::span<const double> weights) const
216{
219
220 const bool predsAreYields = _binw.empty();
221
222 for (std::size_t i = 0; i < preds.size(); ++i) {
223
224 // Calculate log(Poisson(N|mu) for this bin
225 double N = weights[i];
226 double mu = preds[i];
227 if (!predsAreYields) {
228 mu *= _binw[i];
229 }
230
231 if (mu <= 0 && N > 0) {
232 // Catch error condition: data present where zero events are predicted
233 logEvalError(Form("Observed %f events in bin %lu with zero event yield", N, (unsigned long)i));
234 } else {
237 }
238 }
239
241}
242
243void RooNLLVarNew::doEval(RooFit::EvalContext &ctx) const
244{
245 std::span<const double> weights = ctx.at(_weightVar);
246 std::span<const double> weightsSumW2 = ctx.at(_weightSquaredVar);
247
248 if (_binnedL) {
249 return doEvalBinnedL(ctx, ctx.at(&*_pdf), _weightSquared ? weightsSumW2 : weights);
250 }
251
252 auto config = ctx.config(this);
253
254 auto probas = ctx.at(_pdf);
255
256 _sumWeight = weights.size() == 1 ? weights[0] * probas.size()
257 : RooBatchCompute::reduceSum(config, weights.data(), weights.size());
258 if (_expectedEvents && _weightSquared && _sumWeight2 == 0.0) {
259 _sumWeight2 = weights.size() == 1 ? weightsSumW2[0] * probas.size()
260 : RooBatchCompute::reduceSum(config, weightsSumW2.data(), weightsSumW2.size());
261 }
262
263 auto nllOut = RooBatchCompute::reduceNLL(config, probas, _weightSquared ? weightsSumW2 : weights,
264 _doBinOffset ? ctx.at(*_offsetPdf) : std::span<const double>{});
265
266 if (nllOut.nInfiniteValues > 0) {
267 oocoutW(&*_pdf, Eval) << "RooAbsPdf::getLogVal(" << _pdf->GetName()
268 << ") WARNING: top-level pdf has some infinite values" << std::endl;
269 }
270 for (std::size_t i = 0; i < nllOut.nNonPositiveValues; ++i) {
271 _pdf->logEvalError("getLogVal() top-level p.d.f not greater than zero");
272 }
273 for (std::size_t i = 0; i < nllOut.nNaNValues; ++i) {
274 _pdf->logEvalError("getLogVal() top-level p.d.f evaluates to NaN");
275 }
276
277 if (_expectedEvents) {
278 std::span<const double> expected = ctx.at(*_expectedEvents);
279 nllOut.nllSum += _pdf->extendedTerm(_sumWeight, expected[0], _weightSquared ? _sumWeight2 : 0.0, _doBinOffset);
280 }
281
282 finalizeResult(ctx, {nllOut.nllSum, nllOut.nllSumCarry}, _sumWeight);
283}
284
285////////////////////////////////////////////////////////////////////////////////
286/// Sets the prefix for the special variables of this NLL, like weights or bin
287/// volumes.
288/// \param[in] prefix The prefix to add to the observables and weight names.
289void RooNLLVarNew::setPrefix(std::string const &prefix)
290{
291 _prefix = prefix;
292
294}
295
296void RooNLLVarNew::resetWeightVarNames()
297{
298 _weightVar->SetName((_prefix + weightVarName).c_str());
299 _weightSquaredVar->SetName((_prefix + weightVarNameSumW2).c_str());
300 if (_offsetPdf) {
301 (*_offsetPdf)->SetName((_prefix + "_offset_pdf").c_str());
302 }
303}
304
305////////////////////////////////////////////////////////////////////////////////
306/// Toggles the weight square correction.
307void RooNLLVarNew::applyWeightSquared(bool flag)
308{
310}
311
312void RooNLLVarNew::enableOffsetting(bool flag)
313{
314 _doOffset = flag;
316}
317
318void RooNLLVarNew::finalizeResult(RooFit::EvalContext &ctx, ROOT::Math::KahanSum<double> result, double weightSum) const
319{
320 // If part of simultaneous PDF normalize probability over
321 // number of simultaneous PDFs: -sum(log(p/n)) = -sum(log(p)) + N*log(n)
322 // If we do bin-by bin offsetting, we don't do this because it cancels out
323 if (!_doBinOffset && _simCount > 1) {
324 result += weightSum * std::log(static_cast<double>(_simCount));
325 }
326
327 // Check if value offset flag is set.
328 if (_doOffset) {
329
330 // If no offset is stored enable this feature now
331 if (_offset.Sum() == 0 && _offset.Carry() == 0 && (result.Sum() != 0 || result.Carry() != 0)) {
332 _offset = result;
333 }
334 }
335 ctx.setOutputWithOffset(this, result, _offset);
336}
337
338} // namespace Detail
339} // namespace RooFit
340
341/// \endcond
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
#define oocoutW(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define N
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 result
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
The Kahan summation is a compensated summation algorithm, which significantly reduces numerical error...
Definition Util.h:136
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
Storage_t::size_type size() const
RooAbsArg * first() const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const
Returns an object that represents the expected number of events for a given normalization set,...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
virtual std::list< double > * binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
void add(const RooArgSet &row, double wgt=1.0) override
Add wgt to the bin content enclosed by the coordinates passed in row.
Definition RooDataHist.h:72
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
void setOutputWithOffset(RooAbsArg const *arg, ROOT::Math::KahanSum< double > val, ROOT::Math::KahanSum< double > const &offset)
Sets the output value with an offset.
A probability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41
double reduceSum(Config cfg, InputArr input, size_t n)
ReduceNLLOutput reduceNLL(Config cfg, std::span< const double > probas, std::span< const double > weights, std::span< const double > offsetProbas)
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:379
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
OffsetMode
For setting the offset mode with the Offset() command argument to RooAbsPdf::fitTo()
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
void probas(TString dataset, TString fin="TMVA.root", Bool_t useTMVAStyle=kTRUE)
static void output()