Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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
18 RooExtendPdf is a wrapper around an existing PDF that adds a
19 parameteric extended likelihood term to the PDF, optionally divided by a
20 fractional 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]
24 where \f$ N \f$ is supplied as a RooAbsReal to RooExtendPdf.
25 The 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 
31 where \f$ x \f$ is the set of dependents involved in the selection region and \f$ y \f$
32 is the set of remaining dependents.
33 
34 \f$ \mathrm{cutRegion}[x] \f$ is a limited integration range that is contained in
35 the 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 "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 
51 using namespace std;
52 
53 ClassImp(RooExtendPdf);
54 ;
55 
56 
57 RooExtendPdf::RooExtendPdf() : _rangeName(0)
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` if `rangeName` remains empty. If `rangeName` is not empty,
65 /// `norm` will refer to this range, and expectedEvents will return the
66 /// total number of events over the full range of the observables.
67 /// \param[in] name Name of the pdf
68 /// \param[in] title Title of the pdf (for plotting)
69 /// \param[in] pdf The pdf to be extended
70 /// \param[in] norm Expected number of events
71 /// \param[in] rangeName If given, the number of events denoted by `norm` is interpreted as
72 /// the number of events in this range only
73 RooExtendPdf::RooExtendPdf(const char *name, const char *title, const RooAbsPdf& pdf,
74  const RooAbsReal& norm, const char* rangeName) :
75  RooAbsPdf(name,title),
76  _pdf("pdf","PDF",this,(RooAbsReal&)pdf),
77  _n("n","Normalization",this,(RooAbsReal&)norm),
78  _rangeName(RooNameReg::ptr(rangeName))
79 {
80 
81  // Copy various setting from pdf
82  setUnit(_pdf.arg().getUnit()) ;
83  setPlotLabel(_pdf.arg().getPlotLabel()) ;
84 }
85 
86 
87 
88 RooExtendPdf::RooExtendPdf(const RooExtendPdf& other, const char* name) :
89  RooAbsPdf(other,name),
90  _pdf("pdf",this,other._pdf),
91  _n("n",this,other._n),
92  _rangeName(other._rangeName)
93 {
94  // Copy constructor
95 }
96 
97 
98 RooExtendPdf::~RooExtendPdf()
99 {
100  // Destructor
101 
102 }
103 
104 
105 /// Return the number of expected events over the full range of all variables.
106 /// `norm`, the variable set as normalisation constant in the constructor,
107 /// will yield the number of events in the range set in the constructor. That is, the function returns
108 /// \f[
109 /// N = \mathrm{norm} \; \cdot \; \frac{\int_{(x_F,y_F)} \mathrm{pdf}(x,y) }{\int_{(x_C,y_F)} \mathrm{pdf}(x,y)}
110 /// \f]
111 /// Where \f$ x \f$ is the set of dependents with a restricted range (defined by `rangeName` in the constructor),
112 /// and \f$ y \f$ are the other dependents. \f$ x_C \f$ is the integration
113 /// of \f$ x \f$ over the restricted range, and \f$ x_F \f$ is the integration of
114 /// \f$ x \f$ over the full range. `norm` is the number of events given as parameter to the constructor.
115 ///
116 /// If the nested PDF can be extended, \f$ N \f$ is further scaled by its expected number of events.
117 Double_t RooExtendPdf::expectedEvents(const RooArgSet* nset) const
118 {
119  RooAbsPdf& pdf = (RooAbsPdf&)_pdf.arg() ;
120 
121  if (_rangeName && (!nset || nset->getSize()==0)) {
122  coutW(InputArguments) << "RooExtendPdf::expectedEvents(" << GetName() << ") WARNING: RooExtendPdf needs non-null normalization set to calculate fraction in range "
123  << _rangeName << ". Results may be nonsensical" << endl ;
124  }
125 
126  Double_t nExp = _n ;
127 
128  // Optionally multiply with fractional normalization
129  if (_rangeName) {
130 
131  double fracInt;
132  {
133  GlobalSelectComponentRAII globalSelComp(true);
134  fracInt = pdf.getNormObj(nset,nset,_rangeName)->getVal();
135  }
136 
137 
138  if ( fracInt == 0. || _n == 0.) {
139  coutW(Eval) << "RooExtendPdf(" << GetName() << ") WARNING: nExpected = " << _n << " / "
140  << fracInt << " for nset = " << (nset?*nset:RooArgSet()) << endl ;
141  }
142 
143  nExp /= fracInt ;
144 
145 
146  // cout << "RooExtendPdf::expectedEvents(" << GetName() << ") fracInt = " << fracInt << " _n = " << _n << " nExpect = " << nExp << endl ;
147 
148  }
149 
150  // Multiply with original Nexpected, if defined
151  if (pdf.canBeExtended()) nExp *= pdf.expectedEvents(nset) ;
152 
153  return nExp ;
154 }
155 
156 
157