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