Logo ROOT   6.07/09
Reference Guide
rf506_msgservice.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_roofit
3 /// \notebook -nodraw
4 /// 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #506
5 ///
6 /// Tuning and customizing the RooFit message logging facility
7 ///
8 /// \macro_output
9 /// \macro_code
10 /// \author 07/2008 - Wouter Verkerke
11 
12 
13 #include "RooRealVar.h"
14 #include "RooDataSet.h"
15 #include "RooGaussian.h"
16 #include "RooConstVar.h"
17 #include "RooPolynomial.h"
18 #include "RooAddPdf.h"
19 #include "TCanvas.h"
20 #include "TAxis.h"
21 #include "RooPlot.h"
22 #include "RooMsgService.h"
23 
24 using namespace RooFit ;
25 
26 
27 void rf506_msgservice()
28 {
29  // C r e a t e p d f
30  // --------------------
31 
32  // Construct gauss(x,m,s)
33  RooRealVar x("x","x",-10,10) ;
34  RooRealVar m("m","m",0,-10,10) ;
35  RooRealVar s("s","s",1,-10,10) ;
36  RooGaussian gauss("g","g",x,m,s) ;
37 
38  // Construct poly(x,p0)
39  RooRealVar p0("p0","p0",0.01,0.,1.) ;
40  RooPolynomial poly("p","p",x,p0) ;
41 
42  // Construct model = f*gauss(x) + (1-f)*poly(x)
43  RooRealVar f("f","f",0.5,0.,1.) ;
44  RooAddPdf model("model","model",RooArgSet(gauss,poly),f) ;
45 
46  RooDataSet* data = model.generate(x,10) ;
47 
48 
49 
50  // P r i n t c o n f i g u r a t i o n o f m e s s a g e s e r v i c e
51  // ---------------------------------------------------------------------------
52 
53  // Print streams configuration
55  cout << endl ;
56 
57 
58 
59  // A d d i n g I n t e g r a t i o n t o p i c t o e x i s t i n g I N F O s t r e a m
60  // -----------------------------------------------------------------------------------------------
61 
62  // Print streams configuration
64  cout << endl ;
65 
66  // Add Integration topic to existing INFO stream
68 
69  // Construct integral over gauss to demonstrate new message stream
70  RooAbsReal* igauss = gauss.createIntegral(x) ;
71  igauss->Print() ;
72 
73  // Print streams configuration in verbose, which also shows inactive streams
74  cout << endl ;
76  cout << endl ;
77 
78  // Remove stream
80 
81 
82 
83  // E x a m p l e s o f p d f v a l u e t r a c i n g s t r e a m
84  // -----------------------------------------------------------------------
85 
86  // Show DEBUG level message on function tracing, trace RooGaussian only
88 
89  // Perform a fit to generate some tracing messages
90  model.fitTo(*data,Verbose(kTRUE)) ;
91 
92  // Reset message service to default stream configuration
94 
95 
96 
97  // Show DEBUG level message on function tracing on all objects, redirect output to file
99 
100  // Perform a fit to generate some tracing messages
101  model.fitTo(*data,Verbose(kTRUE)) ;
102 
103  // Reset message service to default stream configuration
105 
106 
107 
108  // E x a m p l e o f a n o t h e r d e b u g g i n g s t r e a m
109  // ---------------------------------------------------------------------
110 
111  // Show DEBUG level messages on client/server link state management
114 
115  // Clone composite pdf g to trigger some link state management activity
116  RooAbsArg* gprime = gauss.cloneTree() ;
117  gprime->Print() ;
118 
119  // Reset message service to default stream configuration
121 
122 
123 
124 }
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
void addTopic(RooFit::MsgTopic newTopic)
virtual RooAbsArg * cloneTree(const char *newname=0) const
Clone tree expression of objects.
Definition: RooAbsArg.cxx:2296
StreamConfig & getStream(Int_t id)
void removeTopic(RooFit::MsgTopic oldTopic)
static RooMsgService & instance()
Return reference to singleton instance.
RooCmdArg ClassName(const char *name)
Double_t x[n]
Definition: legend1.C:17
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
Definition: RooAbsArg.h:227
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooCmdArg Topic(Int_t topic)
TMarker * m
Definition: textangle.C:8
Int_t addStream(RooFit::MsgLevel level, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg())
Add a message logging stream for message with given RooFit::MsgLevel or higher (i.e.
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
double f(double x)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooCmdArg Verbose(Bool_t flag=kTRUE)
RooCmdArg OutputFile(const char *fileName)
RooPolynomial implements a polynomial p.d.f of the form By default coefficient a_0 is chosen to be 1...
Definition: RooPolynomial.h:28
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
void Print(Option_t *options=0) const
Print configuration of message service.
const Bool_t kTRUE
Definition: Rtypes.h:91
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:503