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

Detailed Description

View in nbviewer Open in SWAN
Multidimensional models: working with parametrized ranges in a fit.

This an example of a fit with an acceptance that changes per-event

pdf = exp(-t/tau) with t[tmin,5]

where t and tmin are both observables in the dataset

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooGaussian.h"
#include "RooConstVar.h"
#include "RooExponential.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "RooFitResult.h"
using namespace RooFit;
{
// D e f i n e o b s e r v a b l e s a n d d e c a y p d f
// ---------------------------------------------------------------
// Declare observables
RooRealVar t("t", "t", 0, 5);
RooRealVar tmin("tmin", "tmin", 0, 0, 5);
// Make parametrized range in t : [tmin,5]
t.setRange(tmin, RooConst(t.getMax()));
// Make pdf
RooRealVar tau("tau", "tau", -1.54, -10, -0.1);
RooExponential model("model", "model", t, tau);
// C r e a t e i n p u t d a t a
// ------------------------------------
// Generate complete dataset without acceptance cuts (for reference)
std::unique_ptr<RooDataSet> dall{model.generate(t, 10000)};
// Generate a (fake) prototype dataset for acceptance limit values
std::unique_ptr<RooDataSet> tmp{RooGaussian("gmin", "gmin", tmin, 0.0, 0.5).generate(tmin, 5000)};
// Generate dataset with t values that observe (t>tmin)
std::unique_ptr<RooDataSet> dacc{model.generate(t, ProtoData(*tmp))};
// F i t p d f t o d a t a i n a c c e p t a n c e r e g i o n
// -----------------------------------------------------------------------
std::unique_ptr<RooFitResult> r{model.fitTo(*dacc, Save(), PrintLevel(-1))};
// P l o t f i t t e d p d f o n f u l l a n d a c c e p t e d d a t a
// ---------------------------------------------------------------------------------
// Make plot frame, add datasets and overlay model
RooPlot *frame = t.frame(Title("Fit to data with per-event acceptance"));
dall->plotOn(frame, MarkerColor(kRed), LineColor(kRed));
model.plotOn(frame);
dacc->plotOn(frame);
// Print fit results to demonstrate absence of bias
r->Print("v");
new TCanvas("rf314_paramranges", "rf314_paramranges", 600, 600);
gPad->SetLeftMargin(0.15);
frame->GetYaxis()->SetTitleOffset(1.6);
frame->Draw();
return;
}
@ kRed
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
#define gPad
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:57
Exponential PDF.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
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:237
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
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 Save(bool flag=true)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg ProtoData(const RooDataSet &protoData, bool randomizeOrder=false, bool resample=false)
RooCmdArg MarkerColor(Color_t color)
RooCmdArg LineColor(Color_t color)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Title
Definition TXMLSetup.cxx:68
[#1] INFO:Fitting -- RooAbsPdf::fitTo(model_over_model_Int[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_model_over_model_Int[t]_modelData) Summation contains a RooNLLVar, using its error level
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 5000 will supersede previous event count of 10000 for normalization of PDF projections
RooFitResult: minimized FCN value: 2823.97, estimated distance to minimum: 3.17108e-08
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter InitialValue FinalValue +/- Error GblCorr.
-------------------- ------------ -------------------------- --------
tau -1.5400e+00 -1.5335e+00 +/- 2.22e-02 <none>
Date
July 2008
Author
Wouter Verkerke

Definition in file rf314_paramfitrange.C.