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 "RooChebychev.h"
00021 #include "RooExponential.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 using namespace RooFit ;
00026
00027
00028 void rf205_compplot()
00029 {
00030
00031
00032
00033
00034 RooRealVar x("x","x",0,10) ;
00035
00036
00037 RooRealVar mean("mean","mean of gaussians",5) ;
00038 RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
00039 RooRealVar sigma2("sigma2","width of gaussians",1) ;
00040 RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
00041 RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
00042
00043
00044 RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
00045 RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
00046
00047
00048 RooRealVar a0("a0","a0",0.5,0.,1.) ;
00049 RooRealVar a1("a1","a1",-0.2,0.,1.) ;
00050 RooChebychev bkg1("bkg1","Background 1",x,RooArgSet(a0,a1)) ;
00051
00052
00053 RooRealVar alpha("alpha","alpha",-1) ;
00054 RooExponential bkg2("bkg2","Background 2",x,alpha) ;
00055
00056
00057 RooRealVar bkg1frac("sig1frac","fraction of component 1 in background",0.2,0.,1.) ;
00058 RooAddPdf bkg("bkg","Signal",RooArgList(bkg1,bkg2),sig1frac) ;
00059
00060
00061 RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
00062 RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
00063
00064
00065
00066
00067
00068
00069
00070 RooDataSet *data = model.generate(x,1000) ;
00071
00072
00073 RooPlot* xframe = x.frame(Title("Component plotting of pdf=(sig1+sig2)+(bkg1+bkg2)")) ;
00074 data->plotOn(xframe) ;
00075 model.plotOn(xframe) ;
00076
00077
00078 RooPlot* xframe2 = (RooPlot*) xframe->Clone("xframe2") ;
00079
00080
00081
00082
00083
00084
00085 model.plotOn(xframe,Components(bkg),LineColor(kRed)) ;
00086
00087
00088 model.plotOn(xframe,Components(bkg2),LineStyle(kDashed),LineColor(kRed)) ;
00089
00090
00091
00092
00093 model.plotOn(xframe,Components(RooArgSet(bkg,sig2)),LineStyle(kDotted)) ;
00094
00095
00096
00097
00098
00099
00100
00101 model.plotOn(xframe2,Components("bkg"),LineColor(kCyan)) ;
00102
00103
00104 model.plotOn(xframe2,Components("bkg1,sig2"),LineStyle(kDotted),LineColor(kCyan)) ;
00105
00106
00107 model.plotOn(xframe2,Components("sig*"),LineStyle(kDashed),LineColor(kCyan)) ;
00108
00109
00110 model.plotOn(xframe2,Components("bkg1,sig*"),LineStyle(kDashed),LineColor(kYellow),Invisible()) ;
00111
00112
00113
00114 TCanvas* c = new TCanvas("rf205_compplot","rf205_compplot",800,400) ;
00115 c->Divide(2) ;
00116 c->cd(1) ; gPad->SetLeftMargin(0.15) ; xframe->GetYaxis()->SetTitleOffset(1.4) ; xframe->Draw() ;
00117 c->cd(2) ; gPad->SetLeftMargin(0.15) ; xframe2->GetYaxis()->SetTitleOffset(1.4) ; xframe2->Draw() ;
00118
00119
00120 }