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 namespace std;
72
74
75
76////////////////////////////////////////////////////////////////////////////////
77/// Constructor with convolution variable 'x'.
78/// The convolution variable needs to be convertable 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.
81RooResolutionModel::RooResolutionModel(const char *name, const char *title, RooAbsRealLValue& _x) :
82 RooAbsPdf(name,title),
83 x("x","Dependent or convolution variable",this,_x),
84 _basisCode(0), _basis(0),
85 _ownBasis(false)
86{
87
88}
89
90
91
92////////////////////////////////////////////////////////////////////////////////
93/// Copy constructor
94
96 RooAbsPdf(other,name),
97 x("x",this,other.x),
98 _basisCode(other._basisCode), _basis(0),
99 _ownBasis(false)
100{
101 if (other._basis) {
102 _basis = (RooFormulaVar*) other._basis->Clone() ;
103 _ownBasis = true ;
104 //_basis = other._basis ;
105 }
106
107 if (_basis) {
108 for (RooAbsArg * basisServer : _basis->servers()) {
109 addServer(*basisServer,true,false) ;
110 }
111 }
112}
113
114
115
116////////////////////////////////////////////////////////////////////////////////
117/// Destructor
118
120{
121 if (_ownBasis && _basis) {
122 delete _basis ;
123 }
124}
125
126
127
128////////////////////////////////////////////////////////////////////////////////
129/// Return identity formula pointer
130
132{
133 static RooFormulaVar identity("identity","1",RooArgSet(""));
134 return &identity;
135}
136
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// Instantiate a clone of this resolution model representing a convolution with given
141/// basis function. The owners object name is incorporated in the clones name
142/// to avoid multiple convolution objects with the same name in complex PDF structures.
143///
144/// Note: The 'inBasis' formula expression must be a RooFormulaVar that encodes the formula
145/// in the title of the object and this expression must be an exact match against the
146/// implemented basis function strings (see derived class implementation of method basisCode()
147/// for those strings
148
150{
151 // Check that primary variable of basis functions is our convolution variable
152 if (inBasis->getParameter(0) != x.absArg()) {
153 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
154 << ") convolution parameter of basis function and PDF don't match" << endl
155 << "basis->findServer(0) = " << inBasis->findServer(0) << endl
156 << "x.absArg() = " << x.absArg() << endl ;
157 return 0 ;
158 }
159
160 if (basisCode(inBasis->GetTitle())==0) {
161 coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this
162 << ") basis function '" << inBasis->GetTitle() << "' is not supported." << endl ;
163 return 0 ;
164 }
165
166 TString newName(GetName()) ;
167 newName.Append("_conv_") ;
168 newName.Append(inBasis->GetName()) ;
169 newName.Append("_[") ;
170 newName.Append(owner->GetName()) ;
171 newName.Append("]") ;
172
173 RooResolutionModel* conv = (RooResolutionModel*) clone(newName) ;
174
175 TString newTitle(conv->GetTitle()) ;
176 newTitle.Append(" convoluted with basis function ") ;
177 newTitle.Append(inBasis->GetName()) ;
178 conv->SetTitle(newTitle.Data()) ;
179
180 conv->changeBasis(inBasis) ;
181
182 return conv ;
183}
184
185
186
187////////////////////////////////////////////////////////////////////////////////
188/// Change the basis function we convolute with.
189/// For one-time use by convolution() only.
190
192{
193 // Remove client-server link to old basis
194 if (_basis) {
195 for (RooAbsArg* basisServer : _basis->servers()) {
196 removeServer(*basisServer) ;
197 }
198
199 if (_ownBasis) {
200 delete _basis ;
201 }
202 }
203 _ownBasis = false ;
204
205 // Change basis pointer and update client-server link
206 _basis = inBasis ;
207 if (_basis) {
208 for (RooAbsArg* basisServer : _basis->servers()) {
209 addServer(*basisServer,true,false) ;
210 }
211 }
212
213 _basisCode = inBasis?basisCode(inBasis->GetTitle()):0 ;
214}
215
216
217
218////////////////////////////////////////////////////////////////////////////////
219/// Return the convolution variable of the selection basis function.
220/// This is, by definition, the first parameter of the basis function
221
223{
224 // Convolution variable is by definition first server of basis function
225 return *static_cast<RooRealVar const*>(*basis().servers().begin());
226}
227
228
229////////////////////////////////////////////////////////////////////////////////
230/// Modified version of RooAbsPdf::getValF(). If used as regular PDF,
231/// call RooAbsPdf::getValF(), otherwise return unnormalized value
232/// regardless of specified normalization set
233
234double RooResolutionModel::getValV(const RooArgSet* nset) const
235{
236 if (!_basis) return RooAbsPdf::getValV(nset) ;
237
238 // Return value of object. Calculated if dirty, otherwise cached value is returned.
239 if (isValueDirty()) {
240 _value = evaluate() ;
241
242 // WVE insert traceEval traceEval
243 if (_verboseDirty) cxcoutD(Tracing) << "RooResolutionModel(" << GetName() << ") value = " << _value << endl ;
244
247 }
248
249 return _value ;
250}
251
252
253
254////////////////////////////////////////////////////////////////////////////////
255/// Forward redirectServers call to our basis function, which is not connected to either resolution
256/// model or the physics model.
257
258bool RooResolutionModel::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
259{
260 if (!_basis) {
261 _norm = 0 ;
262 return false ;
263 }
264
265 RooFormulaVar* newBasis = (RooFormulaVar*) newServerList.find(_basis->GetName()) ;
266 if (newBasis) {
267
268 if (_ownBasis) {
269 delete _basis ;
270 }
271
272 _basis = newBasis ;
273 _ownBasis = false ;
274 }
275
276 _basis->redirectServers(newServerList,mustReplaceAll,nameChange) ;
277
278 return (mustReplaceAll && !newBasis) || RooAbsPdf::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
279}
280
281
282
283////////////////////////////////////////////////////////////////////////////////
284/// Floating point error checking and tracing for given float value
285
286//bool RooResolutionModel::traceEvalHook(double value) const
287//{
288// // check for a math error or negative value
289// return TMath::IsNaN(value) ;
290//}
291
292
293
294////////////////////////////////////////////////////////////////////////////////
295/// Return the list of servers used by our normalization integral
296
298{
299 _norm->leafNodeServerList(&list) ;
300}
301
302
303
304////////////////////////////////////////////////////////////////////////////////
305/// Return the integral of this PDF over all elements of 'nset'.
306
307double RooResolutionModel::getNorm(const RooArgSet* nset) const
308{
309 if (!nset) {
310 return getVal() ;
311 }
312
313 syncNormalization(nset,false) ;
314 if (_verboseEval>1) cxcoutD(Tracing) << ClassName() << "::getNorm(" << GetName()
315 << "): norm(" << _norm << ") = " << _norm->getVal() << endl ;
316
317 double ret = _norm->getVal() ;
318 return ret ;
319}
320
321
322
323////////////////////////////////////////////////////////////////////////////////
324/// Print info about this object to the specified stream. In addition to the info
325/// from RooAbsArg::printStream() we add:
326///
327/// Shape : value, units, plot range
328/// Verbose : default binning and print label
329
330void RooResolutionModel::printMultiline(ostream& os, Int_t content, bool verbose, TString indent) const
331{
332 RooAbsPdf::printMultiline(os,content,verbose,indent) ;
333
334 if(verbose) {
335 os << indent << "--- RooResolutionModel ---" << endl;
336 os << indent << "basis function = " ;
337 if (_basis) {
339 } else {
340 os << "<none>" << endl ;
341 }
342 }
343}
344
#define cxcoutD(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
@ kName
@ kTitle
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
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:690
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:599
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:204
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:418
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:86
void clearShapeDirty() const
Definition RooAbsArg.h:602
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:208
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 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:376
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
The cache manager.
static Int_t _verboseEval
Definition RooAbsPdf.h:371
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
double _value
Cache for current value of object.
Definition RooAbsReal.h:509
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:47
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
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,...
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
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.
RooFormulaVar * _basis
Basis function convolved with this resolution model.
const RooFormulaVar & basis() const
RooTemplateProxy< RooAbsRealLValue > x
Dependent/convolution variable.
Container_t::const_iterator begin() const
Iterator over contained objects.
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:207
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
TString & Append(const char *cs)
Definition TString.h:576
Double_t x[n]
Definition legend1.C:17