Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 *
38 * A minimal implementation of a RooResolutionModel consists of a
39 * ```
40 * Int_t basisCode(const char* name)
41 * ```
42 * function indicating which basis functions this resolution model supports, and
43 * ```
44 * double evaluate(),
45 * ```
46 * which should implement the resolution model (optionally convoluted with one of the
47 * supported basis functions). RooResolutionModel objects can be used as regular
48 * PDFs (They inherit from RooAbsPdf), or as resolution model convoluted with
49 * a basis function. The implementation of evaluate() can identify the requested
50 * mode using basisCode(). If zero, the regular PDF value
51 * should be calculated. If non-zero, the model's value convoluted with the
52 * basis function identified by the code should be calculated.
53 *
54 * Optionally, analytical integrals can be advertised and implemented, in the
55 * same way as done for regular PDFS (see RooAbsPdf for further details).
56 * Also in getAnalyticalIntegral() / analyticalIntegral(), the implementation
57 * should use basisCode() to determine for which scenario the integral is
58 * requested.
59 *
60 * The choice of basis returned by basisCode() is guaranteed not to change
61 * during the lifetime of a RooResolutionModel object.
62 *
63 */
64
65#include "TClass.h"
66#include "TMath.h"
67#include "Riostream.h"
68#include "RooResolutionModel.h"
69#include "RooMsgService.h"
70
71using std::endl, std::ostream;
72
74
75
76////////////////////////////////////////////////////////////////////////////////
77/// Constructor with convolution variable 'x'.
78/// The convolution variable needs to be convertible to real values, and be able
79/// to give information about its range. This is supported by e.g. RooRealVar or RooLinearVar, which
80/// accepts offsetting and scaling an observable.
82 : RooAbsPdf(name, title), x("x", "Dependent or convolution variable", this, _x), _basisCode(0), _ownBasis(false)
83{
84
85}
86
87
88
89////////////////////////////////////////////////////////////////////////////////
90/// Copy constructor
91
93 : RooAbsPdf(other, name), x("x", this, other.x), _basisCode(other._basisCode), _ownBasis(false)
94{
95 if (other._basis) {
96 _basis = static_cast<RooFormulaVar*>(other._basis->Clone()) ;
97 _ownBasis = true ;
98 //_basis = other._basis ;
99 }
100
101 if (_basis) {
102 for (RooAbsArg * basisServer : _basis->servers()) {
103 addServer(*basisServer,true,false) ;
104 }
105 }
106}
107
108
109
110////////////////////////////////////////////////////////////////////////////////
111/// Destructor
112
114{
115 if (_ownBasis && _basis) {
116 delete _basis ;
117 }
118}
119
120
121
122////////////////////////////////////////////////////////////////////////////////
123/// Return identity formula pointer
124
126{
127 static RooFormulaVar identity("identity","1",RooArgSet(""));
128 return &identity;
129}
130
131
132
133////////////////////////////////////////////////////////////////////////////////
134/// Instantiate a clone of this resolution model representing a convolution with given
135/// basis function. The owners object name is incorporated in the clones name
136/// to avoid multiple convolution objects with the same name in complex PDF structures.
137///
138/// Note: The 'inBasis' formula expression must be a RooFormulaVar that encodes the formula
139/// in the title of the object and this expression must be an exact match against the
140/// implemented basis function strings (see derived class implementation of method basisCode()
141/// for those strings
142
144{
145 // Check that primary variable of basis functions is our convolution variable
146 if (inBasis->getParameter(0) != x.absArg()) {
147 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
148 << ") convolution parameter of basis function and PDF don't match" << endl
149 << "basis->findServer(0) = " << inBasis->findServer(0) << endl
150 << "x.absArg() = " << x.absArg() << endl ;
151 return nullptr ;
152 }
153
154 if (basisCode(inBasis->GetTitle())==0) {
155 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
156 << ") basis function '" << inBasis->GetTitle() << "' is not supported." << endl ;
157 return nullptr ;
158 }
159
160 TString newName(GetName()) ;
161 newName.Append("_conv_") ;
162 newName.Append(inBasis->GetName()) ;
163 newName.Append("_[") ;
164 newName.Append(owner->GetName()) ;
165 newName.Append("]") ;
166
167 RooResolutionModel* conv = static_cast<RooResolutionModel*>(clone(newName)) ;
168
169 TString newTitle(conv->GetTitle()) ;
170 newTitle.Append(" convoluted with basis function ") ;
171 newTitle.Append(inBasis->GetName()) ;
172 conv->SetTitle(newTitle.Data()) ;
173
174 conv->changeBasis(inBasis) ;
175
176 return conv ;
177}
178
179
180
181////////////////////////////////////////////////////////////////////////////////
182/// Change the basis function we convolute with.
183/// For one-time use by convolution() only.
184
186{
187 // Remove client-server link to old basis
188 if (_basis) {
189 for (RooAbsArg* basisServer : _basis->servers()) {
190 removeServer(*basisServer) ;
191 }
192
193 if (_ownBasis) {
194 delete _basis ;
195 }
196 }
197 _ownBasis = false ;
198
199 // Change basis pointer and update client-server link
200 _basis = inBasis ;
201 if (_basis) {
202 for (RooAbsArg* basisServer : _basis->servers()) {
203 addServer(*basisServer,true,false) ;
204 }
205 }
206
207 _basisCode = inBasis?basisCode(inBasis->GetTitle()):0 ;
208}
209
210
211
212////////////////////////////////////////////////////////////////////////////////
213/// Return the convolution variable of the selection basis function.
214/// This is, by definition, the first parameter of the basis function
215
217{
218 // Convolution variable is by definition first server of basis function
219 return *static_cast<RooRealVar const*>(*basis().servers().begin());
220}
221
222
223////////////////////////////////////////////////////////////////////////////////
224/// Modified version of RooAbsPdf::getValF(). If used as regular PDF,
225/// call RooAbsPdf::getValF(), otherwise return unnormalized value
226/// regardless of specified normalization set
227
228double RooResolutionModel::getValV(const RooArgSet* nset) const
229{
230 if (!_basis) return RooAbsPdf::getValV(nset) ;
231
232 // Return value of object. Calculated if dirty, otherwise cached value is returned.
233 if (isValueDirty()) {
234 _value = evaluate() ;
235
236 // WVE insert traceEval traceEval
237 if (_verboseDirty) cxcoutD(Tracing) << "RooResolutionModel(" << GetName() << ") value = " << _value << endl ;
238
241 }
242
243 return _value ;
244}
245
246
247
248////////////////////////////////////////////////////////////////////////////////
249/// Forward redirectServers call to our basis function, which is not connected to either resolution
250/// model or the physics model.
251
252bool RooResolutionModel::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
253{
254 if (!_basis) {
255 _norm = nullptr ;
256 return false ;
257 }
258
259 RooFormulaVar* newBasis = static_cast<RooFormulaVar*>(newServerList.find(_basis->GetName())) ;
260 if (newBasis) {
261
262 if (_ownBasis) {
263 delete _basis ;
264 }
265
266 _basis = newBasis ;
267 _ownBasis = false ;
268 }
269
270 _basis->redirectServers(newServerList,mustReplaceAll,nameChange) ;
271
272 return (mustReplaceAll && !newBasis) || RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
273}
274
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// Floating point error checking and tracing for given float value
279
280//bool RooResolutionModel::traceEvalHook(double value) const
281//{
282// // check for a math error or negative value
283// return TMath::IsNaN(value) ;
284//}
285
286
287
288////////////////////////////////////////////////////////////////////////////////
289/// Return the list of servers used by our normalization integral
290
292{
293 _norm->leafNodeServerList(&list) ;
294}
295
296
297
298////////////////////////////////////////////////////////////////////////////////
299/// Return the integral of this PDF over all elements of 'nset'.
300
301double RooResolutionModel::getNorm(const RooArgSet* nset) const
302{
303 if (!nset) {
304 return getVal() ;
305 }
306
307 syncNormalization(nset,false) ;
308 if (_verboseEval>1) cxcoutD(Tracing) << ClassName() << "::getNorm(" << GetName()
309 << "): norm(" << _norm << ") = " << _norm->getVal() << endl ;
310
311 double ret = _norm->getVal() ;
312 return ret ;
313}
314
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// Print info about this object to the specified stream. In addition to the info
319/// from RooAbsArg::printStream() we add:
320///
321/// Shape : value, units, plot range
322/// Verbose : default binning and print label
323
324void RooResolutionModel::printMultiline(ostream& os, Int_t content, bool verbose, TString indent) const
325{
326 RooAbsPdf::printMultiline(os,content,verbose,indent) ;
327
328 if(verbose) {
329 os << indent << "--- RooResolutionModel ---" << endl;
330 os << indent << "basis function = " ;
331 if (_basis) {
333 } else {
334 os << "<none>" << endl ;
335 }
336 }
337}
338
#define cxcoutD(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:382
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
@ kName
@ kTitle
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
R__DEPRECATED(6, 36, "Use getObservables().") RooFit R__DEPRECATED(6, 36, "Use getObservables().") RooFit R__DEPRECATED(6, 36, "Use getObservables().") RooFit const RooAbsArg &testArg const
Definition RooAbsArg.h:145
void removeServer(RooAbsArg &server, bool force=false)
Unregister another RooAbsArg as a server to us, ie, declare that we no longer depend on its value and...
static bool _verboseDirty
cache of the list of proxies. Avoids type casting.
Definition RooAbsArg.h:667
bool redirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool isRecursionStep=false)
Replace all direct servers of this object with the new servers in newServerList.
void clearValueDirty() const
Definition RooAbsArg.h:576
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:180
void addServer(RooAbsArg &server, bool valueProp=true, bool shapeProp=false, std::size_t refCount=1)
Register another RooAbsArg as a server to us, ie, declare that we depend on it.
bool isValueDirty() const
Definition RooAbsArg.h:393
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:91
void clearShapeDirty() const
Definition RooAbsArg.h:579
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=nullptr, bool recurseNonDerived=false) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition RooAbsArg.h:184
Abstract container object that can hold multiple RooAbsArg objects.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual bool syncNormalization(const RooArgSet *dset, bool adjustProxies=true) const
Verify that the normalization integral cached with this PDF is valid for given set of normalization o...
double getValV(const RooArgSet *set=nullptr) const override
Return current value, normalized by integrating over the observables in nset.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print multi line detailed information of this RooAbsPdf.
RooAbsReal * _norm
Definition RooAbsPdf.h:319
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
static Int_t _verboseEval
Definition RooAbsPdf.h:314
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
double _value
Cache for current value of object.
Definition RooAbsReal.h:536
virtual double evaluate() const =0
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsArg * absArg() const
Return pointer to contained argument.
Definition RooArgProxy.h:46
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooAbsArg * getParameter(const char *name) const
Return pointer to parameter with given name.
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,...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
~RooResolutionModel() override
Destructor.
TObject * clone(const char *newname) const override=0
double getValV(const RooArgSet *nset=nullptr) const override
Modified version of RooAbsPdf::getValF().
virtual void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
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.
double getNorm(const RooArgSet *nset=nullptr) const override
Return the integral of this PDF over all elements of 'nset'.
static RooFormulaVar * identity()
Return identity formula pointer.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override
Forward redirectServers call to our basis function, which is not connected to either resolution model...
const RooRealVar & basisConvVar() const
Return the convolution variable of the selection basis function.
bool _ownBasis
Flag indicating ownership of _basis.
virtual void normLeafServerList(RooArgSet &list) const
Floating point error checking and tracing for given float value.
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Print info about this object to the specified stream.
Int_t _basisCode
Identifier code for selected basis function.
RooResolutionModel()=default
RooFormulaVar * _basis
Basis function convolved with this resolution model.
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:213
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
TString & Append(const char *cs)
Definition TString.h:572
Double_t x[n]
Definition legend1.C:17