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 "RooConstVar.h"
00020 #include "RooAddPdf.h"
00021 #include "RooChebychev.h"
00022 #include "RooFitResult.h"
00023 #include "TCanvas.h"
00024 #include "TAxis.h"
00025 #include "RooPlot.h"
00026 #include "TFile.h"
00027 #include "TStyle.h"
00028 #include "TH2.h"
00029 #include "TH3.h"
00030
00031 using namespace RooFit ;
00032
00033
00034 void rf608_fitresultaspdf()
00035 {
00036
00037
00038
00039
00040 RooRealVar x("x","x",-20,20) ;
00041
00042
00043 RooRealVar mean("mean","mean of g1 and g2",0,-1,1) ;
00044 RooRealVar sigma_g1("sigma_g1","width of g1",2) ;
00045 RooGaussian g1("g1","g1",x,mean,sigma_g1) ;
00046
00047 RooRealVar sigma_g2("sigma_g2","width of g2",4,3.0,5.0) ;
00048 RooGaussian g2("g2","g2",x,mean,sigma_g2) ;
00049
00050 RooRealVar frac("frac","frac",0.5,0.0,1.0) ;
00051 RooAddPdf model("model","model",RooArgList(g1,g2),frac) ;
00052
00053
00054 RooDataSet* data = model.generate(x,1000) ;
00055
00056
00057
00058
00059
00060 RooFitResult* r = model.fitTo(*data,Save()) ;
00061
00062
00063
00064
00065
00066 RooAbsPdf* parabPdf = r->createHessePdf(RooArgSet(frac,mean,sigma_g2)) ;
00067
00068
00069
00070
00071
00072
00073 RooDataSet* d = parabPdf->generate(RooArgSet(mean,sigma_g2,frac),100000) ;
00074
00075
00076
00077 TH3* hh_3d = (TH3*) parabPdf->createHistogram("mean,sigma_g2,frac",25,25,25) ;
00078 hh_3d->SetFillColor(kBlue) ;
00079
00080
00081
00082
00083
00084 RooAbsPdf* pdf_sigmag2_frac = parabPdf->createProjection(mean) ;
00085 RooAbsPdf* pdf_mean_frac = parabPdf->createProjection(sigma_g2) ;
00086 RooAbsPdf* pdf_mean_sigmag2 = parabPdf->createProjection(frac) ;
00087
00088
00089
00090 TH2* hh_sigmag2_frac = (TH2*) pdf_sigmag2_frac->createHistogram("sigma_g2,frac",50,50) ;
00091 TH2* hh_mean_frac = (TH2*) pdf_mean_frac->createHistogram("mean,frac",50,50) ;
00092 TH2* hh_mean_sigmag2 = (TH2*) pdf_mean_sigmag2->createHistogram("mean,sigma_g2",50,50) ;
00093 hh_mean_frac->SetLineColor(kBlue) ;
00094 hh_sigmag2_frac->SetLineColor(kBlue) ;
00095 hh_mean_sigmag2->SetLineColor(kBlue) ;
00096
00097
00098
00099 gStyle->SetCanvasPreferGL(true);
00100 gStyle->SetPalette(1) ;
00101 new TCanvas("rf608_fitresultaspdf_1","rf608_fitresultaspdf_1",600,600) ;
00102 hh_3d->Draw("gliso") ;
00103
00104
00105 TCanvas* c2 = new TCanvas("rf608_fitresultaspdf_2","rf608_fitresultaspdf_2",900,600) ;
00106 c2->Divide(3,2) ;
00107 c2->cd(1) ; gPad->SetLeftMargin(0.15) ; hh_mean_sigmag2->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_sigmag2->Draw("surf3") ;
00108 c2->cd(2) ; gPad->SetLeftMargin(0.15) ; hh_sigmag2_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_sigmag2_frac->Draw("surf3") ;
00109 c2->cd(3) ; gPad->SetLeftMargin(0.15) ; hh_mean_frac->GetZaxis()->SetTitleOffset(1.4) ; hh_mean_frac->Draw("surf3") ;
00110
00111
00112 TH1* tmp1 = d->createHistogram("mean,sigma_g2",50,50) ;
00113 TH1* tmp2 = d->createHistogram("sigma_g2,frac",50,50) ;
00114 TH1* tmp3 = d->createHistogram("mean,frac",50,50) ;
00115
00116 c2->cd(4) ; gPad->SetLeftMargin(0.15) ; tmp1->GetZaxis()->SetTitleOffset(1.4) ; tmp1->Draw("lego3") ;
00117 c2->cd(5) ; gPad->SetLeftMargin(0.15) ; tmp2->GetZaxis()->SetTitleOffset(1.4) ; tmp2->Draw("lego3") ;
00118 c2->cd(6) ; gPad->SetLeftMargin(0.15) ; tmp3->GetZaxis()->SetTitleOffset(1.4) ; tmp3->Draw("lego3") ;
00119
00120 }
00121