Logo ROOT   6.10/00
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
rf211_paramconv.C File Reference

Detailed Description

'ADDITION AND CONVOLUTION' RooFit tutorial macro #211

Working a with a p.d.f. with a convolution operator in terms of a parameter

(require ROOT to be compiled with –enable-fftw3)

pict1_rf211_paramconv.C.png
Processing /builddir/build/BUILD/root-6.10.00/tutorials/roofit/rf211_paramconv.C...
RooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:Eval -- RooRealVar::setRange(mean) new range named 'refrange_fft_model' created with bounds [-3,3]
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x223c100 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean] for nset (x,mean) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x1f012b0 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean] for nset (x,mean) with code 0 from preexisting content.
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#0] WARNING:Minization -- RooMinimizerFcn::synchronize: WARNING: no initial error estimate available for a: using 0.5
[#0] WARNING:Minization -- RooMinimizerFcn::synchronize: WARNING: no initial error estimate available for sigma: using 0.2
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 a 2.00000e+00 5.00000e-01 1.00000e+00 1.00000e+01
2 sigma 5.00000e-01 2.00000e-01 1.00000e-01 1.00000e+01
**********
** 3 **SET ERR 0.5
**********
**********
** 4 **SET PRINT 1
**********
**********
** 5 **SET STR 1
**********
NOW USING STRATEGY 1: TRY TO BALANCE SPEED AGAINST RELIABILITY
**********
** 6 **MIGRAD 1000 1
**********
FIRST CALL TO USER FUNCTION AT NEW START POINT, WITH IFLAG=4.
prevFCN = 2171.275755 START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
a=2.012, sigma=0.5, [#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2.012,
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2, sigma=0.5047,
prevFCN = 2171.286528 sigma=0.4953,
prevFCN = 2171.267762 sigma=0.5029,
prevFCN = 2171.281998 sigma=0.4971,
prevFCN = 2171.270547 a=2.012, sigma=0.5,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.012,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2.012,
prevFCN = 2171.275755 a=1.988,
prevFCN = 2171.275755 a=2.121,
prevFCN = 2171.861215 a=1.886,
prevFCN = 2172.184717 a=2, sigma=0.5029,
prevFCN = 2171.281998 sigma=0.4971,
prevFCN = 2171.270547 FCN=2171.28 FROM MIGRAD STATUS=INITIATE 29 CALLS 30 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a 2.00000e+00 5.00000e-01 0.00000e+00 -3.89301e+00
2 sigma 5.00000e-01 2.00000e-01 0.00000e+00 3.88521e+00
ERR DEF= 0.5
a=2.013, sigma=0.4843,
prevFCN = 2171.259881 a=2.009, sigma=0.4884,
prevFCN = 2171.260992 a=2.025, sigma=0.4843,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.134,
prevFCN = 2171.692149 a=1.898,
prevFCN = 2172.378568 a=2.025,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.013, sigma=0.4871,
prevFCN = 2171.26042 sigma=0.4815,
prevFCN = 2171.260367 MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
sigma=0.4843,
prevFCN = 2171.259881 a=2.025,
prevFCN = 2171.259881 a=2.001,
prevFCN = 2171.259881 a=2.013, sigma=0.4871,
prevFCN = 2171.26042 sigma=0.4815,
prevFCN = 2171.260367 a=2.015, sigma=0.4843,
prevFCN = 2171.259881 a=2.01,
prevFCN = 2171.259881 a=2.013, sigma=0.4848,
prevFCN = 2171.259907 sigma=0.4837,
prevFCN = 2171.259896 a=2.025, sigma=0.4871,
prevFCN = 2171.26042 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=2171.26 FROM MIGRAD STATUS=CONVERGED 49 CALLS 50 TOTAL
EDM=6.76006e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 a 2.01276e+00 6.78267e+00 4.15492e-03 -2.73620e-10
2 sigma 4.84270e-01 8.79333e-02 1.47360e-03 1.78734e-02
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
1.535e+08 0.000e+00
0.000e+00 7.738e-03
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.00000 1.000 0.000
2 0.00000 0.000 1.000
a=2.013, sigma=0.4843, **********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 1000
**********
prevFCN = 2171.259881 a=2.015,
prevFCN = 2171.259881 a=2.01,
prevFCN = 2171.259881 a=2.266,
prevFCN = 2173.490353 a=1.784,
prevFCN = 2175.249474 a=2.114,
prevFCN = 2171.692149 a=1.915,
prevFCN = 2172.378568 a=2.013, sigma=0.4848,
prevFCN = 2171.259907 sigma=0.4837,
prevFCN = 2171.259896 a=2.015, sigma=0.4843,
prevFCN = 2171.259881 a=2.011,
prevFCN = 2171.259881 a=2.013, sigma=0.4844,
prevFCN = 2171.259883 sigma=0.4842,
prevFCN = 2171.25988 a=2.114, sigma=0.4848,
prevFCN = 2171.697792 COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=2171.26 FROM HESSE STATUS=OK 14 CALLS 64 TOTAL
EDM=0.151236 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 a 2.01276e+00 1.12520e-01 3.49974e-02 -8.86625e-01
2 sigma 4.84270e-01 1.23783e-01 2.94719e-04 -1.17417e+00
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
1.267e-02 -9.816e-03
-9.816e-03 1.534e-02
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.70406 1.000 -0.704
2 0.70406 -0.704 1.000
a=2.013, sigma=0.4843, [#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(gx_Int[mean,x]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(model_mean_Int[mean]) using numeric integrator RooIntegrator1D to calculate Int(mean)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x241fb60 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean] for nset (x,mean) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(model) creating new cache 0x2989270 with pdf gx_CONV_model_mean_CACHE_Obs[x,mean] for nset (x) with code 0 from preexisting content.
#include "RooRealVar.h"
#include "RooDataHist.h"
#include "RooGaussian.h"
#include "RooGenericPdf.h"
#include "RooFormulaVar.h"
#include "RooFFTConvPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH2.h"
using namespace RooFit ;
void rf211_paramconv()
{
// S e t u p c o m p o n e n t p d f s
// ---------------------------------------
// Gaussian g(x ; mean,sigma)
RooRealVar x("x","x",-10,10) ;
RooRealVar mean("mean","mean",-3,3) ;
RooRealVar sigma("sigma","sigma",0.5,0.1,10) ;
RooGaussian modelx("gx","gx",x,mean,sigma) ;
// Block function in mean
RooRealVar a("a","a",2,1,10) ;
RooGenericPdf model_mean("model_mean","abs(mean)<a",RooArgList(mean,a)) ;
// Convolution in mean parameter model = g(x,mean,sigma) (x) block(mean)
x.setBins(1000,"cache") ;
mean.setBins(50,"cache") ;
RooFFTConvPdf model("model","model",mean,modelx,model_mean) ;
// Configure convolution to construct a 2-D cache in (x,mean)
// rather than a 1-d cache in mean that needs to be recalculated
// for each value of x
model.setCacheObservables(x) ;
model.setBufferFraction(1.0) ;
// Integrate model over mean projModel = Int model dmean
RooAbsPdf* projModel = model.createProjection(mean) ;
// Generate 1000 toy events
RooDataHist* d = projModel->generateBinned(x,1000) ;
// Fit p.d.f. to toy data
projModel->fitTo(*d,Verbose()) ;
// Plot data and fitted p.d.f.
RooPlot* frame = x.frame(Bins(25)) ;
d->plotOn(frame) ;
projModel->plotOn(frame) ;
// Make 2d histogram of model(x;mean)
TH1* hh = model.createHistogram("hh",x,Binning(50),YVar(mean,Binning(50)),ConditionalObservables(mean)) ;
hh->SetTitle("histogram of model(x|mean)") ;
// Draw frame on canvas
TCanvas* c = new TCanvas("rf211_paramconv","rf211_paramconv",800,400) ;
c->Divide(2) ;
c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
c->cd(2) ; gPad->SetLeftMargin(0.20) ; hh->GetZaxis()->SetTitleOffset(2.5) ; hh->Draw("surf") ;
}
Author
04/2009 - Wouter Verkerke

Definition in file rf211_paramconv.C.