Logo ROOT  
Reference Guide
rf704_amplitudefit.C File Reference

Detailed Description

View in nbviewer Open in SWAN Speecial p.d.f.

's: using a p.d.f defined by a sum of real-valued amplitude components

␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Minization -- RooMinimizer::optimizeConst: activating const optimization
[#1] INFO:Minization -- The following expressions have been identified as constant and will be precalculated and cached: (ampl1,ampl2)
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 f1 1.00000e+00 2.00000e-01 0.00000e+00 2.00000e+00
2 f2 5.00000e-01 2.00000e-01 0.00000e+00 2.00000e+00
**********
** 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.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=24323.4 FROM MIGRAD STATUS=INITIATE 8 CALLS 9 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 f1 1.00000e+00 2.00000e-01 2.01358e-01 2.36561e+01
2 f2 5.00000e-01 2.00000e-01 2.35352e-01 -4.09909e+01
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=24323.1 FROM MIGRAD STATUS=CONVERGED 33 CALLS 34 TOTAL
EDM=1.04489e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 f1 9.85261e-01 4.18474e-01 3.64059e-03 -9.57093e-03
2 f2 5.07031e-01 2.19792e-01 2.15145e-03 1.61567e-02
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
1.864e-01 9.565e-02
9.565e-02 4.937e-02
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.99694 1.000 0.997
2 0.99694 0.997 1.000
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 1000
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=24323.1 FROM HESSE STATUS=OK 10 CALLS 44 TOTAL
EDM=1.04101e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 f1 9.85261e-01 1.40907e+00 1.45624e-04 -1.47397e-02
2 f2 5.07031e-01 1.24234e+00 8.60579e-05 -5.15499e-01
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
4.590e+00 2.362e+00
2.362e+00 1.215e+00
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.99988 1.000 1.000
2 0.99988 1.000 1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#0] WARNING:InputArguments -- RooAbsReal::createHistogram(ampl1) WARNING extended mode requested for a non-pdf object, ignored
[#0] WARNING:InputArguments -- RooAbsReal::createHistogram(ampl2) WARNING extended mode requested for a non-pdf object, ignored
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_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,tm_conv_coshGBasis_[t],coshGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_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,tm_conv_sinhGBasis_[t],sinhGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on t integrates over variables (cosa)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_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(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_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,tm_conv_coshGBasis_[t],coshGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on cosa integrates over variables (t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_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,tm_conv_sinhGBasis_[t],sinhGBasis)
[#1] INFO:Plotting -- RooAbsReal::plotOn(pdf) plot on cosa integrates over variables (t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_coshGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(tm_conv_sinhGBasis_[t]_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooTruthModel.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.);
// Use RooTruthModel to obtain compiled implementation of sinh/cosh modulated decay functions
RooRealVar tau("tau", "#tau", 1.5);
RooRealVar deltaGamma("deltaGamma", "deltaGamma", 0.3);
RooTruthModel truthModel("tm", "tm", t);
RooFormulaVar coshGBasis("coshGBasis", "exp(-@0/ @1)*cosh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
RooFormulaVar sinhGBasis("sinhGBasis", "exp(-@0/ @1)*sinh(@0*@2/2)", RooArgList(t, tau, deltaGamma));
RooAbsReal *coshGConv = truthModel.convolution(&coshGBasis, &t);
RooAbsReal *sinhGConv = truthModel.convolution(&sinhGBasis, &t);
// Construct polynomial amplitudes in cos(a)
RooPolyVar poly1("poly1", "poly1", cosa, RooArgList(RooConst(0.5), RooConst(0.2), RooConst(0.2)), 0);
RooPolyVar poly2("poly2", "poly2", cosa, RooArgList(RooConst(1), RooConst(-0.2), RooConst(3)), 0);
// Construct 2D amplitude as uncorrelated product of amp(t)*amp(cosa)
RooProduct ampl1("ampl1", "amplitude 1", RooArgSet(poly1, *coshGConv));
RooProduct ampl2("ampl2", "amplitude 2", RooArgSet(poly2, *sinhGConv));
// 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
RooDataSet *data = pdf.generate(RooArgSet(t, cosa), 10000);
// Fit pdf to toy data with only amplitude strength floating
pdf.fitTo(*data);
// 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:64
@ kBlue
Definition: Rtypes.h:64
@ kDashed
Definition: TAttLine.h:48
#define gPad
Definition: TVirtualPad.h:286
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
Definition: RooAbsData.cxx:550
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
Definition: RooFormulaVar.h:29
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAxis * GetYaxis() const
Definition: RooPlot.cxx:1277
static RooPlot * frame(const RooAbsRealLValue &var, Double_t xmin, Double_t xmax, Int_t nBins)
Create a new frame for a given variable in x.
Definition: RooPlot.cxx:277
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:712
Class RooPolyVar is a RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficie...
Definition: RooPolyVar.h:28
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition: RooProduct.h:32
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
Definition: RooRealSumPdf.h:24
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:31
The TH1 histogram class.
Definition: TH1.h:56
TAxis * GetZaxis()
Definition: TH1.h:318
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2998
TF1 * f1
Definition: legend1.C:11
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg Components(const RooArgSet &compSet)
RooConstVar & RooConst(Double_t val)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
Author
07/2008 - Wouter Verkerke

Definition in file rf704_amplitudefit.C.