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

Detailed Description

View in nbviewer Open in SWAN
Special pdf's: linear interpolation between pdf shapes using the 'Alex Read' algorithm

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooPolynomial.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "TH1.h"
using namespace RooFit;
{
// C r e a t e e n d p o i n t p d f s h a p e s
// ------------------------------------------------------
// Observable
RooRealVar x("x", "x", -20, 20);
// Lower end point shape: a Gaussian
RooRealVar g1mean("g1mean", "g1mean", -10);
RooGaussian g1("g1", "g1", x, g1mean, RooConst(2));
// Upper end point shape: a Polynomial
RooPolynomial g2("g2", "g2", x, RooArgSet(-0.03, -0.001));
// C r e a t e i n t e r p o l a t i n g p d f
// -----------------------------------------------
// Create interpolation variable
RooRealVar alpha("alpha", "alpha", 0, 1.0);
// Specify sampling density on observable and interpolation variable
x.setBins(1000, "cache");
alpha.setBins(50, "cache");
// Construct interpolating pdf in (x,a) represent g1(x) at a=a_min
// and g2(x) at a=a_max
RooIntegralMorph lmorph("lmorph", "lmorph", g1, g2, x, alpha);
// P l o t i n t e r p o l a t i n g p d f a t v a r i o u s a l p h a
// -----------------------------------------------------------------------------
// Show end points as blue curves
RooPlot *frame1 = x.frame();
g1.plotOn(frame1);
g2.plotOn(frame1);
// Show interpolated shapes in red
alpha.setVal(0.125);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.25);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.375);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.50);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.625);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.75);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.875);
lmorph.plotOn(frame1, LineColor(kRed));
alpha.setVal(0.95);
lmorph.plotOn(frame1, LineColor(kRed));
// S h o w 2 D d i s t r i b u t i o n o f p d f ( x , a l p h a )
// -----------------------------------------------------------------------
// Create 2D histogram
TH1 *hh = lmorph.createHistogram("hh", x, Binning(40), YVar(alpha, Binning(40)));
// F i t p d f t o d a t a s e t w i t h a l p h a = 0 . 8
// -----------------------------------------------------------------
// Generate a toy dataset at alpha = 0.8
alpha = 0.8;
RooDataSet *data = lmorph.generate(x, 1000);
// Fit pdf to toy data
lmorph.setCacheAlpha(true);
lmorph.fitTo(*data, Verbose(true), PrintLevel(-1));
// Plot fitted pdf and data overlaid
RooPlot *frame2 = x.frame(Bins(100));
data->plotOn(frame2);
lmorph.plotOn(frame2);
// S c a n - l o g ( L ) v s a l p h a
// -----------------------------------------
// Show scan -log(L) of dataset w.r.t alpha
RooPlot *frame3 = alpha.frame(Bins(100), Range(0.1, 0.9));
// Make 2D pdf of histogram
std::unique_ptr<RooAbsReal> nll{lmorph.createNLL(*data)};
nll->plotOn(frame3, ShiftToZero());
lmorph.setCacheAlpha(false);
TCanvas *c = new TCanvas("rf705_linearmorph", "rf705_linearmorph", 800, 800);
c->Divide(2, 2);
c->cd(1);
gPad->SetLeftMargin(0.15);
frame1->GetYaxis()->SetTitleOffset(1.6);
frame1->Draw();
c->cd(2);
gPad->SetLeftMargin(0.20);
hh->GetZaxis()->SetTitleOffset(2.5);
hh->Draw("surf");
c->cd(3);
gPad->SetLeftMargin(0.15);
frame3->GetYaxis()->SetTitleOffset(1.4);
frame3->Draw();
c->cd(4);
gPad->SetLeftMargin(0.15);
frame2->GetYaxis()->SetTitleOffset(1.4);
frame2->Draw();
return;
}
#define c(i)
Definition RSha256.hxx:101
PrintLevel
Definition RooMinuit.h:6
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Class RooIntegralMorph is an implementation of the histogram interpolation technique described by Ale...
A RooPlot is a 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:239
TAxis * GetYaxis() const
Definition RooPlot.cxx:1276
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Definition RooPlot.cxx:649
RooPolynomial implements a polynomial p.d.f of the form.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition TAttAxis.cxx:298
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:324
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3060
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
Ta Range(0, 0, 1, 1)
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a671abfa0 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x_alpha for nset (x,alpha) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a69e91370 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a671abfa0 with pdf g1_MORPH_g2_CACHE_Obs[x]_NORM_x for nset (x) with code 0 from preexisting content.
[#0] PROGRESS:Eval -- RooIntegralMorph::fillCacheObject(lmorph) filling multi-dimensional cache..................................................
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a671abfa0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (g1,g2)
[#0] WARNING:Minimization -- RooAbsMinimizerFcn::synchronize: WARNING: no initial error estimate available for alpha: using 0.1
prevFCN = 2639.749221 alpha=0.8026,
prevFCN = 2639.409749 alpha=0.7974,
prevFCN = 2640.17258 alpha=0.8003,
prevFCN = 2639.703144 alpha=0.7997,
prevFCN = 2639.796557 alpha=0.8116,
prevFCN = 2638.919064 alpha=0.8119,
prevFCN = 2638.930985 alpha=0.8113,
prevFCN = 2638.908768 alpha=0.8094,
prevFCN = 2638.901287 alpha=0.8102,
prevFCN = 2638.885143 alpha=0.8103,
prevFCN = 2638.886672 alpha=0.8101,
prevFCN = 2638.883909 alpha=0.8095,
prevFCN = 2638.898442 alpha=0.8098,
prevFCN = 2638.889204 alpha=0.8099,
prevFCN = 2638.884962 alpha=0.8104,
prevFCN = 2638.887231 alpha=0.8098,
prevFCN = 2638.888592 alpha=0.8102,
prevFCN = 2638.885076 alpha=0.81,
prevFCN = 2638.883889 alpha=0.8098,
prevFCN = 2638.888745 alpha=0.81,
prevFCN = 2638.883589 alpha=0.81,
prevFCN = 2638.883345 alpha=0.8095,
prevFCN = 2638.898453 alpha=0.8098,
prevFCN = 2638.89011 alpha=0.8099,
prevFCN = 2638.886243 alpha=0.81,
prevFCN = 2638.884385 alpha=0.81,
prevFCN = 2638.883475 alpha=0.8101,
prevFCN = 2638.884457 alpha=0.8099,
prevFCN = 2638.885715 alpha=0.8101,
prevFCN = 2638.883988 alpha=0.81,
prevFCN = 2638.884446 alpha=0.81,
prevFCN = 2638.883552 alpha=0.81,
prevFCN = 2638.883372 alpha=0.81,
prevFCN = 2638.883349 alpha=0.81,
prevFCN = 2638.883345 alpha=0.81,
prevFCN = 2638.883345 alpha=0.8101,
prevFCN = 2638.883988 alpha=0.81,
prevFCN = 2638.884446 alpha=0.81,
prevFCN = 2638.883467 alpha=0.81,
prevFCN = 2638.883226 alpha=0.81,
prevFCN = 2638.883369 alpha=0.81,
prevFCN = 2638.883321 alpha=0.81,
prevFCN = 2638.883345 alpha=0.81,
prevFCN = 2638.883369 alpha=0.81,
prevFCN = 2638.883321 alpha=0.8103,
prevFCN = 2638.886457 alpha=0.8097,
prevFCN = 2638.89054 alpha=0.8131,
prevFCN = 2638.992315 alpha=0.8069,
prevFCN = 2639.023618 alpha=0.8171,
prevFCN = 2639.367946 alpha=0.8029,
prevFCN = 2639.374671 alpha=0.8102,
prevFCN = 2638.885206 alpha=0.8098,
prevFCN = 2638.887615 alpha=0.8101,
prevFCN = 2638.883674 alpha=0.81,
prevFCN = 2638.883555 alpha=0.81, [#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a671abfa0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a671abfa0 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a68694230 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
[#1] INFO:Caching -- RooAbsCachedPdf::getCache(lmorph) creating new cache 0x555a68694230 with pdf g1_MORPH_g2_CACHE_Obs[alpha,x]_NORM_x for nset (x) with code 0 from preexisting content.
Date
July 2008
Author
Wouter Verkerke

Definition in file rf705_linearmorph.C.