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 <RooRealConstant.h>
17#include <RooRealIntegral.h>
18#include <RooRealVar.h>
19#include <RooRatio.h>
20
21////////////////////////////////////////////////////////////////////////////////
22/// Create a RooAddPdf cache element for a given normalization set and
23/// projection configuration.
24
25AddCacheElem::AddCacheElem(RooAbsPdf const &addPdf, RooArgList const &pdfList, RooArgList const &coefList,
26 const RooArgSet *nset, const RooArgSet *iset, RooArgSet const &refCoefNormSet,
27 std::string const &refCoefNormRange, int verboseEval)
28{
29 // We put the normRange into a std::string to not have to deal with
30 // nullptr vs. "" ambiguities
31 const std::string normRange = addPdf.normRange() ? addPdf.normRange() : "";
32
33 _suppNormList.reserve(pdfList.size());
34
35 // *** PART 1 : Create supplemental normalization list ***
36
37 // Retrieve the combined set of dependents of this PDF ;
38 RooArgSet fullDepList;
39 addPdf.getObservables(nset, fullDepList);
40 if (iset) {
41 fullDepList.remove(*iset, true, true);
42 }
43
44 bool hasPdfWithCustomRange = false;
45
46 // Fill with dummy unit RRVs for now
47 for (std::size_t i = 0; i < pdfList.size(); ++i) {
48 auto pdf = static_cast<const RooAbsPdf *>(pdfList.at(i));
49 auto coef = static_cast<const RooAbsReal *>(coefList.at(i));
50
51 hasPdfWithCustomRange |= pdf->normRange() != nullptr;
52
53 // Start with full list of dependents
54 RooArgSet supNSet(fullDepList);
55
56 // Remove PDF dependents
57 if (auto pdfDeps = std::unique_ptr<RooArgSet>{pdf->getObservables(nset)}) {
58 supNSet.remove(*pdfDeps, true, true);
59 }
60
61 // Remove coef dependents
62 if (auto coefDeps = std::unique_ptr<RooArgSet>{coef ? coef->getObservables(nset) : nullptr}) {
63 supNSet.remove(*coefDeps, true, true);
64 }
65
66 std::unique_ptr<RooAbsReal> snorm;
67 auto name = std::string(addPdf.GetName()) + "_" + pdf->GetName() + "_SupNorm";
68 if (!supNSet.empty()) {
69 snorm = std::make_unique<RooRealIntegral>(name.c_str(), "Supplemental normalization integral",
70 RooRealConstant::value(1.0), supNSet);
71 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << " " << addPdf.GetName()
72 << " making supplemental normalization set " << supNSet << " for pdf component "
73 << pdf->GetName() << std::endl;
74 }
75
76 if (!normRange.empty()) {
77 auto snormTerm = std::unique_ptr<RooAbsReal>(pdf->createIntegral(*nset, *nset, normRange.c_str()));
78 if (snorm) {
79 auto oldSnorm = std::move(snorm);
80 snorm = std::make_unique<RooProduct>("snorm", "snorm", *oldSnorm.get(), *snormTerm.get());
81 snorm->addOwnedComponents(std::move(snormTerm), std::move(oldSnorm));
82 } else {
83 snorm = std::move(snormTerm);
84 }
85 }
86 _suppNormList.emplace_back(std::move(snorm));
87 }
88
89 if (verboseEval > 1) {
90 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << "::syncSuppNormList(" << addPdf.GetName()
91 << ") synching supplemental normalization list for norm"
92 << (nset ? *nset : RooArgSet()) << std::endl;
93 }
94
95 // *** PART 2 : Create projection coefficients ***
96
97 const bool projectCoefsForRangeReasons = !refCoefNormRange.empty() || !normRange.empty() || hasPdfWithCustomRange;
98
99 // If no projections required stop here
100 if (refCoefNormSet.empty() && !projectCoefsForRangeReasons) {
101 return;
102 }
103
104 // Reduce iset/nset to actual dependents of this PDF
105 RooArgSet nset2;
106 if (nset)
107 addPdf.getObservables(nset, nset2);
108 oocxcoutD(&addPdf, Caching) << addPdf.ClassName() << "(" << addPdf.GetName()
109 << ")::getPC nset = " << (nset ? *nset : RooArgSet()) << " nset2 = " << nset2
110 << std::endl;
111
112 if (nset2.empty() && !refCoefNormSet.empty()) {
113 // cout << "WVE: evaluating RooAddPdf without normalization, but have reference normalization for coefficient
114 // definition" << std::endl ;
115
116 nset2.add(refCoefNormSet);
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 printVector(_suppNormList, "_suppNormList");
183 printVector(_projList, "_projList");
184 printVector(_suppProjList, "_suppProjList");
185 printVector(_rangeProjList, "_rangeProjList");
186}
187
188////////////////////////////////////////////////////////////////////////////////
189/// List all RooAbsArg derived contents in this cache element
190
192{
193 RooArgList allNodes;
194 // need to iterate manually because _suppProjList can contain nullptr
195 for (auto const &arg : _projList) {
196 if (arg)
197 allNodes.add(*arg);
198 }
199 for (auto const &arg : _suppProjList) {
200 if (arg)
201 allNodes.add(*arg);
202 }
203 for (auto const &arg : _rangeProjList) {
204 if (arg)
205 allNodes.add(*arg);
206 }
207
208 return allNodes;
209}
210
211////////////////////////////////////////////////////////////////////////////////
212/// Update the RooAddPdf coefficients for a given normalization set and
213/// projection configuration. The `coefCache` argument should have the same
214/// size as `pdfList`. It needs to be initialized with the raw values of the
215/// coefficients, as obtained from the `_coefList` proxy in the RooAddPdf. If
216/// the last coefficient is not given, the initial value of the last element of
217/// `_coefCache` does not matter. After this function, the `_coefCache` will be
218/// filled with the correctly scaled coefficients for each pdf.
219
220void RooAddHelpers::updateCoefficients(RooAbsPdf const &addPdf, std::vector<double> &coefCache,
221 RooArgList const &pdfList, bool haveLastCoef, AddCacheElem &cache,
222 const RooArgSet *nset, RooArgSet const &refCoefNormSet, bool allExtendable,
223 int &coefErrCount)
224{
225 // Straight coefficients
226 if (allExtendable) {
227
228 // coef[i] = expectedEvents[i] / SUM(expectedEvents)
229 double coefSum(0);
230 std::size_t i = 0;
231 for (auto arg : pdfList) {
232 auto pdf = static_cast<RooAbsPdf *>(arg);
233 coefCache[i] = pdf->expectedEvents(!refCoefNormSet.empty() ? &refCoefNormSet : nset);
234 coefSum += coefCache[i];
235 i++;
236 }
237
238 if (coefSum == 0.) {
239 oocoutW(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
240 << ") WARNING: total number of expected events is 0" << std::endl;
241 } else {
242 for (std::size_t j = 0; j < pdfList.size(); j++) {
243 coefCache[j] /= coefSum;
244 }
245 }
246
247 } else {
248 if (haveLastCoef) {
249
250 // coef[i] = coef[i] / SUM(coef)
251 double coefSum = std::accumulate(coefCache.begin(), coefCache.end(), 0.0);
252 if (coefSum == 0.) {
253 oocoutW(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
254 << ") WARNING: sum of coefficients is zero 0" << std::endl;
255 } else {
256 const double invCoefSum = 1. / coefSum;
257 for (std::size_t j = 0; j < coefCache.size(); j++) {
258 coefCache[j] *= invCoefSum;
259 }
260 }
261 } else {
262
263 // coef[i] = coef[i] ; coef[n] = 1-SUM(coef[0...n-1])
264 double lastCoef = 1.0 - std::accumulate(coefCache.begin(), coefCache.end() - 1, 0.0);
265 coefCache.back() = lastCoef;
266
267 // Treat coefficient degeneration
268 const float coefDegen = lastCoef < 0. ? -lastCoef : (lastCoef > 1. ? lastCoef - 1. : 0.);
269 if (coefDegen > 1.E-5) {
270 coefCache.back() = RooNaNPacker::packFloatIntoNaN(100.f * coefDegen);
271
272 std::stringstream msg;
273 if (coefErrCount-- > 0) {
274 msg << "RooAddPdf::updateCoefCache(" << addPdf.GetName()
275 << " WARNING: sum of PDF coefficients not in range [0-1], value=" << 1 - lastCoef;
276 if (coefErrCount == 0) {
277 msg << " (no more will be printed)";
278 }
279 oocoutW(&addPdf, Eval) << msg.str() << std::endl;
280 }
281 }
282 }
283 }
284
285 // Stop here if not projection is required or needed
286 if (!cache.doProjection()) {
287 return;
288 }
289
290 // Adjust coefficients for given projection
291 double coefSum(0);
292 {
294
295 for (std::size_t i = 0; i < pdfList.size(); i++) {
296 coefCache[i] *= cache.projVal(i) / cache.projSuppNormVal(i) * cache.rangeProjScaleFactor(i);
297 coefSum += coefCache[i];
298 }
299 }
300
301 if ((RooMsgService::_debugCount > 0) &&
303 for (std::size_t i = 0; i < pdfList.size(); ++i) {
304 ooccoutD(&addPdf, Caching) << " ALEX: POST-SYNC coef[" << i << "] = " << coefCache[i]
305 << " ( _coefCache[i]/coefSum = " << coefCache[i] * coefSum << "/" << coefSum
306 << " ) " << std::endl;
307 }
308 }
309
310 if (coefSum == 0.) {
311 oocoutE(&addPdf, Eval) << addPdf.ClassName() << "::updateCoefCache(" << addPdf.GetName()
312 << ") sum of coefficients is zero." << std::endl;
313 }
314
315 for (std::size_t i = 0; i < pdfList.size(); i++) {
316 coefCache[i] /= coefSum;
317 }
318}
#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:481
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
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:310
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
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.