Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAddHelpers.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*
4 * Project: RooFit
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 "RooAddHelpers.h"
14
15#include <RooAbsPdf.h>
16#include <RooArgSet.h>
17#include <RooNaNPacker.h>
18#include <RooProduct.h>
19#include <RooRatio.h>
20#include <RooRealConstant.h>
21#include <RooRealIntegral.h>
22#include <RooRealVar.h>
23
24////////////////////////////////////////////////////////////////////////////////
25/// Create a RooAddPdf cache element for a given normalization set and
26/// projection configuration.
27
28AddCacheElem::AddCacheElem(RooAbsPdf const &addPdf, RooArgList const &pdfList, RooArgList const &coefList,
29 const RooArgSet *nset, const RooArgSet *iset, RooArgSet const &refCoefNormSet,
30 std::string const &refCoefNormRange, int verboseEval)
31{
32 // Projection integrals are always over all pdf components. Overriding the
33 // global component selection temporarily makes all RooRealIntegrals created
34 // during that time always include all components.
36
37 // We put the normRange into a std::string to not have to deal with
38 // nullptr vs. "" ambiguities
39 const std::string normRange = addPdf.normRange() ? addPdf.normRange() : "";
40
41 _suppNormList.reserve(pdfList.size());
42
43 // *** PART 1 : Create supplemental normalization list ***
44
45 // Retrieve the combined set of dependents of this PDF ;
47 addPdf.getObservables(nset, fullDepList);
48 if (iset) {
49 fullDepList.remove(*iset, true, true);
50 }
51
52 // Reduce iset/nset to actual dependents of this PDF
54 addPdf.getObservables(nset, nset2);
55
56 if (nset2.empty() && !refCoefNormSet.empty()) {
57 // Evaluating RooAddPdf without normalization, but have reference normalization for coefficient definition
59 }
60
61 bool hasPdfWithDifferentRange = false;
62
63 // Fill with dummy unit RRVs for now
64 for (std::size_t i = 0; i < pdfList.size(); ++i) {
65 auto pdf = static_cast<const RooAbsPdf *>(pdfList.at(i));
66 auto coef = static_cast<const RooAbsReal *>(coefList.at(i));
67
68 const std::string normRangeComponent = pdf->normRange() ? pdf->normRange() : "";
69 const bool componentHasDifferentNormRange = normRangeComponent != normRange;
70
72
73 // Start with full list of dependents
75
76 // Remove PDF dependents
77 if (auto pdfDeps = std::unique_ptr<RooArgSet>{pdf->getObservables(nset)}) {
78 supNSet.remove(*pdfDeps, true, true);
79 }
80
81 // Remove coef dependents
82 if (auto coefDeps = std::unique_ptr<RooArgSet>{coef ? coef->getObservables(nset) : nullptr}) {
83 supNSet.remove(*coefDeps, true, true);
84 }
85
86 std::unique_ptr<RooAbsReal> snorm;
87 auto name = std::string(addPdf.GetName()) + "_" + pdf->GetName() + "_SupNorm";
88 if (!supNSet.empty()) {
89 snorm = std::make_unique<RooRealIntegral>(name.c_str(), "Supplemental normalization integral",
91 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << " " << addPdf.GetName()
92 << " making supplemental normalization set " << supNSet << " for pdf component "
93 << pdf->GetName() << std::endl;
94 }
95
97 auto snormTerm = std::unique_ptr<RooAbsReal>(pdf->createIntegral(nset2, nset2, normRange.c_str()));
98 if (snorm) {
99 auto oldSnorm = std::move(snorm);
100 snorm = std::make_unique<RooProduct>("snorm", "snorm", *oldSnorm.get(), *snormTerm);
101 snorm->addOwnedComponents(std::move(snormTerm), std::move(oldSnorm));
102 } else {
103 snorm = std::move(snormTerm);
104 }
105 }
106 _suppNormList.emplace_back(std::move(snorm));
107 }
108
109 if (verboseEval > 1) {
110 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << "::syncSuppNormList(" << addPdf.GetName()
111 << ") synching supplemental normalization list for norm"
112 << (nset ? *nset : RooArgSet()) << std::endl;
113 }
114
115 // *** PART 2 : Create projection coefficients ***
116
117 const bool projectCoefsForRangeReasons = !refCoefNormRange.empty() || !normRange.empty() || hasPdfWithDifferentRange;
118
119 // If no projections required stop here
121 return;
122 }
123
125
126 // Recalculate projection integrals of PDFs
127 for (auto *pdf : static_range_cast<const RooAbsPdf *>(pdfList)) {
128
129 // Calculate projection integral
130 std::unique_ptr<RooAbsReal> pdfProj;
131 if (!refCoefNormSet.empty() && !nset2.equals(refCoefNormSet)) {
132 pdfProj = std::unique_ptr<RooAbsReal>{pdf->createIntegral(nset2, refCoefNormSet, normRange.c_str())};
133 pdfProj->setOperMode(addPdf.operMode());
134 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << "(" << addPdf.GetName() << ")::getPC nset2(" << nset2
135 << ")!=_refCoefNormSet(" << refCoefNormSet
136 << ") --> pdfProj = " << pdfProj->GetName() << std::endl;
137 oocxcoutD(&addPdf, Caching) << " " << addPdf.ClassName() << "::syncCoefProjList(" << addPdf.GetName()
138 << ") PP = " << pdfProj->GetName() << std::endl;
139 }
140
141 _projList.emplace_back(std::move(pdfProj));
142
143 // Calculation optional supplemental normalization term
145 auto deps = std::unique_ptr<RooArgSet>{pdf->getParameters(RooArgSet())};
146 supNormSet.remove(*deps, true, true);
147
148 std::unique_ptr<RooAbsReal> snorm;
149 auto name = std::string(addPdf.GetName()) + "_" + pdf->GetName() + "_ProjSupNorm";
150 if (!supNormSet.empty() && !nset2.equals(refCoefNormSet)) {
151 snorm = std::make_unique<RooRealIntegral>(name.c_str(), "Projection Supplemental normalization integral",
153 oocxcoutD(&addPdf, Caching) << " " << addPdf.ClassName() << "::syncCoefProjList(" << addPdf.GetName()
154 << ") SN = " << snorm->GetName() << std::endl;
155 }
156 _suppProjList.emplace_back(std::move(snorm));
157
158 // Calculate range adjusted projection integral
159 std::unique_ptr<RooAbsReal> rangeProj2;
160 if (normRange != refCoefNormRange) {
161 tmp;
162 pdf->getObservables(tmp);
163 auto int1 = std::unique_ptr<RooAbsReal>{pdf->createIntegral(tmp, tmp, normRange.c_str())};
164 auto int2 = std::unique_ptr<RooAbsReal>{pdf->createIntegral(tmp, tmp, refCoefNormRange.c_str())};
165 rangeProj2 = std::make_unique<RooRatio>("rangeProj", "rangeProj", *int1, *int2);
166 rangeProj2->addOwnedComponents(std::move(int1), std::move(int2));
167 }
168
169 _rangeProjList.emplace_back(std::move(rangeProj2));
170 }
171 }
172}
173
174void AddCacheElem::print() const
175{
176 auto printVector = [](auto const &vec, const char *name) {
177 std::cout << "+++ " << name << ":" << std::endl;
178 for (auto const &arg : vec) {
179 std::cout << " ";
180 if (arg) {
181 arg->Print();
182 } else {
183 std::cout << "nullptr" << std::endl;
184 }
185 }
186 };
187
188 printVector(_suppNormList, "_suppNormList");
189 printVector(_projList, "_projList");
190 printVector(_suppProjList, "_suppProjList");
191 printVector(_rangeProjList, "_rangeProjList");
192}
193
194////////////////////////////////////////////////////////////////////////////////
195/// List all RooAbsArg derived contents in this cache element
196
197RooArgList AddCacheElem::containedArgs(Action)
198{
200 // need to iterate manually because _suppProjList can contain nullptr
201 for (auto const &arg : _projList) {
202 if (arg)
203 allNodes.add(*arg);
204 }
205 for (auto const &arg : _suppProjList) {
206 if (arg)
207 allNodes.add(*arg);
208 }
209 for (auto const &arg : _rangeProjList) {
210 if (arg)
211 allNodes.add(*arg);
212 }
213
214 return allNodes;
215}
216
217////////////////////////////////////////////////////////////////////////////////
218/// Update the RooAddPdf coefficients for a given normalization set and
219/// projection configuration. The `coefCache` argument should have the same
220/// size as `pdfList`. It needs to be initialized with the raw values of the
221/// coefficients, as obtained from the `_coefList` proxy in the RooAddPdf. If
222/// the last coefficient is not given, the initial value of the last element of
223/// `_coefCache` does not matter. After this function, the `_coefCache` will be
224/// filled with the correctly scaled coefficients for each pdf.
225
226void RooAddHelpers::updateCoefficients(RooAbsPdf const &addPdf, std::vector<double> &coefCache,
227 RooArgList const &pdfList, bool haveLastCoef, AddCacheElem &cache,
228 const RooArgSet *nset, RooArgSet const &refCoefNormSet, bool allExtendable,
229 int &coefErrCount)
230{
231 // Straight coefficients
232 if (allExtendable) {
233
234 // coef[i] = expectedEvents[i] / SUM(expectedEvents)
235 double coefSum(0);
236 std::size_t i = 0;
237 for (auto arg : pdfList) {
238 auto pdf = static_cast<RooAbsPdf *>(arg);
239 coefCache[i] = pdf->expectedEvents(!refCoefNormSet.empty() ? &refCoefNormSet : nset);
240 coefSum += coefCache[i];
241 i++;
242 }
243
244 if (coefSum == 0.) {
245 oocoutW(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
246 << ") WARNING: total number of expected events is 0" << std::endl;
247 } else {
248 for (std::size_t j = 0; j < pdfList.size(); j++) {
249 coefCache[j] /= coefSum;
250 }
251 }
252
253 } else {
254 if (haveLastCoef) {
255
256 // coef[i] = coef[i] / SUM(coef)
257 double coefSum = std::accumulate(coefCache.begin(), coefCache.end(), 0.0);
258 if (coefSum == 0.) {
259 oocoutW(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
260 << ") WARNING: sum of coefficients is zero 0" << std::endl;
261 } else {
262 const double invCoefSum = 1. / coefSum;
263 for (std::size_t j = 0; j < coefCache.size(); j++) {
265 }
266 }
267 } else {
268
269 // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
270 double lastCoef = 1.0 - std::accumulate(coefCache.begin(), coefCache.end() - 1, 0.0);
271 coefCache.back() = lastCoef;
272
273 // Treat coefficient degeneration
274 const float coefDegen = lastCoef < 0. ? -lastCoef : (lastCoef > 1. ? lastCoef - 1. : 0.);
275 if (coefDegen > 1.E-5) {
277
278 std::stringstream msg;
279 if (coefErrCount-- > 0) {
280 msg << "RooAddPdf::updateCoefCache(" << addPdf.GetName()
281 << " WARNING: sum of PDF coefficients not in range [0-1], value=" << 1 - lastCoef;
282 if (coefErrCount == 0) {
283 msg << " (no more will be printed)";
284 }
285 oocoutW(&addPdf, Eval) << msg.str() << std::endl;
286 }
287 }
288 }
289 }
290
291 // Stop here if not projection is required or needed
292 if (!cache.doProjection()) {
293 return;
294 }
295
296 // Adjust coefficients for given projection
297 double coefSum(0);
298 for (std::size_t i = 0; i < pdfList.size(); i++) {
299 coefCache[i] *= cache.projVal(i) / cache.projSuppNormVal(i) * cache.rangeProjScaleFactor(i);
300 coefSum += coefCache[i];
301 }
302
303 if ((RooMsgService::_debugCount > 0) &&
305 for (std::size_t i = 0; i < pdfList.size(); ++i) {
306 ooccoutD(&addPdf, Caching) << " ALEX: POST-SYNC coef[" << i << "] = " << coefCache[i]
307 << " ( _coefCache[i]/coefSum = " << coefCache[i] * coefSum << "/" << coefSum
308 << " ) " << std::endl;
309 }
310 }
311
312 if (coefSum == 0.) {
313 oocoutE(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
314 << ") sum of coefficients is zero." << std::endl;
315 }
316
317 for (std::size_t i = 0; i < pdfList.size(); i++) {
318 coefCache[i] /= coefSum;
319 }
320}
321
322/// \endcond
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
#define oocoutW(o, a)
#define oocxcoutD(o, a)
#define oocoutE(o, a)
#define ooccoutD(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char name[80]
Definition TGX11.cxx:110
const_iterator begin() const
const_iterator end() const
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.
OperMode operMode() const
Query the operation mode of this node.
Definition RooAbsArg.h:425
Storage_t::size_type size() const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
const char * normRange() const
Definition RooAbsPdf.h:250
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
static RooMsgService & instance()
Return reference to singleton instance.
static Int_t _debugCount
static RooConstVar & value(double value)
Return a constant value object with given value.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
static double packFloatIntoNaN(float payload)
Pack float into mantissa of a NaN.