Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooExtendPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/** \class RooExtendPdf
18RooExtendPdf is a wrapper around an existing PDF that adds a
19parameteric extended likelihood term to the PDF, optionally divided by a
20fractional term from a partial normalization of the PDF:
21\f[
22 n_\mathrm{Expected} = N \quad \text{or} \quad n_\mathrm{Expected} = N / \mathrm{frac}
23\f]
24where \f$ N \f$ is supplied as a RooAbsReal to RooExtendPdf.
25The fractional term is defined as
26\f[
27 \mathrm{frac} = \frac{\int_\mathrm{cutRegion[x]} \mathrm{pdf}(x,y) \; \mathrm{d}x \mathrm{d}y}{
28 \int_\mathrm{normRegion[x]} \mathrm{pdf}(x,y) \; \mathrm{d}x \mathrm{d}y}
29\f]
30
31where \f$ x \f$ is the set of dependents involved in the selection region and \f$ y \f$
32is the set of remaining dependents.
33
34\f$ \mathrm{cutRegion}[x] \f$ is a limited integration range that is contained in
35the nominal integration range \f$ \mathrm{normRegion}[x] \f$.
36*/
37
38#include "Riostream.h"
39
40#include "RooExtendPdf.h"
41#include "RooArgList.h"
42#include "RooRealVar.h"
43#include "RooFormulaVar.h"
44#include "RooNameReg.h"
45#include "RooConstVar.h"
46#include "RooProduct.h"
47#include "RooRatio.h"
48#include "RooMsgService.h"
49
50
51
52using namespace std;
53
55
56/// Constructor. The ExtendPdf behaves identical to the supplied input pdf,
57/// but adds an extended likelihood term. expectedEvents() will return
58/// `norm` if `rangeName` remains empty. If `rangeName` is not empty,
59/// `norm` will refer to this range, and expectedEvents will return the
60/// total number of events over the full range of the observables.
61/// \param[in] name Name of the pdf
62/// \param[in] title Title of the pdf (for plotting)
63/// \param[in] pdf The pdf to be extended
64/// \param[in] norm Expected number of events
65/// \param[in] rangeName If given, the number of events denoted by `norm` is interpreted as
66/// the number of events in this range only
67RooExtendPdf::RooExtendPdf(const char *name, const char *title, RooAbsPdf& pdf,
68 RooAbsReal& norm, const char* rangeName) :
69 RooAbsPdf(name,title),
70 _pdf("pdf", "PDF", this, pdf),
71 _n("n","Normalization",this,norm),
72 _rangeName(RooNameReg::ptr(rangeName))
73{
74
75 // Copy various setting from pdf
76 setUnit(_pdf->getUnit()) ;
78}
79
80
81
82RooExtendPdf::RooExtendPdf(const RooExtendPdf& other, const char* name) :
83 RooAbsPdf(other,name),
84 _pdf("pdf",this,other._pdf),
85 _n("n",this,other._n),
86 _rangeName(other._rangeName)
87{
88 // Copy constructor
89}
90
91
92/// Return the number of expected events over the full range of all variables.
93/// `norm`, the variable set as normalisation constant in the constructor,
94/// will yield the number of events in the range set in the constructor. That is, the function returns
95/// \f[
96/// N = \mathrm{norm} \; \cdot \; \frac{\int_{(x_F,y_F)} \mathrm{pdf}(x,y) }{\int_{(x_C,y_F)} \mathrm{pdf}(x,y)}
97/// \f]
98/// Where \f$ x \f$ is the set of dependents with a restricted range (defined by `rangeName` in the constructor),
99/// and \f$ y \f$ are the other dependents. \f$ x_C \f$ is the integration
100/// of \f$ x \f$ over the restricted range, and \f$ x_F \f$ is the integration of
101/// \f$ x \f$ over the full range. `norm` is the number of events given as parameter to the constructor.
102///
103/// If the nested PDF can be extended, \f$ N \f$ is further scaled by its expected number of events.
105{
106 const RooAbsPdf& pdf = *_pdf;
107
108 if (_rangeName && (!nset || nset->empty())) {
109 coutW(InputArguments) << "RooExtendPdf::expectedEvents(" << GetName() << ") WARNING: RooExtendPdf needs non-null normalization set to calculate fraction in range "
110 << _rangeName << ". Results may be nonsensical" << endl ;
111 }
112
113 double nExp = _n ;
114
115 // Optionally multiply with fractional normalization
116 if (_rangeName) {
117
118 double fracInt;
119 {
120 GlobalSelectComponentRAII globalSelComp(true);
121 fracInt = pdf.getNormObj(nset,nset,_rangeName)->getVal();
122 }
123
124
125 if ( fracInt == 0. || _n == 0.) {
126 coutW(Eval) << "RooExtendPdf(" << GetName() << ") WARNING: nExpected = " << _n << " / "
127 << fracInt << " for nset = " << (nset?*nset:RooArgSet()) << endl ;
128 }
129
130 nExp /= fracInt ;
131 }
132
133 // Multiply with original Nexpected, if defined
134 if (pdf.canBeExtended()) nExp *= pdf.expectedEvents(nset) ;
135
136 return nExp ;
137}
138
139
140std::unique_ptr<RooAbsReal> RooExtendPdf::createExpectedEventsFunc(const RooArgSet *nset) const
141{
142 const RooAbsPdf& pdf = *_pdf;
143
144 RooArgList prodList;
145 prodList.add(*_n);
146
147 // Optionally multiply with fractional normalization
148 std::unique_ptr<RooAbsReal> rangeFactor;
149 if (_rangeName) {
150 std::unique_ptr<RooAbsReal> fracInteg{pdf.createIntegral(*nset, *nset, RooNameReg::str(_rangeName))};
151 // Create one over integral term
152 auto rangeFactorName = std::string("one_over_") + fracInteg->GetName();
153 rangeFactor = std::make_unique<RooRatio>(rangeFactorName.c_str(), rangeFactorName.c_str(), RooFit::RooConst(1.0), *fracInteg);
154 rangeFactor->addOwnedComponents(std::move(fracInteg));
155 prodList.add(*rangeFactor);
156 }
157
158 // Multiply with original Nexpected, if defined
159 std::unique_ptr<RooAbsReal> pdfExpectedEvents;
160 if (pdf.canBeExtended()) {
161 pdfExpectedEvents = pdf.createExpectedEventsFunc(nset);
162 prodList.add(*pdfExpectedEvents);
163 }
164
165 auto name = std::string(GetName()) + "_expectedEvents";
166 auto out = std::make_unique<RooProduct>(name.c_str(), name.c_str(), prodList);
167 if(rangeFactor) {
168 out->addOwnedComponents(std::move(rangeFactor));
169 }
170 if(pdfExpectedEvents) {
171 out->addOwnedComponents(std::move(pdfExpectedEvents));
172 }
173 return out;
174}
#define coutW(a)
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:278
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,...
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=nullptr) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables std::liste...
void setUnit(const char *unit)
Definition RooAbsReal.h:151
const char * getPlotLabel() const
Get the label associated with the variable.
const Text_t * getUnit() const
Definition RooAbsReal.h:147
void setPlotLabel(const char *label)
Set the label associated with this variable.
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:55
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
std::unique_ptr< RooAbsReal > createExpectedEventsFunc(const RooArgSet *nset) const override
Returns an object that represents the expected number of events for a given normalization set,...
RooTemplateProxy< RooAbsPdf > _pdf
Input p.d.f.
double expectedEvents(const RooArgSet *nset) const override
Return the number of expected events over the full range of all variables.
const TNamed * _rangeName
Name of subset range.
RooExtendPdf()=default
RooTemplateProxy< RooAbsReal > _n
Number of expected events.
RooNameReg is a registry for const char* names.
Definition RooNameReg.h:25
static const char * str(const TNamed *ptr)
Return C++ string corresponding to given TNamed pointer.
Definition RooNameReg.h:37
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
RooConstVar & RooConst(double val)