Logo ROOT   6.16/01
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 x is the set of dependents involved in the selection region and y
32is the set of remaining dependents.
33
34cutRegion[x] is a limited integration range that is contained in
35the nominal integration range normRegion[x]
36*/
37
38#include "RooFit.h"
39#include "Riostream.h"
40
41#include "RooExtendPdf.h"
42#include "RooExtendPdf.h"
43#include "RooArgList.h"
44#include "RooRealVar.h"
45#include "RooFormulaVar.h"
46#include "RooNameReg.h"
47#include "RooMsgService.h"
48
49
50
51using namespace std;
52
54;
55
56
58{
59 // Default constructor
60}
61
62/// Constructor. The ExtendPdf behaves identical to the supplied input pdf,
63/// but adds an extended likelihood term. expectedEvents() will return
64/// `norm` (or `norm` scaled to the full range of the variable if `rangeName` is used)
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 is interpreted as
70/// the number of events in this range only
71RooExtendPdf::RooExtendPdf(const char *name, const char *title, const RooAbsPdf& pdf,
72 const RooAbsReal& norm, const char* rangeName) :
73 RooAbsPdf(name,title),
74 _pdf("pdf","PDF",this,(RooAbsReal&)pdf),
75 _n("n","Normalization",this,(RooAbsReal&)norm),
76 _rangeName(RooNameReg::ptr(rangeName))
77{
78
79 // Copy various setting from pdf
80 setUnit(_pdf.arg().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. That is
104 /// \f[
105 /// N = \mathrm{norm} \; / \; \frac{\int_{(x_C,y_F)} \mathrm{pdf}(x,y)}{\int_{(x_F,y_F)} \mathrm{pdf}(x,y) }
106 /// \f]
107 /// Where \f$ x \f$ is the set of dependents with cuts defined (`rangeName` in the constructor)
108 /// and \f$ y \f$ are the other dependents. \f$ x_C \f$ is the integration
109 /// of \f$ x \f$ over the cut range, \f$ x_F \f$ is the integration of
110 /// \f$ x \f$ over the full range, `norm` is the number of events from the constructor.
111 /// If the nested PDF can be extended, \f$ N \f$ is scaled by its expected number of events.
113{
114 RooAbsPdf& pdf = (RooAbsPdf&)_pdf.arg() ;
115
116 if (_rangeName && (!nset || nset->getSize()==0)) {
117 coutW(InputArguments) << "RooExtendPdf::expectedEvents(" << GetName() << ") WARNING: RooExtendPdf needs non-null normalization set to calculate fraction in range "
118 << _rangeName << ". Results may be nonsensical" << endl ;
119 }
120
121 Double_t nExp = _n ;
122
123 // Optionally multiply with fractional normalization
124 if (_rangeName) {
125
127 Double_t fracInt = pdf.getNormObj(nset,nset,_rangeName)->getVal() ;
129
130
131 if ( fracInt == 0. || _n == 0.) {
132 coutW(Eval) << "RooExtendPdf(" << GetName() << ") WARNING: nExpected = " << _n << " / "
133 << fracInt << " for nset = " << (nset?*nset:RooArgSet()) << endl ;
134 }
135
136 nExp /= fracInt ;
137
138
139 // cout << "RooExtendPdf::expectedEvents(" << GetName() << ") fracInt = " << fracInt << " _n = " << _n << " nExpect = " << nExp << endl ;
140
141 }
142
143 // Multiply with original Nexpected, if defined
144 if (pdf.canBeExtended()) nExp *= pdf.expectedEvents(nset) ;
145
146 return nExp ;
147}
148
149
150
#define coutW(a)
Definition: RooMsgService.h:33
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
friend class RooArgSet
Definition: RooAbsArg.h:471
Int_t getSize() const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:230
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:401
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
Definition: RooAbsPdf.cxx:2855
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
void setUnit(const char *unit)
Definition: RooAbsReal.h:88
static void globalSelectComp(Bool_t flag)
Global switch controlling the activation of the selectComp() functionality.
Double_t getVal(const RooArgSet *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:374
const Text_t * getUnit() const
Definition: RooAbsReal.h:84
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:384
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Definition: RooExtendPdf.h:22
RooRealProxy _pdf
Definition: RooExtendPdf.h:52
const TNamed * _rangeName
Definition: RooExtendPdf.h:54
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return the number of expected events.
RooRealProxy _n
Definition: RooExtendPdf.h:53
virtual ~RooExtendPdf()
RooNameReg is a registry for 'const char*' name.
Definition: RooNameReg.h:23
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
@ InputArguments
Definition: RooGlobalFunc.h:58
STL namespace.