Logo ROOT   6.14/05
Reference Guide
rf107_plotstyles.C File Reference

Detailed Description

View in nbviewer Open in SWAN 'BASIC FUNCTIONALITY' RooFit tutorial macro #107

Demonstration of various plotting styles of data, functions in a RooPlot

pict1_rf107_plotstyles.C.png
Processing /mnt/build/workspace/root-makedoc-v614/rootspi/rdoc/src/v6-14-00-patches/tutorials/roofit/rf107_plotstyles.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:Minization -- RooMinimizer::optimizeConst: activating const optimization
**********
** 1 **SET PRINT 1
**********
**********
** 2 **SET NOGRAD
**********
PARAMETER DEFINITIONS:
NO. NAME VALUE STEP SIZE LIMITS
1 mean -3.00000e+00 2.00000e+00 -1.00000e+01 1.00000e+01
2 sigma 3.00000e+00 9.90000e-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.
START MIGRAD MINIMIZATION. STRATEGY 1. CONVERGENCE WHEN EDM .LT. 1.00e-03
FCN=244.778 FROM MIGRAD STATUS=INITIATE 6 CALLS 7 TOTAL
EDM= unknown STRATEGY= 1 NO ERROR MATRIX
EXT PARAMETER CURRENT GUESS STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 mean -3.00000e+00 2.00000e+00 2.11716e-01 7.88402e+00
2 sigma 3.00000e+00 9.90000e-01 2.22742e-01 8.68850e+00
ERR DEF= 0.5
MIGRAD MINIMIZATION HAS CONVERGED.
MIGRAD WILL VERIFY CONVERGENCE AND ERROR MATRIX.
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=244.648 FROM MIGRAD STATUS=CONVERGED 27 CALLS 28 TOTAL
EDM=6.12289e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER STEP FIRST
NO. NAME VALUE ERROR SIZE DERIVATIVE
1 mean -3.06106e+00 3.00167e-01 3.38614e-04 -1.01280e-02
2 sigma 2.89572e+00 2.28664e-01 5.51106e-04 1.31676e-02
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
9.013e-02 -8.498e-03
-8.498e-03 5.233e-02
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.12374 1.000 -0.124
2 0.12374 -0.124 1.000
**********
** 7 **SET ERR 0.5
**********
**********
** 8 **SET PRINT 1
**********
**********
** 9 **HESSE 1000
**********
COVARIANCE MATRIX CALCULATED SUCCESSFULLY
FCN=244.648 FROM HESSE STATUS=OK 10 CALLS 38 TOTAL
EDM=6.13161e-07 STRATEGY= 1 ERROR MATRIX ACCURATE
EXT PARAMETER INTERNAL INTERNAL
NO. NAME VALUE ERROR STEP SIZE VALUE
1 mean -3.06106e+00 3.00196e-01 6.77227e-05 -3.11100e-01
2 sigma 2.89572e+00 2.28685e-01 1.10221e-04 -4.50268e-01
ERR DEF= 0.5
EXTERNAL ERROR MATRIX. NDIM= 25 NPAR= 2 ERR DEF=0.5
9.015e-02 -8.552e-03
-8.552e-03 5.234e-02
PARAMETER CORRELATION COEFFICIENTS
NO. GLOBAL 1 2
1 0.12449 1.000 -0.124
2 0.12449 -0.124 1.000
[#1] INFO:Minization -- RooMinimizer::optimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(gauss) only plotting range [-8,3], curve is normalized to data in given given range
[#1] INFO:Eval -- RooRealVar::setRange(x) new range named 'plotRange' created with bounds [-8,3]
#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
using namespace RooFit ;
void rf107_plotstyles()
{
// S e t u p m o d e l
// ---------------------
// Create observables
RooRealVar x("x","x",-10,10) ;
// Create Gaussian
RooRealVar sigma("sigma","sigma",3,0.1,10) ;
RooRealVar mean("mean","mean",-3,-10,10) ;
RooGaussian gauss("gauss","gauss",x,mean,sigma) ;
// Generate a sample of 100 events with sigma=3
RooDataSet* data = gauss.generate(x,100) ;
// Fit pdf to data
gauss.fitTo(*data) ;
// M a k e p l o t f r a m e s
// -------------------------------
// Make four plot frames to demonstrate various plotting features
RooPlot* frame1 = x.frame(Name("xframe"),Title("Red Curve / SumW2 Histo errors"),Bins(20)) ;
RooPlot* frame2 = x.frame(Name("xframe"),Title("Dashed Curve / No XError bars"),Bins(20)) ;
RooPlot* frame3 = x.frame(Name("xframe"),Title("Filled Curve / Blue Histo"),Bins(20)) ;
RooPlot* frame4 = x.frame(Name("xframe"),Title("Partial Range / Filled Bar chart"),Bins(20)) ;
// D a t a p l o t t i n g s t y l e s
// ---------------------------------------
// Use sqrt(sum(weights^2)) error instead of Poisson errors
// Remove horizontal error bars
data->plotOn(frame2,XErrorSize(0)) ;
// Blue markers and error bors
data->plotOn(frame3,MarkerColor(kBlue),LineColor(kBlue)) ;
// Filled bar chart
// F u n c t i o n p l o t t i n g s t y l e s
// -----------------------------------------------
// Change line color to red
gauss.plotOn(frame1,LineColor(kRed)) ;
// Change line style to dashed
gauss.plotOn(frame2,LineStyle(kDashed)) ;
// Filled shapes in green color
gauss.plotOn(frame3,DrawOption("F"),FillColor(kOrange),MoveToBack()) ;
//
gauss.plotOn(frame4,Range(-8,3),LineColor(kMagenta)) ;
TCanvas* c = new TCanvas("rf107_plotstyles","rf107_plotstyles",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.15) ; frame2->GetYaxis()->SetTitleOffset(1.6) ; frame2->Draw() ;
c->cd(3) ; gPad->SetLeftMargin(0.15) ; frame3->GetYaxis()->SetTitleOffset(1.6) ; frame3->Draw() ;
c->cd(4) ; gPad->SetLeftMargin(0.15) ; frame4->GetYaxis()->SetTitleOffset(1.6) ; frame4->Draw() ;
}
Author
07/2008 - Wouter Verkerke

Definition in file rf107_plotstyles.C.