Logo ROOT   6.18/05
Reference Guide
RooNumConvolution.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $Id: RooNumConvolution.h,v 1.4 2007/05/11 09:11:30 verkerke Exp $
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#ifndef ROO_NUM_CONVOLUTION
17#define ROO_NUM_CONVOLUTION
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21#include "RooSetProxy.h"
22#include "RooListProxy.h"
23#include "RooNumIntConfig.h"
24
26class RooAbsIntegrator ;
27class TH2 ;
28
30public:
31
33
34 RooNumConvolution(const char *name, const char *title,
35 RooRealVar& convVar, RooAbsReal& pdf, RooAbsReal& resmodel, const RooNumConvolution* proto=0) ;
36
37 RooNumConvolution(const RooNumConvolution& other, const char* name=0) ;
38
39 virtual TObject* clone(const char* newname) const { return new RooNumConvolution(*this,newname) ; }
40 virtual ~RooNumConvolution() ;
41
42 Double_t evaluate() const ;
43
46
48 void setConvolutionWindow(RooAbsReal& centerParam, RooAbsReal& widthParam, Double_t widthScaleFactor=1) ;
49
50 void setCallWarning(Int_t threshold=2000) ;
51 void setCallProfiling(Bool_t flag, Int_t nbinX = 40, Int_t nbinCall = 40, Int_t nCallHigh=1000) ;
52 const TH2* profileData() const { return _doProf ? _callHist : 0 ; }
53
54 // Access components
55 RooRealVar& var() const { return (RooRealVar&) _origVar.arg() ; }
56 RooAbsReal& pdf() const { return (RooAbsReal&) _origPdf.arg() ; }
57 RooAbsReal& model() const { return (RooAbsReal&) _origModel.arg() ; }
58
59protected:
60
61 friend class RooNumConvPdf ;
62
63 mutable Bool_t _init ;
64 void initialize() const ;
65 Bool_t redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive) ;
66
67 virtual void printCompactTreeHook(std::ostream& os, const char* indent="") ;
68
69 RooNumIntConfig _convIntConfig ; // Configuration of numeric convolution integral ;
70 mutable RooConvIntegrandBinding* _integrand ; //! Binding of Convolution Integrand function
71 mutable RooAbsIntegrator* _integrator ; //! Numeric integrator of convolution integrand
72
73 RooRealProxy _origVar ; // Original convolution variable
74 RooRealProxy _origPdf ; // Original input PDF
75 RooRealProxy _origModel ; // Original resolution model
76
77 mutable RooArgSet _ownedClonedPdfSet ; // Owning set of cloned PDF components
78 mutable RooArgSet _ownedClonedModelSet ; // Owning set of cloned model components
79
80 mutable RooAbsReal* _cloneVar ; // Pointer to cloned convolution variable
81 mutable RooAbsReal* _clonePdf ; // Pointer to cloned PDF
82 mutable RooAbsReal* _cloneModel ; // Pointer to cloned model
83
84 friend class RooConvGenContext ;
85 RooRealVar& cloneVar() const { if (!_init) initialize() ; return (RooRealVar&) *_cloneVar ; }
86 RooAbsReal& clonePdf() const { if (!_init) initialize() ; return (RooAbsReal&) *_clonePdf ; }
87 RooAbsReal& cloneModel() const { if (!_init) initialize() ; return (RooAbsReal&) *_cloneModel ; }
88
89 Bool_t _useWindow ; // Switch to activate window convolution
90 Double_t _windowScale ; // Scale factor for window parameter
91 RooListProxy _windowParam ; // Holder for optional convolution integration window scaling parameter
92
93 Int_t _verboseThresh ; // Call count threshold for verbose printing
94 Bool_t _doProf ; // Switch to activate profiling option
95 TH2* _callHist ; //! Histogram recording number of calls per convolution integral calculation
96
97 ClassDef(RooNumConvolution,1) // Operator PDF implementing numeric convolution of 2 input functions
98};
99
100#endif
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
#define ClassDef(name, id)
Definition: Rtypes.h:326
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition: TGX11.cxx:109
const char * proto
Definition: civetweb.c:16604
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
RooAbsIntegrator is the abstract interface for integrators of real-valued functions that implement th...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooConvGenContext is an efficient implementation of the generator context specific for RooAbsAnaConvP...
Implementation of RooAbsFunc that represent the the integrand of a generic (numeric) convolution A (x...
RooListProxy is the concrete proxy for RooArgList objects.
Definition: RooListProxy.h:25
Numeric 1-dimensional convolution operator PDF.
Definition: RooNumConvPdf.h:26
Numeric 1-dimensional convolution operator PDF.
void clearConvolutionWindow()
Removes previously defined convolution window, reverting to convolution from -inf to +inf.
RooArgSet _ownedClonedModelSet
const TH2 * profileData() const
RooAbsReal * _cloneModel
Bool_t redirectServersHook(const RooAbsCollection &newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
Intercept server redirects. Throw away cache, as figuring out redirections on the cache is an unsolva...
void setConvolutionWindow(RooAbsReal &centerParam, RooAbsReal &widthParam, Double_t widthScaleFactor=1)
Restrict convolution integral to finite range [ x - C - S*W, x - C + S*W ] where x is current value o...
Double_t evaluate() const
Calculate convolution integral.
void initialize() const
One-time initialization of object.
RooRealProxy _origVar
Numeric integrator of convolution integrand.
virtual void printCompactTreeHook(std::ostream &os, const char *indent="")
Hook function to intercept printCompactTree() calls so that it can print out the content of its priva...
RooConvIntegrandBinding * _integrand
virtual TObject * clone(const char *newname) const
RooAbsIntegrator * _integrator
Binding of Convolution Integrand function.
RooRealProxy _origPdf
RooAbsReal & cloneModel() const
RooRealVar & cloneVar() const
void setCallProfiling(Bool_t flag, Int_t nbinX=40, Int_t nbinCall=40, Int_t nCallHigh=1000)
Activate call profile if flag is set to true.
RooListProxy _windowParam
const RooNumIntConfig & convIntConfig() const
RooArgSet _ownedClonedPdfSet
RooRealVar & var() const
RooAbsReal * _cloneVar
RooRealProxy _origModel
RooAbsReal * _clonePdf
RooNumIntConfig & convIntConfig()
virtual ~RooNumConvolution()
Destructor.
void setCallWarning(Int_t threshold=2000)
Activate warning messages if number of function calls needed for evaluation of convolution integral e...
RooNumIntConfig _convIntConfig
RooAbsReal & model() const
RooAbsReal & pdf() const
RooAbsReal & clonePdf() const
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
Service class for 2-Dim histogram classes.
Definition: TH2.h:30
Mother of all ROOT objects.
Definition: TObject.h:37