00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooGaussian.h"
00019 #include "RooAddPdf.h"
00020 #include "RooNLLVar.h"
00021 #include "RooProfileLL.h"
00022 #include "RooMinuit.h"
00023 #include "TCanvas.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026
00027
00028 class TestBasic605 : public RooFitTestUnit
00029 {
00030 public:
00031 TestBasic605(TFile* refFile, Bool_t writeRef, Int_t verbose) : RooFitTestUnit("Profile Likelihood operator",refFile,writeRef,verbose) {} ;
00032 Bool_t testCode() {
00033
00034
00035
00036
00037
00038 RooRealVar x("x","x",-20,20) ;
00039
00040
00041 RooRealVar mean("mean","mean of g1 and g2",0,-10,10) ;
00042 RooRealVar sigma_g1("sigma_g1","width of g1",3) ;
00043 RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
00044
00045 RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,6.0) ;
00046 RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
00047
00048 RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
00049 RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
00050
00051
00052 RooDataSet* data = model.generate(x,1000) ;
00053
00054
00055
00056
00057
00058
00059
00060 RooNLLVar nll("nll","nll",model,*data) ;
00061
00062
00063 RooMinuit(nll).migrad() ;
00064
00065
00066 RooPlot* frame1 = frac.frame(Bins(10),Range(0.01,0.95),Title("LL and profileLL in frac")) ;
00067 nll.plotOn(frame1,ShiftToZero()) ;
00068
00069
00070 RooPlot* frame2 = sigma_g2.frame(Bins(10),Range(3.3,5.0),Title("LL and profileLL in sigma_g2")) ;
00071 nll.plotOn(frame2,ShiftToZero()) ;
00072
00073
00074
00075
00076
00077
00078
00079
00080 RooProfileLL pll_frac("pll_frac","pll_frac",nll,frac) ;
00081
00082
00083 pll_frac.plotOn(frame1,LineColor(kRed)) ;
00084
00085
00086 frame1->SetMinimum(0) ;
00087 frame1->SetMaximum(3) ;
00088
00089
00090
00091
00092
00093
00094
00095
00096 RooProfileLL pll_sigmag2("pll_sigmag2","pll_sigmag2",nll,sigma_g2) ;
00097
00098
00099 pll_sigmag2.plotOn(frame2,LineColor(kRed)) ;
00100
00101
00102 frame2->SetMinimum(0) ;
00103 frame2->SetMaximum(3) ;
00104
00105
00106 regPlot(frame1,"rf605_plot1") ;
00107 regPlot(frame2,"rf605_plot2") ;
00108
00109 delete data ;
00110
00111 return kTRUE ;
00112 }
00113 } ;
00114