Logo ROOT  
Reference Guide
RooResolutionModel.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/**
19 * \class RooResolutionModel
20 * RooResolutionModel is the base class for PDFs that represent a
21 * resolution model that can be convoluted with a physics model of the form
22 * \f[
23 * \mathrm{Phys}(x,a,b) = \sum_k \mathrm{coef}_k(a) * \mathrm{basis}_k(x,b)
24 * \f]
25 * where basis_k are a limited number of functions in terms of the variable
26 * to be convoluted and coef_k are coefficients independent of the convolution
27 * variable.
28 *
29 * Classes derived from RooResolutionModel implement
30 * \f[
31 * R_k(x,\bar{b},\bar{c}) = \int \mathrm{basis}_k(x',\bar{b}) * \mathrm{resModel}(x-x',\bar{c}) \; \mathrm{d} x',
32 * \f]
33 * which RooAbsAnaConvPdf uses to construct the pdf for [ Phys (x) R ] :
34 * \f[
35 * \mathrm{PDF}(x,\bar a, \bar b, \bar c) = \sum_k \mathrm{coef}_k(\bar a) * R_k(x, \bar b, \bar c)
36 * \f]
37 * A minimal implementation of a RooResolutionModel consists of a
38 * ```
39 * Int_t basisCode(const char* name)
40 * ```
41 * function indication which basis functions this resolution model supports, and
42 * ```
43 * Double_t evaluate()
44 * ```
45 * Implementing the resolution model, optionally convoluted with one of the
46 * supported basis functions. RooResolutionModel objects can be used as regular
47 * PDFs (They inherit from RooAbsPdf), or as resolution model convoluted with
48 * a basis function. The implementation of evaluate() can identify the requested
49 * from of use from the basisCode() function. If zero, the regular PDF value
50 * should be calculated. If non-zero, the models value convoluted with the
51 * basis function identified by the code should be calculated.
52 *
53 * Optionally, analytical integrals can be advertised and implemented, in the
54 * same way as done for regular PDFS (see RooAbsPdf for further details).
55 * Also in getAnalyticalIntegral()/analyticalIntegral() the implementation
56 * should use basisCode() to determine for which scenario the integral is
57 * requested.
58 *
59 * The choice of basis returned by basisCode() is guaranteed not to change
60 * of the lifetime of a RooResolutionModel object.
61 *
62 */
63
64#include "RooFit.h"
65
66#include "TClass.h"
67#include "TMath.h"
68#include "Riostream.h"
69#include "RooResolutionModel.h"
70#include "RooMsgService.h"
71
72using namespace std;
73
75
76
77////////////////////////////////////////////////////////////////////////////////
78/// Constructor with convolution variable 'x'
79
80RooResolutionModel::RooResolutionModel(const char *name, const char *title, RooRealVar& _x) :
81 RooAbsPdf(name,title),
82 x("x","Dependent or convolution variable",this,_x),
83 _basisCode(0), _basis(0),
84 _ownBasis(kFALSE)
85{
86
87}
88
89
90
91////////////////////////////////////////////////////////////////////////////////
92/// Copy constructor
93
95 RooAbsPdf(other,name),
96 x("x",this,other.x),
97 _basisCode(other._basisCode), _basis(0),
98 _ownBasis(kFALSE)
99{
100 if (other._basis) {
101 _basis = (RooFormulaVar*) other._basis->Clone() ;
102 _ownBasis = kTRUE ;
103 //_basis = other._basis ;
104 }
105
106 if (_basis) {
107 TIterator* bsIter = _basis->serverIterator() ;
108 RooAbsArg* basisServer ;
109 while((basisServer = (RooAbsArg*)bsIter->Next())) {
110 addServer(*basisServer,kTRUE,kFALSE) ;
111 }
112 delete bsIter ;
113 }
114}
115
116
117
118////////////////////////////////////////////////////////////////////////////////
119/// Destructor
120
122{
123 if (_ownBasis && _basis) {
124 delete _basis ;
125 }
126}
127
128
129
130////////////////////////////////////////////////////////////////////////////////
131/// Return identity formula pointer
132
134{
135 static RooFormulaVar identity("identity","1",RooArgSet(""));
136 return &identity;
137}
138
139
140
141////////////////////////////////////////////////////////////////////////////////
142/// Instantiate a clone of this resolution model representing a convolution with given
143/// basis function. The owners object name is incorporated in the clones name
144/// to avoid multiple convolution objects with the same name in complex PDF structures.
145///
146/// Note: The 'inBasis' formula expression must be a RooFormulaVar that encodes the formula
147/// in the title of the object and this expression must be an exact match against the
148/// implemented basis function strings (see derived class implementation of method basisCode()
149/// for those strings
150
152{
153 // Check that primary variable of basis functions is our convolution variable
154 if (inBasis->getParameter(0) != x.absArg()) {
155 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
156 << ") convolution parameter of basis function and PDF don't match" << endl
157 << "basis->findServer(0) = " << inBasis->findServer(0) << endl
158 << "x.absArg() = " << x.absArg() << endl ;
159 return 0 ;
160 }
161
162 if (basisCode(inBasis->GetTitle())==0) {
163 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
164 << ") basis function '" << inBasis->GetTitle() << "' is not supported." << endl ;
165 return 0 ;
166 }
167
168 TString newName(GetName()) ;
169 newName.Append("_conv_") ;
170 newName.Append(inBasis->GetName()) ;
171 newName.Append("_[") ;
172 newName.Append(owner->GetName()) ;
173 newName.Append("]") ;
174
175 RooResolutionModel* conv = (RooResolutionModel*) clone(newName) ;
176
177 TString newTitle(conv->GetTitle()) ;
178 newTitle.Append(" convoluted with basis function ") ;
179 newTitle.Append(inBasis->GetName()) ;
180 conv->SetTitle(newTitle.Data()) ;
181
182 conv->changeBasis(inBasis) ;
183
184 return conv ;
185}
186
187
188
189////////////////////////////////////////////////////////////////////////////////
190/// Change the basis function we convolute with.
191/// For one-time use by convolution() only.
192
194{
195 // Remove client-server link to old basis
196 if (_basis) {
197 TIterator* bsIter = _basis->serverIterator() ;
198 RooAbsArg* basisServer ;
199 while((basisServer = (RooAbsArg*)bsIter->Next())) {
200 removeServer(*basisServer) ;
201 }
202 delete bsIter ;
203
204 if (_ownBasis) {
205 delete _basis ;
206 }
207 }
208 _ownBasis = kFALSE ;
209
210 // Change basis pointer and update client-server link
211 _basis = inBasis ;
212 if (_basis) {
213 TIterator* bsIter = _basis->serverIterator() ;
214 RooAbsArg* basisServer ;
215 while((basisServer = (RooAbsArg*)bsIter->Next())) {
216 addServer(*basisServer,kTRUE,kFALSE) ;
217 }
218 delete bsIter ;
219 }
220
221 _basisCode = inBasis?basisCode(inBasis->GetTitle()):0 ;
222}
223
224
225
226////////////////////////////////////////////////////////////////////////////////
227/// Return the convolution variable of the selection basis function.
228/// This is, by definition, the first parameter of the basis function
229
231{
232 // Convolution variable is by definition first server of basis function
233 TIterator* sIter = basis().serverIterator() ;
234 RooRealVar* var = (RooRealVar*) sIter->Next() ;
235 delete sIter ;
236
237 return *var ;
238}
239
240
241
242////////////////////////////////////////////////////////////////////////////////
243/// Return the convolution variable of the resolution model
244
246{
247 return (RooRealVar&) x.arg() ;
248}
249
250
251
252////////////////////////////////////////////////////////////////////////////////
253/// Modified version of RooAbsPdf::getValF(). If used as regular PDF,
254/// call RooAbsPdf::getValF(), otherwise return unnormalized value
255/// regardless of specified normalization set
256
258{
259 if (!_basis) return RooAbsPdf::getValV(nset) ;
260
261 // Return value of object. Calculated if dirty, otherwise cached value is returned.
262 if (isValueDirty()) {
263 _value = evaluate() ;
264
265 // WVE insert traceEval traceEval
266 if (_verboseDirty) cxcoutD(Tracing) << "RooResolutionModel(" << GetName() << ") value = " << _value << endl ;
267
268 clearValueDirty() ;
269 clearShapeDirty() ;
270 }
271
272 return _value ;
273}
274
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// Forward redirectServers call to our basis function, which is not connected to either resolution
279/// model or the physics model.
280
281Bool_t RooResolutionModel::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t /*isRecursive*/)
282{
283 if (!_basis) {
284 _norm = 0 ;
285 return kFALSE ;
286 }
287
288 RooFormulaVar* newBasis = (RooFormulaVar*) newServerList.find(_basis->GetName()) ;
289 if (newBasis) {
290
291 if (_ownBasis) {
292 delete _basis ;
293 }
294
295 _basis = newBasis ;
296 _ownBasis = kFALSE ;
297 }
298
299 _basis->redirectServers(newServerList,mustReplaceAll,nameChange) ;
300
301 return (mustReplaceAll && !newBasis) ;
302}
303
304
305
306////////////////////////////////////////////////////////////////////////////////
307/// Floating point error checking and tracing for given float value
308
309//Bool_t RooResolutionModel::traceEvalHook(Double_t value) const
310//{
311// // check for a math error or negative value
312// return TMath::IsNaN(value) ;
313//}
314
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// Return the list of servers used by our normalization integral
319
321{
322 _norm->leafNodeServerList(&list) ;
323}
324
325
326
327////////////////////////////////////////////////////////////////////////////////
328/// Return the integral of this PDF over all elements of 'nset'.
329
331{
332 if (!nset) {
333 return getVal() ;
334 }
335
337 if (_verboseEval>1) cxcoutD(Tracing) << IsA()->GetName() << "::getNorm(" << GetName()
338 << "): norm(" << _norm << ") = " << _norm->getVal() << endl ;
339
340 Double_t ret = _norm->getVal() ;
341 return ret ;
342}
343
344
345
346////////////////////////////////////////////////////////////////////////////////
347/// Print info about this object to the specified stream. In addition to the info
348/// from RooAbsArg::printStream() we add:
349///
350/// Shape : value, units, plot range
351/// Verbose : default binning and print label
352
354{
356
357 if(verbose) {
358 os << indent << "--- RooResolutionModel ---" << endl;
359 os << indent << "basis function = " ;
360 if (_basis) {
362 } else {
363 os << "<none>" << endl ;
364 }
365 }
366}
367
#define cxcoutD(a)
Definition: RooMsgService.h:82
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition: TGX11.cxx:109
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:71
Bool_t redirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t isRecursionStep=kFALSE)
Substitute our servers with those listed in newSet.
Definition: RooAbsArg.cxx:913
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
Definition: RooAbsArg.cxx:479
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:83
friend class RooArgSet
Definition: RooAbsArg.h:551
void addServer(RooAbsArg &server, Bool_t valueProp=kTRUE, Bool_t shapeProp=kFALSE, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
Definition: RooAbsArg.cxx:358
Bool_t isValueDirty() const
Definition: RooAbsArg.h:390
void clearValueDirty() const
Definition: RooAbsArg.h:526
static Bool_t _verboseDirty
Definition: RooAbsArg.h:619
void removeServer(RooAbsArg &server, Bool_t force=kFALSE)
Unregister another RooAbsArg as a server to us, ie, declare that we no longer depend on its value and...
Definition: RooAbsArg.cxx:409
void clearShapeDirty() const
Definition: RooAbsArg.h:529
TIterator * serverIterator() const R__SUGGEST_ALTERNATIVE("Use servers() and begin()
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition: RooAbsArg.h:167
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print multi line detailed information of this RooAbsPdf.
Definition: RooAbsPdf.cxx:2006
RooAbsReal * _norm
Definition: RooAbsPdf.h:321
virtual Bool_t syncNormalization(const RooArgSet *dset, Bool_t adjustProxies=kTRUE) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
Definition: RooAbsPdf.cxx:530
virtual Double_t getValV(const RooArgSet *set=0) const
Return current value, normalized by integrating over the observables in nset.
Definition: RooAbsPdf.cxx:281
static Int_t _verboseEval
Definition: RooAbsPdf.h:315
virtual Double_t evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Double_t _value
Definition: RooAbsReal.h:446
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
RooAbsArg * absArg() const
Definition: RooArgProxy.h:37
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition: RooFormulaVar.h:29
RooAbsArg * getParameter(const char *name) const
Definition: RooFormulaVar.h:40
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
virtual TObject * clone(const char *newname) const =0
virtual void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
Double_t getValV(const RooArgSet *nset=0) const
Modified version of RooAbsPdf::getValF().
virtual Int_t basisCode(const char *name) const =0
virtual RooResolutionModel * convolution(RooFormulaVar *basis, RooAbsArg *owner) const
Instantiate a clone of this resolution model representing a convolution with given basis function.
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Print info about this object to the specified stream.
static RooFormulaVar * identity()
Return identity formula pointer.
const RooRealVar & basisConvVar() const
Return the convolution variable of the selection basis function.
virtual Bool_t redirectServersHook(const RooAbsCollection &newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
Forward redirectServers call to our basis function, which is not connected to either resolution model...
virtual ~RooResolutionModel()
Destructor.
RooRealVar & convVar() const
Return the convolution variable of the resolution model.
virtual void normLeafServerList(RooArgSet &list) const
Floating point error checking and tracing for given float value.
Double_t getNorm(const RooArgSet *nset=0) const
Return the integral of this PDF over all elements of 'nset'.
RooFormulaVar * _basis
const RooFormulaVar & basis() const
const T & arg() const
Return reference to object held in proxy.
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
TString & Append(const char *cs)
Definition: TString.h:559
Double_t x[n]
Definition: legend1.C:17
@ InputArguments
Definition: RooGlobalFunc.h:68