Logo ROOT  
Reference Guide
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 "RooMsgService.h"
46
47
48
49using namespace std;
50
52;
53
54
56{
57 // Default constructor
58}
59
60/// Constructor. The ExtendPdf behaves identical to the supplied input pdf,
61/// but adds an extended likelihood term. expectedEvents() will return
62/// `norm` if `rangeName` remains empty. If `rangeName` is not empty,
63/// `norm` will refer to this range, and expectedEvents will return the
64/// total number of events over the full range of the observables.
65/// \param[in] name Name of the pdf
66/// \param[in] title Title of the pdf (for plotting)
67/// \param[in] pdf The pdf to be extended
68/// \param[in] norm Expected number of events
69/// \param[in] rangeName If given, the number of events denoted by `norm` is interpreted as
70/// the number of events in this range only
71RooExtendPdf::RooExtendPdf(const char *name, const char *title, RooAbsPdf& pdf,
72 RooAbsReal& norm, const char* rangeName) :
73 RooAbsPdf(name,title),
74 _pdf("pdf", "PDF", this, pdf),
75 _n("n","Normalization",this,norm),
76 _rangeName(RooNameReg::ptr(rangeName))
77{
78
79 // Copy various setting from pdf
80 setUnit(_pdf->getUnit()) ;
82}
83
84
85
86RooExtendPdf::RooExtendPdf(const RooExtendPdf& other, const char* name) :
87 RooAbsPdf(other,name),
88 _pdf("pdf",this,other._pdf),
89 _n("n",this,other._n),
90 _rangeName(other._rangeName)
91{
92 // Copy constructor
93}
94
95
97{
98 // Destructor
99
100}
101
102
103/// Return the number of expected events over the full range of all variables.
104/// `norm`, the variable set as normalisation constant in the constructor,
105/// will yield the number of events in the range set in the constructor. That is, the function returns
106/// \f[
107/// N = \mathrm{norm} \; \cdot \; \frac{\int_{(x_F,y_F)} \mathrm{pdf}(x,y) }{\int_{(x_C,y_F)} \mathrm{pdf}(x,y)}
108/// \f]
109/// Where \f$ x \f$ is the set of dependents with a restricted range (defined by `rangeName` in the constructor),
110/// and \f$ y \f$ are the other dependents. \f$ x_C \f$ is the integration
111/// of \f$ x \f$ over the restricted range, and \f$ x_F \f$ is the integration of
112/// \f$ x \f$ over the full range. `norm` is the number of events given as parameter to the constructor.
113///
114/// If the nested PDF can be extended, \f$ N \f$ is further scaled by its expected number of events.
116{
117 const RooAbsPdf& pdf = *_pdf;
118
119 if (_rangeName && (!nset || nset->getSize()==0)) {
120 coutW(InputArguments) << "RooExtendPdf::expectedEvents(" << GetName() << ") WARNING: RooExtendPdf needs non-null normalization set to calculate fraction in range "
121 << _rangeName << ". Results may be nonsensical" << endl ;
122 }
123
124 Double_t nExp = _n ;
125
126 // Optionally multiply with fractional normalization
127 if (_rangeName) {
128
129 double fracInt;
130 {
131 GlobalSelectComponentRAII globalSelComp(true);
132 fracInt = pdf.getNormObj(nset,nset,_rangeName)->getVal();
133 }
134
135
136 if ( fracInt == 0. || _n == 0.) {
137 coutW(Eval) << "RooExtendPdf(" << GetName() << ") WARNING: nExpected = " << _n << " / "
138 << fracInt << " for nset = " << (nset?*nset:RooArgSet()) << endl ;
139 }
140
141 nExp /= fracInt ;
142
143
144 // cout << "RooExtendPdf::expectedEvents(" << GetName() << ") fracInt = " << fracInt << " _n = " << _n << " nExpect = " << nExp << endl ;
145
146 }
147
148 // Multiply with original Nexpected, if defined
149 if (pdf.canBeExtended()) nExp *= pdf.expectedEvents(nset) ;
150
151 return nExp ;
152}
153
154
155
#define coutW(a)
Definition: RooMsgService.h:32
#define ClassImp(name)
Definition: Rtypes.h:375
char name[80]
Definition: TGX11.cxx:110
friend class RooArgSet
Definition: RooAbsArg.h:651
Int_t getSize() const
Return the number of elements in the collection.
Bool_t canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition: RooAbsPdf.h:263
virtual const RooAbsReal * getNormObj(const RooArgSet *set, const RooArgSet *iset, const TNamed *rangeName=0) const
Return pointer to RooAbsReal object that implements calculation of integral over observables iset in ...
Definition: RooAbsPdf.cxx:476
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
Definition: RooAbsPdf.cxx:3231
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:63
void setUnit(const char *unit)
Definition: RooAbsReal.h:153
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:93
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:446
const Text_t * getUnit() const
Definition: RooAbsReal.h:149
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:456
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:57
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Definition: RooExtendPdf.h:22
RooTemplateProxy< RooAbsPdf > _pdf
Input p.d.f.
Definition: RooExtendPdf.h:50
const TNamed * _rangeName
Name of subset range.
Definition: RooExtendPdf.h:52
~RooExtendPdf() override
RooTemplateProxy< RooAbsReal > _n
Number of expected events.
Definition: RooExtendPdf.h:51
Double_t expectedEvents(const RooArgSet *nset) const override
Return the number of expected events over the full range of all variables.
RooNameReg is a registry for const char* names.
Definition: RooNameReg.h:25
const char * GetName() const override
Returns name of object.
Definition: TNamed.h:47
@ InputArguments
Definition: RooGlobalFunc.h:64