rf501_simultaneouspdf.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'ORGANIZATION AND SIMULTANEOUS FITS' RooFit tutorial macro #501
00004 // 
00005 // Using simultaneous p.d.f.s to describe simultaneous fits to multiple
00006 // datasets
00007 //
00008 //
00009 //
00010 // 07/2008 - Wouter Verkerke 
00011 // 
00012 /////////////////////////////////////////////////////////////////////////
00013 
00014 #ifndef __CINT__
00015 #include "RooGlobalFunc.h"
00016 #endif
00017 #include "RooRealVar.h"
00018 #include "RooDataSet.h"
00019 #include "RooGaussian.h"
00020 #include "RooConstVar.h"
00021 #include "RooChebychev.h"
00022 #include "RooAddPdf.h"
00023 #include "RooSimultaneous.h"
00024 #include "RooCategory.h"
00025 #include "TCanvas.h"
00026 #include "TAxis.h"
00027 #include "RooPlot.h"
00028 using namespace RooFit ;
00029 
00030 
00031 void rf501_simultaneouspdf()
00032 {
00033   // C r e a t e   m o d e l   f o r   p h y s i c s   s a m p l e
00034   // -------------------------------------------------------------
00035 
00036   // Create observables
00037   RooRealVar x("x","x",-8,8) ;
00038 
00039   // Construct signal pdf
00040   RooRealVar mean("mean","mean",0,-8,8) ;
00041   RooRealVar sigma("sigma","sigma",0.3,0.1,10) ;
00042   RooGaussian gx("gx","gx",x,mean,sigma) ;
00043 
00044   // Construct background pdf
00045   RooRealVar a0("a0","a0",-0.1,-1,1) ;
00046   RooRealVar a1("a1","a1",0.004,-1,1) ;
00047   RooChebychev px("px","px",x,RooArgSet(a0,a1)) ;
00048 
00049   // Construct composite pdf
00050   RooRealVar f("f","f",0.2,0.,1.) ;
00051   RooAddPdf model("model","model",RooArgList(gx,px),f) ;
00052 
00053 
00054 
00055   // C r e a t e   m o d e l   f o r   c o n t r o l   s a m p l e
00056   // --------------------------------------------------------------
00057 
00058   // Construct signal pdf. 
00059   // NOTE that sigma is shared with the signal sample model
00060   RooRealVar mean_ctl("mean_ctl","mean_ctl",-3,-8,8) ;
00061   RooGaussian gx_ctl("gx_ctl","gx_ctl",x,mean_ctl,sigma) ;
00062 
00063   // Construct the background pdf
00064   RooRealVar a0_ctl("a0_ctl","a0_ctl",-0.1,-1,1) ;
00065   RooRealVar a1_ctl("a1_ctl","a1_ctl",0.5,-0.1,1) ;
00066   RooChebychev px_ctl("px_ctl","px_ctl",x,RooArgSet(a0_ctl,a1_ctl)) ;
00067 
00068   // Construct the composite model
00069   RooRealVar f_ctl("f_ctl","f_ctl",0.5,0.,1.) ;
00070   RooAddPdf model_ctl("model_ctl","model_ctl",RooArgList(gx_ctl,px_ctl),f_ctl) ;
00071 
00072 
00073 
00074   // G e n e r a t e   e v e n t s   f o r   b o t h   s a m p l e s 
00075   // ---------------------------------------------------------------
00076 
00077   // Generate 1000 events in x and y from model
00078   RooDataSet *data = model.generate(RooArgSet(x),100) ;
00079   RooDataSet *data_ctl = model_ctl.generate(RooArgSet(x),2000) ;
00080 
00081 
00082 
00083   // C r e a t e   i n d e x   c a t e g o r y   a n d   j o i n   s a m p l e s 
00084   // ---------------------------------------------------------------------------
00085 
00086   // Define category to distinguish physics and control samples events
00087   RooCategory sample("sample","sample") ;
00088   sample.defineType("physics") ;
00089   sample.defineType("control") ;
00090 
00091   // Construct combined dataset in (x,sample)
00092   RooDataSet combData("combData","combined data",x,Index(sample),Import("physics",*data),Import("control",*data_ctl)) ;
00093 
00094 
00095 
00096   // C o n s t r u c t   a   s i m u l t a n e o u s   p d f   i n   ( x , s a m p l e )
00097   // -----------------------------------------------------------------------------------
00098 
00099   // Construct a simultaneous pdf using category sample as index
00100   RooSimultaneous simPdf("simPdf","simultaneous pdf",sample) ;
00101 
00102   // Associate model with the physics state and model_ctl with the control state
00103   simPdf.addPdf(model,"physics") ;
00104   simPdf.addPdf(model_ctl,"control") ;
00105 
00106 
00107 
00108   // P e r f o r m   a   s i m u l t a n e o u s   f i t
00109   // ---------------------------------------------------
00110 
00111   // Perform simultaneous fit of model to data and model_ctl to data_ctl
00112   simPdf.fitTo(combData) ;
00113 
00114 
00115 
00116   // P l o t   m o d e l   s l i c e s   o n   d a t a    s l i c e s 
00117   // ----------------------------------------------------------------
00118 
00119   // Make a frame for the physics sample
00120   RooPlot* frame1 = x.frame(Bins(30),Title("Physics sample")) ;
00121 
00122   // Plot all data tagged as physics sample
00123   combData.plotOn(frame1,Cut("sample==sample::physics")) ;
00124 
00125   // Plot "physics" slice of simultaneous pdf. 
00126   // NBL You _must_ project the sample index category with data using ProjWData 
00127   // as a RooSimultaneous makes no prediction on the shape in the index category 
00128   // and can thus not be integrated
00129   simPdf.plotOn(frame1,Slice(sample,"physics"),ProjWData(sample,combData)) ;
00130   simPdf.plotOn(frame1,Slice(sample,"physics"),Components("px"),ProjWData(sample,combData),LineStyle(kDashed)) ;
00131 
00132   // The same plot for the control sample slice
00133   RooPlot* frame2 = x.frame(Bins(30),Title("Control sample")) ;
00134   combData.plotOn(frame2,Cut("sample==sample::control")) ;
00135   simPdf.plotOn(frame2,Slice(sample,"control"),ProjWData(sample,combData)) ;
00136   simPdf.plotOn(frame2,Slice(sample,"control"),Components("px_ctl"),ProjWData(sample,combData),LineStyle(kDashed)) ;
00137 
00138 
00139 
00140   TCanvas* c = new TCanvas("rf501_simultaneouspdf","rf403_simultaneouspdf",800,400) ;
00141   c->Divide(2) ;
00142   c->cd(1) ; gPad->SetLeftMargin(0.15) ; frame1->GetYaxis()->SetTitleOffset(1.4) ; frame1->Draw() ;
00143   c->cd(2) ; gPad->SetLeftMargin(0.15) ; frame2->GetYaxis()->SetTitleOffset(1.4) ; frame2->Draw() ;
00144 
00145 
00146 }

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