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

Detailed Description

View in nbviewer Open in SWAN
Performing a simple fit with RooLagrangianMorphFunc.

a morphing function is setup as a function of three variables and a fit is performed on a pseudo-dataset.

#include <RooDataHist.h>
#include <RooFitResult.h>
#include <RooPlot.h>
#include <RooRealVar.h>
#include <TAxis.h>
#include <TCanvas.h>
#include <TH2.h>
#include <TStyle.h>
using namespace RooFit;
{
// C r e a t e v a r i a b l e s f o r
// m o r p h i n g f u n c t i o n
// ---------------------------------------------
std::string observablename = "pTV";
RooRealVar obsvar(observablename.c_str(), "observable of pTV", 10, 600);
RooRealVar kSM("kSM", "sm modifier", 1.0);
RooRealVar cHq3("cHq3", "EFT modifier", -10.0, 10.0);
cHq3.setAttribute("NewPhysics", true);
RooRealVar cHl3("cHl3", "EFT modifier", -10.0, 10.0);
cHl3.setAttribute("NewPhysics", true);
RooRealVar cHDD("cHDD", "EFT modifier", -10.0, 10.0);
cHDD.setAttribute("NewPhysics", true);
// I n p u t s n e e d e d f o r c o n f i g
// ---------------------------------------------
std::string infilename = std::string(gROOT->GetTutorialDir()) + "/roofit/input_histos_rf_lagrangianmorph.root";
std::vector<std::string> samplelist = {"SM_NPsq0", "cHq3_NPsq1", "cHq3_NPsq2", "cHl3_NPsq1",
"cHl3_NPsq2", "cHDD_NPsq1", "cHDD_NPsq2", "cHl3_cHDD_NPsq2",
"cHq3_cHDD_NPsq2", "cHl3_cHq3_NPsq2"};
// S e t u p C o n f i g
// ---------------------------------------------
config.fileName = infilename;
config.observableName = observablename;
config.folderNames = samplelist;
config.couplings.add(cHq3);
config.couplings.add(cHl3);
config.couplings.add(cHDD);
config.couplings.add(kSM);
// C r e a t e m o r p h i n g f u n c t i o n
// ---------------------------------------------
RooLagrangianMorphFunc morphfunc("morphfunc", "morphed dist. of pTV", config);
// C r e a t e p s e u d o d a t a h i s t o g r a m
// f o r f i t
// ---------------------------------------------
morphfunc.setParameter("cHq3", 0.01);
morphfunc.setParameter("cHl3", 1.0);
morphfunc.setParameter("cHDD", 0.2);
auto pseudo_hist = morphfunc.createTH1("pseudo_hist");
auto pseudo_dh = new RooDataHist("pseudo_dh", "pseudo_dh", RooArgList(obsvar), pseudo_hist);
// reset parameters to zeros before fit
morphfunc.setParameter("cHq3", 0.0);
morphfunc.setParameter("cHl3", 0.0);
morphfunc.setParameter("cHDD", 0.0);
// error set used as initial step size
cHq3.setError(0.1);
cHl3.setError(0.1);
cHDD.setError(0.1);
// W r a p p d f o n m o r p h f u n c a n d
// f i t t o d a t a h i s t o g r a m
// ---------------------------------------------
// wrapper pdf to normalise morphing function to a morphing pdf
RooWrapperPdf model("wrap_pdf", "wrap_pdf", morphfunc);
auto fitres = model.fitTo(*pseudo_dh, SumW2Error(true), Optimize(false), Save(), PrintLevel(-1));
auto hcorr = fitres->correlationHist();
// E x t r a c t p o s t f i t d i s t r i b u t i o n
// a n d p l o t w i t h i n i t i a l
// h i s t o g r a m
// ---------------------------------------------
auto postfit_hist = morphfunc.createTH1("morphing_postfit_hist");
RooDataHist postfit_dh("morphing_postfit_dh", "morphing_postfit_dh", RooArgList(obsvar), postfit_hist);
auto frame0 = obsvar.frame(Title("Fitted histogram of p_{T}^{V}"));
postfit_dh.plotOn(frame0, Name("postfit_dist"), DrawOption("C"), LineColor(kBlue), DataError(RooAbsData::None),
XErrorSize(0));
pseudo_dh->plotOn(frame0, Name("input"));
// D r a w p l o t s o n c a n v a s
// ---------------------------------------------
TCanvas *c1 = new TCanvas("fig3", "fig3", 800, 400);
c1->Divide(2, 1);
c1->cd(1);
gPad->SetLeftMargin(0.15);
gPad->SetRightMargin(0.05);
model.paramOn(frame0, Layout(0.50, 0.75, 0.9), Parameters(config.couplings));
frame0->GetXaxis()->SetTitle("p_{T}^{V}");
frame0->Draw();
c1->cd(2);
gPad->SetLeftMargin(0.15);
gPad->SetRightMargin(0.15);
hcorr->SetMarkerSize(3.);
hcorr->SetTitle("correlation matrix");
hcorr->GetYaxis()->SetTitleOffset(1.4);
hcorr->GetYaxis()->SetLabelSize(0.1);
hcorr->GetXaxis()->SetLabelSize(0.1);
hcorr->GetYaxis()->SetBinLabel(1, "c_{HDD}");
hcorr->GetYaxis()->SetBinLabel(2, "c_{Hl^{(3)}}");
hcorr->GetYaxis()->SetBinLabel(3, "c_{Hq^{(3)}}");
hcorr->GetXaxis()->SetBinLabel(3, "c_{HDD}");
hcorr->GetXaxis()->SetBinLabel(2, "c_{Hl^{(3)}}");
hcorr->GetXaxis()->SetBinLabel(1, "c_{Hq^{(3)}}");
hcorr->Draw("colz text");
c1->SaveAs("rf712_lagrangianmorphfit.png");
}
@ kBlue
Definition Rtypes.h:66
#define gROOT
Definition TROOT.h:406
R__EXTERN TStyle * gStyle
Definition TStyle.h:436
#define gPad
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Class RooLagrangianMorphing is a implementation of the method of Effective Lagrangian Morphing,...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The RooWrapperPdf is a class that can be used to convert a function into a PDF.
The Canvas class.
Definition TCanvas.h:23
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1640
void SetPaintTextFormat(const char *format="g")
Definition TStyle.h:386
return c1
Definition legend1.C:41
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
std::vector< std::string > folderNames
[#0] PROGRESS:InputArguments -- initializing physics inputs from file /home/sftnight/build/workspace/root-makedoc-master/rootspi/rdoc/src/master.build/tutorials/roofit/input_histos_rf_lagrangianmorph.root with object name(s) 'pTV'
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x7ffd6ad98448
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
[#1] INFO:DataHandling -- RooDataHist::adjustBinning(pseudo_dh): fit range of variable pTV expanded to nearest bin boundaries: [10,600] --> [0,600]
[#0] PROGRESS:Caching -- creating cache from getCache function for 0x55a2b827b690
[#0] PROGRESS:Caching -- current storage has size 10
[#0] PROGRESS:ObjectHandling -- observable: pTV
[#0] PROGRESS:ObjectHandling -- binWidth: binWidth_pTV
[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf_over_wrap_pdf_Int[pTV]) 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_wrap_pdf_over_wrap_pdf_Int[pTV]_pseudo_dh) Summation contains a RooNLLVar, using its error level
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0 cHl3=0 cHq3=-0.0202918
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0376648, denominator=wrap_pdf_Int[pTV]=10358.6
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.213672 cHl3=1.97898 cHq3=0.00773174
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=0.0882577, denominator=wrap_pdf_Int[pTV]=4536.67
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.0861399 cHl3=-9.50561 cHq3=0.0801661
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-0.675078, denominator=wrap_pdf_Int[pTV]=86190.2
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.186765 cHl3=8.8591 cHq3=-0.971282
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-6.32705, denominator=wrap_pdf_Int[pTV]=46316
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.218731 cHl3=0.37397 cHq3=-2.08166
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=38.6705, denominator=wrap_pdf_Int[pTV]=90131.3
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=3.22694 cHl3=-7.04051 cHq3=0.54016
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=-2.21831, denominator=wrap_pdf_Int[pTV]=89722.5
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=3.44258 cHl3=4.96668 cHq3=0.0273884
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=0.403042, denominator=wrap_pdf_Int[pTV]=242.998
... (remaining 14 messages suppressed)
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=-3.35967 cHl3=-9.58704 cHq3=-6.27461
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=486.968, denominator=wrap_pdf_Int[pTV]=829530
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.154263 cHl3=2.95902 cHq3=-2.78828
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=59.1285, denominator=wrap_pdf_Int[pTV]=200921
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.338546 cHl3=0.879879 cHq3=-1.35856
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
p.d.f value is less than zero, trying to recover @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=14.3535, denominator=wrap_pdf_Int[pTV]=34082.1
RooAbsMinimizerFcn: Minimized function has error status.
Returning maximum FCN so far (86348.3) to force MIGRAD to back out of this region. Error log follows.
Parameter values: cHDD=0.372361 cHl3=0.491134 cHq3=-0.886807
RooFit::Detail::RooNormalizedPdf::wrap_pdf_over_wrap_pdf_Int[pTV][ numerator=wrap_pdf denominator=wrap_pdf_Int[pTV] ]
getLogVal() top-level p.d.f evaluates to NaN @ numerator=wrap_pdf=5.8312, denominator=wrap_pdf_Int[pTV]=12183.6
[#1] INFO:Fitting -- RooAbsPdf::fitTo(wrap_pdf) Calculating sum-of-weights-squared correction matrix for covariance matrix
[#1] INFO:InputArguments -- RooAbsData::plotOn(pseudo_dh) INFO: dataset has non-integer weights, auto-selecting SumW2 errors instead of Poisson errors
[#1] INFO:Plotting -- RooPlot::updateFitRangeNorm: New event count of 7389.24 will supersede previous event count of 9313.81 for normalization of PDF projections
Date
January 2022
Author
Rahul Balasubramanian

Definition in file rf712_lagrangianmorphfit.C.