Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf712_lagrangianmorphfit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roofit
3/// \notebook -js
4/// Performing a simple fit with RooLagrangianMorphFunc.
5/// a morphing function is setup as a function of three variables and
6/// a fit is performed on a pseudo-dataset.
7///
8/// \macro_image
9/// \macro_code
10/// \macro_output
11///
12/// \date January 2022
13/// \author Rahul Balasubramanian
14
15#include <RooDataHist.h>
16#include <RooFitResult.h>
18#include <RooPlot.h>
19#include <RooRealVar.h>
20
21#include <TAxis.h>
22#include <TCanvas.h>
23#include <TH2.h>
24#include <TStyle.h>
25
26using namespace RooFit;
27
29{
30 // C r e a t e v a r i a b l e s f o r
31 // m o r p h i n g f u n c t i o n
32 // ---------------------------------------------
33
34 std::string observablename = "pTV";
35 RooRealVar obsvar(observablename.c_str(), "observable of pTV", 10, 600);
36 RooRealVar kSM("kSM", "sm modifier", 1.0);
37 RooRealVar cHq3("cHq3", "EFT modifier", -10.0, 10.0);
38 cHq3.setAttribute("NewPhysics", true);
39 RooRealVar cHl3("cHl3", "EFT modifier", -10.0, 10.0);
40 cHl3.setAttribute("NewPhysics", true);
41 RooRealVar cHDD("cHDD", "EFT modifier", -10.0, 10.0);
42 cHDD.setAttribute("NewPhysics", true);
43
44 // I n p u t s n e e d e d f o r c o n f i g
45 // ---------------------------------------------
46 std::string infilename = std::string(gROOT->GetTutorialDir()) + "/roofit/input_histos_rf_lagrangianmorph.root";
47 std::vector<std::string> samplelist = {"SM_NPsq0", "cHq3_NPsq1", "cHq3_NPsq2", "cHl3_NPsq1",
48 "cHl3_NPsq2", "cHDD_NPsq1", "cHDD_NPsq2", "cHl3_cHDD_NPsq2",
49 "cHq3_cHDD_NPsq2", "cHl3_cHq3_NPsq2"};
50
51 // S e t u p C o n f i g
52 // ---------------------------------------------
54 config.fileName = infilename;
55 config.observableName = observablename;
56 config.folderNames = samplelist;
57 config.couplings.add(cHq3);
58 config.couplings.add(cHl3);
59 config.couplings.add(cHDD);
60 config.couplings.add(kSM);
61
62 // C r e a t e m o r p h i n g f u n c t i o n
63 // ---------------------------------------------
64 RooLagrangianMorphFunc morphfunc("morphfunc", "morphed dist. of pTV", config);
65
66 // 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
67 // f o r f i t
68 // ---------------------------------------------
69 morphfunc.setParameter("cHq3", 0.01);
70 morphfunc.setParameter("cHl3", 1.0);
71 morphfunc.setParameter("cHDD", 0.2);
72
73 auto pseudo_hist = morphfunc.createTH1("pseudo_hist");
74 auto pseudo_dh = new RooDataHist("pseudo_dh", "pseudo_dh", RooArgList(obsvar), pseudo_hist);
75
76 // reset parameters to zeros before fit
77 morphfunc.setParameter("cHq3", 0.0);
78 morphfunc.setParameter("cHl3", 0.0);
79 morphfunc.setParameter("cHDD", 0.0);
80
81 // error set used as initial step size
82 cHq3.setError(0.1);
83 cHl3.setError(0.1);
84 cHDD.setError(0.1);
85
86 // W r a p p d f o n m o r p h f u n c a n d
87 // f i t t o d a t a h i s t o g r a m
88 // ---------------------------------------------
89 // wrapper pdf to normalise morphing function to a morphing pdf
90 RooWrapperPdf model("wrap_pdf", "wrap_pdf", morphfunc);
91 auto fitres = model.fitTo(*pseudo_dh, SumW2Error(true), Optimize(false), Save(), PrintLevel(-1));
92 auto hcorr = fitres->correlationHist();
93
94 // 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
95 // a n d p l o t w i t h i n i t i a l
96 // h i s t o g r a m
97 // ---------------------------------------------
98 auto postfit_hist = morphfunc.createTH1("morphing_postfit_hist");
99 RooDataHist postfit_dh("morphing_postfit_dh", "morphing_postfit_dh", RooArgList(obsvar), postfit_hist);
100
101 auto frame0 = obsvar.frame(Title("Fitted histogram of p_{T}^{V}"));
102 postfit_dh.plotOn(frame0, Name("postfit_dist"), DrawOption("C"), LineColor(kBlue), DataError(RooAbsData::None),
103 XErrorSize(0));
104 pseudo_dh->plotOn(frame0, Name("input"));
105
106 // D r a w p l o t s o n c a n v a s
107 // ---------------------------------------------
108 TCanvas *c1 = new TCanvas("fig3", "fig3", 800, 400);
109 c1->Divide(2, 1);
110
111 c1->cd(1);
112 gPad->SetLeftMargin(0.15);
113 gPad->SetRightMargin(0.05);
114
115 model.paramOn(frame0, Layout(0.50, 0.75, 0.9), Parameters(config.couplings));
116 frame0->GetXaxis()->SetTitle("p_{T}^{V}");
117 frame0->Draw();
118
119 c1->cd(2);
120 gPad->SetLeftMargin(0.15);
121 gPad->SetRightMargin(0.15);
122 gStyle->SetPaintTextFormat("4.1f");
123 gStyle->SetOptStat(0);
124 hcorr->SetMarkerSize(3.);
125 hcorr->SetTitle("correlation matrix");
126 hcorr->GetYaxis()->SetTitleOffset(1.4);
127 hcorr->GetYaxis()->SetLabelSize(0.1);
128 hcorr->GetXaxis()->SetLabelSize(0.1);
129 hcorr->GetYaxis()->SetBinLabel(1, "c_{HDD}");
130 hcorr->GetYaxis()->SetBinLabel(2, "c_{Hl^{(3)}}");
131 hcorr->GetYaxis()->SetBinLabel(3, "c_{Hq^{(3)}}");
132 hcorr->GetXaxis()->SetBinLabel(3, "c_{HDD}");
133 hcorr->GetXaxis()->SetBinLabel(2, "c_{Hl^{(3)}}");
134 hcorr->GetXaxis()->SetBinLabel(1, "c_{Hq^{(3)}}");
135 hcorr->Draw("colz text");
136 c1->SaveAs("rf712_lagrangianmorphfit.png");
137}
@ kBlue
Definition Rtypes.h:66
#define gROOT
Definition TROOT.h:407
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
#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
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
Class RooLagrangianMorphing is a implementation of the method of Effective Lagrangian Morphing,...
RooRealVar represents a 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:1636
void SetPaintTextFormat(const char *format="g")
Definition TStyle.h:383
RooCmdArg Parameters(const RooArgSet &params)
RooCmdArg Layout(double xmin, double xmax=0.99, double ymin=0.95)
RooCmdArg Save(bool flag=true)
RooCmdArg SumW2Error(bool flag)
RooCmdArg Optimize(Int_t flag=2)
RooCmdArg DataError(Int_t)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg DrawOption(const char *opt)
RooCmdArg XErrorSize(double width)
RooCmdArg LineColor(Color_t color)
return c1
Definition legend1.C:41
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
const char * Name
Definition TXMLSetup.cxx:67
const char * Title
Definition TXMLSetup.cxx:68
std::vector< std::string > folderNames