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

Detailed Description

View in nbviewer Open in SWAN
Special pdf's: using a pdf defined by a sum of real-valued amplitude components

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooFormulaVar.h"
#include "RooRealSumPdf.h"
#include "RooPolyVar.h"
#include "RooProduct.h"
#include "TH1.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit;
{
// S e t u p 2 D a m p l i t u d e f u n c t i o n s
// -------------------------------------------------------
// Observables
RooRealVar t("t", "time", -1., 15.);
RooRealVar cosa("cosa", "cos(alpha)", -1., 1.);
RooRealVar tau("tau", "#tau", 1.5);
RooRealVar deltaGamma("deltaGamma", "deltaGamma", 0.3);
RooFormulaVar coshG("coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", {t, tau, deltaGamma});
RooFormulaVar sinhG("sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", {t, tau, deltaGamma});
// Construct polynomial amplitudes in cos(a)
RooPolyVar poly1("poly1", "poly1", cosa, RooArgList{0.5, 0.2, 0.2}, 0);
RooPolyVar poly2("poly2", "poly2", cosa, RooArgList{1.0, -0.2, 3.0}, 0);
// Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
RooProduct ampl1("ampl1", "amplitude 1", {poly1, coshG});
RooProduct ampl2("ampl2", "amplitude 2", {poly2, sinhG});
// C o n s t r u c t a m p l i t u d e s u m p d f
// -----------------------------------------------------
// Amplitude strengths
RooRealVar f1("f1", "f1", 1, 0, 2);
RooRealVar f2("f2", "f2", 0.5, 0, 2);
// Construct pdf
RooRealSumPdf pdf("pdf", "pdf", RooArgList(ampl1, ampl2), RooArgList(f1, f2));
// Generate some toy data from pdf
std::unique_ptr<RooDataSet> data{pdf.generate({t, cosa}, 10000)};
// Fit pdf to toy data with only amplitude strength floating
pdf.fitTo(*data, PrintLevel(-1));
// P l o t a m p l i t u d e s u m p d f
// -------------------------------------------
// Make 2D plots of amplitudes
TH1 *hh_cos = ampl1.createHistogram("hh_cos", t, Binning(50), YVar(cosa, Binning(50)));
TH1 *hh_sin = ampl2.createHistogram("hh_sin", t, Binning(50), YVar(cosa, Binning(50)));
hh_cos->SetLineColor(kBlue);
hh_sin->SetLineColor(kRed);
// Make projection on t, plot data, pdf and its components
// Note component projections may be larger than sum because amplitudes can be negative
RooPlot *frame1 = t.frame();
data->plotOn(frame1);
pdf.plotOn(frame1);
pdf.plotOn(frame1, Components(ampl1), LineStyle(kDashed));
pdf.plotOn(frame1, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
// Make projection on cosa, plot data, pdf and its components
// Note that components projection may be larger than sum because amplitudes can be negative
RooPlot *frame2 = cosa.frame();
data->plotOn(frame2);
pdf.plotOn(frame2);
pdf.plotOn(frame2, Components(ampl1), LineStyle(kDashed));
pdf.plotOn(frame2, Components(ampl2), LineStyle(kDashed), LineColor(kRed));
TCanvas *c = new TCanvas("rf704_amplitudefit", "rf704_amplitudefit", 800, 800);
c->Divide(2, 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();
c->cd(3);
gPad->SetLeftMargin(0.20);
hh_cos->GetZaxis()->SetTitleOffset(2.3);
hh_cos->Draw("surf");
c->cd(4);
gPad->SetLeftMargin(0.20);
hh_sin->GetZaxis()->SetTitleOffset(2.3);
hh_sin->Draw("surf");
}
#define c(i)
Definition RSha256.hxx:101
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kDashed
Definition TAttLine.h:48
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
#define gPad
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
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
A RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficients.
Definition RooPolyVar.h:25
Represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
Implements a PDF constructed from a sum of functions:
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
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:59
TAxis * GetZaxis()
Definition TH1.h:327
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3068
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Components(Args_t &&... argsOrArgSet)
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
TF1 * f1
Definition legend1.C:11
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Fitting -- RooAbsPdf::fitTo(pdf_over_pdf_Int[cosa,t]) 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_pdf_over_pdf_Int[cosa,t]_pdfData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) directly selected PDF components: (ampl1)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) indirectly selected PDF components: (poly1,coshGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) directly selected PDF components: (ampl2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) indirectly selected PDF components: (poly2,sinhGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on cosa integrates over variables (t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) directly selected PDF components: (ampl1)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) indirectly selected PDF components: (poly1,coshGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on cosa integrates over variables (t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) directly selected PDF components: (ampl2)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(pdf) indirectly selected PDF components: (poly2,sinhGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on cosa integrates over variables (t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(coshGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(sinhGBasis_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
Date
July 2008
Author
Wouter Verkerke

Definition in file rf704_amplitudefit.C.