Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf506_msgservice.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -nodraw
4## Organization and simultaneous fits: tuning and customizing the ROOT.RooFit message logging facility
5##
6## \macro_code
7##
8## \date February 2018
9## \authors Clemens Lange, Wouter Verkerke (C++ version)
10
11import ROOT
12
13# Create pdf
14# --------------------
15
16# Construct gauss(x,m,s)
17x = ROOT.RooRealVar("x", "x", -10, 10)
18m = ROOT.RooRealVar("m", "m", 0, -10, 10)
19s = ROOT.RooRealVar("s", "s", 1, -10, 10)
20gauss = ROOT.RooGaussian("g", "g", x, m, s)
21
22# Construct poly(x,p0)
23p0 = ROOT.RooRealVar("p0", "p0", 0.01, 0., 1.)
24poly = ROOT.RooPolynomial("p", "p", x, ROOT.RooArgList(p0))
25
26# model = f*gauss(x) + (1-f)*poly(x)
27f = ROOT.RooRealVar("f", "f", 0.5, 0., 1.)
28model = ROOT.RooAddPdf("model", "model", ROOT.RooArgList(
29 gauss, poly), ROOT.RooArgList(f))
30
31data = model.generate(ROOT.RooArgSet(x), 10)
32
33# Print configuration of message service
34# ------------------------------------------
35
36# Print streams configuration
37ROOT.RooMsgService.instance().Print()
38
39# Adding integration topic to existing INFO stream
40# ---------------------------------------------------
41
42# Print streams configuration
43ROOT.RooMsgService.instance().Print()
44
45# Add Integration topic to existing INFO stream
46ROOT.RooMsgService.instance().getStream(1).addTopic(ROOT.RooFit.Integration)
47
48# Construct integral over gauss to demonstrate message stream
49igauss = gauss.createIntegral(ROOT.RooArgSet(x))
50igauss.Print()
51
52# Print streams configuration in verbose, also shows inactive streams
53ROOT.RooMsgService.instance().Print()
54
55# Remove stream
56ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.Integration)
57
58# Examples of pdf value tracing
59# -----------------------------------------------------------------------
60
61# Show DEBUG level message on function tracing, ROOT.RooGaussian only
62ROOT.RooMsgService.instance().addStream(
63 ROOT.RooFit.DEBUG,
64 ROOT.RooFit.Topic(
65 ROOT.RooFit.Tracing),
66 ROOT.RooFit.ClassName("RooGaussian"))
67
68# Perform a fit to generate some tracing messages
69model.fitTo(data, ROOT.RooFit.Verbose(ROOT.kTRUE))
70
71# Reset message service to default stream configuration
72ROOT.RooMsgService.instance().reset()
73
74# Show DEBUG level message on function tracing on all objects, output to
75# file
76ROOT.RooMsgService.instance().addStream(
77 ROOT.RooFit.DEBUG,
78 ROOT.RooFit.Topic(
79 ROOT.RooFit.Tracing),
80 ROOT.RooFit.OutputFile("rf506_debug.log"))
81
82# Perform a fit to generate some tracing messages
83model.fitTo(data, ROOT.RooFit.Verbose(ROOT.kTRUE))
84
85# Reset message service to default stream configuration
86ROOT.RooMsgService.instance().reset()
87
88# Example of another debugging stream
89# ---------------------------------------------------------------------
90
91# Show DEBUG level messages on client/server link state management
92ROOT.RooMsgService.instance().addStream(
93 ROOT.RooFit.DEBUG, ROOT.RooFit.Topic(ROOT.RooFit.LinkStateMgmt))
94ROOT.RooMsgService.instance().Print("v")
95
96# Clone composite pdf g to trigger some link state management activity
97gprime = gauss.cloneTree()
98gprime.Print()
99
100# Reset message service to default stream configuration
101ROOT.RooMsgService.instance().reset()