00001
00002
00003
00004
00005
00006
00007
00008
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 "RooChebychev.h"
00020 #include "RooAddPdf.h"
00021 #include "RooWorkspace.h"
00022 #include "RooPlot.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "TFile.h"
00026 #include "TH1.h"
00027 using namespace RooFit ;
00028
00029
00030 void fillWorkspace(RooWorkspace& w) ;
00031
00032 void rf510_wsnamedsets()
00033 {
00034
00035
00036
00037 RooWorkspace* w = new RooWorkspace("w") ;
00038 fillWorkspace(*w) ;
00039
00040
00041
00042 RooAbsPdf* model = w->pdf("model") ;
00043
00044
00045 RooDataSet* data = model->generate(*w->set("observables"),1000) ;
00046
00047
00048 model->fitTo(*data) ;
00049
00050
00051 RooPlot* frame = ((RooRealVar*)w->set("observables")->first())->frame() ;
00052 data->plotOn(frame) ;
00053 model->plotOn(frame) ;
00054
00055
00056 w->loadSnapshot("reference_fit") ;
00057 model->plotOn(frame,LineColor(kRed)) ;
00058 w->loadSnapshot("reference_fit_bkgonly") ;
00059 model->plotOn(frame,LineColor(kRed),LineStyle(kDashed)) ;
00060
00061
00062
00063 new TCanvas("rf510_wsnamedsets","rf503_wsnamedsets",600,600) ;
00064 gPad->SetLeftMargin(0.15) ; frame->GetYaxis()->SetTitleOffset(1.4) ; frame->Draw() ;
00065
00066
00067
00068 w->Print() ;
00069
00070
00071
00072 gDirectory->Add(w) ;
00073
00074 }
00075
00076
00077
00078 void fillWorkspace(RooWorkspace& w)
00079 {
00080
00081
00082
00083
00084 RooRealVar x("x","x",0,10) ;
00085
00086
00087 RooRealVar mean("mean","mean of gaussians",5,0,10) ;
00088 RooRealVar sigma1("sigma1","width of gaussians",0.5) ;
00089 RooRealVar sigma2("sigma2","width of gaussians",1) ;
00090
00091 RooGaussian sig1("sig1","Signal component 1",x,mean,sigma1) ;
00092 RooGaussian sig2("sig2","Signal component 2",x,mean,sigma2) ;
00093
00094
00095 RooRealVar a0("a0","a0",0.5,0.,1.) ;
00096 RooRealVar a1("a1","a1",-0.2,0.,1.) ;
00097 RooChebychev bkg("bkg","Background",x,RooArgSet(a0,a1)) ;
00098
00099
00100 RooRealVar sig1frac("sig1frac","fraction of component 1 in signal",0.8,0.,1.) ;
00101 RooAddPdf sig("sig","Signal",RooArgList(sig1,sig2),sig1frac) ;
00102
00103
00104 RooRealVar bkgfrac("bkgfrac","fraction of background",0.5,0.,1.) ;
00105 RooAddPdf model("model","g1+g2+a",RooArgList(bkg,sig),bkgfrac) ;
00106
00107
00108 w.import(model) ;
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122 RooArgSet* params = (RooArgSet*) model.getParameters(x) ;
00123 w.defineSet("parameters",*params) ;
00124 w.defineSet("observables",x) ;
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141 RooDataSet* refData = model.generate(x,10000) ;
00142 model.fitTo(*refData,PrintLevel(-1)) ;
00143
00144
00145
00146 w.saveSnapshot("reference_fit",*params,kTRUE) ;
00147
00148
00149
00150
00151 bkgfrac.setVal(1) ;
00152 bkgfrac.setConstant(kTRUE) ;
00153 bkgfrac.removeError() ;
00154 model.fitTo(*refData,PrintLevel(-1)) ;
00155
00156 w.saveSnapshot("reference_fit_bkgonly",*params,kTRUE) ;
00157
00158
00159 }