Logo ROOT   6.14/05
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 ///////////////////////////////////////////////////////////////////////////
18 // RooExtendPdf is a wrappper 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 //
22 // nExpected = N _or Expected = N / frac
23 //
24 // where N is supplied as a RooAbsReal to RooExtendPdf.
25 // The fractional term is defined as
26 // _ _ _ _ _
27 // Int(cutRegion[x]) pdf(x,y) dx dy
28 // frac = ---------------_-------_-_---_--_
29 // Int(normRegion[x]) pdf(x,y) dx dy
30 //
31 // _ _
32 // where x is the set of dependents involved in the selection region and y
33 // is the set of remaining dependents.
34 // _
35 // cutRegion[x] is an limited integration range that is contained in
36 // the nominal integration range normRegion[x[]
37 //
38 
39 #include "RooFit.h"
40 #include "Riostream.h"
41 
42 #include "RooExtendPdf.h"
43 #include "RooExtendPdf.h"
44 #include "RooArgList.h"
45 #include "RooRealVar.h"
46 #include "RooFormulaVar.h"
47 #include "RooNameReg.h"
48 #include "RooMsgService.h"
49 
50 
51 
52 using namespace std;
53 
55 ;
56 
57 
58 RooExtendPdf::RooExtendPdf() : _rangeName(0)
59 {
60  // Default constructor
61 }
62 
63 RooExtendPdf::RooExtendPdf(const char *name, const char *title, const RooAbsPdf& pdf,
64  const RooAbsReal& norm, const char* rangeName) :
65  RooAbsPdf(name,title),
66  _pdf("pdf","PDF",this,(RooAbsReal&)pdf),
67  _n("n","Normalization",this,(RooAbsReal&)norm),
68  _rangeName(RooNameReg::ptr(rangeName))
69 {
70  // Constructor. The ExtendedPdf behaves identical to the supplied input pdf,
71  // but adds an extended likelihood term. The expected number of events return
72  // is 'norm'. If a rangename is given, the number of events is interpreted as
73 # // the number of events in the given range
74 
75  // Copy various setting from pdf
76  setUnit(_pdf.arg().getUnit()) ;
78 }
79 
80 
81 
82 RooExtendPdf::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 
93 {
94  // Destructor
95 
96 }
97 
98 
99 
101 {
102  // Return the number of expected events, which is
103  //
104  // n / [ Int(xC,yF) pdf(x,y) / Int(xF,yF) pdf(x,y) ]
105  //
106  // Where x is the set of dependents with cuts defined
107  // and y are the other dependents. xC is the integration
108  // of x over the cut range, xF is the integration of
109  // x over the full range.
110 
111  RooAbsPdf& pdf = (RooAbsPdf&)_pdf.arg() ;
112 
113  if (_rangeName && (!nset || nset->getSize()==0)) {
114  coutW(InputArguments) << "RooExtendPdf::expectedEvents(" << GetName() << ") WARNING: RooExtendPdf needs non-null normalization set to calculate fraction in range "
115  << _rangeName << ". Results may be nonsensical" << endl ;
116  }
117 
118  Double_t nExp = _n ;
119 
120  // Optionally multiply with fractional normalization
121  if (_rangeName) {
122 
124  Double_t fracInt = pdf.getNormObj(nset,nset,_rangeName)->getVal() ;
126 
127 
128  if ( fracInt == 0. || _n == 0.) {
129  coutW(Eval) << "RooExtendPdf(" << GetName() << ") WARNING: nExpected = " << _n << " / "
130  << fracInt << " for nset = " << (nset?*nset:RooArgSet()) << endl ;
131  }
132 
133  nExp /= fracInt ;
134 
135 
136  // cout << "RooExtendPdf::expectedEvents(" << GetName() << ") fracInt = " << fracInt << " _n = " << _n << " nExpect = " << nExp << endl ;
137 
138  }
139 
140  // Multiply with original Nexpected, if defined
141  if (pdf.canBeExtended()) nExp *= pdf.expectedEvents(nset) ;
142 
143  return nExp ;
144 }
145 
146 
147 
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
static void globalSelectComp(Bool_t flag)
Global switch controlling the activation of the selectComp() functionality.
const Text_t * getUnit() const
Definition: RooAbsReal.h:83
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:383
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual ~RooExtendPdf()
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:33
friend class RooArgSet
Definition: RooAbsArg.h:471
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:2953
const TNamed * _rangeName
Definition: RooExtendPdf.h:53
Int_t getSize() const
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:373
Bool_t canBeExtended() const
Definition: RooAbsPdf.h:216
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
const Bool_t kFALSE
Definition: RtypesCore.h:88
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooRealProxy _n
Definition: RooExtendPdf.h:52
RooRealProxy _pdf
Definition: RooExtendPdf.h:51
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooNameReg is a registry for &#39;const char*&#39; name.
Definition: RooNameReg.h:23
void setUnit(const char *unit)
Definition: RooAbsReal.h:87
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
const Bool_t kTRUE
Definition: RtypesCore.h:87
char name[80]
Definition: TGX11.cxx:109