Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf210_angularconv.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Addition and convolution: convolution in cyclical angular observables theta

and construction of pdf in terms of transformed angular coordinates, e.g. cos(theta), where the convolution is performed in theta rather than cos(theta)

pdf(theta) = T(theta) (x) gauss(theta)
pdf(cosTheta) = T(acos(cosTheta)) (x) gauss(acos(cosTheta))
Double_t x[n]
Definition legend1.C:17

This tutorial requires FFT3 to be enabled.

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooGenericPdf.h"
#include "RooFormulaVar.h"
#include "RooFFTConvPdf.h"
#include "RooPlot.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "TH1.h"
using namespace RooFit;
{
// S e t u p c o m p o n e n t p d f s
// ---------------------------------------
// Define angle psi
RooRealVar psi("psi", "psi", 0, 3.14159268);
// Define physics pdf T(psi)
RooGenericPdf Tpsi("Tpsi", "1+sin(2*@0)", psi);
// Define resolution R(psi)
RooRealVar gbias("gbias", "gbias", 0.2, 0., 1);
RooRealVar greso("greso", "greso", 0.3, 0.1, 1.0);
RooGaussian Rpsi("Rpsi", "Rpsi", psi, gbias, greso);
// Define cos(psi) and function psif that calculates psi from cos(psi)
RooRealVar cpsi("cpsi", "cos(psi)", -1, 1);
RooFormulaVar psif("psif", "acos(cpsi)", cpsi);
// Define physics pdf also as function of cos(psi): T(psif(cpsi)) = T(cpsi) ;
RooGenericPdf Tcpsi("T", "1+sin(2*@0)", psif);
// C o n s t r u c t c o n v o l u t i o n p d f i n p s i
// --------------------------------------------------------------
// Define convoluted pdf as function of psi: M=[T(x)R](psi) = M(psi)
RooFFTConvPdf Mpsi("Mf", "Mf", psi, Tpsi, Rpsi);
// Set the buffer fraction to zero to obtain a true cyclical convolution
Mpsi.setBufferFraction(0);
// S a m p l e , f i t a n d p l o t c o n v o l u t e d p d f ( p s i )
// --------------------------------------------------------------------------------
// Generate some events in observable psi
std::unique_ptr<RooDataSet> data_psi{Mpsi.generate(psi, 10000)};
// Fit convoluted model as function of angle psi
Mpsi.fitTo(*data_psi, PrintLevel(-1));
// Plot cos(psi) frame with Mf(cpsi)
RooPlot *frame1 = psi.frame(Title("Cyclical convolution in angle psi"));
data_psi->plotOn(frame1);
Mpsi.plotOn(frame1);
// Overlay comparison to unsmeared physics pdf T(psi)
Tpsi.plotOn(frame1, LineColor(kRed));
// C o n s t r u c t c o n v o l u t i o n p d f i n c o s ( p s i )
// --------------------------------------------------------------------------
// Define convoluted pdf as function of cos(psi): M=[T(x)R](psif(cpsi)) = M(cpsi)
//
// Need to give both observable psi here (for definition of convolution)
// and function psif here (for definition of observables, ultimately in cpsi)
RooFFTConvPdf Mcpsi("Mf", "Mf", psif, psi, Tpsi, Rpsi);
// Set the buffer fraction to zero to obtain a true cyclical convolution
Mcpsi.setBufferFraction(0);
// S a m p l e , f i t a n d p l o t c o n v o l u t e d p d f ( c o s p s i )
// --------------------------------------------------------------------------------
// Generate some events
std::unique_ptr<RooDataSet> data_cpsi{Mcpsi.generate(cpsi, 10000)};
// set psi constant to exclude to be a parameter of the fit
psi.setConstant(true);
// Fit convoluted model as function of cos(psi)
Mcpsi.fitTo(*data_cpsi, PrintLevel(-1));
// Plot cos(psi) frame with Mf(cpsi)
RooPlot *frame2 = cpsi.frame(Title("Same convolution in psi, expressed in cos(psi)"));
data_cpsi->plotOn(frame2);
Mcpsi.plotOn(frame2);
// Overlay comparison to unsmeared physics pdf Tf(cpsi)
Tcpsi.plotOn(frame2, LineColor(kRed));
// Draw frame on canvas
TCanvas *c = new TCanvas("rf210_angularconv", "rf210_angularconv", 800, 400);
c->Divide(2);
c->cd(1);
gPad->SetLeftMargin(0.15);
frame1->GetYaxis()->SetTitleOffset(1.4);
frame1->Draw();
c->cd(2);
gPad->SetLeftMargin(0.15);
frame2->GetYaxis()->SetTitleOffset(1.4);
frame2->Draw();
}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
#define gPad
PDF for the numerical (FFT) convolution of two PDFs.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Implementation of a probability density function that takes a RooArgList of servers and a C++ express...
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
static RooPlot * frame(const RooAbsRealLValue &var, double xmin, double xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition RooPlot.cxx:185
TAxis * GetYaxis() const
Definition RooPlot.cxx:1224
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:597
Variable that can be changed from the outside.
Definition RooRealVar.h:37
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
The Canvas class.
Definition TCanvas.h:23
RooCmdArg PrintLevel(Int_t code)
RooCmdArg LineColor(Color_t color)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
const char * Title
Definition TXMLSetup.cxx:68
[#1] INFO:Caching -- Changing internal binning of variable 'psi' in FFT 'Mf' from 100 to 930 to improve the precision of the numerical FFT. This can be done manually by setting an additional binning named 'cache'.
[#1] INFO:Eval -- RooRealVar::setRange(psi) new range named 'refrange_fft_Mf' created with bounds [0,3.14159]
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564411f0b790 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564413930510 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564413930510 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0 from preexisting content.
[#1] INFO:Fitting -- RooAbsPdf::fitTo(Mf_over_Mf_Int[psi]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- using CPU computation library compiled with -mavx2
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_Mf_over_Mf_Int[psi]_MfData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564413966520 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi] for nset () with code 1 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564411ef6dc0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[psi]_NORM_psi for nset (psi) with code 0
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564411ef6dc0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x5644120fdb80 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x5644120fdb80 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0 from preexisting content.
[#1] INFO:Fitting -- RooAbsPdf::fitTo(Mf_over_Mf_Int[cpsi]) fixing normalization set for coefficient determination to observables in data
[#1] INFO:Fitting -- RooAddition::defaultErrorLevel(nll_Mf_over_Mf_Int[cpsi]_MfData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564411d564f0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi] for nset () with code 1 from preexisting content.
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Mf_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_Int[psi]) using numeric integrator RooIntegrator1D to calculate Int(psi)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(Mf) creating new cache 0x564411d564f0 with pdf Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi for nset (cpsi) with code 0
[#1] INFO:NumericIntegration -- RooRealIntegral::init(Tpsi_CONV_Rpsi_CACHE_Obs[cpsi]_NORM_cpsi_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(T_Int[cpsi]) using numeric integrator RooIntegrator1D to calculate Int(cpsi)
Date
April 2009
Author
Wouter Verkerke

Definition in file rf210_angularconv.C.