rf604_constraints.C

Go to the documentation of this file.
00001 /////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'LIKELIHOOD AND MINIMIZATION' RooFit tutorial macro #604
00004 // 
00005 // Fitting with constraints
00006 //
00007 //
00008 // 07/2008 - Wouter Verkerke 
00009 //
00010 /////////////////////////////////////////////////////////////////////////
00011 
00012 #ifndef __CINT__
00013 #include "RooGlobalFunc.h"
00014 #endif
00015 #include "RooRealVar.h"
00016 #include "RooDataSet.h"
00017 #include "RooGaussian.h"
00018 #include "RooConstVar.h"
00019 #include "RooPolynomial.h"
00020 #include "RooAddPdf.h"
00021 #include "RooProdPdf.h"
00022 #include "RooFitResult.h"
00023 #include "RooPlot.h"
00024 #include "TCanvas.h"
00025 #include "TAxis.h"
00026 #include "TH1.h"
00027 using namespace RooFit ;
00028 
00029 
00030 void rf604_constraints()
00031 {
00032 
00033   // C r e a t e   m o d e l  a n d   d a t a s e t 
00034   // ----------------------------------------------
00035 
00036   // Construct a Gaussian p.d.f
00037   RooRealVar x("x","x",-10,10) ;
00038 
00039   RooRealVar m("m","m",0,-10,10) ;
00040   RooRealVar s("s","s",2,0.1,10) ;
00041   RooGaussian gauss("gauss","gauss(x,m,s)",x,m,s) ;
00042 
00043   // Construct a flat p.d.f (polynomial of 0th order)
00044   RooPolynomial poly("poly","poly(x)",x) ;
00045 
00046   // Construct model = f*gauss + (1-f)*poly
00047   RooRealVar f("f","f",0.5,0.,1.) ;
00048   RooAddPdf model("model","model",RooArgSet(gauss,poly),f) ;
00049 
00050   // Generate small dataset for use in fitting below
00051   RooDataSet* d = model.generate(x,50) ;
00052 
00053 
00054 
00055   // C r e a t e   c o n s t r a i n t   p d f 
00056   // -----------------------------------------
00057 
00058   // Construct Gaussian constraint p.d.f on parameter f at 0.8 with resolution of 0.1
00059   RooGaussian fconstraint("fconstraint","fconstraint",f,RooConst(0.8),RooConst(0.1)) ;
00060 
00061 
00062 
00063   // M E T H O D   1   -   A d d   i n t e r n a l   c o n s t r a i n t   t o   m o d e l 
00064   // -------------------------------------------------------------------------------------
00065 
00066   // Multiply constraint term with regular p.d.f using RooProdPdf
00067   // Specify in fitTo() that internal constraints on parameter f should be used
00068 
00069   // Multiply constraint with p.d.f
00070   RooProdPdf modelc("modelc","model with constraint",RooArgSet(model,fconstraint)) ;
00071   
00072   // Fit model (without use of constraint term)
00073   RooFitResult* r1 = model.fitTo(*d,Save()) ;
00074 
00075   // Fit modelc with constraint term on parameter f
00076   RooFitResult* r2 = modelc.fitTo(*d,Constrain(f),Save()) ;
00077 
00078 
00079 
00080   // M E T H O D   2   -     S p e c i f y   e x t e r n a l   c o n s t r a i n t   w h e n   f i t t i n g
00081   // -------------------------------------------------------------------------------------------------------
00082 
00083   // Construct another Gaussian constraint p.d.f on parameter f at 0.8 with resolution of 0.1
00084   RooGaussian fconstext("fconstext","fconstext",f,RooConst(0.2),RooConst(0.1)) ;
00085 
00086   // Fit with external constraint
00087   RooFitResult* r3 = model.fitTo(*d,ExternalConstraints(fconstext),Save()) ;
00088 
00089 
00090 
00091   // Print the fit results
00092   cout << "fit result without constraint (data generated at f=0.5)" << endl ;
00093   r1->Print("v") ;
00094   cout << "fit result with internal constraint (data generated at f=0.5, constraint is f=0.8+/-0.2)" << endl ;
00095   r2->Print("v") ;
00096   cout << "fit result with (another) external constraint (data generated at f=0.5, constraint is f=0.2+/-0.1)" << endl ;
00097   r3->Print("v") ;
00098   
00099 }

Generated on Tue Jul 5 15:45:08 2011 for ROOT_528-00b_version by  doxygen 1.5.1