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